deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
|
#include <deal.II/lac/petsc_solver.h>
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) |
SolverControl & | control () 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 |
ObserverPointer< SolverControl, SolverBase > | solver_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 |
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.
Definition at line 91 of file petsc_solver.h.
SolverBase< VectorType >::SolverBase | ( | ) |
Default constructor without SolverControl. The resulting solver will use PETSc's default convergence testing routines.
Definition at line 41 of file petsc_solver.cc.
SolverBase< VectorType >::SolverBase | ( | SolverControl & | cn | ) |
Constructor with a SolverControl instance.
Definition at line 47 of file petsc_solver.cc.
|
virtual |
Destructor.
Definition at line 60 of file petsc_solver.cc.
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 83 of file petsc_solver.cc.
|
virtual |
Resets the contained preconditioner and solver object. See class description for more details.
Definition at line 148 of file petsc_solver.cc.
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 141 of file petsc_solver.cc.
SolverControl & SolverBase< VectorType >::control | ( | ) | const |
Access to object that controls convergence.
Definition at line 155 of file petsc_solver.cc.
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 231 of file petsc_solver.cc.
KSP SolverBase< VectorType >::petsc_ksp | ( | ) |
Return the PETSc KSP object.
Definition at line 68 of file petsc_solver.cc.
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 75 of file petsc_solver.cc.
|
protected |
Utility to create the KSP object and attach convergence test.
Definition at line 206 of file petsc_solver.cc.
|
protectedvirtual |
Function that takes a Krylov Subspace Solver context object, and sets the type of solver that is requested by the derived class.
Reimplemented in PETScWrappers::SolverRichardson, PETScWrappers::SolverChebychev, PETScWrappers::SolverCG, PETScWrappers::SolverBiCG, PETScWrappers::SolverGMRES, PETScWrappers::SolverBicgstab, PETScWrappers::SolverCGS, PETScWrappers::SolverTFQMR, PETScWrappers::SolverTCQMR, PETScWrappers::SolverCR, PETScWrappers::SolverLSQR, PETScWrappers::SolverPreOnly, and PETScWrappers::SparseDirectMUMPS.
Definition at line 55 of file petsc_solver.cc.
|
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 222 of file petsc_solver.cc.
|
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 166 of file petsc_solver.cc.
|
friend |
Definition at line 225 of file petsc_solver.h.
|
protected |
The PETSc KSP object.
Definition at line 167 of file petsc_solver.h.
|
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 175 of file petsc_solver.h.
|
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 205 of file petsc_solver.h.