IBAMR  IBAMR version 0.19.
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
IBTK::CCPoissonHypreLevelSolver Class Reference

Class CCPoissonHypreLevelSolver is a concrete LinearSolver for solving elliptic equations of the form \( \mbox{$L u$} = \mbox{$(C I + \nabla \cdot D \nabla) u$} = f \) on a single SAMRAI::hier::PatchLevel using hypre. More...

#include <ibtk/CCPoissonHypreLevelSolver.h>

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

Public Member Functions

 CCPoissonHypreLevelSolver (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, const std::string &default_options_prefix)
 Constructor. More...
 
 ~CCPoissonHypreLevelSolver ()
 Destructor. More...
 
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. More...
 
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< PoissonSolverallocate_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. More...
 

Protected Member Functions

void initSpecialized (const std::string &object_name, bool homogeneous_bc) override
 

Protected Attributes

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
 

Private Attributes

Problem specification.
bool d_grid_aligned_anisotropy = true
 
hypre objects.
unsigned int d_depth = 0
 
HYPRE_StructGrid d_grid = nullptr
 
HYPRE_StructStencil d_stencil = nullptr
 
std::vector< HYPRE_StructMatrix > d_matrices
 
std::vector< HYPRE_StructVector > d_rhs_vecs
 
std::vector< HYPRE_StructVector > d_sol_vecs
 
std::vector< HYPRE_StructSolver > d_solvers
 
std::vector< HYPRE_StructSolver > d_preconds
 
std::vector< SAMRAI::hier::Index< NDIM > > d_stencil_offsets
 
std::string d_solver_type = "PFMG"
 
std::string d_precond_type = "none"
 
int d_rel_change = 0
 
int d_num_pre_relax_steps = 1
 
int d_num_post_relax_steps = 1
 
int d_memory_use = 0
 
int d_rap_type
 
int d_relax_type
 
int d_skip_relax = 1
 
int d_two_norm = 1
 

Linear solver functionality.

SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > d_hierarchy
 Associated hierarchy. More...
 
int d_level_num = IBTK::invalid_level_number
 Associated patch level and C-F boundary (for level numbers > 0). More...
 
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchLevel< NDIM > > d_level
 
SAMRAI::tbox::Pointer< SAMRAI::hier::CoarseFineBoundary< NDIM > > d_cf_boundary
 
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...
 
 CCPoissonHypreLevelSolver ()=delete
 Default constructor. More...
 
 CCPoissonHypreLevelSolver (const CCPoissonHypreLevelSolver &from)=delete
 Copy constructor. More...
 
CCPoissonHypreLevelSolveroperator= (const CCPoissonHypreLevelSolver &that)=delete
 Assignment operator. More...
 
void allocateHypreData ()
 Functions to allocate, initialize, access, and deallocate hypre data structures. More...
 
void setMatrixCoefficients_aligned ()
 
void setMatrixCoefficients_nonaligned ()
 
void setupHypreSolver ()
 
bool solveSystem (int x_idx, int b_idx)
 
void destroyHypreSolver ()
 
void deallocateHypreData ()
 

Functions to access solver parameters.

virtual void setMaxIterations (int max_iterations)
 Set the maximum number of nonlinear iterations to use per solve. More...
 
virtual int getMaxIterations () const
 Get the maximum number of nonlinear iterations to use per solve. More...
 
virtual void setAbsoluteTolerance (double abs_residual_tol)
 Set the absolute residual tolerance for convergence. More...
 
virtual double getAbsoluteTolerance () const
 Get the absolute residual tolerance for convergence. More...
 
virtual void setRelativeTolerance (double rel_residual_tol)
 Set the relative residual tolerance for convergence. More...
 
virtual double getRelativeTolerance () const
 Get the relative residual tolerance for convergence. More...
 

Logging functions.

virtual void setLoggingEnabled (bool enable_logging=true)
 Enable or disable logging. More...
 
virtual bool getLoggingEnabled () const
 Determine whether logging is enabled or disabled. More...
 
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
 
void init (const std::string &object_name, bool homogeneous_bc)
 
virtual void initSpecialized (const std::string &object_name, bool homogeneous_bc)
 

Functions to access solver parameters.

virtual void setInitialGuessNonzero (bool initial_guess_nonzero=true)
 Set whether the initial guess is non-zero. More...
 
