Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | Static Private Member Functions | Private Attributes | List of all members
PETScWrappers::SolverTFQMR Class Reference

#include <deal.II/lac/petsc_solver.h>

Inheritance diagram for PETScWrappers::SolverTFQMR:
[legend]

Classes

struct  AdditionalData
 

Public Member Functions

 SolverTFQMR (SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
 
void solve (const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionBase &preconditioner)
 
virtual void reset ()
 
void set_prefix (const std::string &prefix)
 
SolverControlcontrol () const
 
void initialize (const PreconditionBase &preconditioner)
 

Protected Member Functions

virtual void set_solver_type (KSP &ksp) const override
 

Protected Attributes

const AdditionalData additional_data
 
SolverControlsolver_control
 
const MPI_Comm mpi_communicator
 
std::string prefix_name
 

Static Private Member Functions

static PetscErrorCode convergence_test (KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
 

Private Attributes

std::unique_ptr< SolverDatasolver_data
 

Detailed Description

An implementation of the solver interface using the PETSc TFQMR solver.

Definition at line 646 of file petsc_solver.h.

Constructor & Destructor Documentation

◆ SolverTFQMR()

PETScWrappers::SolverTFQMR::SolverTFQMR ( SolverControl cn,
const MPI_Comm mpi_communicator = PETSC_COMM_SELF,
const AdditionalData data = AdditionalData() 
)

Constructor. In contrast to deal.II's own solvers, there is no need to give a vector memory object. However, PETSc solvers want to have an MPI communicator context over which computations are parallelized. By default, PETSC_COMM_SELF is used here, but you can change this. Note that for single processor (non-MPI) versions, this parameter does not have any effect.

The last argument takes a structure with additional, solver dependent flags for tuning.

Note that the communicator used here must match the communicator used in the system matrix, solution, and right hand side object of the solve to be done with this solver. Otherwise, PETSc will generate hard to track down errors, see the documentation of the SolverBase class.

Definition at line 496 of file petsc_solver.cc.

Member Function Documentation

◆ set_solver_type()

void PETScWrappers::SolverTFQMR::set_solver_type ( KSP &  ksp) const
overrideprotectedvirtual

Function that takes a Krylov Subspace Solver context object, and sets the type of solver that is appropriate for this class.

Implements PETScWrappers::SolverBase.

Definition at line 505 of file petsc_solver.cc.

◆ solve()

void SolverBase< VectorType >::solve ( const MatrixBase A,
VectorBase x,
const VectorBase b,
const PreconditionBase preconditioner 
)
inherited

Solve the linear system Ax=b. Depending on the information provided by derived classes and the object passed as a preconditioner, one of the linear solvers and preconditioners of PETSc is chosen. Repeated calls to solve() do not reconstruct the preconditioner for performance reasons. See class Documentation.

Definition at line 53 of file petsc_solver.cc.

◆ reset()

void SolverBase< VectorType >::reset ( )
virtualinherited

Resets the contained preconditioner and solver object. See class description for more details.

Definition at line 161 of file petsc_solver.cc.

◆ set_prefix()

void SolverBase< VectorType >::set_prefix ( const std::string &  prefix)
inherited

Sets a prefix name for the solver object. Useful when customizing the PETSc KSP object with command-line options.

Definition at line 154 of file petsc_solver.cc.

◆ control()

SolverControl & SolverBase< VectorType >::control ( ) const
inherited

Access to object that controls convergence.

Definition at line 168 of file petsc_solver.cc.

◆ initialize()

void SolverBase< VectorType >::initialize ( const PreconditionBase preconditioner)
inherited

initialize the solver with the preconditioner. This function is intended for use with SLEPc spectral transformation class.

Definition at line 213 of file petsc_solver.cc.

◆ convergence_test()

int SolverBase< VectorType >::convergence_test ( KSP  ksp,
const PetscInt  iteration,
const PetscReal  residual_norm,
KSPConvergedReason *  reason,
void *  solver_control 
)
staticprivateinherited

A function that is used in PETSc as a callback to check on convergence. It takes the information provided from PETSc and checks it against deal.II's own SolverControl objects to see if convergence has been reached.

Definition at line 175 of file petsc_solver.cc.

Member Data Documentation

◆ additional_data

const AdditionalData PETScWrappers::SolverTFQMR::additional_data
protected

Store a copy of the flags for this particular solver.

Definition at line 679 of file petsc_solver.h.

◆ solver_control

SolverControl& PETScWrappers::SolverBase::solver_control
protectedinherited

Reference to the object that controls convergence of the iterative solver. In fact, for these PETSc wrappers, PETSc does so itself, but we copy the data from this object before starting the solution process, and copy the data back into it afterwards.

Definition at line 181 of file petsc_solver.h.

◆ mpi_communicator

const MPI_Comm PETScWrappers::SolverBase::mpi_communicator
protectedinherited

Copy of the MPI communicator object to be used for the solver.

Definition at line 186 of file petsc_solver.h.

◆ prefix_name

std::string PETScWrappers::SolverBase::prefix_name
protectedinherited

Solver prefix name to qualify options specific to the PETSc KSP object in the current context. Note: A hyphen (-) must NOT be given at the beginning of the prefix name. The first character of all runtime options is AUTOMATICALLY the hyphen.

Definition at line 201 of file petsc_solver.h.

◆ solver_data

std::unique_ptr<SolverData> PETScWrappers::SolverBase::solver_data
privateinherited

Pointer to an object that stores the solver context. This is recreated in the main solver routine if necessary.

Definition at line 251 of file petsc_solver.h.


The documentation for this class was generated from the following files: