Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
amrex::HypreMLABecLap Class Reference

Multilevel HYPRE solver that mirrors MLMG's ABecLaplacian structure. More...

#include <AMReX_HypreMLABecLap.H>

Public Member Functions

 HypreMLABecLap (HypreMLABecLap const &)=delete
 
 HypreMLABecLap (HypreMLABecLap &&)=delete
 
HypreMLABecLapoperator= (HypreMLABecLap const &)=delete
 
HypreMLABecLapoperator= (HypreMLABecLap &&)=delete
 
 HypreMLABecLap (Vector< Geometry > a_geom, Vector< BoxArray > a_grids, Vector< DistributionMapping > a_dmap, HypreSolverID a_hypre_solver_id, std::string a_parmparse_prefix="hypre_mlabeclap")
 Construct the multilevel HYPRE wrapper and copy the hierarchy metadata.
 
 ~HypreMLABecLap ()
 
void setVerbose (int v)
 Set HYPRE's logging verbosity (v).
 
void setMaxIter (int v)
 Override the maximum number of HYPRE iterations (v).
 
void setIsSingular (bool v)
 Flag whether the assembled operator is singular (v).
 
void setup (Real a_ascalar, Real a_bscalar, Vector< MultiFab const * > const &a_acoefs, Vector< Array< MultiFab const *, 3 > > const &a_bcoefs, Array< LinOpBCType, 3 > const &a_lobc, Array< LinOpBCType, 3 > const &a_hibc, Vector< MultiFab const * > const &a_levelbcdata, std::pair< MultiFab const *, IntVect > const &a_coarse_bc={nullptr, IntVect(0)})
 Assemble the HYPRE operators/storage for the provided coefficients.
 
void solve (Vector< MultiFab * > const &a_sol, Vector< MultiFab const * > const &a_rhs, Real a_reltol, Real a_abstol)
 Solve the assembled system with HYPRE (single-component).
 
void commBCoefs (int flev, Array< MultiFab const *, 3 > const &a_bcoefs)
 Communicate b coefficients needed by level flev at coarse/fine interfaces.
 
void commBCoefs_local (int flev, Array< MultiFab const *, 3 > const &a_bcoefs, Vector< FabArrayBase::CopyComTag > const &tags)
 Local portion of commBCoefs used to issue device-friendly copies.
 

Detailed Description

Multilevel HYPRE solver that mirrors MLMG's ABecLaplacian structure.

Constructor & Destructor Documentation

◆ HypreMLABecLap() [1/3]

amrex::HypreMLABecLap::HypreMLABecLap ( HypreMLABecLap const &  )
delete

◆ HypreMLABecLap() [2/3]

amrex::HypreMLABecLap::HypreMLABecLap ( HypreMLABecLap &&  )
delete

◆ HypreMLABecLap() [3/3]

amrex::HypreMLABecLap::HypreMLABecLap ( Vector< Geometry a_geom,
Vector< BoxArray a_grids,
Vector< DistributionMapping a_dmap,
HypreSolverID  a_hypre_solver_id,
std::string  a_parmparse_prefix = "hypre_mlabeclap" 
)

Construct the multilevel HYPRE wrapper and copy the hierarchy metadata.

Parameters
a_geomPer-level geometries (fine to coarse).
a_gridsPer-level BoxArrays.
a_dmapPer-level DistributionMappings.
a_hypre_solver_idSolver type (BoomerAMG or SSAMG).
a_parmparse_prefixParmParse prefix for runtime options.

◆ ~HypreMLABecLap()

amrex::HypreMLABecLap::~HypreMLABecLap ( )

Member Function Documentation

◆ commBCoefs()

void amrex::HypreMLABecLap::commBCoefs ( int  flev,
Array< MultiFab const *, 3 > const &  a_bcoefs 
)

Communicate b coefficients needed by level flev at coarse/fine interfaces.

Parameters
flevFinest AMR level whose coarse/fine stencils are being prepared.
a_bcoefsFace-centered coefficients for that level (per direction).

◆ commBCoefs_local()

void amrex::HypreMLABecLap::commBCoefs_local ( int  flev,
Array< MultiFab const *, 3 > const &  a_bcoefs,
Vector< FabArrayBase::CopyComTag > const &  tags 
)

Local portion of commBCoefs used to issue device-friendly copies.

Parameters
flevAMR level whose coarse/fine data are exchanged.
a_bcoefsFace arrays for that level.
tagsCached copy tags extracted from FabArray communications.

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ setIsSingular()

void amrex::HypreMLABecLap::setIsSingular ( bool  v)
inline

Flag whether the assembled operator is singular (v).

◆ setMaxIter()

void amrex::HypreMLABecLap::setMaxIter ( int  v)
inline

Override the maximum number of HYPRE iterations (v).

◆ setup()

void amrex::HypreMLABecLap::setup ( Real  a_ascalar,
Real  a_bscalar,
Vector< MultiFab const * > const &  a_acoefs,
Vector< Array< MultiFab const *, 3 > > const &  a_bcoefs,
Array< LinOpBCType, 3 > const &  a_lobc,
Array< LinOpBCType, 3 > const &  a_hibc,
Vector< MultiFab const * > const &  a_levelbcdata,
std::pair< MultiFab const *, IntVect > const &  a_coarse_bc = {nullptr, IntVect(0)} 
)

Assemble the HYPRE operators/storage for the provided coefficients.

Parameters
a_ascalarScalar applied to the diagonal term.
a_bscalarScalar applied to the diffusive term.
a_acoefsLevel-by-level a coefficients.
a_bcoefsLevel-by-level arrays of face-centered b coefficients.
a_lobcLow-side boundary-condition types.
a_hibcHigh-side boundary-condition types.
a_levelbcdataBC data per level.
a_coarse_bcOptional BC data from even coarser AMR level.

◆ setVerbose()

void amrex::HypreMLABecLap::setVerbose ( int  v)
inline

Set HYPRE's logging verbosity (v).

◆ solve()

void amrex::HypreMLABecLap::solve ( Vector< MultiFab * > const &  a_sol,
Vector< MultiFab const * > const &  a_rhs,
Real  a_reltol,
Real  a_abstol 
)

Solve the assembled system with HYPRE (single-component).

Parameters
a_solVector of level solution MultiFabs (filled on return).
a_rhsVector of level right-hand sides.
a_reltolRelative tolerance.
a_abstolAbsolute tolerance.

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