Block-Structured AMR Software Framework
amrex::MLLinOpT< MF > Class Template Referenceabstract

#include <AMReX_MLLinOp.H>

Inheritance diagram for amrex::MLLinOpT< MF >:
amrex::MLCellLinOpT< MF > amrex::MLNodeLinOp amrex::MLCellABecLapT< MF > amrex::MLEBNodeFDLaplacian amrex::MLNodeABecLaplacian amrex::MLNodeLaplacian amrex::MLNodeTensorLaplacian amrex::MLABecLaplacianT< MF > amrex::MLALaplacianT< MF > amrex::MLEBABecLap amrex::MLPoissonT< MF > amrex::MLTensorOp amrex::MLEBTensorOp

Classes

struct  CommContainer
 

Public Types

using MFType = MF
 
using FAB = typename FabDataType< MF >::fab_type
 
using RT = typename FabDataType< MF >::value_type
 
using BCType = LinOpBCType
 
using BCMode = LinOpEnumType::BCMode
 
using StateMode = LinOpEnumType::StateMode
 
using Location = LinOpEnumType::Location
 

Public Member Functions

 MLLinOpT ()=default
 
virtual ~MLLinOpT ()=default
 
 MLLinOpT (const MLLinOpT< MF > &)=delete
 
 MLLinOpT (MLLinOpT< MF > &&)=delete
 
MLLinOpT< MF > & operator= (const MLLinOpT< MF > &)=delete
 
MLLinOpT< MF > & operator= (MLLinOpT< MF > &&)=delete
 
void define (const Vector< Geometry > &a_geom, const Vector< BoxArray > &a_grids, const Vector< DistributionMapping > &a_dmap, const LPInfo &a_info, const Vector< FabFactory< FAB > const * > &a_factory, bool eb_limit_coarsening=true)
 
virtual std::string name () const
 
void setDomainBC (const Array< BCType, AMREX_SPACEDIM > &lobc, const Array< BCType, AMREX_SPACEDIM > &hibc) noexcept
 Boundary of the whole domain. More...
 
void setDomainBC (const Vector< Array< BCType, AMREX_SPACEDIM > > &lobc, const Vector< Array< BCType, AMREX_SPACEDIM > > &hibc) noexcept
 Boundary of the whole domain. More...
 
void setDomainBCLoc (const Array< Real, AMREX_SPACEDIM > &lo_bcloc, const Array< Real, AMREX_SPACEDIM > &hi_bcloc) noexcept
 Set location of domain boundaries. More...
 
bool needsCoarseDataForBC () const noexcept
 Needs coarse data for bc? More...
 