virtual bool getInitialGuessNonzero () const
 Get whether the initial guess is non-zero. More...
 

Logging functions.

virtual void printClassData (std::ostream &stream) override
 Print class data to stream. More...
 
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
 

Detailed Description

This solver class uses the hypre library to solve linear equations of the form \( (C I + \nabla \cdot D \nabla ) u = f \), where \(C\) is a cell-centered array, \(D\) is a side-centered array, and \(u\) and \(f\) 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:

  1. Create a CCPoissonHypreLevelSolver object.
  2. Set the problem specification via setPoissonSpecifications(), setPhysicalBcCoef(), and setHomogeneousBc().
  3. Initialize CCPoissonHypreLevelSolver object using the function initializeSolverState().
  4. Solve the linear system using the member function solveSystem(), passing in SAMRAI::solv::SAMRAIVectorReal objects corresponding to \(u\) and \(f\).

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.

Constructor & Destructor Documentation

◆ CCPoissonHypreLevelSolver() [1/3]

IBTK::CCPoissonHypreLevelSolver::CCPoissonHypreLevelSolver ( const std::string &  object_name,
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database input_db,
const std::string &  default_options_prefix 
)

◆ ~CCPoissonHypreLevelSolver()

IBTK::CCPoissonHypreLevelSolver::~CCPoissonHypreLevelSolver ( )

◆ CCPoissonHypreLevelSolver() [2/3]

IBTK::CCPoissonHypreLevelSolver::CCPoissonHypreLevelSolver ( )
privatedelete
Note
This constructor is not implemented and should not be used.

◆ CCPoissonHypreLevelSolver() [3/3]

IBTK::CCPoissonHypreLevelSolver::CCPoissonHypreLevelSolver ( const CCPoissonHypreLevelSolver from)
privatedelete
Note
This constructor is not implemented and should not be used.
Parameters
fromThe value to copy to this object.

Member Function Documentation

◆ allocate_solver()

static SAMRAI::tbox::Pointer<PoissonSolver> IBTK::CCPoissonHypreLevelSolver::allocate_solver ( const std::string &  object_name,
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database input_db,
const std::string &  default_options_prefix 
)
inlinestatic

◆ solveSystem() [1/2]

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

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.

◆ initializeSolverState()

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

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.
See also
deallocateSolverState

Reimplemented from IBTK::GeneralSolver.

◆ deallocateSolverState()

void IBTK::CCPoissonHypreLevelSolver::deallocateSolverState ( )
overridevirtual
Note
It is safe to call deallocateSolverState() when the solver state is already deallocated.
See also
initializeSolverState

Reimplemented from IBTK::GeneralSolver.

◆ operator=()

CCPoissonHypreLevelSolver& IBTK::CCPoissonHypreLevelSolver::operator= ( const CCPoissonHypreLevelSolver that)
privatedelete
Note
This operator is not implemented and should not be used.
Parameters
thatThe value to assign to this object.
Returns
A reference to this object.

◆ allocateHypreData()

void IBTK::CCPoissonHypreLevelSolver::allocateHypreData ( )
private

◆ setMatrixCoefficients_aligned()

void IBTK::CCPoissonHypreLevelSolver::setMatrixCoefficients_aligned ( )
private

◆ setMatrixCoefficients_nonaligned()

void IBTK::CCPoissonHypreLevelSolver::setMatrixCoefficients_nonaligned ( )
private

◆ setupHypreSolver()

void IBTK::CCPoissonHypreLevelSolver::setupHypreSolver ( )
private

◆ solveSystem() [2/2]

bool IBTK::CCPoissonHypreLevelSolver::solveSystem ( int  x_idx,
int  b_idx 
)
private

◆ destroyHypreSolver()

void IBTK::CCPoissonHypreLevelSolver::destroyHypreSolver ( )
private

◆ deallocateHypreData()

void IBTK::CCPoissonHypreLevelSolver::deallocateHypreData ( )
private

◆ setNullSpace()

virtual void IBTK::LinearSolver::setNullSpace ( bool  nullspace_contains_constant_vec,
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, double > > > &  nullspace_basis_vecs = std::vector< SAMRAI::tbox::PointerSAMRAI::solv::SAMRAIVectorReal< NDIM, double > > >() 
)
virtualinherited

