|
IBAMR
IBAMR version 0.19.
|
Class SCPoissonHypreLevelSolver 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/SCPoissonHypreLevelSolver.h>

Public Member Functions | |
| SCPoissonHypreLevelSolver (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, const std::string &default_options_prefix) | |
| Constructor. More... | |
| ~SCPoissonHypreLevelSolver () | |
| 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< 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 SCPoissonHypreLevelSolver. 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 |
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< HierarchyMathOps > | d_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 |
hypre objects. | |
| HYPRE_SStructGrid | d_grid = nullptr |
| HYPRE_SStructStencil | d_stencil [NVARS] |
| HYPRE_SStructGraph | d_graph = nullptr |
| HYPRE_SStructMatrix | d_matrix = nullptr |
| HYPRE_SStructVector | d_rhs_vec = nullptr |
| HYPRE_SStructVector | d_sol_vec = nullptr |
| HYPRE_SStructSolver | d_solver = nullptr |
| HYPRE_SStructSolver | d_precond = nullptr |
| std::vector< SAMRAI::hier::Index< NDIM > > | d_stencil_offsets |
| std::string | d_solver_type = "Split" |
| std::string | d_precond_type = "none" |
| std::string | d_split_solver_type = "PFMG" |
| int | d_rel_change = 0 |
| int | d_num_pre_relax_steps = 1 |
| int | d_num_post_relax_steps = 1 |
| int | d_relax_type |
| int | d_skip_relax = 1 |
| int | d_two_norm = 1 |
| static const int | PART = 0 |
| static const int | NPARTS = 1 |
| static const int | NVARS = NDIM |
| static const int | X_VAR = 0 |
| static const int | Y_VAR = 1 |
| static const int | Z_VAR = 2 |
This solver class uses the hypre library to solve linear equations of the form \( (C I + \nabla \cdot D \nabla ) u = f \), where \(C\) and \(D\) are constants, and \(u\) and \(f\) are side-centered arrays. The discretization is second-order accurate.
Robin boundary conditions may be specified through the interface class SAMRAI::solv::RobinBcCoefStrategy, but boundary conditions must be of either Dirichlet or Neumann type. In particular, mixed boundary conditions are not presently supported.
The user must perform the following steps to use class SCPoissonHypreLevelSolver:
Sample parameters for initialization from database (and their default values):
enable_logging = FALSE // see setLoggingEnabled() solver_type = "Split" // choices are: "Split", "SysPFMG", "PCG", "GMRES", "FlexGMRES" , "LGMRES", "BiCGSTAB" precond_type = "none" // choices are: "Split", "SysPFMG" split_solver_type = "PFMG" // choices are: "PFMG", "SMG", "Jacobi" 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 (only used by SysPFMG or PCG solver) num_pre_relax_steps = 1 // number of pre-sweeps (only used by SysPFMG solver) num_post_relax_steps = 1 // number of post-sweeps (only used by SysPFMG solver) relax_type = 1 // see hypre User's Manual (only used by SysPFMG solver or preconditioner) skip_relax = 1 // see hypre User's Manual (only used by SysPFMG 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.
| IBTK::SCPoissonHypreLevelSolver::SCPoissonHypreLevelSolver | ( | const std::string & | object_name, |
| SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | input_db, | ||
| const std::string & | default_options_prefix | ||
| ) |
| IBTK::SCPoissonHypreLevelSolver::~SCPoissonHypreLevelSolver | ( | ) |
|
privatedelete |
|
privatedelete |
| from | The value to copy to this object. |
|
inlinestatic |
|
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.
| 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.
|
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:
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 |
Reimplemented from IBTK::GeneralSolver.
|
privatedelete |
| that | The value to assign to this object. |
|
private |
|
private |
|
private |
|
private |
|
private |
|
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.
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
Reimplemented in IBTK::PETScPCLSWrapper.
|
overridevirtualinherited |
Reimplemented from IBTK::GeneralSolver.
|
inherited |
|
virtualinherited |
|
virtualinherited |
Reimplemented in IBTK::FACPreconditioner, IBTK::NewtonKrylovSolver, and IBTK::KrylovLinearSolver.
|
virtualinherited |
|
virtualinherited |
Reimplemented in IBAMR::CIBStaggeredStokesSolver, IBTK::FACPreconditioner, IBTK::NewtonKrylovSolver, and IBTK::KrylovLinearSolver.
|
virtualinherited |
|
virtualinherited |
Reimplemented in IBAMR::CIBStaggeredStokesSolver, IBTK::FACPreconditioner, IBTK::NewtonKrylovSolver, and IBTK::KrylovLinearSolver.
|
virtualinherited |
|
virtualinherited |
Reimplemented in IBTK::NewtonKrylovSolver, and IBTK::KrylovLinearSolver.
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
Reimplemented in IBTK::PETScPCLSWrapper.
|
virtualinherited |
Reimplemented in IBTK::PETScPCLSWrapper.
|
virtualinherited |
Reimplemented in IBTK::PETScPCLSWrapper.
|
virtualinherited |
Reimplemented in IBTK::PETScPCLSWrapper.
|
virtualinherited |
Reimplemented in IBTK::PETScPCLSWrapper.
|
virtualinherited |
Reimplemented in IBTK::BGaussSeidelPreconditioner, IBTK::PETScPCLSWrapper, and IBTK::BJacobiPreconditioner.
|
virtualinherited |
Reimplemented in IBTK::BGaussSeidelPreconditioner, IBTK::PETScPCLSWrapper, and IBTK::BJacobiPreconditioner.
|
virtualinherited |
|
virtualinherited |
|
protectedinherited |
|
protectedvirtualinherited |
Reimplemented in IBTK::PoissonSolver.
|
virtualinherited |
Reimplemented in IBTK::PoissonFACPreconditioner, and IBTK::KrylovLinearSolverPoissonSolverInterface.
|
virtualinherited |
| bc_coef | Pointer to an object that can set the Robin boundary condition coefficients |
Reimplemented in IBTK::PoissonFACPreconditioner, and IBTK::KrylovLinearSolverPoissonSolverInterface.
|
virtualinherited |
| bc_coefs | Vector of pointers to objects that can set the Robin boundary condition coefficients |
Reimplemented in IBTK::PoissonFACPreconditioner, and IBTK::KrylovLinearSolverPoissonSolverInterface.
|
overrideprotectedvirtualinherited |
Reimplemented from IBTK::GeneralSolver.
|
private |
|
private |
|
private |
|
private |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
1.8.17