void setCoarseFineBC (const MF *crse, int crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
 Set coarse/fine boundary conditions. For cell-centered solves only. More...
 
void setCoarseFineBC (const MF *crse, IntVect const &crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
 
template<typename AMF , std::enable_if_t<!std::is_same_v< MF, AMF >, int > = 0>
void setCoarseFineBC (const AMF *crse, int crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
 
template<typename AMF , std::enable_if_t<!std::is_same_v< MF, AMF >, int > = 0>
void setCoarseFineBC (const AMF *crse, IntVect const &crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
 
virtual void setLevelBC (int, const MF *, const MF *=nullptr, const MF *=nullptr, const MF *=nullptr)=0
 Set boundary conditions for given level. For cell-centered solves only. More...
 
template<typename AMF , std::enable_if_t<!std::is_same_v< MF, AMF >, int > = 0>
void setLevelBC (int amrlev, const AMF *levelbcdata, const AMF *robinbc_a=nullptr, const AMF *robinbc_b=nullptr, const AMF *robinbc_f=nullptr)
 
void setVerbose (int v) noexcept
 Set verbosity. More...
 
void setMaxOrder (int o) noexcept
 Set order of interpolation at coarse/fine boundary. More...
 
int getMaxOrder () const noexcept
 Get order of interpolation at coarse/fine boundary. More...
 
void setEnforceSingularSolvable (bool o) noexcept
 
bool getEnforceSingularSolvable () const noexcept
 
virtual BottomSolver getDefaultBottomSolver () const
 
virtual int getNComp () const
 Return number of components. More...
 
virtual int getNGrow (int=0, int=0) const
 
virtual bool needsUpdate () const
 Does it need update if it's reused? More...
 
virtual void update ()
 Update for reuse. More...
 
virtual void restriction (int amrlev, int cmglev, MF &crse, MF &fine) const =0
 Restriction onto coarse MG level. More...
 
virtual void interpolation (int amrlev, int fmglev, MF &fine, const MF &crse) const =0
 Add interpolated coarse MG level data to fine MG level data. More...
 
virtual void interpAssign (int amrlev, int fmglev, MF &fine, MF &crse) const
 Overwrite fine MG level data with interpolated coarse data. More...
 
virtual void interpolationAmr (int famrlev, MF &fine, const MF &crse, IntVect const &nghost) const
 Interpolation between AMR levels. More...
 
virtual void averageDownSolutionRHS (int camrlev, MF &crse_sol, MF &crse_rhs, const MF &fine_sol, const MF &fine_rhs)
 Average-down data from fine AMR level to coarse AMR level. More...
 
virtual void apply (int amrlev, int mglev, MF &out, MF &in, BCMode bc_mode, StateMode s_mode, const MLMGBndryT< MF > *bndry=nullptr) const =0
 Apply the linear operator, out = L(in) More...
 
virtual void smooth (int amrlev, int mglev, MF &sol, const MF &rhs, bool skip_fillboundary=false) const =0
 Smooth. More...
 
virtual void normalize (int amrlev, int mglev, MF &mf) const
 Divide mf by the diagonal component of the operator. Used by bicgstab. More...
 
virtual void solutionResidual (int amrlev, MF &resid, MF &x, const MF &b, const MF *crse_bcdata=nullptr)=0
 Compute residual for solution. More...
 
virtual void prepareForFluxes (int, const MF *=nullptr)
 
virtual void correctionResidual (int amrlev, int mglev, MF &resid, MF &x, const MF &b, BCMode bc_mode, const MF *crse_bcdata=nullptr)=0
 Compute residual for the residual-correction form, resid = b - L(x) More...
 
virtual void reflux (int crse_amrlev, MF &res, const MF &crse_sol, const MF &crse_rhs, MF &fine_res, MF &fine_sol, const MF &fine_rhs) const
 Reflux at AMR coarse/fine boundary. More...
 
virtual void compFlux (int, const Array< MF *, AMREX_SPACEDIM > &, MF &, Location) const
 Compute fluxes. More...
 
virtual void compGrad (int, const Array< MF *, AMREX_SPACEDIM > &, MF &, Location) const
 Compute gradients of the solution. More...
 
virtual void applyMetricTerm (int, int, MF &) const
 apply metric terms if there are any More...
 
virtual void unapplyMetricTerm (int, int, MF &) const
 unapply metric terms if there are any More...
 
virtual void unimposeNeumannBC (int, MF &) const
 This is needed for our nodal projection solver. More...
 
virtual void applyInhomogNeumannTerm (int, MF &) const
 Extra terms introduced when we treat inhomogeneous Nuemann BC as homogeneous. More...
 
virtual void applyOverset (int, MF &) const
 for overset solver only More...
 
virtual void scaleRHS (int, MF &) const
 scale RHS to fix solvability More...
 
virtual Vector< RTgetSolvabilityOffset (int, int, MF const &) const
 get offset for fixing solvability More...
 
virtual void fixSolvabilityByOffset (int, int, MF &, Vector< RT > const &) const
 fix solvability by subtracting offset from RHS More...
 
virtual void prepareForSolve ()=0
 
virtual void prepareForGMRES ()
 
virtual void setDirichletNodesToZero (int, int, MF &) const
 
virtual bool isSingular (int amrlev) const =0
 Is it singular on given AMR level? More...
 
virtual bool isBottomSingular () const =0
 Is the bottom of MG singular? More...
 
virtual RT xdoty (int amrlev, int mglev, const MF &x, const MF &y, bool local) const =0
 x dot y, used by the bottom solver More...
 
virtual std::unique_ptr< MLLinOpT< MF > > makeNLinOp (int) const
 
virtual void getFluxes (const Vector< Array< MF *, AMREX_SPACEDIM > > &, const Vector< MF * > &, Location) const
 
virtual void getFluxes (const Vector< MF * > &, const Vector< MF * > &) const
 
virtual bool supportNSolve () const
 
virtual void copyNSolveSolution (MF &, MF const &) const
 
virtual void postSolve (Vector< MF > &) const
 
virtual RT normInf (int amrlev, MF const &mf, bool local) const =0
 
virtual void averageDownAndSync (Vector< MF > &sol) const =0
 
virtual void avgDownResAmr (int clev, MF &cres, MF const &fres) const
 
virtual void avgDownResMG (int clev, MF &cres, MF const &fres) const
 
bool isMFIterSafe (int amrlev, int mglev1, int mglev2) const
 
int NAMRLevels () const noexcept
 Return the number of AMR levels. More...
 
int NMGLevels (int amrlev) const noexcept
 Return the number of MG levels at given AMR level. More...
 
const GeometryGeom (int amr_lev, int mglev=0) const noexcept
 

Public Attributes

Vector< Array< BCType, AMREX_SPACEDIM > > m_lobc
 
Vector< Array< BCType, AMREX_SPACEDIM > > m_hibc
 
Vector< Array< BCType, AMREX_SPACEDIM > > m_lobc_orig
 
Vector< Array< BCType, AMREX_SPACEDIM > > m_hibc_orig
 

Protected Member Functions

const Vector< int > & AMRRefRatio () const noexcept
 Return AMR refinement ratios. More...
 
int AMRRefRatio (int amr_lev) const noexcept
 Return AMR refinement ratio at given AMR level. More...
 
FabFactory< FAB > const * Factory (int amr_lev, int mglev=0) const noexcept
 
GpuArray< BCType, AMREX_SPACEDIM > LoBC (int icomp=0) const noexcept
 
GpuArray< BCType, AMREX_SPACEDIM > HiBC (int icomp=0) const noexcept
 
bool hasBC (BCType bct) const noexcept
 
bool hasInhomogNeumannBC () const noexcept
 
bool hasRobinBC () const noexcept
 
virtual bool supportRobinBC () const noexcept
 
virtual bool supportInhomogNeumannBC () const noexcept
 
bool isBottomActive () const noexcept
 
MPI_Comm BottomCommunicator () const noexcept
 
MPI_Comm Communicator () const noexcept
 
void setCoarseFineBCLocation (const RealVect &cloc) noexcept
 
bool doAgglomeration () const noexcept
 
bool doConsolidation () const noexcept
 
bool doSemicoarsening () const noexcept
 
bool isCellCentered () const noexcept
 
virtual IntVect getNGrowVectRestriction () const
 
virtual void make (Vector< Vector< MF > > &mf, IntVect const &ng) const
 
virtual MF make (int amrlev, int mglev, IntVect const &ng) const
 
virtual MF makeAlias (MF const &mf) const
 
virtual MF makeCoarseMG (int amrlev, int mglev, IntVect const &ng) const
 
virtual MF makeCoarseAmr (int famrlev, IntVect const &ng) const
 
virtual std::unique_ptr< FabFactory< FAB > > makeFactory (int, int) const
 
virtual void resizeMultiGrid (int new_size)
 
bool hasHiddenDimension () const noexcept
 
int hiddenDirection () const noexcept
 
Box compactify (Box const &b) const noexcept
 
template<typename T >
Array4< T > compactify (Array4< T > const &a) const noexcept
 
template<typename T >
get_d0 (T const &dx, T const &dy, T const &) const noexcept
 
template<typename T >
get_d1 (T const &, T const &dy, T const &dz) const noexcept
 

Protected Attributes

int mg_domain_min_width = 2
 
LPInfo info
 
int verbose = 0
 
int maxorder = 3
 
bool enforceSingularSolvable = true
 
int m_num_amr_levels = 0
 
Vector< intm_amr_ref_ratio
 
Vector< intm_num_mg_levels
 
const MLLinOpT< MF > * m_parent = nullptr
 
IntVect m_ixtype
 
bool m_do_agglomeration = false
 
bool m_do_consolidation = false
 
bool m_do_semicoarsening = false
 
Vector< IntVectmg_coarsen_ratio_vec
 
Vector< Vector< Geometry > > m_geom
 first Vector is for amr level and second is mg level More...
 
Vector< Vector< BoxArray > > m_grids
 
Vector< Vector< DistributionMapping > > m_dmap
 
Vector< Vector< std::unique_ptr< FabFactory< FAB > > > > m_factory
 
Vector< intm_domain_covered
 
MPI_Comm m_default_comm = MPI_COMM_NULL
 
MPI_Comm m_bottom_comm = MPI_COMM_NULL
 
std::unique_ptr< CommContainerm_raii_comm
 
Array< Real, AMREX_SPACEDIM > m_domain_bloc_lo {{AMREX_D_DECL(0._rt,0._rt,0._rt)}}
 
Array< Real, AMREX_SPACEDIM > m_domain_bloc_hi {{AMREX_D_DECL(0._rt,0._rt,0._rt)}}
 
bool m_needs_coarse_data_for_bc = false
 
LinOpBCType m_coarse_fine_bc_type = LinOpBCType::Dirichlet
 
IntVect m_coarse_data_crse_ratio = IntVect(-1)
 
RealVect m_coarse_bc_loc
 
const MF * m_coarse_data_for_bc = nullptr
 
MF m_coarse_data_for_bc_raii
 

Static Protected Attributes

static constexpr int mg_coarsen_ratio = 2
 
static constexpr int mg_box_min_width = 2
 

Private Member Functions

void defineGrids (const Vector< Geometry > &a_geom, const Vector< BoxArray > &a_grids, const Vector< DistributionMapping > &a_dmap, const Vector< FabFactory< FAB > const * > &a_factory)
 
void defineBC ()
 
MPI_Comm makeSubCommunicator (const DistributionMapping &dm)
 
virtual void checkPoint (std::string const &) const
 

Static Private Member Functions

static void makeAgglomeratedDMap (const Vector< BoxArray > &ba, Vector< DistributionMapping > &dm)
 
static void makeConsolidatedDMap (const Vector< BoxArray > &ba, Vector< DistributionMapping > &dm, int ratio, int strategy)
 

Private Attributes

Vector< std::unique_ptr< MF > > levelbc_raii
 
Vector< std::unique_ptr< MF > > robin_a_raii
 
Vector< std::unique_ptr< MF > > robin_b_raii
 
Vector< std::unique_ptr< MF > > robin_f_raii
 

Friends

template<typename T >
class MLMGT
 
template<typename T >
class MLCGSolverT
 
template<typename T >
class MLPoissonT
 
template<typename T >
class MLABecLaplacianT
 
template<typename T >
class GMRESMLMGT
 

Member Typedef Documentation

◆ BCMode

template<typename MF >
using amrex::MLLinOpT< MF >::BCMode = LinOpEnumType::BCMode

◆ BCType

template<typename MF >
using amrex::MLLinOpT< MF >::BCType = LinOpBCType

◆ FAB

template<typename MF >
using amrex::MLLinOpT< MF >::FAB = typename FabDataType<MF>::fab_type

◆ Location

template<typename MF >
using amrex::MLLinOpT< MF >::Location = LinOpEnumType::Location

◆ MFType

template<typename MF >
using amrex::MLLinOpT< MF >::MFType = MF

◆ RT

template<typename MF >
using amrex::MLLinOpT< MF >::RT = typename FabDataType<MF>::value_type

◆ StateMode

template<typename MF >
using amrex::MLLinOpT< MF >::StateMode = LinOpEnumType::StateMode

Constructor & Destructor Documentation

◆ MLLinOpT() [1/3]

template<typename MF >
amrex::MLLinOpT< MF >::MLLinOpT ( )
default

◆ ~MLLinOpT()

template<typename MF >
virtual amrex::MLLinOpT< MF >::~MLLinOpT ( )
virtualdefault

◆ MLLinOpT() [2/3]

template<typename MF >
amrex::MLLinOpT< MF >::MLLinOpT ( const MLLinOpT< MF > &  )
delete

◆ MLLinOpT() [3/3]

template<typename MF >
amrex::MLLinOpT< MF >::MLLinOpT ( MLLinOpT< MF > &&  )
delete

Member Function Documentation

◆ AMRRefRatio() [1/2]

template<typename MF >
const Vector<int>& amrex::MLLinOpT< MF >::AMRRefRatio ( ) const
inlineprotectednoexcept

Return AMR refinement ratios.

◆ AMRRefRatio() [2/2]

template<typename MF >
int amrex::MLLinOpT< MF >::AMRRefRatio ( int  amr_lev) const
inlineprotectednoexcept

Return AMR refinement ratio at given AMR level.

◆ apply()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::apply ( int  amrlev,
int  mglev,
MF &  out,
MF &  in,
BCMode  bc_mode,
StateMode  s_mode,
const MLMGBndryT< MF > *  bndry = nullptr 
) const
pure virtual

Apply the linear operator, out = L(in)

Parameters
amrlevAMR level
mglevMG level
outoutput
ininput
bc_modeIs the BC homogeneous or inhomogeneous?
s_modeAre data data solution or correction?
bndryobject for handling coarse/fine and physical boundaries

Implemented in amrex::MLCellLinOpT< MF >.

◆ applyInhomogNeumannTerm()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::applyInhomogNeumannTerm ( int  ,
MF &   
) const
inlinevirtual

Extra terms introduced when we treat inhomogeneous Nuemann BC as homogeneous.

Reimplemented in amrex::MLCellABecLapT< MF >.

◆ applyMetricTerm()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::applyMetricTerm ( int  ,
int  ,
MF &   
) const
inlinevirtual

apply metric terms if there are any

Reimplemented in amrex::MLCellLinOpT< MF >.

◆ applyOverset()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::applyOverset ( int  ,
MF &   
) const
inlinevirtual

for overset solver only

Reimplemented in amrex::MLCellABecLapT< MF >.

◆ averageDownAndSync()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::averageDownAndSync ( Vector< MF > &  sol) const
pure virtual

Implemented in amrex::MLCellLinOpT< MF >.

◆ averageDownSolutionRHS()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::averageDownSolutionRHS ( int  camrlev,
MF &  crse_sol,
MF &  crse_rhs,
const MF &  fine_sol,
const MF &  fine_rhs 
)
inlinevirtual

Average-down data from fine AMR level to coarse AMR level.

Parameters
camrlevcoarse AMR level
crse_solsolutoin on coarse AMR level
crse_rhsRHS on coarse AMR level
fine_solsolution on fine AMR level
fine_rhsRHS on fine AMR level

Reimplemented in amrex::MLCellLinOpT< MF >.

◆ avgDownResAmr()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::avgDownResAmr ( int  clev,
MF &  cres,
MF const &  fres 
) const
inlinevirtual

Reimplemented in amrex::MLCellLinOpT< MF >.

◆ avgDownResMG()

template<typename MF >
void amrex::MLLinOpT< MF >::avgDownResMG ( int  clev,
MF &  cres,
MF const &  fres 
) const
virtual

◆ BottomCommunicator()

template<typename MF >
MPI_Comm amrex::MLLinOpT< MF >::BottomCommunicator ( ) const
inlineprotectednoexcept

◆ checkPoint()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::checkPoint ( std::string const &  ) const
inlineprivatevirtual

Reimplemented in amrex::MLNodeLaplacian.

◆ Communicator()

template<typename MF >
MPI_Comm amrex::MLLinOpT< MF >::Communicator ( ) const
inlineprotectednoexcept

◆ compactify() [1/2]

template<typename MF >
template<typename T >
Array4<T> amrex::MLLinOpT< MF >::compactify ( Array4< T > const &  a) const
inlineprotectednoexcept

◆ compactify() [2/2]

template<typename MF >
Box amrex::MLLinOpT< MF >::compactify ( Box const &  b) const
protectednoexcept

◆ compFlux()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::compFlux ( int  ,
const Array< MF *, AMREX_SPACEDIM > &  ,
MF &  ,
Location   
) const
inlinevirtual

Compute fluxes.

Parameters
amrlevAMR level
fluxesfluxes
solsolution
loclocation of the fluxes

Reimplemented in amrex::MLCellLinOpT< MF >.

◆ compGrad()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::compGrad ( int  ,
const Array< MF *, AMREX_SPACEDIM > &  ,
MF &  ,
Location   
) const
inlinevirtual

Compute gradients of the solution.

Parameters
amrlevAMR level
gradgrad(sol)
solsolution
loclocation of the gradients

Reimplemented in amrex::MLCellLinOpT< MF >.

◆ copyNSolveSolution()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::copyNSolveSolution ( MF &  ,
MF const &   
) const
inlinevirtual

◆ correctionResidual()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::correctionResidual ( int  amrlev,
int  mglev,
MF &  resid,
MF &  x,
const MF &  b,
BCMode  bc_mode,
const MF *  crse_bcdata = nullptr 
)
pure virtual

Compute residual for the residual-correction form, resid = b - L(x)

Parameters
amrlevAMR level
mglevMG level
residresidual
xunknown in the residual-correction form
bRHS in the residual-correction form
bc_modeIs the BC homogeneous or inhomogeneous?
crse_bc_dataoptional argument providing BC at coarse/fine boundary.

Implemented in amrex::MLCellLinOpT< MF >.

◆ define()

template<typename MF >
void amrex::MLLinOpT< MF >::define ( const Vector< Geometry > &  a_geom,
const Vector< BoxArray > &  a_grids,
const Vector< DistributionMapping > &  a_dmap,
const LPInfo a_info,
const Vector< FabFactory< FAB > const * > &  a_factory,
bool  eb_limit_coarsening = true 
)

◆ defineBC()

template<typename MF >
void amrex::MLLinOpT< MF >::defineBC
private

◆ defineGrids()

template<typename MF >
void amrex::MLLinOpT< MF >::defineGrids ( const Vector< Geometry > &  a_geom,
const Vector< BoxArray > &  a_grids,
const Vector< DistributionMapping > &  a_dmap,
const Vector< FabFactory< FAB > const * > &  a_factory 
)
private

◆ doAgglomeration()

template<typename MF >
bool amrex::MLLinOpT< MF >::doAgglomeration ( ) const
inlineprotectednoexcept

◆ doConsolidation()

template<typename MF >
bool amrex::MLLinOpT< MF >::doConsolidation ( ) const
inlineprotectednoexcept

◆ doSemicoarsening()

template<typename MF >
bool amrex::MLLinOpT< MF >::doSemicoarsening ( ) const
inlineprotectednoexcept

◆ Factory()

template<typename MF >
FabFactory<FAB> const* amrex::MLLinOpT< MF >::Factory ( int  amr_lev,
int  mglev = 0 
) const
inlineprotectednoexcept

◆ fixSolvabilityByOffset()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::fixSolvabilityByOffset ( int  ,
int  ,
MF &  ,
Vector< RT > const &   
) const
inlinevirtual

fix solvability by subtracting offset from RHS

◆ Geom()

template<typename MF >
const Geometry& amrex::MLLinOpT< MF >::Geom ( int  amr_lev,
int  mglev = 0 
) const
inlinenoexcept

◆ get_d0()

template<typename MF >
template<typename T >
T amrex::MLLinOpT< MF >::get_d0 ( T const &  dx,
T const &  dy,
T const &   
) const
inlineprotectednoexcept

◆ get_d1()

template<typename MF >
template<typename T >
T amrex::MLLinOpT< MF >::get_d1 ( T const &  ,
T const &  dy,
T const &  dz 
) const
inlineprotectednoexcept

◆ getDefaultBottomSolver()

template<typename MF >
virtual BottomSolver amrex::MLLinOpT< MF >::getDefaultBottomSolver ( ) const
inlinevirtual

Reimplemented in amrex::MLNodeLaplacian.

◆ getEnforceSingularSolvable()

template<typename MF >
bool amrex::MLLinOpT< MF >::getEnforceSingularSolvable ( ) const
inlinenoexcept

Get the flag for whether the solver should try to make singular problem solvable.

◆ getFluxes() [1/2]

template<typename MF >
virtual void amrex::MLLinOpT< MF >::getFluxes ( const Vector< Array< MF *, AMREX_SPACEDIM > > &  ,
const Vector< MF * > &  ,
Location   
) const
inlinevirtual

Reimplemented in amrex::MLCellABecLapT< MF >.

◆ getFluxes() [2/2]

template<typename MF >
virtual void amrex::MLLinOpT< MF >::getFluxes ( const Vector< MF * > &  ,
const Vector< MF * > &   
) const
inlinevirtual

Reimplemented in amrex::MLCellABecLapT< MF >.

◆ getMaxOrder()

template<typename MF >
int amrex::MLLinOpT< MF >::getMaxOrder ( ) const
inlinenoexcept

Get order of interpolation at coarse/fine boundary.

◆ getNComp()

template<typename MF >
virtual int amrex::MLLinOpT< MF >::getNComp ( ) const
inlinevirtual

◆ getNGrow()

template<typename MF >
virtual int amrex::MLLinOpT< MF >::getNGrow ( int  = 0,
int  = 0 
) const
inlinevirtual

◆ getNGrowVectRestriction()

template<typename MF >
virtual IntVect amrex::MLLinOpT< MF >::getNGrowVectRestriction ( ) const
inlineprotectedvirtual

Reimplemented in amrex::MLCurlCurl.

◆ getSolvabilityOffset()

template<typename MF >
virtual Vector<RT> amrex::MLLinOpT< MF >::getSolvabilityOffset ( int  ,
int  ,
MF const &   
) const
inlinevirtual

get offset for fixing solvability

Reimplemented in amrex::MLCellLinOpT< MF >.

◆ hasBC()

template<typename MF >
bool amrex::MLLinOpT< MF >::hasBC ( BCType  bct) const
protectednoexcept

◆ hasHiddenDimension()

template<typename MF >
bool amrex::MLLinOpT< MF >::hasHiddenDimension ( ) const
inlineprotectednoexcept

◆ hasInhomogNeumannBC()

template<typename MF >
bool amrex::MLLinOpT< MF >::hasInhomogNeumannBC
protectednoexcept

◆ hasRobinBC()

template<typename MF >
bool amrex::MLLinOpT< MF >::hasRobinBC
protectednoexcept

◆ HiBC()

template<typename MF >
GpuArray<BCType,AMREX_SPACEDIM> amrex::MLLinOpT< MF >::HiBC ( int  icomp = 0) const
inlineprotectednoexcept

◆ hiddenDirection()

template<typename MF >
int amrex::MLLinOpT< MF >::hiddenDirection ( ) const
inlineprotectednoexcept

◆ interpAssign()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::interpAssign ( int  amrlev,
int  fmglev,
MF &  fine,
MF &  crse 
) const
inlinevirtual

Overwrite fine MG level data with interpolated coarse data.

Parameters
amrlevAMR level
fmglevfine MG level
finefine MG level data
crsecoarse MG level data

Reimplemented in amrex::MLCellLinOpT< MF >.

◆ interpolation()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::interpolation ( int  amrlev,
int  fmglev,
MF &  fine,
const MF &  crse 
) const
pure virtual

Add interpolated coarse MG level data to fine MG level data.

Parameters
amrlevAMR level
fmglevfine MG level
crsecoarse data.
finefine data.

Implemented in amrex::MLCellLinOpT< MF >.

◆ interpolationAmr()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::interpolationAmr ( int  famrlev,
MF &  fine,
const MF &  crse,
IntVect const &  nghost 
) const
inlinevirtual

Interpolation between AMR levels.

Parameters
famrlevfine AMR level
finefine level data
crsecoarse level data
nghostnumber of ghost cells

Reimplemented in amrex::MLCellLinOpT< MF >.

◆ isBottomActive()

template<typename MF >
bool amrex::MLLinOpT< MF >::isBottomActive ( ) const
inlineprotectednoexcept

◆ isBottomSingular()

◆ isCellCentered()

template<typename MF >
bool amrex::MLLinOpT< MF >::isCellCentered ( ) const
inlineprotectednoexcept

◆ isMFIterSafe()

template<typename MF >
bool amrex::MLLinOpT< MF >::isMFIterSafe ( int  amrlev,
int  mglev1,
int  mglev2 
) const

◆ isSingular()

template<typename MF >
virtual bool amrex::MLLinOpT< MF >::isSingular ( int  amrlev) const
pure virtual

◆ LoBC()

template<typename MF >
GpuArray<BCType,AMREX_SPACEDIM> amrex::MLLinOpT< MF >::LoBC ( int  icomp = 0) const
inlineprotectednoexcept

◆ make() [1/2]

template<typename MF >
MF amrex::MLLinOpT< MF >::make ( int  amrlev,
int  mglev,
IntVect const &  ng 
) const
protectedvirtual

Reimplemented in amrex::MLCurlCurl.

◆ make() [2/2]

template<typename MF >
void amrex::MLLinOpT< MF >::make ( Vector< Vector< MF > > &  mf,
IntVect const &  ng 
) const
protectedvirtual

◆ makeAgglomeratedDMap()

template<typename MF >
void amrex::MLLinOpT< MF >::makeAgglomeratedDMap ( const Vector< BoxArray > &  ba,
Vector< DistributionMapping > &  dm 
)
staticprivate

◆ makeAlias()

template<typename MF >
MF amrex::MLLinOpT< MF >::makeAlias ( MF const &  mf) const
protectedvirtual

◆ makeCoarseAmr()

template<typename MF >
MF amrex::MLLinOpT< MF >::makeCoarseAmr ( int  famrlev,
IntVect const &  ng 
) const
protectedvirtual

Reimplemented in amrex::MLCurlCurl.

◆ makeCoarseMG()

template<typename MF >
MF amrex::MLLinOpT< MF >::makeCoarseMG ( int  amrlev,
int  mglev,
IntVect const &  ng 
) const
protectedvirtual

Reimplemented in amrex::MLCurlCurl.

◆ makeConsolidatedDMap()

template<typename MF >
void amrex::MLLinOpT< MF >::makeConsolidatedDMap ( const Vector< BoxArray > &  ba,
Vector< DistributionMapping > &  dm,
int  ratio,
int  strategy 
)
staticprivate

◆ makeFactory()

template<typename MF >
virtual std::unique_ptr<FabFactory<FAB> > amrex::MLLinOpT< MF >::makeFactory ( int  ,
int   
) const
inlineprotectedvirtual

Reimplemented in amrex::MLEBABecLap.

◆ makeNLinOp()

template<typename MF >
virtual std::unique_ptr<MLLinOpT<MF> > amrex::MLLinOpT< MF >::makeNLinOp ( int  ) const
inlinevirtual

◆ makeSubCommunicator()

template<typename MF >
MPI_Comm amrex::MLLinOpT< MF >::makeSubCommunicator ( const DistributionMapping dm)
private

◆ name()

template<typename MF >
virtual std::string amrex::MLLinOpT< MF >::name ( ) const
inlinevirtual

◆ NAMRLevels()

template<typename MF >
int amrex::MLLinOpT< MF >::NAMRLevels ( ) const
inlinenoexcept

Return the number of AMR levels.

◆ needsCoarseDataForBC()

template<typename MF >
bool amrex::MLLinOpT< MF >::needsCoarseDataForBC ( ) const
inlinenoexcept

Needs coarse data for bc?

If the lowest level grids does not cover the entire domain, coarse level data are needed for supplying Dirichlet bc at coarse/fine boundary, even when the domain bc is not Dirichlet.

◆ needsUpdate()

template<typename MF >
virtual bool amrex::MLLinOpT< MF >::needsUpdate ( ) const
inlinevirtual

◆ NMGLevels()

template<typename MF >
int amrex::MLLinOpT< MF >::NMGLevels ( int  amrlev) const
inlinenoexcept

Return the number of MG levels at given AMR level.

◆ normalize()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::normalize ( int  amrlev,
int  mglev,
MF &  mf 
) const
inlinevirtual

Divide mf by the diagonal component of the operator. Used by bicgstab.

Reimplemented in amrex::MLPoissonT< MF >, amrex::MLALaplacianT< MF >, and amrex::MLABecLaplacianT< MF >.

◆ normInf()

template<typename MF >
virtual RT amrex::MLLinOpT< MF >::normInf ( int  amrlev,
MF const &  mf,
bool  local 
) const
pure virtual

Implemented in amrex::MLCellLinOpT< MF >.

◆ operator=() [1/2]

template<typename MF >
MLLinOpT<MF>& amrex::MLLinOpT< MF >::operator= ( const MLLinOpT< MF > &  )
delete

◆ operator=() [2/2]

template<typename MF >
MLLinOpT<MF>& amrex::MLLinOpT< MF >::operator= ( MLLinOpT< MF > &&  )
delete

◆ postSolve()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::postSolve ( Vector< MF > &  ) const
inlinevirtual

◆ prepareForFluxes()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::prepareForFluxes ( int  ,
const MF *  = nullptr 
)
inlinevirtual

Reimplemented in amrex::MLCellLinOpT< MF >.

◆ prepareForGMRES()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::prepareForGMRES ( )
inlinevirtual

Reimplemented in amrex::MLNodeLinOp.

◆ prepareForSolve()

◆ reflux()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::reflux ( int  crse_amrlev,
MF &  res,
const MF &  crse_sol,
const MF &  crse_rhs,
MF &  fine_res,
MF &  fine_sol,
const MF &  fine_rhs 
) const
inlinevirtual

Reflux at AMR coarse/fine boundary.

Parameters
crse_amrlevcoarse AMR level
rescoarse level residual
crse_solcoarse level solution
crse_rhscoarse level RHS
fine_resfine level residual
fine_solfine level solution
fine_rhsfine level RHS

Reimplemented in amrex::MLCellLinOpT< MF >.

◆ resizeMultiGrid()

template<typename MF >
void amrex::MLLinOpT< MF >::resizeMultiGrid ( int  new_size)
protectedvirtual

◆ restriction()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::restriction ( int  amrlev,
int  cmglev,
MF &  crse,
MF &  fine 
) const
pure virtual

Restriction onto coarse MG level.

Parameters
amrlevAMR level
cmglevcoarse MG level
crsecoarse data. This is the output.
finefine data. This is the input. Some operators might need to fill ghost cells. This is why it's not a const reference.

Implemented in amrex::MLCellLinOpT< MF >.

◆ scaleRHS()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::scaleRHS ( int  ,
MF &   
) const
inlinevirtual

scale RHS to fix solvability

◆ setCoarseFineBC() [1/4]

template<typename MF >
template<typename AMF , std::enable_if_t<!std::is_same_v< MF, AMF >, int > >
void amrex::MLLinOpT< MF >::setCoarseFineBC ( const AMF *  crse,
int  crse_ratio,
LinOpBCType  bc_type = LinOpBCType::Dirichlet 
)
noexcept

◆ setCoarseFineBC() [2/4]

template<typename MF >
template<typename AMF , std::enable_if_t<!std::is_same_v< MF, AMF >, int > >
void amrex::MLLinOpT< MF >::setCoarseFineBC ( const AMF *  crse,
IntVect const &  crse_ratio,
LinOpBCType  bc_type = LinOpBCType::Dirichlet 
)
noexcept

◆ setCoarseFineBC() [3/4]

template<typename MF >
void amrex::MLLinOpT< MF >::setCoarseFineBC ( const MF *  crse,
int  crse_ratio,
LinOpBCType  bc_type = LinOpBCType::Dirichlet 
)
noexcept

Set coarse/fine boundary conditions. For cell-centered solves only.

If we want to do a linear solve where the boundary conditions on the coarsest AMR level of the solve come from a coarser level (e.g. the base AMR level of the solve is > 0 and does not cover the entire domain), we must explicitly provide the coarser data. Boundary conditions from a coarser level are Dirichlet by default. The MultiFab crse does not need to have ghost cells and is at a coarser resolution than the coarsest AMR level of the solve; it is used to supply (interpolated) boundary conditions for the solve. NOTE: If this is called, it must be called before setLevelBC. If crse is nullptr, then the bc values are assumed to be zero. The coarse/fine BC type can be changed to homogeneous Neumann by the bc_type argument. In that case, use nullptr for the crse argument.

Parameters
crsethe coarse AMR level data
crse_ratiothe coarsening ratio between fine and coarse AMR levels.
bc_typeoptional. It's Dirichlet by default, and can be Neumann.

◆ setCoarseFineBC() [4/4]

template<typename MF >
void amrex::MLLinOpT< MF >::setCoarseFineBC ( const MF *  crse,
IntVect const &  crse_ratio,
LinOpBCType  bc_type = LinOpBCType::Dirichlet 
)
noexcept

◆ setCoarseFineBCLocation()

template<typename MF >
void amrex::MLLinOpT< MF >::setCoarseFineBCLocation ( const RealVect cloc)
inlineprotectednoexcept

◆ setDirichletNodesToZero()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::setDirichletNodesToZero ( int  ,
int  ,
MF &   
) const
inlinevirtual

For GMRES to work, this might need to be implemented to mask out Dirichlet nodes or cells (e.g., EB covered cells or overset cells)

Reimplemented in amrex::MLCellABecLapT< MF >.

◆ setDomainBC() [1/2]

template<typename MF >
void amrex::MLLinOpT< MF >::setDomainBC ( const Array< BCType, AMREX_SPACEDIM > &  lobc,
const Array< BCType, AMREX_SPACEDIM > &  hibc 
)
noexcept

Boundary of the whole domain.

This functions must be called, and must be called before other bc functions. This version is for single-component solve or when all the components have the same BC types.

Parameters
lobclower boundaries
hibcupper boundaries

◆ setDomainBC() [2/2]

template<typename MF >
void amrex::MLLinOpT< MF >::setDomainBC ( const Vector< Array< BCType, AMREX_SPACEDIM > > &  lobc,
const Vector< Array< BCType, AMREX_SPACEDIM > > &  hibc 
)
noexcept

Boundary of the whole domain.

This functions must be called, and must be called before other bc functions. This version is for multi-component solve.

Parameters
lobclower boundaries
hibcupper boundaries

◆ setDomainBCLoc()

template<typename MF >
void amrex::MLLinOpT< MF >::setDomainBCLoc ( const Array< Real, AMREX_SPACEDIM > &  lo_bcloc,
const Array< Real, AMREX_SPACEDIM > &  hi_bcloc 
)
noexcept

Set location of domain boundaries.

By default, domain BC is on the domain face. If that's the case, this function doesn't need to be called. However, one could use this function to set non-zero domain BC locations. Note all values should be >= 0. If this function is called, it MUST be called before setLevelBC.

◆ setEnforceSingularSolvable()

template<typename MF >
void amrex::MLLinOpT< MF >::setEnforceSingularSolvable ( bool  o)
inlinenoexcept

Set the flag for whether the solver should try to make singular problem solvable, which is on by default.

◆ setLevelBC() [1/2]

template<typename MF >
template<typename AMF , std::enable_if_t<!std::is_same_v< MF, AMF >, int > >
void amrex::MLLinOpT< MF >::setLevelBC ( int  amrlev,
const AMF *  levelbcdata,
const AMF *  robinbc_a = nullptr,
const AMF *  robinbc_b = nullptr,
const AMF *  robinbc_f = nullptr 
)

◆ setLevelBC() [2/2]

template<typename MF >
virtual void amrex::MLLinOpT< MF >::setLevelBC ( int  ,
const MF *  ,
const MF *  = nullptr,
const MF *  = nullptr,
const MF *  = nullptr 
)
pure virtual

Set boundary conditions for given level. For cell-centered solves only.

This must be called for each level. Argument levelbcdata is used to supply Dirichlet or Neumann bc at the physical domain; if those data are homogeneous we can pass nullptr instead of levelbcdata. Regardless, this function must be called. If used, the MultiFab levelbcdata must have one ghost cell. Only the data outside the physical domain will be used. It is assumed that the data in those ghost cells outside the domain live exactly on the face of the physical domain. Argument amrlev is relative level such that the lowest to the solver is always 0. The optional arguments robinbc_[a|b|f] provide Robin boundary condition a*phi + b*dphi/dn = f. Note that d./dn is d./dx at the upper boundary and -d./dx at the lower boundary, for Robin BC. However, for inhomogeneous Neumann BC, the value in leveldata is assumed to be d./dx.

Implemented in amrex::MLCellLinOpT< MF >.

◆ setMaxOrder()

template<typename MF >
void amrex::MLLinOpT< MF >::setMaxOrder ( int  o)
inlinenoexcept

Set order of interpolation at coarse/fine boundary.

◆ setVerbose()

template<typename MF >
void amrex::MLLinOpT< MF >::setVerbose ( int  v)
inlinenoexcept

Set verbosity.

◆ smooth()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::smooth ( int  amrlev,
int  mglev,
MF &  sol,
const MF &  rhs,
bool  skip_fillboundary = false 
) const
pure virtual

Smooth.

Parameters
amrlevAMR level
mglevMG level
solunknowns
rhsRHS
skip_fillboundaryflag controlling whether ghost cell filling can be skipped.

Implemented in amrex::MLCellLinOpT< MF >.

◆ solutionResidual()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::solutionResidual ( int  amrlev,
MF &  resid,
MF &  x,
const MF &  b,
const MF *  crse_bcdata = nullptr 
)
pure virtual

Compute residual for solution.

Parameters
amrlevAMR level
residresidual
xsolution
bRHS
crse_bc_dataoptional argument providing BC at coarse/fine boundary.

Implemented in amrex::MLCellLinOpT< MF >.

◆ supportInhomogNeumannBC()

template<typename MF >
virtual bool amrex::MLLinOpT< MF >::supportInhomogNeumannBC ( ) const
inlineprotectedvirtualnoexcept

Reimplemented in amrex::MLCellABecLapT< MF >.

◆ supportNSolve()

template<typename MF >
virtual bool amrex::MLLinOpT< MF >::supportNSolve ( ) const
inlinevirtual

◆ supportRobinBC()

template<typename MF >
virtual bool amrex::MLLinOpT< MF >::supportRobinBC ( ) const
inlineprotectedvirtualnoexcept

◆ unapplyMetricTerm()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::unapplyMetricTerm ( int  ,
int  ,
MF &   
) const
inlinevirtual

unapply metric terms if there are any

Reimplemented in amrex::MLCellLinOpT< MF >.

◆ unimposeNeumannBC()

template<typename MF >
virtual void amrex::MLLinOpT< MF >::unimposeNeumannBC ( int  ,
MF &   
) const
inlinevirtual

This is needed for our nodal projection solver.

◆ update()

◆ xdoty()

template<typename MF >
virtual RT amrex::MLLinOpT< MF >::xdoty ( int  amrlev,
int  mglev,
const MF &  x,
const MF &  y,
bool  local 
) const
pure virtual

x dot y, used by the bottom solver

Implemented in amrex::MLCellLinOpT< MF >.

Friends And Related Function Documentation

◆ GMRESMLMGT

template<typename MF >
template<typename T >
friend class GMRESMLMGT
friend

◆ MLABecLaplacianT

template<typename MF >
template<typename T >
friend class MLABecLaplacianT
friend

◆ MLCGSolverT

template<typename MF >
template<typename T >
friend class MLCGSolverT
friend

◆ MLMGT

template<typename MF >
template<typename T >
friend class MLMGT
friend

◆ MLPoissonT

template<typename MF >
template<typename T >
friend class MLPoissonT
friend

Member Data Documentation

◆ enforceSingularSolvable

template<typename MF >
bool amrex::MLLinOpT< MF >::enforceSingularSolvable = true
protected

◆ info

template<typename MF >
LPInfo amrex::MLLinOpT< MF >::info
protected

◆ levelbc_raii

template<typename MF >
Vector<std::unique_ptr<MF> > amrex::MLLinOpT< MF >::levelbc_raii
private

◆ m_amr_ref_ratio

template<typename MF >
Vector<int> amrex::MLLinOpT< MF >::m_amr_ref_ratio
protected

◆ m_bottom_comm

template<typename MF >
MPI_Comm amrex::MLLinOpT< MF >::m_bottom_comm = MPI_COMM_NULL
protected

◆ m_coarse_bc_loc

template<typename MF >
RealVect amrex::MLLinOpT< MF >::m_coarse_bc_loc
protected

◆ m_coarse_data_crse_ratio

template<typename MF >
IntVect amrex::MLLinOpT< MF >::m_coarse_data_crse_ratio = IntVect(-1)
protected

◆ m_coarse_data_for_bc

template<typename MF >
const MF* amrex::MLLinOpT< MF >::m_coarse_data_for_bc = nullptr
protected

◆ m_coarse_data_for_bc_raii

template<typename MF >
MF amrex::MLLinOpT< MF >::m_coarse_data_for_bc_raii
protected

◆ m_coarse_fine_bc_type

template<typename MF >
LinOpBCType amrex::MLLinOpT< MF >::m_coarse_fine_bc_type = LinOpBCType::Dirichlet
protected

◆ m_default_comm

template<typename MF >
MPI_Comm amrex::MLLinOpT< MF >::m_default_comm = MPI_COMM_NULL
protected

◆ m_dmap

template<typename MF >
Vector<Vector<DistributionMapping> > amrex::MLLinOpT< MF >::m_dmap
protected

◆ m_do_agglomeration

template<typename MF >
bool amrex::MLLinOpT< MF >::m_do_agglomeration = false
protected

◆ m_do_consolidation

template<typename MF >
bool amrex::MLLinOpT< MF >::m_do_consolidation = false
protected

◆ m_do_semicoarsening

template<typename MF >
bool amrex::MLLinOpT< MF >::m_do_semicoarsening = false
protected

◆ m_domain_bloc_hi

template<typename MF >
Array<Real, AMREX_SPACEDIM> amrex::MLLinOpT< MF >::m_domain_bloc_hi {{AMREX_D_DECL(0._rt,0._rt,0._rt)}}
protected

◆ m_domain_bloc_lo

template<typename MF >
Array<Real, AMREX_SPACEDIM> amrex::MLLinOpT< MF >::m_domain_bloc_lo {{AMREX_D_DECL(0._rt,0._rt,0._rt)}}
protected

◆ m_domain_covered

template<typename MF >
Vector<int> amrex::MLLinOpT< MF >::m_domain_covered
protected

◆ m_factory

template<typename MF >
Vector<Vector<std::unique_ptr<FabFactory<FAB> > > > amrex::MLLinOpT< MF >::m_factory
protected

◆ m_geom

template<typename MF >
Vector<Vector<Geometry> > amrex::MLLinOpT< MF >::m_geom
protected

first Vector is for amr level and second is mg level

◆ m_grids

template<typename MF >
Vector<Vector<BoxArray> > amrex::MLLinOpT< MF >::m_grids
protected

◆ m_hibc

template<typename MF >
Vector<Array<BCType, AMREX_SPACEDIM> > amrex::MLLinOpT< MF >::m_hibc

◆ m_hibc_orig

template<typename MF >
Vector<Array<BCType, AMREX_SPACEDIM> > amrex::MLLinOpT< MF >::m_hibc_orig

◆ m_ixtype

template<typename MF >
IntVect amrex::MLLinOpT< MF >::m_ixtype
protected

◆ m_lobc

template<typename MF >
Vector<Array<BCType, AMREX_SPACEDIM> > amrex::MLLinOpT< MF >::m_lobc

◆ m_lobc_orig

template<typename MF >
Vector<Array<BCType, AMREX_SPACEDIM> > amrex::MLLinOpT< MF >::m_lobc_orig

◆ m_needs_coarse_data_for_bc

template<typename MF >
bool amrex::MLLinOpT< MF >::m_needs_coarse_data_for_bc = false
protected

◆ m_num_amr_levels

template<typename MF >
int amrex::MLLinOpT< MF >::m_num_amr_levels = 0
protected

◆ m_num_mg_levels

template<typename MF >
Vector<int> amrex::MLLinOpT< MF >::m_num_mg_levels
protected

◆ m_parent

template<typename MF >
const MLLinOpT<MF>* amrex::MLLinOpT< MF >::m_parent = nullptr
protected

◆ m_raii_comm

template<typename MF >
std::unique_ptr<CommContainer> amrex::MLLinOpT< MF >::m_raii_comm
protected

◆ maxorder

template<typename MF >
int amrex::MLLinOpT< MF >::maxorder = 3
protected

◆ mg_box_min_width

template<typename MF >
constexpr int amrex::MLLinOpT< MF >::mg_box_min_width = 2
staticconstexprprotected

◆ mg_coarsen_ratio

template<typename MF >
constexpr int amrex::MLLinOpT< MF >::mg_coarsen_ratio = 2
staticconstexprprotected

◆ mg_coarsen_ratio_vec

template<typename MF >
Vector<IntVect> amrex::MLLinOpT< MF >::mg_coarsen_ratio_vec
protected

◆ mg_domain_min_width

template<typename MF >
int amrex::MLLinOpT< MF >::mg_domain_min_width = 2
protected

◆ robin_a_raii

template<typename MF >
Vector<std::unique_ptr<MF> > amrex::MLLinOpT< MF >::robin_a_raii
private

◆ robin_b_raii

template<typename MF >
Vector<std::unique_ptr<MF> > amrex::MLLinOpT< MF >::robin_b_raii
private

◆ robin_f_raii

template<typename MF >
Vector<std::unique_ptr<MF> > amrex::MLLinOpT< MF >::robin_f_raii
private

◆ verbose

template<typename MF >
int amrex::MLLinOpT< MF >::verbose = 0
protected

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