IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
|
Class BGaussSeidelPreconditioner is a block Gauss-Seidel preconditioner which implements the abstract LinearSolver interface. More...
#include </home/runner/work/IBAMR/IBAMR/ibtk/include/ibtk/BGaussSeidelPreconditioner.h>
Public Member Functions | |
BGaussSeidelPreconditioner (std::string object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, const std::string &default_options_prefix) | |
Constructor. | |
~BGaussSeidelPreconditioner () | |
Destructor. | |
void | setComponentPreconditioner (SAMRAI::tbox::Pointer< LinearSolver > preconditioner, unsigned int component) |
Set the preconditioner to be employed on the specified vector component. | |
void | setComponentOperators (const std::vector< SAMRAI::tbox::Pointer< LinearOperator > > &linear_ops, unsigned int component) |
Set the linear operators to be employed on the specified vector component. | |
void | setSymmetricPreconditioner (bool symmetric_preconditioner) |
Indicate whether to apply the component preconditioners symmetrically. | |
void | setReversedOrder (bool reverse_order) |
Indicate whether to apply the component preconditioners in reversed order (i.e., starting with the last component and ending with the first component). | |
Linear solver functionality. | |
bool | solveSystem (SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &b) override |
Solve the linear system of equations for . 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 . More... | |
void | deallocateSolverState () override |
Remove all hierarchy dependent data allocated by initializeSolverState(). More... | |
Functions to access solver parameters. | |
void | setInitialGuessNonzero (bool initial_guess_nonzero=true) override |
Set whether the initial guess is non-zero. | |
void | setMaxIterations (int max_iterations) override |
Set the maximum number of iterations to use per solve. | |
Public Member Functions inherited from IBTK::LinearSolver | |
LinearSolver () | |
Constructor. | |
virtual | ~LinearSolver () |
Empty virtual destructor. | |
virtual void | setNullspace (bool nullspace_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 > > >()) |
Set the nullspace of the linear system. More... | |
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 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::string & | getName () 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< HierarchyMathOps > | getHierarchyMathOps () const |
Get the HierarchyMathOps object used by the solver. | |
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 void | setLoggingEnabled (bool enable_logging=true) |
Enable or disable logging. | |
virtual bool | getLoggingEnabled () const |
Determine whether logging is enabled or disabled. | |
Functions to access data on the most recent solve. | |
int | getNumIterations () const override |
Return the iteration count from the most recent linear solve. | |
double | getResidualNorm () const override |
Return the residual norm from the most recent iteration. | |
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::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< HierarchyMathOps > | d_hier_math_ops |
bool | d_hier_math_ops_external = false |
bool | d_enable_logging = false |
Class BGaussSeidelPreconditioner is a block Gauss-Seidel preconditioner which implements the abstract LinearSolver interface.
This solver class performs a single block Gauss-Seidel sweep, applying specified component LinearSolver and LinearOperator objects to the components of a supplied SAMRAI::solv::SAMRAIVectorReal vector, and doing so in a multiplicative fashion. Note that the block Gauss-Seidel algorithm is not generally convergent, but can be used as a preconditioner for a KrylovLinearSolver.
Note that the default block Gauss-Seidel algorithm is not a symmetric linear operator, even if the individual component linear operators and solvers are symmetric. Instead, the algorithm applies the component preconditioners to the vector components starting with the first vector component and ending with the last vector component. The algorithm can be symmetrized via the setSymmetricPreconditioner() member function, and the order in which vector components are visited can be reversed via the setReversedOrder() member function.
Sample parameters for initialization from database (and their default values):
symmetric_preconditioner = FALSE // see setSymmetricPreconditioner() reverse_order = FALSE // see setReversedOrder() initial_guess_nonzero = FALSE // see setInitialGuessNonzero() rel_residual_tol = 1.0e-6 // see setRelativeTolerance() abs_residual_tol = 1.0e-30 // see setAbsoluteTolerance() max_iterations = 1 // see setMaxIterations()
|
overridevirtual |
Remove all hierarchy dependent data allocated by initializeSolverState().
Reimplemented from IBTK::GeneralSolver.
|
overridevirtual |
Compute hierarchy dependent data required for solving .
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:
Note that it is generally necessary to reinitialize the solver state when the hierarchy configuration changes.
x | solution vector |
b | right-hand-side vector |
Conditions on Parameters:
Reimplemented from IBTK::GeneralSolver.
|
overridevirtual |
Solve the linear system of equations for .
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.
x | solution vector |
b | right-hand-side vector |
Conditions on Parameters:
true
if the solver converged to the specified tolerances, false
otherwise Implements IBTK::GeneralSolver.