IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
Public Member Functions | Static Public Member Functions | Friends | List of all members
IBTK::PETScKrylovLinearSolver Class Reference

Class PETScKrylovLinearSolver provides a KrylovLinearSolver interface for a PETSc Krylov subspace iterative linear solver (KSP). More...

#include </home/runner/work/IBAMR/IBAMR/ibtk/include/ibtk/PETScKrylovLinearSolver.h>

Inheritance diagram for IBTK::PETScKrylovLinearSolver:
Inheritance graph
[legend]

Public Member Functions

 PETScKrylovLinearSolver (std::string object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, std::string default_options_prefix, MPI_Comm petsc_comm=PETSC_COMM_WORLD)
 Constructor for a concrete KrylovLinearSolver that employs the PETSc KSP solver framework.
 
 PETScKrylovLinearSolver (std::string object_name, const KSP &petsc_ksp)
 Constructor for a concrete KrylovLinearSolver that acts as a "wrapper" for a supplied PETSc KSP object. More...
 
 ~PETScKrylovLinearSolver ()
 Destructor.
 
void setKSPType (const std::string &ksp_type)
 Set the KSP type.
 
void setOptionsPrefix (const std::string &options_prefix)
 Set the options prefix used by this PETSc solver object.
 
Functions to access the underlying PETSc objects.
const KSP & getPETScKSP () const
 Get the PETSc KSP object.
 
- Public Member Functions inherited from IBTK::KrylovLinearSolver
 KrylovLinearSolver ()=default
 Default constructor.
 
 ~KrylovLinearSolver ()=default
 Empty destructor.
 
void setHierarchyMathOps (SAMRAI::tbox::Pointer< HierarchyMathOps > hier_math_ops) override
 Set the HierarchyMathOps object used by the solver.
 
void setHomogeneousBc (bool homogeneous_bc) override
 Set whether the solver should use homogeneous boundary conditions.
 
void setSolutionTime (double solution_time) override
 Set the time at which the solution is to be evaluated.
 
void setTimeInterval (double current_time, double new_time) override
 Set the current time interval.
 
virtual SAMRAI::tbox::Pointer< LinearOperatorgetOperator () const
 Retrieve the linear operator used when solving $Ax=b$.
 
virtual SAMRAI::tbox::Pointer< LinearSolvergetPreconditioner () const
 Retrieve the preconditioner used by the Krylov subspace method when solving $Ax=b$.
 
- Public Member Functions inherited from IBTK::LinearSolver
 LinearSolver ()
 Constructor.
 
virtual ~LinearSolver ()
 Empty virtual destructor.
 
virtual bool getNullspaceContainsConstantVector () const
 Get whether the nullspace of the linear system contains th constant vector.
 
virtual const std::vector< SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > > & getNullspaceBasisVectors () const
 Get the basis vectors for the nullspace of the linear system.
 
virtual void setInitialGuessNonzero (bool initial_guess_nonzero=true)
 Set whether the initial guess is non-zero.
 
virtual bool getInitialGuessNonzero () const
 Get whether the initial guess is non-zero.
 
virtual void printClassData (std::ostream &stream) override
 Print class data to stream.
 
- Public Member Functions inherited from IBTK::GeneralSolver
 GeneralSolver ()=default
 Constructor.
 
virtual ~GeneralSolver ()=default
 Empty virtual destructor.
 
const std::stringgetName () const
 Return the object name.
 
virtual bool getIsInitialized () const
 Return whether the operator is initialized.
 
virtual bool getHomogeneousBc () const
 Return whether the solver is using homogeneous boundary conditions.
 
virtual double getSolutionTime () const
 Get the time at which the solution is being evaluated.
 
virtual std::pair< double, double > getTimeInterval () const
 Get the current time interval.
 
virtual double getDt () const
 Get the current time step size.
 
virtual SAMRAI::tbox::Pointer< HierarchyMathOpsgetHierarchyMathOps () const
 Get the HierarchyMathOps object used by the solver.
 
virtual void setMaxIterations (int max_iterations)
 Set the maximum number of nonlinear iterations to use per solve.
 
virtual int getMaxIterations () const
 Get the maximum number of nonlinear iterations to use per solve.
 
virtual void setAbsoluteTolerance (double abs_residual_tol)
 Set the absolute residual tolerance for convergence.
 
virtual double getAbsoluteTolerance () const
 Get the absolute residual tolerance for convergence.
 
virtual void setRelativeTolerance (double rel_residual_tol)
 Set the relative residual tolerance for convergence.
 
virtual double getRelativeTolerance () const
 Get the relative residual tolerance for convergence.
 
virtual int getNumIterations () const
 Return the iteration count from the most recent solve.
 
virtual double getResidualNorm () const
 Return the residual norm from the most recent iteration.
 
