Reference documentation for deal.II version 9.1.0-pre
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | Static Private Member Functions | Private Attributes | Friends | List of all members
PETScWrappers::SolverBase Class Referenceabstract

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

Inheritance diagram for PETScWrappers::SolverBase:
[legend]

Classes

struct  SolverData
 

Public Member Functions

 SolverBase (SolverControl &cn, const MPI_Comm &mpi_communicator)
 
virtual ~SolverBase ()=default
 
void solve (const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionerBase &preconditioner)
 
virtual void reset ()
 
void set_prefix (const std::string &prefix)
 
SolverControlcontrol () const
 
void initialize (const PreconditionerBase &preconditioner)
 

Protected Member Functions

virtual void set_solver_type (KSP &ksp) const =0
 

Protected Attributes

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
 

Friends

class SLEPcWrappers::TransformationBase
 

Detailed Description

Base class for solver classes using the PETSc solvers. Since solvers in PETSc are selected based on flags passed to a generic solver object, basically all the actual solver calls happen in this class, and derived classes simply set the right flags to select one solver or another, or to set certain parameters for individual solvers.

Optionally, the user can create a solver derived from the SolverBase class and can set the default arguments necessary to solve the linear system of equations with SolverControl. These default options can be overridden by specifying command line arguments of the form -ksp_*. For example, -ksp_monitor_true_residual prints out true residual norm (unpreconditioned) at each iteration and -ksp_view provides information about the linear solver and the preconditioner used in the current context. The type of the solver can also be changed during runtime by specifying -ksp_type {richardson, cg, gmres, fgmres, ..} to dynamically test the optimal solver along with a suitable preconditioner set using -pc_type {jacobi, bjacobi, ilu, lu, ..}. There are several other command line options available to modify the behavior of the PETSc linear solver and can be obtained from the documentation and manual pages.

Note
Repeated calls to solve() on a solver object with a Preconditioner must be used with care. The preconditioner is initialized in the first call to solve() and subsequent calls reuse the solver and preconditioner object. This is done for performance reasons. The solver and preconditioner can be reset by calling reset().

One of the gotchas of PETSc is that – in particular in MPI mode – it often does not produce very helpful error messages. In order to save other users some time in searching a hard to track down error, here is one situation and the error message one gets there: when you don't specify an MPI communicator to your solver's constructor. In this case, you will get an error of the following form from each of your parallel processes:

*   [1]PETSC ERROR: PCSetVector() line 1173 in src/ksp/pc/interface/precon.c
*   [1]PETSC ERROR:   Arguments must have same communicators!
*   [1]PETSC ERROR:   Different communicators in the two objects: Argument #
* 1 and 2! [1]PETSC ERROR: KSPSetUp() line 195 in
* src/ksp/ksp/interface/itfunc.c
* 

This error, on which one can spend a very long time figuring out what exactly goes wrong, results from not specifying an MPI communicator. Note that the communicator must match that of the matrix and all vectors in the linear system which we want to solve. Aggravating the situation is the fact that the default argument to the solver classes, PETSC_COMM_SELF, is the appropriate argument for the sequential case (which is why it is the default argument), so this error only shows up in parallel mode.

Author
Wolfgang Bangerth, 2004

Definition at line 108 of file petsc_solver.h.

Constructor & Destructor Documentation

PETScWrappers::SolverBase::SolverBase ( SolverControl cn,
const MPI_Comm &  mpi_communicator 
)

Constructor. Takes the solver control object and the MPI communicator over which parallel computations are to happen.

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 45 of file petsc_solver.cc.

virtual PETScWrappers::SolverBase::~SolverBase ( )
virtualdefault

Destructor.

Member Function Documentation

void PETScWrappers::SolverBase::solve ( const MatrixBase A,
VectorBase x,
const VectorBase b,
const PreconditionerBase preconditioner 
)

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.

void PETScWrappers::SolverBase::reset ( )
virtual

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

Definition at line 161 of file petsc_solver.cc.

void PETScWrappers::SolverBase::set_prefix ( const std::string &  prefix)

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.

SolverControl & PETScWrappers::SolverBase::control ( ) const

Access to object that controls convergence.

Definition at line 168 of file petsc_solver.cc.

void PETScWrappers::SolverBase::initialize ( const PreconditionerBase preconditioner)

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.

virtual void PETScWrappers::SolverBase::set_solver_type ( KSP &  ksp) const
protectedpure virtual
int PETScWrappers::SolverBase::convergence_test ( KSP  ksp,
const PetscInt  iteration,
const PetscReal  residual_norm,
KSPConvergedReason *  reason,
void *  solver_control 
)
staticprivate

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.

Friends And Related Function Documentation

friend class SLEPcWrappers::TransformationBase
friend

Make the transformation class a friend, since it needs to set the KSP solver.

Definition at line 254 of file petsc_solver.h.

Member Data Documentation

SolverControl& PETScWrappers::SolverBase::solver_control
protected

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 177 of file petsc_solver.h.

const MPI_Comm PETScWrappers::SolverBase::mpi_communicator
protected

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

Definition at line 182 of file petsc_solver.h.

std::string PETScWrappers::SolverBase::prefix_name
protected

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 197 of file petsc_solver.h.

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

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

Definition at line 247 of file petsc_solver.h.


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