![]() |
Block-Structured AMR Software Framework
|
Multilevel HYPRE solver that mirrors MLMG's ABecLaplacian structure. More...
#include <AMReX_HypreMLABecLap.H>
Public Member Functions | |
| HypreMLABecLap (HypreMLABecLap const &)=delete | |
| HypreMLABecLap (HypreMLABecLap &&)=delete | |
| HypreMLABecLap & | operator= (HypreMLABecLap const &)=delete |
| HypreMLABecLap & | operator= (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. | |
Multilevel HYPRE solver that mirrors MLMG's ABecLaplacian structure.
|
delete |
|
delete |
| 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.
| a_geom | Per-level geometries (fine to coarse). |
| a_grids | Per-level BoxArrays. |
| a_dmap | Per-level DistributionMappings. |
| a_hypre_solver_id | Solver type (BoomerAMG or SSAMG). |
| a_parmparse_prefix | ParmParse prefix for runtime options. |
| amrex::HypreMLABecLap::~HypreMLABecLap | ( | ) |
Communicate b coefficients needed by level flev at coarse/fine interfaces.
| flev | Finest AMR level whose coarse/fine stencils are being prepared. |
| a_bcoefs | Face-centered coefficients for that level (per direction). |
| 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.
| flev | AMR level whose coarse/fine data are exchanged. |
| a_bcoefs | Face arrays for that level. |
| tags | Cached copy tags extracted from FabArray communications. |
|
delete |
|
delete |
|
inline |
Flag whether the assembled operator is singular (v).
|
inline |
Override the maximum number of HYPRE iterations (v).
| 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.
| a_ascalar | Scalar applied to the diagonal term. |
| a_bscalar | Scalar applied to the diffusive term. |
| a_acoefs | Level-by-level a coefficients. |
| a_bcoefs | Level-by-level arrays of face-centered b coefficients. |
| a_lobc | Low-side boundary-condition types. |
| a_hibc | High-side boundary-condition types. |
| a_levelbcdata | BC data per level. |
| a_coarse_bc | Optional BC data from even coarser AMR level. |
|
inline |
Set HYPRE's logging verbosity (v).