virtual void setLoggingEnabled (bool enable_logging=true)
 Enable or disable logging.
 
virtual bool getLoggingEnabled () const
 Determine whether logging is enabled or disabled.
 

Static Public Member Functions

static SAMRAI::tbox::Pointer< KrylovLinearSolverallocate_solver (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, const std::string &default_options_prefix)
 Static function to construct a PETScKrylovLinearSolver.
 

Friends

class PETScNewtonKrylovSolver
 

Krylov solver functionality.

void setOperator (SAMRAI::tbox::Pointer< LinearOperator > A) override
 Set the linear operator used when solving $Ax=b$.
 
void setPreconditioner (SAMRAI::tbox::Pointer< LinearSolver > pc_solver=NULL) override
 Set the preconditioner used by the Krylov subspace method when solving $Ax=b$. More...
 
void setNullspace (bool contains_constant_vec, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > > &nullspace_basis_vecs=std::vector< SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > >()) override
 Set the nullspace of the linear system. More...
 
bool solveSystem (SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &b) override
 Solve the linear system of equations $Ax=b$ for $x$. More...
 
void initializeSolverState (const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &b) override
 Compute hierarchy dependent data required for solving $Ax=b$. More...
 
void deallocateSolverState () override
 Remove all hierarchy dependent data allocated by initializeSolverState(). More...
 

Additional Inherited Members

- Protected Member Functions inherited from IBTK::GeneralSolver
void init (const std::string &object_name, bool homogeneous_bc)
 
virtual void initSpecialized (const std::string &object_name, bool homogeneous_bc)
 
- Protected Attributes inherited from IBTK::KrylovLinearSolver
SAMRAI::tbox::Pointer< LinearOperatord_A
 
SAMRAI::tbox::Pointer< LinearSolverd_pc_solver
 
SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > d_x
 
SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > d_b
 
- Protected Attributes inherited from IBTK::LinearSolver
bool d_initial_guess_nonzero = true
 
bool d_nullspace_contains_constant_vec = false
 
std::vector< SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > > d_nullspace_basis_vecs
 
- Protected Attributes inherited from IBTK::GeneralSolver
std::string d_object_name = "unitialized"
 
bool d_is_initialized = false
 
bool d_homogeneous_bc = false
 
double d_solution_time = std::numeric_limits<double>::quiet_NaN()
 
double d_current_time = std::numeric_limits<double>::quiet_NaN()
 
double d_new_time = std::numeric_limits<double>::quiet_NaN()
 
double d_rel_residual_tol = 0.0
 
double d_abs_residual_tol = 0.0
 
int d_max_iterations = 100
 
int d_current_iterations = 0
 
double d_current_residual_norm = std::numeric_limits<double>::quiet_NaN()
 
SAMRAI::tbox::Pointer< HierarchyMathOpsd_hier_math_ops
 
bool d_hier_math_ops_external = false
 
bool d_enable_logging = false
 

Detailed Description

Class PETScKrylovLinearSolver provides a KrylovLinearSolver interface for a PETSc Krylov subspace iterative linear solver (KSP).

This solver class provides access to a large number of Krylov subspace solvers for linear problems of the form $Ax=b$ using the PETSc KSP linear solver interface. See http://www.mcs.anl.gov/petsc/documentation/linearsolvertable.html for a complete list of the Krylov solvers provided by this class. Note that solver configuration is typically done at runtime via command line options.

Note
  • Users have the option of supplying the PETScKrylovLinearSolver class constructor with a KSP object. It is important to emphasize that doing so is optional. When a KSP object is provided to the class constructor, the PETScKrylovLinearSolver class acts as a wrapper for the supplied KSP object.
  • The functionality of the PETScKrylovLinearSolver class is essentially identical regardless as to whether a PETSc KSP object is provided to the class constructor. The main exception is that when an external KSP object is provided to the class constructor, memory management of that object is NOT handled by the PETScKrylovLinearSolver. In particular, it is the caller's responsibility to ensure that the supplied KSP object is properly destroyed via KSPDestroy().

Sample parameters for initialization from database (and their default values):

options_prefix = ""           // see setOptionsPrefix()
ksp_type = "gmres"            // see setKSPType()
initial_guess_nonzero = TRUE  // see setInitialGuessNonzero()
rel_residual_tol = 1.0e-5     // see setRelativeTolerance()
abs_residual_tol = 1.0e-50    // see setAbsoluteTolerance()
max_iterations = 10000        // see setMaxIterations()
enable_logging = FALSE        // see setLoggingEnabled()

PETSc is developed in the Mathematics and Computer Science (MCS) Division at Argonne National Laboratory (ANL). For more information about PETSc, see http://www.mcs.anl.gov/petsc.

Constructor & Destructor Documentation

◆ PETScKrylovLinearSolver()

