Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
amrex::HypreABecLap Class Referencefinal

Struct-interface HYPRE solver for (a - b div beta grad) phi = rhs. More...

#include <AMReX_HypreABecLap.H>

Inheritance diagram for amrex::HypreABecLap:
amrex::Hypre

Public Member Functions

 HypreABecLap (const BoxArray &grids, const DistributionMapping &dmap, const Geometry &geom_, MPI_Comm comm_)
 Bind the solver to a single AMR level described by grids/dmap.
 
 ~HypreABecLap () override
 
 HypreABecLap (HypreABecLap const &)=delete
 
 HypreABecLap (HypreABecLap &&)=delete
 
HypreABecLapoperator= (HypreABecLap const &)=delete
 
HypreABecLapoperator= (HypreABecLap &&)=delete
 
void solve (MultiFab &soln, const MultiFab &rhs, Real reltol, Real abstol, int maxiter, const BndryData &bndry, int max_bndry_order) override
 Solve the struct-interface linear system with the configured coefficients.
 
void prepareSolver ()
 Build the HYPRE structures and apply runtime parameters.
 
void loadVectors (MultiFab &soln, const MultiFab &rhs)
 Copy AMReX solution guess soln and RHS rhs into HYPRE vectors prior to a solve.
 
- 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
 
Hypreoperator= (Hypre &&) noexcept=default
 
 Hypre (Hypre const &)=delete
 
Hypreoperator= (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 }
 

Detailed Description

Struct-interface HYPRE solver for (a - b div beta grad) phi = rhs.

Constructor & Destructor Documentation

◆ HypreABecLap() [1/3]

amrex::HypreABecLap::HypreABecLap ( const BoxArray grids,
const DistributionMapping dmap,
const Geometry geom_,
MPI_Comm  comm_ 
)

Bind the solver to a single AMR level described by grids/dmap.

Parameters
gridsGrid layout for the level.
dmapDistributionMapping for grids.
geom_Geometry describing the level.
comm_MPI communicator handed to HYPRE.

◆ ~HypreABecLap()

amrex::HypreABecLap::~HypreABecLap ( )
override

◆ HypreABecLap() [2/3]

amrex::HypreABecLap::HypreABecLap ( HypreABecLap const &  )
delete

◆ HypreABecLap() [3/3]

amrex::HypreABecLap::HypreABecLap ( HypreABecLap &&  )
delete

Member Function Documentation

◆ loadVectors()

void amrex::HypreABecLap::loadVectors ( MultiFab soln,
const MultiFab rhs 
)

Copy AMReX solution guess soln and RHS rhs into HYPRE vectors prior to a solve.

◆ operator=() [1/2]

HypreABecLap & amrex::HypreABecLap::operator= ( HypreABecLap &&  )
delete

◆ operator=() [2/2]

HypreABecLap & amrex::HypreABecLap::operator= ( HypreABecLap const &  )
delete

◆ prepareSolver()

void amrex::HypreABecLap::prepareSolver ( )

Build the HYPRE structures and apply runtime parameters.

◆ solve()

void amrex::HypreABecLap::solve ( MultiFab soln,
const MultiFab rhs,
Real  reltol,
Real  abstol,
int  maxiter,
const BndryData bndry,
int  max_bndry_order 
)
overridevirtual

Solve the struct-interface linear system with the configured coefficients.

Parameters
solnSolution MultiFab (filled on return).
rhsRight-hand side MultiFab.
reltolRelative tolerance.
abstolAbsolute tolerance.
maxiterMaximum number of iterations.
bndryBoundary data.
max_bndry_orderMaximum order at the boundary.

Implements amrex::Hypre.


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