Implementations can require the nullspace basis vectors to be orthogonal but should not assume the basis vectors to be orthonormal. If the basis vectors are not orthonormal, the solver may normalize them in place.

Reimplemented in IBTK::PETScKrylovLinearSolver, and IBTK::PETScLevelSolver.

◆ getNullSpaceContainsConstantVector()

virtual bool IBTK::LinearSolver::getNullSpaceContainsConstantVector ( ) const
virtualinherited

◆ getNullSpaceBasisVectors()

virtual const std::vector<SAMRAI::tbox::Pointer<SAMRAI::solv::SAMRAIVectorReal<NDIM, double> > >& IBTK::LinearSolver::getNullSpaceBasisVectors ( ) const
virtualinherited

◆ setInitialGuessNonzero()

virtual void IBTK::LinearSolver::setInitialGuessNonzero ( bool  initial_guess_nonzero = true)
virtualinherited

◆ getInitialGuessNonzero()

virtual bool IBTK::LinearSolver::getInitialGuessNonzero ( ) const
virtualinherited

Reimplemented in IBTK::PETScPCLSWrapper.

◆ printClassData()

virtual void IBTK::LinearSolver::printClassData ( std::ostream &  stream)
overridevirtualinherited

Reimplemented from IBTK::GeneralSolver.

◆ getName()

const std::string& IBTK::GeneralSolver::getName ( ) const
inherited

◆ getIsInitialized()

virtual bool IBTK::GeneralSolver::getIsInitialized ( ) const
virtualinherited

◆ setHomogeneousBc()

virtual void IBTK::GeneralSolver::setHomogeneousBc ( bool  homogeneous_bc)
virtualinherited

◆ getHomogeneousBc()

virtual bool IBTK::GeneralSolver::getHomogeneousBc ( ) const
virtualinherited

◆ setSolutionTime()

virtual void IBTK::GeneralSolver::setSolutionTime ( double  solution_time)
virtualinherited

◆ getSolutionTime()

virtual double IBTK::GeneralSolver::getSolutionTime ( ) const
virtualinherited

◆ setTimeInterval()

virtual void IBTK::GeneralSolver::setTimeInterval ( double  current_time,
double  new_time 
)
virtualinherited

◆ getTimeInterval()

virtual std::pair<double, double> IBTK::GeneralSolver::getTimeInterval ( ) const
virtualinherited

◆ getDt()

virtual double IBTK::GeneralSolver::getDt ( ) const
virtualinherited

◆ setHierarchyMathOps()

virtual void IBTK::GeneralSolver::setHierarchyMathOps ( SAMRAI::tbox::Pointer< HierarchyMathOps hier_math_ops)
virtualinherited

◆ getHierarchyMathOps()

virtual SAMRAI::tbox::Pointer<HierarchyMathOps> IBTK::GeneralSolver::getHierarchyMathOps ( ) const
virtualinherited

◆ setMaxIterations()

virtual void IBTK::GeneralSolver::setMaxIterations ( int  max_iterations)
virtualinherited

◆ getMaxIterations()

virtual int IBTK::GeneralSolver::getMaxIterations ( ) const
virtualinherited

Reimplemented in IBTK::PETScPCLSWrapper.

◆ setAbsoluteTolerance()

virtual void IBTK::GeneralSolver::setAbsoluteTolerance ( double  abs_residual_tol)
virtualinherited

Reimplemented in IBTK::PETScPCLSWrapper.

◆ getAbsoluteTolerance()

virtual double IBTK::GeneralSolver::getAbsoluteTolerance ( ) const
virtualinherited

Reimplemented in IBTK::PETScPCLSWrapper.

◆ setRelativeTolerance()

virtual void IBTK::GeneralSolver::setRelativeTolerance ( double  rel_residual_tol)
virtualinherited

Reimplemented in IBTK::PETScPCLSWrapper.

◆ getRelativeTolerance()

virtual double IBTK::GeneralSolver::getRelativeTolerance ( ) const
virtualinherited

Reimplemented in IBTK::PETScPCLSWrapper.

◆ getNumIterations()

virtual int IBTK::GeneralSolver::getNumIterations ( ) const
virtualinherited

◆ getResidualNorm()

virtual double IBTK::GeneralSolver::getResidualNorm ( ) const
virtualinherited