IBTK::PETScKrylovLinearSolver::PETScKrylovLinearSolver ( std::string  object_name,
const KSP &  petsc_ksp 
)

Constructor for a concrete KrylovLinearSolver that acts as a "wrapper" for a supplied PETSc KSP object.

Parameters
object_nameName of the solver
petsc_kspPETSc KSP object
Note
This constructor initializes a PETScKrylovLinearSolver object that acts as a "wrapper" for the provided KSP object. Note that memory management of the provided KSP object is NOT handled by this class.

Member Function Documentation

◆ deallocateSolverState()

void IBTK::PETScKrylovLinearSolver::deallocateSolverState ( )
overridevirtual

Remove all hierarchy dependent data allocated by initializeSolverState().

When linear operator or preconditioner objects have been registered with this class via setOperator() and setPreconditioner(), they are also deallocated by this member function.

Note
It is safe to call deallocateSolverState() when the solver state is already deallocated.
See also
initializeSolverState

Reimplemented from IBTK::GeneralSolver.

◆ initializeSolverState()

void IBTK::PETScKrylovLinearSolver::initializeSolverState ( const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &  x,
const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &  b 
)
overridevirtual

Compute hierarchy dependent data required for solving $Ax=b$.

By default, the solveSystem() method computes some required hierarchy dependent data before solving and removes that data after the solve. For multiple solves that use the same hierarchy configuration, it is more efficient to:

  1. initialize the hierarchy-dependent data required by the solver via initializeSolverState(),
  2. solve the system one or more times via solveSystem(), and
  3. remove the hierarchy-dependent data via deallocateSolverState().

Note that it is generally necessary to reinitialize the solver state when the hierarchy configuration changes.

When linear operator or preconditioner objects have been registered with this class via setOperator() and setPreconditioner(), they are also initialized by this member function.

Parameters
xsolution vector
bright-hand-side vector

Conditions on Parameters:

  • vectors x and b must have same patch hierarchy
  • vectors x and b must have same structure, depth, etc.
Note
The vector arguments for solveSystem() need not match those for initializeSolverState(). However, there must be a certain degree of similarity, including:
  • hierarchy configuration (hierarchy pointer and range of levels)
  • number, type and alignment of vector component data
  • ghost cell widths of data in the solution x and right-hand-side b vectors
Note
It is safe to call initializeSolverState() when the state is already initialized. In this case, the solver state is first deallocated and then reinitialized.
See also
deallocateSolverState

Reimplemented from IBTK::GeneralSolver.

◆ setNullspace()

void IBTK::PETScKrylovLinearSolver::setNullspace ( bool  contains_constant_vec,
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > > &  nullspace_basis_vecs = std::vector<SAMRAI::tbox::Pointer<SAMRAI::solv::SAMRAIVectorReal<NDIM, double> > >() 
)
overridevirtual

Set the nullspace of the linear system.

Basis vectors must be orthogonal but are not required to be orthonormal. Basis vectors will be normalized automatically.

Reimplemented from IBTK::LinearSolver.

◆ setPreconditioner()

void IBTK::PETScKrylovLinearSolver::setPreconditioner ( SAMRAI::tbox::Pointer< LinearSolver pc_solver = NULL)
overridevirtual

Set the preconditioner used by the Krylov subspace method when solving $Ax=b$.

Note
If the preconditioner is NULL, no preconditioning is performed.

Reimplemented from IBTK::KrylovLinearSolver.

◆ solveSystem()

bool IBTK::PETScKrylovLinearSolver::solveSystem ( SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &  x,
SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &  b 
)
overridevirtual

Solve the linear system of equations $Ax=b$ for $x$.

Before calling solveSystem(), the form of the solution x and right-hand-side b vectors must be set properly by the user on all patch interiors on the specified range of levels in the patch hierarchy. The user is responsible for all data management for the quantities associated with the solution and right-hand-side vectors. In particular, patch data in these vectors must be allocated prior to calling this method.

Parameters
xsolution vector
bright-hand-side vector

Conditions on Parameters:

  • vectors x and b must have same patch hierarchy
  • vectors x and b must have same structure, depth, etc.
Note
The vector arguments for solveSystem() need not match those for initializeSolverState(). However, there must be a certain degree of similarity, including:
  • hierarchy configuration (hierarchy pointer and range of levels)
  • number, type and alignment of vector component data
  • ghost cell widths of data in the solution x and right-hand-side b vectors
Note
The solver need not be initialized prior to calling solveSystem(); however, see initializeSolverState() and deallocateSolverState() for opportunities to save overhead when performing multiple consecutive solves.
See also
initializeSolverState
deallocateSolverState
Returns
true if the solver converged to the specified tolerances, false otherwise

Implements IBTK::GeneralSolver.


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