![]() |
Block-Structured AMR Software Framework
|
HYPRE IJ solver for ABecLaplacian (optionally with overset masks). More...
#include <AMReX_HypreABecLap3.H>
Public Member Functions | |
| HypreABecLap3 (const BoxArray &grids, const DistributionMapping &dmap, const Geometry &geom_, MPI_Comm comm_, const iMultiFab *overset_mask=nullptr) | |
| Construct the IJ-based solver and capture overset/EB context. | |
| ~HypreABecLap3 () override | |
| HypreABecLap3 (HypreABecLap3 const &)=delete | |
| HypreABecLap3 (HypreABecLap3 &&)=delete | |
| HypreABecLap3 & | operator= (HypreABecLap3 const &)=delete |
| HypreABecLap3 & | operator= (HypreABecLap3 &&)=delete |
| void | solve (MultiFab &soln, const MultiFab &rhs, Real rel_tol, Real abs_tol, int max_iter, const BndryData &bndry, int max_bndry_order) override |
| Solve with the IJ wrapper (supports overset/EB contributions). | |
| void | setEBDirichlet (MultiFab const *beb) |
Provide b coefficients stored in beb for embedded-boundary faces. | |
| void | prepareSolver () |
| Build IJ structures and configure the HYPRE solver/preconditioner. | |
| void | getSolution (MultiFab &soln) |
Copy HYPRE solution back into AMReX storage (soln). | |
| void | loadVectors (MultiFab &soln, const MultiFab &rhs) |
Load RHS rhs and initial guess soln into HYPRE vectors. | |
Public Member Functions inherited from amrex::Hypre | |
| 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). | |
| 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. | |
Additional Inherited Members | |
Public Types inherited from amrex::Hypre | |
| enum class | Interface : int { structed , semi_structed , ij } |
| HYPRE interface modes supported. More... | |
Static Public Member Functions inherited from amrex::Hypre | |
| 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 Public Attributes inherited from amrex::Hypre | |
| 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). | |
Protected Attributes inherited from amrex::Hypre | |
| MPI_Comm | comm = MPI_COMM_NULL |
| Geometry | geom |
| int | verbose = 0 |
| bool | old_default = true |
| int | relax_type = 6 |
| int | relax_order = 1 |
| int | num_sweeps = 2 |
| Real | strong_threshold = Real(0.25) |
| std::string | options_namespace {"hypre"} |
| MultiFab | acoefs |
| Array< MultiFab, 3 > | bcoefs |
| Real | scalar_a |
| Real | scalar_b |
| MultiFab | diaginv |
| FabFactory< FArrayBox > const * | m_factory = nullptr |
| BndryData const * | m_bndry = nullptr |
| int | m_maxorder = -1 |
| bool | is_matrix_singular { false } |
HYPRE IJ solver for ABecLaplacian (optionally with overset masks).
| amrex::HypreABecLap3::HypreABecLap3 | ( | const BoxArray & | grids, |
| const DistributionMapping & | dmap, | ||
| const Geometry & | geom_, | ||
| MPI_Comm | comm_, | ||
| const iMultiFab * | overset_mask = nullptr |
||
| ) |
Construct the IJ-based solver and capture overset/EB context.
| grids | Grid layout for the AMR level. |
| dmap | Distribution map for grids. |
| geom_ | Geometry describing the domain. |
| comm_ | MPI communicator forwarded to HYPRE. |
| overset_mask | Optional overset mask used to drop covered cells. |
|
overridedefault |
|
delete |
|
delete |
| void amrex::HypreABecLap3::getSolution | ( | MultiFab & | soln | ) |
Copy HYPRE solution back into AMReX storage (soln).
Load RHS rhs and initial guess soln into HYPRE vectors.
|
delete |
|
delete |
| void amrex::HypreABecLap3::prepareSolver | ( | ) |
Build IJ structures and configure the HYPRE solver/preconditioner.
|
inline |
Provide b coefficients stored in beb for embedded-boundary faces.
|
overridevirtual |
Solve with the IJ wrapper (supports overset/EB contributions).
| soln | Solution MultiFab (updated in place). |
| rhs | Right-hand side data. |
| rel_tol | Relative tolerance. |
| abs_tol | Absolute tolerance. |
| max_iter | Maximum iterations. |
| bndry | Boundary data. |
| max_bndry_order | Maximum order at the boundary. |
Implements amrex::Hypre.