◆ setLoggingEnabled()

virtual void IBTK::GeneralSolver::setLoggingEnabled ( bool  enable_logging = true)
virtualinherited

◆ getLoggingEnabled()

virtual bool IBTK::GeneralSolver::getLoggingEnabled ( ) const
virtualinherited

◆ init()

void IBTK::GeneralSolver::init ( const std::string &  object_name,
bool  homogeneous_bc 
)
protectedinherited

◆ initSpecialized() [1/2]

virtual void IBTK::GeneralSolver::initSpecialized ( const std::string &  object_name,
bool  homogeneous_bc 
)
protectedvirtualinherited

Reimplemented in IBTK::PoissonSolver.

◆ setPoissonSpecifications()

virtual void IBTK::PoissonSolver::setPoissonSpecifications ( const SAMRAI::solv::PoissonSpecifications poisson_spec)
virtualinherited

◆ setPhysicalBcCoef()

virtual void IBTK::PoissonSolver::setPhysicalBcCoef ( SAMRAI::solv::RobinBcCoefStrategy< NDIM > *  bc_coef)
virtualinherited
Note
bc_coef may be nullptr. In this case, default boundary conditions (as supplied to the class constructor) are employed.
Parameters
bc_coefPointer to an object that can set the Robin boundary condition coefficients

Reimplemented in IBTK::PoissonFACPreconditioner, and IBTK::KrylovLinearSolverPoissonSolverInterface.

◆ setPhysicalBcCoefs()

virtual void IBTK::PoissonSolver::setPhysicalBcCoefs ( const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &  bc_coefs)
virtualinherited
Note
Any of the elements of bc_coefs may be nullptr. In this case, default boundary conditions (as supplied to the class constructor) are employed for that data depth.
Parameters
bc_coefsVector of pointers to objects that can set the Robin boundary condition coefficients

Reimplemented in IBTK::PoissonFACPreconditioner, and IBTK::KrylovLinearSolverPoissonSolverInterface.

◆ initSpecialized() [2/2]

void IBTK::PoissonSolver::initSpecialized ( const std::string &  object_name,
bool  homogeneous_bc 
)
overrideprotectedvirtualinherited

Reimplemented from IBTK::GeneralSolver.

Member Data Documentation

◆ d_hierarchy

SAMRAI::tbox::Pointer<SAMRAI::hier::PatchHierarchy<NDIM> > IBTK::CCPoissonHypreLevelSolver::d_hierarchy
private

◆ d_level_num

int IBTK::CCPoissonHypreLevelSolver::d_level_num = IBTK::invalid_level_number
private

◆ d_level

SAMRAI::tbox::Pointer<SAMRAI::hier::PatchLevel<NDIM> > IBTK::CCPoissonHypreLevelSolver::d_level
private

◆ d_cf_boundary

SAMRAI::tbox::Pointer<SAMRAI::hier::CoarseFineBoundary<NDIM> > IBTK::CCPoissonHypreLevelSolver::d_cf_boundary
private

◆ d_grid_aligned_anisotropy

bool IBTK::CCPoissonHypreLevelSolver::d_grid_aligned_anisotropy = true
private

◆ d_depth

unsigned int IBTK::CCPoissonHypreLevelSolver::d_depth = 0
private

◆ d_grid

HYPRE_StructGrid IBTK::CCPoissonHypreLevelSolver::d_grid = nullptr
private

◆ d_stencil

HYPRE_StructStencil IBTK::CCPoissonHypreLevelSolver::d_stencil = nullptr
private

◆ d_matrices

std::vector<HYPRE_StructMatrix> IBTK::CCPoissonHypreLevelSolver::d_matrices
private

◆ d_rhs_vecs

std::vector<HYPRE_StructVector> IBTK::CCPoissonHypreLevelSolver::d_rhs_vecs
private

◆ d_sol_vecs

std::vector<HYPRE_StructVector> IBTK::CCPoissonHypreLevelSolver::d_sol_vecs
private

◆ d_solvers

std::vector<HYPRE_StructSolver> IBTK::CCPoissonHypreLevelSolver::d_solvers
private

◆ d_preconds

std::vector<HYPRE_StructSolver> IBTK::CCPoissonHypreLevelSolver::d_preconds
private

◆ d_stencil_offsets

