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

Class VCSCViscousPETScLevelSolver is a subclass of SCPoissonPETScLevelSolver class which solves vector-valued elliptic equation of the form $ \mbox{$L u$} = C u + \nabla \cdot \mu (\nabla u + (\nabla u)^T) = f $ on a single SAMRAI::hier::PatchLevel using PETSc. More...

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

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

Public Member Functions

 VCSCViscousPETScLevelSolver (std::string object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, std::string default_options_prefix)
 Constructor.
 
 ~VCSCViscousPETScLevelSolver ()
 Destructor.
 
void setViscosityInterpolationType (IBTK::VCInterpType mu_interp_type)
 Set the interpolation type to be used in interpolating the viscosity.
 
- Public Member Functions inherited from IBTK::SCPoissonPETScLevelSolver
 SCPoissonPETScLevelSolver (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, std::string default_options_prefix)
 Constructor.
 
 ~SCPoissonPETScLevelSolver ()
 Destructor.
 
- Public Member Functions inherited from IBTK::PETScLevelSolver
 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.
 
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...
 
- 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.
 
- 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< 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 VCSCViscousPETScLevelSolver.
 
- Static Public Member Functions inherited from IBTK::SCPoissonPETScLevelSolver
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 SCPoissonPETScLevelSolver.
 

Protected Member Functions

void initializeSolverStateSpecialized (const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &b) override
 Compute hierarchy dependent data required for solving $Ax=b$.
 
void setupKSPVecs (Vec &petsc_x, Vec &petsc_b, SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x, SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &b) override
 Copy solution and right-hand-side data to the PETSc representation, including any modifications to account for boundary conditions.
 
- Protected Member Functions inherited from IBTK::SCPoissonPETScLevelSolver
void generateASMSubdomains (std::vector< std::set< int > > &overlap_is, std::vector< std::set< int > > &nonoverlap_is) override
 Generate IS/subdomains for Schwartz type preconditioners.
 
void deallocateSolverStateSpecialized () override
 Remove all hierarchy dependent data allocated by initializeSolverStateSpecialized().
 
void copyToPETScVec (Vec &petsc_x, SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x) override
 Copy a generic vector to the PETSc representation.
 
void copyFromPETScVec (Vec &petsc_x, SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &x) override
 Copy a generic vector from the PETSc representation.
 
- Protected Member Functions inherited from IBTK::PETScLevelSolver
void init (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, const std::string &default_options_prefix)
 Basic initialization.
 
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 setupNullspace ()
 Setup the solver nullspace (if any).
 
- 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
 

Additional Inherited Members

- Protected Attributes inherited from IBTK::SCPoissonPETScLevelSolver
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContextd_context
 
std::vector< int > d_num_dofs_per_proc
 
int d_dof_index_idx = IBTK::invalid_index
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::SideVariable< NDIM, int > > d_dof_index_var
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > d_data_synch_sched
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > d_ghost_fill_sched
 
- Protected Attributes inherited from IBTK::PETScLevelSolver
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
 
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
 
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.
 
std::vector< std::stringd_field_name
 
std::vector< IS > d_field_is
 
- 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
 
- 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
 

Detailed Description

Class VCSCViscousPETScLevelSolver is a subclass of SCPoissonPETScLevelSolver class which solves vector-valued elliptic equation of the form $ \mbox{$L u$} = C u + \nabla \cdot \mu (\nabla u + (\nabla u)^T) = f $ on a single SAMRAI::hier::PatchLevel using PETSc.

This solver class uses the PETSc library to solve linear equations of the form $ \beta C u + \alpha \nabla \cdot \mu (\nabla u + (\nabla u)^T) = f $, in which $ C $ and $ \mu $ are spacially varying coefficients, and $u$ and $f$ are side-centered arrays. The scaling factors of $ C $ and $ \mu $ are stored separately in the class and are denoted by $ \beta $ and $ \alpha $, respectively. 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 VCSCViscousPETScLevelSolver:

  1. Create a VCSCViscousPETScLevelSolver object.
  2. Set the problem specification via setPoissonSpecifications(), setPhysicalBcCoef(), and setHomogeneousBc().
  3. Initialize VCSCViscousPETScLevelSolver 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):

max_iterations = 10            // see setMaxIterations()
absolute_residual_tol = 0.0    // see setAbsoluteTolerance() (only used by hypre Krylov
solvers)
rel_residual_tol = 1.0e-6      // see setRelativeTolerance()
enable_logging = FALSE         // see setLoggingEnabled()
options_prefix = ""            // see setOptionsPrefix()

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.


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