IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
|
Class CCPoissonHypreLevelSolver is a concrete LinearSolver for solving elliptic equations of the form on a single SAMRAI::hier::PatchLevel using hypre. More...
#include </home/runner/work/IBAMR/IBAMR/ibtk/include/ibtk/CCPoissonHypreLevelSolver.h>
Public Member Functions | |
CCPoissonHypreLevelSolver (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, const std::string &default_options_prefix) | |
Constructor. | |
~CCPoissonHypreLevelSolver () | |
Destructor. | |
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 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::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 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. | |
Public Member Functions inherited from IBTK::PoissonSolver | |
PoissonSolver ()=default | |
Default constructor. | |
~PoissonSolver ()=default | |
Default destructor. | |
virtual void | setPoissonSpecifications (const SAMRAI::solv::PoissonSpecifications &poisson_spec) |
Set the SAMRAI::solv::PoissonSpecifications object used to specify the coefficients for the scalar-valued or vector-valued Laplace operator. | |
virtual void | setPhysicalBcCoef (SAMRAI::solv::RobinBcCoefStrategy< NDIM > *bc_coef) |
Set the SAMRAI::solv::RobinBcCoefStrategy object used to specify physical boundary conditions. More... | |
virtual void | setPhysicalBcCoefs (const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &bc_coefs) |
Set the SAMRAI::solv::RobinBcCoefStrategy objects used to specify physical boundary conditions. More... | |
Static Public Member Functions | |
static SAMRAI::tbox::Pointer< PoissonSolver > | allocate_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 CCPoissonHypreLevelSolver. | |
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... | |
Additional Inherited Members | |
Protected Member Functions inherited from IBTK::GeneralSolver | |
void | init (const std::string &object_name, bool homogeneous_bc) |
Protected Member Functions inherited from IBTK::PoissonSolver | |
void | initSpecialized (const std::string &object_name, bool homogeneous_bc) override |
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 |
Protected Attributes inherited from IBTK::PoissonSolver | |
SAMRAI::solv::PoissonSpecifications | d_poisson_spec = SAMRAI::solv::PoissonSpecifications("") |
std::unique_ptr< SAMRAI::solv::RobinBcCoefStrategy< NDIM > > | d_default_bc_coef |
std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > | d_bc_coefs |
Class CCPoissonHypreLevelSolver is a concrete LinearSolver for solving elliptic equations of the form on a single SAMRAI::hier::PatchLevel using hypre.
This solver class uses the hypre library to solve linear equations of the form , where is a cell-centered array, is a side-centered array, and and are cell-centered arrays. The discretization is second-order accurate.
Robin boundary conditions may be specified through the interface class SAMRAI::solv::RobinBcCoefStrategy.
The user must perform the following steps to use class CCPoissonHypreLevelSolver:
Sample parameters for initialization from database (and their default values):
enable_logging = FALSE // see setLoggingEnabled() solver_type = "PFMG" // choices are: "PFMG", "SMG", "PCG", "GMRES", "FlexGMRES" , "LGMRES", "BiCGSTAB" precond_type = "none" // choices are: "PFMG", "SMG", "Jacobi", "none" max_iterations = 25 // see setMaxIterations() abs_residual_tol = 1.e-50 // see setAbsoluteTolerance() (only used by hypre Krylov solvers) rel_residual_tol = 1.0e-5 // see setRelativeTolerance() initial_guess_nonzero = FALSE // see setInitialGuessNonzero() rel_change = 0 // see hypre User's Manual num_pre_relax_steps = 1 // number of pre-sweeps (only used by SMG or PFMG solver or preconditioner) num_post_relax_steps = 1 // number of post-sweeps (only used by SMG or PFMG solver or preconditioner) memory_use = 0 // see hypre User's Manual (only used by SMG solver or preconditioner) rap_type = 0 // see hypre User's Manual (only used by PFMG solver or preconditioner) relax_type = 1 // see hypre User's Manual (only used by PFMG solver or preconditioner) skip_relax = 1 // see hypre User's Manual (only used by PFMG solver or preconditioner) two_norm = 1 // see hypre User's Manual (only used by PCG solver)
hypre is developed in the Center for Applied Scientific Computing (CASC) at Lawrence Livermore National Laboratory (LLNL). For more information about hypre, see https://computation.llnl.gov/casc/linear_solvers/sls_hypre.html.
|
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.