std::vector<SAMRAI::hier::Index<NDIM> > IBTK::CCPoissonHypreLevelSolver::d_stencil_offsets
private

◆ d_solver_type

std::string IBTK::CCPoissonHypreLevelSolver::d_solver_type = "PFMG"
private

◆ d_precond_type

std::string IBTK::CCPoissonHypreLevelSolver::d_precond_type = "none"
private

◆ d_rel_change

int IBTK::CCPoissonHypreLevelSolver::d_rel_change = 0
private

◆ d_num_pre_relax_steps

int IBTK::CCPoissonHypreLevelSolver::d_num_pre_relax_steps = 1
private

◆ d_num_post_relax_steps

int IBTK::CCPoissonHypreLevelSolver::d_num_post_relax_steps = 1
private

◆ d_memory_use

int IBTK::CCPoissonHypreLevelSolver::d_memory_use = 0
private

◆ d_rap_type

int IBTK::CCPoissonHypreLevelSolver::d_rap_type
private

◆ d_relax_type

int IBTK::CCPoissonHypreLevelSolver::d_relax_type
private

◆ d_skip_relax

int IBTK::CCPoissonHypreLevelSolver::d_skip_relax = 1
private

◆ d_two_norm

int IBTK::CCPoissonHypreLevelSolver::d_two_norm = 1
private

◆ d_initial_guess_nonzero

bool IBTK::LinearSolver::d_initial_guess_nonzero = true
protectedinherited

◆ d_nullspace_contains_constant_vec

bool IBTK::LinearSolver::d_nullspace_contains_constant_vec = false
protectedinherited

◆ d_nullspace_basis_vecs

std::vector<SAMRAI::tbox::Pointer<SAMRAI::solv::SAMRAIVectorReal<NDIM, double> > > IBTK::LinearSolver::d_nullspace_basis_vecs
protectedinherited

◆ d_object_name

std::string IBTK::GeneralSolver::d_object_name = "unitialized"
protectedinherited

◆ d_is_initialized

bool IBTK::GeneralSolver::d_is_initialized = false
protectedinherited

◆ d_homogeneous_bc

bool IBTK::GeneralSolver::d_homogeneous_bc = false
protectedinherited

◆ d_solution_time

double IBTK::GeneralSolver::d_solution_time = std::numeric_limits<double>::quiet_NaN()
protectedinherited

◆ d_current_time

double IBTK::GeneralSolver::d_current_time = std::numeric_limits<double>::quiet_NaN()
protectedinherited

◆ d_new_time

double IBTK::GeneralSolver::d_new_time = std::numeric_limits<double>::quiet_NaN()
protectedinherited

◆ d_rel_residual_tol

double IBTK::GeneralSolver::d_rel_residual_tol = 0.0
protectedinherited

◆ d_abs_residual_tol

double IBTK::GeneralSolver::d_abs_residual_tol = 0.0
protectedinherited

◆ d_max_iterations

int IBTK::GeneralSolver::d_max_iterations = 100
protectedinherited

◆ d_current_iterations

int IBTK::GeneralSolver::d_current_iterations = 0
protectedinherited

◆ d_current_residual_norm

double IBTK::GeneralSolver::d_current_residual_norm = std::numeric_limits<double>::quiet_NaN()
protectedinherited

◆ d_hier_math_ops

SAMRAI::tbox::Pointer<HierarchyMathOps> IBTK::GeneralSolver::d_hier_math_ops
protectedinherited

◆ d_hier_math_ops_external

bool IBTK::GeneralSolver::d_hier_math_ops_external = false
protectedinherited

◆ d_enable_logging

bool IBTK::GeneralSolver::d_enable_logging = false
protectedinherited

◆ d_poisson_spec

SAMRAI::solv::PoissonSpecifications IBTK::PoissonSolver::d_poisson_spec = SAMRAI::solv::PoissonSpecifications("")
protectedinherited

◆ d_default_bc_coef

std::unique_ptr<SAMRAI::solv::RobinBcCoefStrategy<NDIM> > IBTK::PoissonSolver::d_default_bc_coef
protectedinherited

◆ d_bc_coefs

std::vector<SAMRAI::solv::RobinBcCoefStrategy<NDIM>*> IBTK::PoissonSolver::d_bc_coefs
protectedinherited

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