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

Class PETScLevelSolver is an abstract LinearSolver for solving systems of linear equations on a single SAMRAI::hier::PatchLevel using PETSc. More...

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

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

Public Member Functions

 PETScLevelSolver ()
 Default constructor.
 
 ~PETScLevelSolver ()
 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.
 
const KSP & getPETScKSP () const
 Get the PETSc KSP object.
 
void getASMSubdomains (std::vector< IS > **nonoverlapping_subdomains, std::vector< IS > **overlapping_subdomains)
 Get ASM subdomains.
 
- 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 void setHomogeneousBc (bool homogeneous_bc)
 Set whether the solver should use homogeneous boundary conditions.
 
virtual bool getHomogeneousBc () const
 Return whether the solver is using homogeneous boundary conditions.
 
virtual void setSolutionTime (double solution_time)
 Set the time at which the solution is to be evaluated.
 
virtual double getSolutionTime () const
 Get the time at which the solution is being evaluated.
 
virtual void setTimeInterval (double current_time, double new_time)
 Set the current time interval.
 
virtual std::pair< double, double > getTimeInterval () const
 Get the current time interval.
 
virtual double getDt () const
 Get the current time step size.
 
virtual void setHierarchyMathOps (SAMRAI::tbox::Pointer< HierarchyMathOps > hier_math_ops)
 Set the HierarchyMathOps object used by the solver.
 
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.
 

Protected Attributes

PETSc objects.
std::string d_ksp_type = KSPGMRES
 
std::string d_pc_type = PCILU
 
std::string d_shell_pc_type
 
std::string d_options_prefix
 
KSP d_petsc_ksp = nullptr
 
Mat d_petsc_mat = nullptr
 
Mat d_petsc_pc = nullptr
 
MatNullSpace d_petsc_nullsp
 
Vec d_petsc_x = nullptr
 
Vec d_petsc_b = nullptr
 
Support for additive and multiplicative Schwarz preconditioners.
Vec d_local_x
 
Vec d_local_y
 
SAMRAI::hier::IntVector< NDIM > d_box_size
 
SAMRAI::hier::IntVector< NDIM > d_overlap_size
 
int d_n_local_subdomains
 
int d_n_subdomains_max
 
std::vector< IS > d_overlap_is
 
std::vector< IS > d_nonoverlap_is
 
std::vector< IS > d_local_overlap_is
 
std::vector< IS > d_local_nonoverlap_is
 
std::vector< VecScatter > d_restriction
 
std::vector< VecScatter > d_prolongation
 
std::vector< KSP > d_sub_ksp
 
Mat * d_sub_mat
 
Mat * d_sub_bc_mat
 
std::vector< Vec > d_sub_x
 
std::vector< Vec > d_sub_y
 
- 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
 

Linear solver functionality.

SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > d_hierarchy
 Associated hierarchy.
 
int d_level_num = IBTK::invalid_level_number
 Associated patch level and C-F boundary (for level numbers > 0).
 
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > > d_level
 
SAMRAI::tbox::Pointer< SAMRAI::hier::CoarseFineBoundary< NDIM > > d_cf_boundary
 
SAMRAIDataCache d_cached_eulerian_data
 Scratch data.
 
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.
 
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...
 
void init (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, const std::string &default_options_prefix)
 Basic initialization.
 
virtual void generateASMSubdomains (std::vector< std::set< int > > &overlap_is, std::vector< std::set< int > > &nonoverlap_is)
 Generate IS/subdomains for Schwartz type preconditioners.
 
virtual void generateFieldSplitSubdomains (std::vector< std::string > &field_names, std::vector< std::set< int > > &field_is)
 Generate IS/subdomains for fieldsplit type preconditioners.
 
virtual void initializeSolverStateSpecialized (const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &b)=0
 Compute hierarchy dependent data required for solving $Ax=b$.
 
virtual void deallocateSolverStateSpecialized ()=0
 Remove all hierarchy dependent data allocated by initializeSolverStateSpecialized().
 
virtual void copyToPETScVec (Vec &petsc_x, SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x)=0
 Copy a generic vector to the PETSc representation.
 
virtual void copyFromPETScVec (Vec &petsc_x, SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x)=0
 Copy a generic vector from the PETSc representation.
 
virtual void setupKSPVecs (Vec &petsc_x, Vec &petsc_b, SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &b)=0
 Copy solution and right-hand-side data to the PETSc representation, including any modifications to account for boundary conditions.
 
virtual void setupNullspace ()
 Setup the solver nullspace (if any).
 

Field split preconditioner.

std::vector< std::stringd_field_name
 
std::vector< IS > d_field_is
 

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)
 

Detailed Description

Class PETScLevelSolver is an abstract LinearSolver for solving systems of linear equations on a single SAMRAI::hier::PatchLevel using PETSc.

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 at the Argonne National Laboratory Mathematics and Computer Science Division. For more information about PETSc, see http://www.mcs.anl.gov/petsc.

Member Function Documentation

◆ deallocateSolverState()

void IBTK::PETScLevelSolver::deallocateSolverState ( )
overridevirtual

Remove all hierarchy dependent data allocated by initializeSolverState().

Note
It is safe to call deallocateSolverState() when the solver state is already deallocated.
Subclasses of class PETScLevelSolver should not override this method. Instead, they should override the protected method deallocatedSolverStateSpecialized().
See also
initializeSolverState

Reimplemented from IBTK::GeneralSolver.

◆ initializeSolverState()

void IBTK::PETScLevelSolver::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.

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.
Subclasses of class PETScLevelSolver should not override this method. Instead, they should override the protected method initializeSolverStateSpecialized().
See also
deallocateSolverState

Reimplemented from IBTK::GeneralSolver.

◆ solveSystem()

bool IBTK::PETScLevelSolver::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: