Reference documentation for deal.II version 9.5.0
\(\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\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Protected Attributes | Static Private Member Functions | Friends | List of all members
PETScWrappers::SolverBase Class Reference

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

Inheritance diagram for PETScWrappers::SolverBase:
[legend]

Public Member Functions

 SolverBase ()
 
 SolverBase (SolverControl &cn)
 
virtual ~SolverBase ()
 
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)
 
KSP petsc_ksp ()
 
 operator KSP () const
 

Protected Member Functions

void initialize_ksp_with_comm (const MPI_Comm comm)
 
virtual void set_solver_type (KSP &ksp) const
 
void perhaps_set_convergence_test () const
 

Protected Attributes

KSP ksp
 
SmartPointer< SolverControl, SolverBasesolver_control
 
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)
 

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
Command line options relative to convergence tolerances are ignored when the solver is instantiated with a SolverControl.
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() or by using the command line option "-ksp_reuse_preconditioner false".

Definition at line 92 of file petsc_solver.h.

Constructor & Destructor Documentation

◆ SolverBase() [1/2]

SolverBase< VectorType >::SolverBase ( )

Default constructor without SolverControl. The resulting solver will use PETSc's default convergence testing routines.

Definition at line 42 of file petsc_solver.cc.

◆ SolverBase() [2/2]

SolverBase< VectorType >::SolverBase ( SolverControl cn)

Constructor with a SolverControl instance.

Definition at line 48 of file petsc_solver.cc.

◆ ~SolverBase()

SolverBase< VectorType >::~SolverBase ( )
virtual

Destructor.

Definition at line 61 of file petsc_solver.cc.

Member Function Documentation

◆ solve()

void SolverBase< VectorType >::solve ( const MatrixBase A,
VectorBase x,
const VectorBase b,
const PreconditionBase 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 84 of file petsc_solver.cc.

◆ reset()

void SolverBase< VectorType >::reset ( )
virtual

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

Definition at line 149 of file petsc_solver.cc.

◆ set_prefix()

void SolverBase< VectorType >::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 142 of file petsc_solver.cc.

◆ control()

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

Access to object that controls convergence.

Definition at line 156 of file petsc_solver.cc.

◆ initialize()

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

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

Definition at line 232 of file petsc_solver.cc.

◆ petsc_ksp()

KSP SolverBase< VectorType >::petsc_ksp ( )

Return the PETSc KSP object.

Definition at line 69 of file petsc_solver.cc.

◆ operator KSP()

SolverBase< VectorType >::operator KSP ( ) const

Conversion operator to gain access to the underlying PETSc type. If you do this, you cut this class off some information it may need, so this conversion operator should only be used if you know what you do.

Definition at line 76 of file petsc_solver.cc.

◆ initialize_ksp_with_comm()

void SolverBase< VectorType >::initialize_ksp_with_comm ( const MPI_Comm  comm)
protected

Utility to create the KSP object and attach convergence test.

Definition at line 207 of file petsc_solver.cc.

◆ set_solver_type()

void SolverBase< VectorType >::set_solver_type ( KSP &  ksp) const
protectedvirtual

◆ perhaps_set_convergence_test()

void SolverBase< VectorType >::perhaps_set_convergence_test ( ) const
protected

Utility to use deal.II convergence testing.

This call changes the convergence criterion when the instance of the class has a SolverControl object associated.

Definition at line 223 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 
)
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 167 of file petsc_solver.cc.

Friends And Related Symbol Documentation

◆ SLEPcWrappers::TransformationBase

friend class SLEPcWrappers::TransformationBase
friend

Definition at line 226 of file petsc_solver.h.

Member Data Documentation

◆ ksp

KSP PETScWrappers::SolverBase::ksp
protected

The PETSc KSP object.

Definition at line 168 of file petsc_solver.h.

◆ solver_control

SmartPointer<SolverControl, SolverBase> 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 176 of file petsc_solver.h.

◆ prefix_name

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


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