Lightweight base class shared by HYPRE-backed linear solvers, HypreABecLap, HypreABecLap2, and HypreABecLap3.
More...
#include <AMReX_Hypre.H>
|
| | Hypre (const BoxArray &grids, const DistributionMapping &dmap, const Geometry &geom, MPI_Comm comm_) |
| | Construct a HYPRE wrapper bound to the supplied grids.
|
| |
| virtual | ~Hypre () |
| |
| | Hypre (Hypre &&) noexcept=default |
| |
| Hypre & | operator= (Hypre &&) noexcept=default |
| |
| | Hypre (Hypre const &)=delete |
| |
| Hypre & | operator= (Hypre const &)=delete |
| |
| void | setScalars (Real sa, Real sb) |
| | Store constant scalars applied to the a (diagonal) and b (Laplacian) terms.
|
| |
| void | setACoeffs (const MultiFab &alpha) |
| | Install per-cell a coefficients.
|
| |
| void | setBCoeffs (const Array< const MultiFab *, 3 > &beta) |
| | Install per-face b coefficients for each spatial direction.
|
| |
| void | setVerbose (int _verbose) |
| | Set HYPRE verbosity level (_verbose).
|
| |
| void | setIsMatrixSingular (bool flag) |
| | Flag a singular operator (enables null-space handling in derived classes).
|
| |
| virtual void | solve (MultiFab &soln, const MultiFab &rhs, Real rel_tol, Real abs_tol, int max_iter, const BndryData &bndry, int max_bndry_order)=0 |
| | Solve the linear system using the concrete interface implementation.
|
| |
| void | setHypreOptionsNamespace (const std::string &ns) noexcept |
| | Override the ParmParse namespace consulted for HYPRE runtime options (ns).
|
| |
| void | setHypreOldDefault (bool l) noexcept |
| | Toggle HYPRE's "old default" AMG setup (Falgout + modified classical interpolation) via l.
|
| |
| void | setHypreRelaxType (int n) noexcept |
| | Choose the relaxation type handed to HYPRE (see HYPRE docs for valid codes) via n.
|
| |
| void | setHypreRelaxOrder (int n) noexcept |
| | Choose the relaxation ordering (C/F, lexicographic, etc.) via n.
|
| |
| void | setHypreNumSweeps (int n) noexcept |
| | Set the number of pre/post smoothing sweeps HYPRE executes on each level via n.
|
| |
| void | setHypreStrongThreshold (Real t) noexcept |
| | Adjust HYPRE's AMG strong-connection threshold used during coarsening via t.
|
| |
|
| static HYPRE_Int | ispow2 (HYPRE_Int i) |
| | Return 1 if i is a power of two, otherwise 0.
|
| |
| static Array< HYPRE_Int, 3 > | loV (const Box &b) |
| | Convenience helper that casts a box's low corner to HYPRE_Int.
|
| |
| static Array< HYPRE_Int, 3 > | hiV (const Box &b) |
| | Convenience helper that casts a box's high corner to HYPRE_Int.
|
| |
|
| static constexpr HYPRE_Int | regular_stencil_size = 2*3 + 1 |
| | Stencil size used on regular cells.
|
| |
| static constexpr HYPRE_Int | eb_stencil_size = 3 *3 *3 |
| | Stencil size used on EB cells (3^D).
|
| |
Lightweight base class shared by HYPRE-backed linear solvers, HypreABecLap, HypreABecLap2, and HypreABecLap3.
Stores geometry, coefficients, and solver-wide options that concrete implementations (Struct, SStruct, IJ) consume in their respective solve() routines.
◆ Interface
HYPRE interface modes supported.
| Enumerator |
|---|
| structed | |
| semi_structed | |
| ij | |
◆ Hypre() [1/3]
Construct a HYPRE wrapper bound to the supplied grids.
- Parameters
-
| grids | Grid layout on the AMR level being solved. |
| dmap | Distribution map for grids. |
| geom | Geometry describing the domain. |
| comm_ | MPI communicator handed to HYPRE. |
◆ ~Hypre()
◆ Hypre() [2/3]
| amrex::Hypre::Hypre |
( |
Hypre && |
| ) |
|
|
defaultnoexcept |
◆ Hypre() [3/3]
| amrex::Hypre::Hypre |
( |
Hypre const & |
| ) |
|
|
delete |
◆ hiV()
| static Array< HYPRE_Int, 3 > amrex::Hypre::hiV |
( |
const Box & |
b | ) |
|
|
inlinestatic |
Convenience helper that casts a box's high corner to HYPRE_Int.
◆ ispow2()
| static HYPRE_Int amrex::Hypre::ispow2 |
( |
HYPRE_Int |
i | ) |
|
|
inlinestatic |
Return 1 if i is a power of two, otherwise 0.
◆ loV()
| static Array< HYPRE_Int, 3 > amrex::Hypre::loV |
( |
const Box & |
b | ) |
|
|
inlinestatic |
Convenience helper that casts a box's low corner to HYPRE_Int.
◆ operator=() [1/2]
◆ operator=() [2/2]
◆ setACoeffs()
| void amrex::Hypre::setACoeffs |
( |
const MultiFab & |
alpha | ) |
|
Install per-cell a coefficients.
◆ setBCoeffs()
| void amrex::Hypre::setBCoeffs |
( |
const Array< const MultiFab *, 3 > & |
beta | ) |
|
Install per-face b coefficients for each spatial direction.
◆ setHypreNumSweeps()
| void amrex::Hypre::setHypreNumSweeps |
( |
int |
n | ) |
|
|
inlinenoexcept |
Set the number of pre/post smoothing sweeps HYPRE executes on each level via n.
◆ setHypreOldDefault()
| void amrex::Hypre::setHypreOldDefault |
( |
bool |
l | ) |
|
|
inlinenoexcept |
Toggle HYPRE's "old default" AMG setup (Falgout + modified classical interpolation) via l.
◆ setHypreOptionsNamespace()
| void amrex::Hypre::setHypreOptionsNamespace |
( |
const std::string & |
ns | ) |
|
|
inlinenoexcept |
Override the ParmParse namespace consulted for HYPRE runtime options (ns).
◆ setHypreRelaxOrder()
| void amrex::Hypre::setHypreRelaxOrder |
( |
int |
n | ) |
|
|
inlinenoexcept |
Choose the relaxation ordering (C/F, lexicographic, etc.) via n.
◆ setHypreRelaxType()
| void amrex::Hypre::setHypreRelaxType |
( |
int |
n | ) |
|
|
inlinenoexcept |
Choose the relaxation type handed to HYPRE (see HYPRE docs for valid codes) via n.
◆ setHypreStrongThreshold()
| void amrex::Hypre::setHypreStrongThreshold |
( |
Real |
t | ) |
|
|
inlinenoexcept |
Adjust HYPRE's AMG strong-connection threshold used during coarsening via t.
◆ setIsMatrixSingular()
| void amrex::Hypre::setIsMatrixSingular |
( |
bool |
flag | ) |
|
|
inline |
Flag a singular operator (enables null-space handling in derived classes).
◆ setScalars()
| void amrex::Hypre::setScalars |
( |
Real |
sa, |
|
|
Real |
sb |
|
) |
| |
Store constant scalars applied to the a (diagonal) and b (Laplacian) terms.
◆ setVerbose()
| void amrex::Hypre::setVerbose |
( |
int |
_verbose | ) |
|
Set HYPRE verbosity level (_verbose).
◆ solve()
Solve the linear system using the concrete interface implementation.
- Parameters
-
| soln | Solution MultiFab updated in place. |
| rhs | Right-hand side. |
| rel_tol | Relative convergence tolerance. |
| abs_tol | Absolute convergence tolerance. |
| max_iter | Maximum iterations allowed. |
| bndry | Boundary condition data. |
| max_bndry_order | Maximum order at the boundary. |
Implemented in amrex::HypreABecLap3, amrex::HypreABecLap, and amrex::HypreABecLap2.
◆ acoefs
◆ bcoefs
◆ comm
◆ diaginv
◆ eb_stencil_size
| constexpr HYPRE_Int amrex::Hypre::eb_stencil_size = 3 *3 *3 |
|
staticconstexpr |
Stencil size used on EB cells (3^D).
◆ geom
◆ is_matrix_singular
| bool amrex::Hypre::is_matrix_singular { false } |
|
protected |
◆ m_bndry
| BndryData const* amrex::Hypre::m_bndry = nullptr |
|
protected |
◆ m_factory
◆ m_maxorder
| int amrex::Hypre::m_maxorder = -1 |
|
protected |
◆ num_sweeps
| int amrex::Hypre::num_sweeps = 2 |
|
protected |
◆ old_default
| bool amrex::Hypre::old_default = true |
|
protected |
◆ options_namespace
| std::string amrex::Hypre::options_namespace {"hypre"} |
|
protected |
◆ regular_stencil_size
| constexpr HYPRE_Int amrex::Hypre::regular_stencil_size = 2*3 + 1 |
|
staticconstexpr |
Stencil size used on regular cells.
◆ relax_order
| int amrex::Hypre::relax_order = 1 |
|
protected |
◆ relax_type
| int amrex::Hypre::relax_type = 6 |
|
protected |
◆ scalar_a
| Real amrex::Hypre::scalar_a |
|
protected |
◆ scalar_b
| Real amrex::Hypre::scalar_b |
|
protected |
◆ strong_threshold
| Real amrex::Hypre::strong_threshold = Real(0.25) |
|
protected |
◆ verbose
| int amrex::Hypre::verbose = 0 |
|
protected |
The documentation for this class was generated from the following files: