Block-Structured AMR Software Framework
amrex::MLLinOp Class Referenceabstract

#include <AMReX_MLLinOp.H>

Inheritance diagram for amrex::MLLinOp:
amrex::MLCellLinOp amrex::MLNodeLinOp amrex_temp::MLLinOpTemp amrex::MLCellABecLap amrex::MLEBNodeFDLaplacian amrex::MLNodeLaplacian amrex::MLNodeTensorLaplacian amrex::MLABecLaplacian amrex::MLALaplacian amrex::MLEBABecLap amrex::MLPoisson amrex::MLTensorOp amrex::MLEBTensorOp

Classes

struct  CommContainer
 

Public Types

enum class  BCMode { Homogeneous , Inhomogeneous }
 
enum class  StateMode { Solution , Correction }
 
enum class  Location { FaceCenter , FaceCentroid , CellCenter , CellCentroid }
 
using BCType = LinOpBCType
 

Public Member Functions

 MLLinOp ()
 
virtual ~MLLinOp ()
 
 MLLinOp (const MLLinOp &)=delete
 
 MLLinOp (MLLinOp &&)=delete
 
MLLinOpoperator= (const MLLinOp &)=delete
 
MLLinOpoperator= (MLLinOp &&)=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< FArrayBox > 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. This functions must be called, and must be called before other bc functions. More...
 
void setDomainBC (const Vector< Array< BCType, AMREX_SPACEDIM > > &lobc, const Vector< Array< BCType, AMREX_SPACEDIM > > &hibc) noexcept
 
void setDomainBCLoc (const Array< Real, AMREX_SPACEDIM > &lo_bcloc, const Array< Real, AMREX_SPACEDIM > &hi_bcloc) noexcept
 
bool needsCoarseDataForBC () const noexcept
 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. More...
 
void setCoarseFineBC (const MultiFab *crse, int crse_ratio) noexcept
 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 always Dirichlet. 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. More...
 
virtual void setLevelBC (int, const MultiFab *, const MultiFab *=nullptr, const MultiFab *=nullptr, const MultiFab *=nullptr)
 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. More...
 
void setVerbose (int v) noexcept
 
void setMaxOrder (int o) noexcept
 
int getMaxOrder () const noexcept
 
void setEnforceSingularSolvable (bool o) noexcept
 
bool getEnforceSingularSolvable () const noexcept
 
virtual BottomSolver getDefaultBottomSolver () const
 
virtual int getNComp () const
 
virtual int getNGrow (int=0, int=0) const
 
virtual bool needsUpdate () const
 
virtual void update ()
 
virtual void restriction (int, int, MultiFab &, MultiFab &) const
 
virtual void interpolation (int, int, MultiFab &, const MultiFab &) const
 
virtual void interpAssign (int, int, MultiFab &, MultiFab &) const
 
virtual void averageDownSolutionRHS (int, MultiFab &, MultiFab &, const MultiFab &, const MultiFab &)
 
virtual void apply (int, int, MultiFab &, MultiFab &, BCMode, StateMode, const MLMGBndry *=nullptr) const
 
virtual void smooth (int, int, MultiFab &, const MultiFab &, bool=false) const
 
virtual void normalize (int, int, MultiFab &) const
 
virtual void solutionResidual (int, MultiFab &, MultiFab &, const MultiFab &, const MultiFab *=nullptr)
 
virtual void correctionResidual (int, int, MultiFab &, MultiFab &, const MultiFab &, BCMode, const MultiFab *=nullptr)
 
virtual void reflux (int, MultiFab &, const MultiFab &, const MultiFab &, MultiFab &, MultiFab &, const MultiFab &) const
 
virtual void compFlux (int, const Array< MultiFab *, AMREX_SPACEDIM > &, MultiFab &, Location) const
 
virtual void compGrad (int, const Array< MultiFab *, AMREX_SPACEDIM > &, MultiFab &, Location) const
 
virtual void applyMetricTerm (int, int, Any &) const
 
virtual void unapplyMetricTerm (int, int, MultiFab &) const
 
virtual void unimposeNeumannBC (int, Any &) const
 
virtual void applyInhomogNeumannTerm (int, Any &) const
 
virtual void applyOverset (int, Any &) const
 
virtual void scaleRHS (int, Any &) const
 
virtual Vector< Real > getSolvabilityOffset (int, int, Any const &) const
 
virtual void fixSolvabilityByOffset (int, int, Any &, Vector< Real > const &) const
 
virtual void prepareForSolve ()=0
 
virtual bool isSingular (int) const
 
virtual bool isBottomSingular () const
 
virtual Real xdoty (int, int, const MultiFab &, const MultiFab &, bool) const
 
virtual std::unique_ptr< MLLinOpmakeNLinOp (int) const
 
virtual void getFluxes (const Vector< Array< MultiFab *, AMREX_SPACEDIM > > &, const Vector< MultiFab * > &, Location) const
 
virtual void getFluxes (const Vector< MultiFab * > &, const Vector< MultiFab * > &) const
 
virtual bool supportNSolve () const
 
virtual void copyNSolveSolution (MultiFab &, MultiFab const &) const
 
virtual Any AnyMake (int amrlev, int mglev, IntVect const &ng) const
 
virtual Any AnyMakeCoarseMG (int amrlev, int mglev, IntVect const &ng) const
 
virtual Any AnyMakeCoarseAmr (int famrlev, IntVect const &ng) const
 
virtual Any AnyMakeAlias (Any const &a) const
 
virtual IntVect AnyGrowVect (Any const &a) const
 
virtual void AnyCopy (Any &dst, Any const &src, IntVect const &ng) const
 
virtual void AnyAdd (Any &dst, Any const &src, IntVect const &ng) const
 
virtual void AnySetToZero (Any &a) const
 
virtual void AnySetBndryToZero (Any &a) const
 
virtual void AnyParallelCopy (Any &dst, Any const &src, IntVect const &src_nghost, IntVect const &dst_nghost, Periodicity const &period=Periodicity::NonPeriodic()) const
 
virtual Real AnyNormInf (Any &a) const
 
virtual Real AnyNormInfMask (int amrlev, Any const &a, bool local) const =0
 
virtual void AnySolutionResidual (int amrlev, Any &resid, Any &x, Any const &b, Any const *crse_bcdata=nullptr)
 
virtual void AnyCorrectionResidual (int amrlev, int mglev, Any &resid, Any &x, const Any &b, BCMode bc_mode, const Any *crse_bcdata=nullptr)
 
virtual void AnyReflux (int crse_amrlev, Any &res, const Any &crse_sol, const Any &crse_rhs, Any &fine_res, Any &fine_sol, const Any &fine_rhs)
 
virtual void AnyAvgDownResAmr (int clev, Any &cres, Any const &fres) const =0
 
virtual void AnyAvgDownResMG (int clev, Any &cres, Any const &fres) const
 
virtual void AnySmooth (int amrlev, int mglev, Any &sol, const Any &rhs, bool skip_fillboundary=false) const
 
virtual void AnyRestriction (int amrlev, int cmglev, Any &crse, Any &fine) const
 
virtual void AnyInterpolationMG (int amrlev, int fmglev, Any &fine, const Any &crse) const
 
virtual void AnyInterpAssignMG (int amrlev, int fmglev, Any &fine, Any &crse) const
 
virtual void AnyInterpolationAmr (int famrlev, Any &fine, const Any &crse, IntVect const &) const =0
 
virtual void AnyAverageDownSolutionRHS (int camrlev, Any &crse_sol, Any &crse_rhs, const Any &fine_sol, const Any &fine_rhs)
 
virtual void AnyAverageDownAndSync (Vector< Any > &sol) const =0
 
virtual void postSolve (Vector< Any > &sol) const
 
Real MFNormInf (MultiFab const &mf, iMultiFab const *fine_mask, bool local) const
 
bool isMFIterSafe (int amrlev, int mglev1, int mglev2) const
 

Static Public Member Functions

static void Initialize ()
 
static void Finalize ()
 

Protected Member Functions

int NAMRLevels () const noexcept
 functions More...
 
int NMGLevels (int amrlev) const noexcept
 
const Vector< int > & AMRRefRatio () const noexcept
 
int AMRRefRatio (int amr_lev) const noexcept
 
const GeometryGeom (int amr_lev, int mglev=0) const noexcept
 
FabFactory< FArrayBox > 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 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
 
void make (Vector< Vector< Any > > &mf, IntVect const &ng) const
 
virtual std::unique_ptr< FabFactory< FArrayBox > > 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

LPInfo info
 
int verbose = 0
 
int maxorder = 3
 
bool enforceSingularSolvable = true
 
int m_num_amr_levels
 
Vector< intm_amr_ref_ratio
 
Vector< intm_num_mg_levels
 
const MLLinOpm_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< FArrayBox > > > > 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
 
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
 
Array< Real, AMREX_SPACEDIM > m_domain_bloc_lo {{AMREX_D_DECL(0.,0.,0.)}}
 
Array< Real, AMREX_SPACEDIM > m_domain_bloc_hi {{AMREX_D_DECL(0.,0.,0.)}}
 
bool m_needs_coarse_data_for_bc
 
int m_coarse_data_crse_ratio = -1
 
RealVect m_coarse_bc_loc
 
const MultiFabm_coarse_data_for_bc = nullptr
 

Static Protected Attributes

static constexpr int mg_coarsen_ratio = 2
 
static constexpr int mg_box_min_width = 2
 
static constexpr int mg_domain_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< FArrayBox > const * > &a_factory)
 
void defineAuxData ()
 
void defineBC ()
 
MPI_Comm makeSubCommunicator (const DistributionMapping &dm)
 
void remapNeighborhoods (Vector< DistributionMapping > &dms)
 
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)
 

Friends

class MLMG
 
class MLCGSolver
 
class MLPoisson
 
class MLABecLaplacian
 

Member Typedef Documentation

◆ BCType

using amrex::MLLinOp::BCType = LinOpBCType

Member Enumeration Documentation

◆ BCMode

Enumerator
Homogeneous 
Inhomogeneous 

◆ Location

Enumerator
FaceCenter 
FaceCentroid 
CellCenter 
CellCentroid 

◆ StateMode

Enumerator
Solution 
Correction 

Constructor & Destructor Documentation

◆ MLLinOp() [1/3]

amrex::MLLinOp::MLLinOp ( )

◆ ~MLLinOp()

amrex::MLLinOp::~MLLinOp ( )
virtual

◆ MLLinOp() [2/3]

amrex::MLLinOp::MLLinOp ( const MLLinOp )
delete

◆ MLLinOp() [3/3]

amrex::MLLinOp::MLLinOp ( MLLinOp &&  )
delete

Member Function Documentation

◆ AMRRefRatio() [1/2]

const Vector<int>& amrex::MLLinOp::AMRRefRatio ( ) const
inlineprotectednoexcept

◆ AMRRefRatio() [2/2]

int amrex::MLLinOp::AMRRefRatio ( int  amr_lev) const
inlineprotectednoexcept

◆ AnyAdd()

void amrex::MLLinOp::AnyAdd ( Any dst,
Any const &  src,
IntVect const &  ng 
) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyAverageDownAndSync()

virtual void amrex::MLLinOp::AnyAverageDownAndSync ( Vector< Any > &  sol) const
pure virtual

◆ AnyAverageDownSolutionRHS()

void amrex::MLLinOp::AnyAverageDownSolutionRHS ( int  camrlev,
Any crse_sol,
Any crse_rhs,
const Any fine_sol,
const Any fine_rhs 
)
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyAvgDownResAmr()

virtual void amrex::MLLinOp::AnyAvgDownResAmr ( int  clev,
Any cres,
Any const &  fres 
) const
pure virtual

◆ AnyAvgDownResMG()

void amrex::MLLinOp::AnyAvgDownResMG ( int  clev,
Any cres,
Any const &  fres 
) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyCopy()

void amrex::MLLinOp::AnyCopy ( Any dst,
Any const &  src,
IntVect const &  ng 
) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyCorrectionResidual()

void amrex::MLLinOp::AnyCorrectionResidual ( int  amrlev,
int  mglev,
Any resid,
Any x,
const Any b,
BCMode  bc_mode,
const Any crse_bcdata = nullptr 
)
virtual

◆ AnyGrowVect()

IntVect amrex::MLLinOp::AnyGrowVect ( Any const &  a) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyInterpAssignMG()

void amrex::MLLinOp::AnyInterpAssignMG ( int  amrlev,
int  fmglev,
Any fine,
Any crse 
) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyInterpolationAmr()

virtual void amrex::MLLinOp::AnyInterpolationAmr ( int  famrlev,
Any fine,
const Any crse,
IntVect const &   
) const
pure virtual

◆ AnyInterpolationMG()

void amrex::MLLinOp::AnyInterpolationMG ( int  amrlev,
int  fmglev,
Any fine,
const Any crse 
) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyMake()

Any amrex::MLLinOp::AnyMake ( int  amrlev,
int  mglev,
IntVect const &  ng 
) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyMakeAlias()

Any amrex::MLLinOp::AnyMakeAlias ( Any const &  a) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyMakeCoarseAmr()

Any amrex::MLLinOp::AnyMakeCoarseAmr ( int  famrlev,
IntVect const &  ng 
) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyMakeCoarseMG()

Any amrex::MLLinOp::AnyMakeCoarseMG ( int  amrlev,
int  mglev,
IntVect const &  ng 
) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyNormInf()

Real amrex::MLLinOp::AnyNormInf ( Any a) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyNormInfMask()

virtual Real amrex::MLLinOp::AnyNormInfMask ( int  amrlev,
Any const &  a,
bool  local 
) const
pure virtual

◆ AnyParallelCopy()

void amrex::MLLinOp::AnyParallelCopy ( Any dst,
Any const &  src,
IntVect const &  src_nghost,
IntVect const &  dst_nghost,
Periodicity const &  period = Periodicity::NonPeriodic() 
) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyReflux()

void amrex::MLLinOp::AnyReflux ( int  crse_amrlev,
Any res,
const Any crse_sol,
const Any crse_rhs,
Any fine_res,
Any fine_sol,
const Any fine_rhs 
)
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnyRestriction()

void amrex::MLLinOp::AnyRestriction ( int  amrlev,
int  cmglev,
Any crse,
Any fine 
) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnySetBndryToZero()

void amrex::MLLinOp::AnySetBndryToZero ( Any a) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnySetToZero()

void amrex::MLLinOp::AnySetToZero ( Any a) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnySmooth()

void amrex::MLLinOp::AnySmooth ( int  amrlev,
int  mglev,
Any sol,
const Any rhs,
bool  skip_fillboundary = false 
) const
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ AnySolutionResidual()

void amrex::MLLinOp::AnySolutionResidual ( int  amrlev,
Any resid,
Any x,
Any const &  b,
Any const *  crse_bcdata = nullptr 
)
virtual

Reimplemented in amrex_temp::MLLinOpTemp.

◆ apply()

virtual void amrex::MLLinOp::apply ( int  ,
int  ,
MultiFab ,
MultiFab ,
BCMode  ,
StateMode  ,
const MLMGBndry = nullptr 
) const
inlinevirtual

◆ applyInhomogNeumannTerm()

virtual void amrex::MLLinOp::applyInhomogNeumannTerm ( int  ,
Any  
) const
inlinevirtual

Reimplemented in amrex::MLCellABecLap.

◆ applyMetricTerm()

virtual void amrex::MLLinOp::applyMetricTerm ( int  ,
int  ,
Any  
) const
inlinevirtual

Reimplemented in amrex::MLNodeLinOp, and amrex::MLCellLinOp.

◆ applyOverset()

virtual void amrex::MLLinOp::applyOverset ( int  ,
Any  
) const
inlinevirtual

Reimplemented in amrex::MLCellABecLap.

◆ averageDownSolutionRHS()

virtual void amrex::MLLinOp::averageDownSolutionRHS ( int  ,
MultiFab ,
MultiFab ,
const MultiFab ,
const MultiFab  
)
inlinevirtual

◆ BottomCommunicator()

MPI_Comm amrex::MLLinOp::BottomCommunicator ( ) const
inlineprotectednoexcept

◆ checkPoint()

virtual void amrex::MLLinOp::checkPoint ( std::string const &  ) const
inlineprivatevirtual

Reimplemented in amrex::MLNodeLaplacian.

◆ Communicator()

MPI_Comm amrex::MLLinOp::Communicator ( ) const
inlineprotectednoexcept

◆ compactify() [1/2]

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

◆ compactify() [2/2]

Box amrex::MLLinOp::compactify ( Box const &  b) const
protectednoexcept

◆ compFlux()

virtual void amrex::MLLinOp::compFlux ( int  ,
const Array< MultiFab *, AMREX_SPACEDIM > &  ,
MultiFab ,
Location   
) const
inlinevirtual

◆ compGrad()

virtual void amrex::MLLinOp::compGrad ( int  ,
const Array< MultiFab *, AMREX_SPACEDIM > &  ,
MultiFab ,
Location   
) const
inlinevirtual

◆ copyNSolveSolution()

virtual void amrex::MLLinOp::copyNSolveSolution ( MultiFab ,
MultiFab const &   
) const
inlinevirtual

Reimplemented in amrex::MLPoisson, and amrex::MLABecLaplacian.

◆ correctionResidual()

virtual void amrex::MLLinOp::correctionResidual ( int  ,
int  ,
MultiFab ,
MultiFab ,
const MultiFab ,
BCMode  ,
const MultiFab = nullptr 
)
inlinevirtual

Reimplemented in amrex::MLNodeLinOp, and amrex::MLCellLinOp.

◆ define()

void amrex::MLLinOp::define ( const Vector< Geometry > &  a_geom,
const Vector< BoxArray > &  a_grids,
const Vector< DistributionMapping > &  a_dmap,
const LPInfo a_info,
const Vector< FabFactory< FArrayBox > const * > &  a_factory,
bool  eb_limit_coarsening = true 
)

◆ defineAuxData()

void amrex::MLLinOp::defineAuxData ( )
private

◆ defineBC()

void amrex::MLLinOp::defineBC ( )
private

◆ defineGrids()

void amrex::MLLinOp::defineGrids ( const Vector< Geometry > &  a_geom,
const Vector< BoxArray > &  a_grids,
const Vector< DistributionMapping > &  a_dmap,
const Vector< FabFactory< FArrayBox > const * > &  a_factory 
)
private

◆ doAgglomeration()

bool amrex::MLLinOp::doAgglomeration ( ) const
inlineprotectednoexcept

◆ doConsolidation()

bool amrex::MLLinOp::doConsolidation ( ) const
inlineprotectednoexcept

◆ doSemicoarsening()

bool amrex::MLLinOp::doSemicoarsening ( ) const
inlineprotectednoexcept

◆ Factory()

FabFactory<FArrayBox> const* amrex::MLLinOp::Factory ( int  amr_lev,
int  mglev = 0 
) const
inlineprotectednoexcept

◆ Finalize()

void amrex::MLLinOp::Finalize ( )
static

◆ fixSolvabilityByOffset()

virtual void amrex::MLLinOp::fixSolvabilityByOffset ( int  ,
int  ,
Any ,
Vector< Real > const &   
) const
inlinevirtual

◆ Geom()

const Geometry& amrex::MLLinOp::Geom ( int  amr_lev,
int  mglev = 0 
) const
inlineprotectednoexcept

◆ get_d0()

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

◆ get_d1()

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

◆ getDefaultBottomSolver()

virtual BottomSolver amrex::MLLinOp::getDefaultBottomSolver ( ) const
inlinevirtual

◆ getEnforceSingularSolvable()

bool amrex::MLLinOp::getEnforceSingularSolvable ( ) const
inlinenoexcept

◆ getFluxes() [1/2]

virtual void amrex::MLLinOp::getFluxes ( const Vector< Array< MultiFab *, AMREX_SPACEDIM > > &  ,
const Vector< MultiFab * > &  ,
Location   
) const
inlinevirtual

◆ getFluxes() [2/2]

virtual void amrex::MLLinOp::getFluxes ( const Vector< MultiFab * > &  ,
const Vector< MultiFab * > &   
) const
inlinevirtual

◆ getMaxOrder()

int amrex::MLLinOp::getMaxOrder ( ) const
inlinenoexcept

◆ getNComp()

virtual int amrex::MLLinOp::getNComp ( ) const
inlinevirtual

◆ getNGrow()

virtual int amrex::MLLinOp::getNGrow ( int  = 0,
int  = 0 
) const
inlinevirtual

◆ getSolvabilityOffset()

virtual Vector<Real> amrex::MLLinOp::getSolvabilityOffset ( int  ,
int  ,
Any const &   
) const
inlinevirtual

◆ hasHiddenDimension()

bool amrex::MLLinOp::hasHiddenDimension ( ) const
inlineprotectednoexcept

◆ hasInhomogNeumannBC()

bool amrex::MLLinOp::hasInhomogNeumannBC ( ) const
protectednoexcept

◆ hasRobinBC()

bool amrex::MLLinOp::hasRobinBC ( ) const
protectednoexcept

◆ HiBC()

GpuArray<BCType,AMREX_SPACEDIM> amrex::MLLinOp::HiBC ( int  icomp = 0) const
inlineprotectednoexcept

◆ hiddenDirection()

int amrex::MLLinOp::hiddenDirection ( ) const
inlineprotectednoexcept

◆ Initialize()

void amrex::MLLinOp::Initialize ( )
static

◆ interpAssign()

virtual void amrex::MLLinOp::interpAssign ( int  ,
int  ,
MultiFab ,
MultiFab  
) const
inlinevirtual

Reimplemented in amrex::MLNodeLinOp, and amrex::MLCellLinOp.

◆ interpolation()

virtual void amrex::MLLinOp::interpolation ( int  ,
int  ,
MultiFab ,
const MultiFab  
) const
inlinevirtual

◆ isBottomActive()

bool amrex::MLLinOp::isBottomActive ( ) const
inlineprotectednoexcept

◆ isBottomSingular()

virtual bool amrex::MLLinOp::isBottomSingular ( ) const
inlinevirtual

◆ isCellCentered()

bool amrex::MLLinOp::isCellCentered ( ) const
inlineprotectednoexcept

◆ isMFIterSafe()

bool amrex::MLLinOp::isMFIterSafe ( int  amrlev,
int  mglev1,
int  mglev2 
) const

◆ isSingular()

◆ LoBC()

GpuArray<BCType,AMREX_SPACEDIM> amrex::MLLinOp::LoBC ( int  icomp = 0) const
inlineprotectednoexcept

◆ make()

void amrex::MLLinOp::make ( Vector< Vector< Any > > &  mf,
IntVect const &  ng 
) const
protected

◆ makeAgglomeratedDMap()

void amrex::MLLinOp::makeAgglomeratedDMap ( const Vector< BoxArray > &  ba,
Vector< DistributionMapping > &  dm 
)
staticprivate

◆ makeConsolidatedDMap()

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

◆ makeFactory()

virtual std::unique_ptr<FabFactory<FArrayBox> > amrex::MLLinOp::makeFactory ( int  ,
int   
) const
inlineprotectedvirtual

Reimplemented in amrex::MLEBABecLap.

◆ makeNLinOp()

virtual std::unique_ptr<MLLinOp> amrex::MLLinOp::makeNLinOp ( int  ) const
inlinevirtual

◆ makeSubCommunicator()

MPI_Comm amrex::MLLinOp::makeSubCommunicator ( const DistributionMapping dm)
private

◆ MFNormInf()

Real amrex::MLLinOp::MFNormInf ( MultiFab const &  mf,
iMultiFab const *  fine_mask,
bool  local 
) const

◆ name()

virtual std::string amrex::MLLinOp::name ( ) const
inlinevirtual

◆ NAMRLevels()

int amrex::MLLinOp::NAMRLevels ( ) const
inlineprotectednoexcept

functions

◆ needsCoarseDataForBC()

bool amrex::MLLinOp::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()

virtual bool amrex::MLLinOp::needsUpdate ( ) const
inlinevirtual

◆ NMGLevels()

int amrex::MLLinOp::NMGLevels ( int  amrlev) const
inlineprotectednoexcept

◆ normalize()

virtual void amrex::MLLinOp::normalize ( int  ,
int  ,
MultiFab  
) const
inlinevirtual

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ postSolve()

void amrex::MLLinOp::postSolve ( Vector< Any > &  sol) const
virtual

Reimplemented in amrex::MLEBNodeFDLaplacian.

◆ prepareForSolve()

◆ reflux()

virtual void amrex::MLLinOp::reflux ( int  ,
MultiFab ,
const MultiFab ,
const MultiFab ,
MultiFab ,
MultiFab ,
const MultiFab  
) const
inlinevirtual

◆ remapNeighborhoods()

void amrex::MLLinOp::remapNeighborhoods ( Vector< DistributionMapping > &  dms)
private

◆ resizeMultiGrid()

void amrex::MLLinOp::resizeMultiGrid ( int  new_size)
protectedvirtual

◆ restriction()

virtual void amrex::MLLinOp::restriction ( int  ,
int  ,
MultiFab ,
MultiFab  
) const
inlinevirtual

◆ scaleRHS()

virtual void amrex::MLLinOp::scaleRHS ( int  ,
Any  
) const
inlinevirtual

◆ setCoarseFineBC()

void amrex::MLLinOp::setCoarseFineBC ( const MultiFab crse,
int  crse_ratio 
)
noexcept

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 always Dirichlet. 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.

Parameters
crse
crse_ratio

◆ setCoarseFineBCLocation()

void amrex::MLLinOp::setCoarseFineBCLocation ( const RealVect cloc)
inlineprotectednoexcept

◆ setDomainBC() [1/2]

void amrex::MLLinOp::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.

Parameters
lobc
hibc

◆ setDomainBC() [2/2]

void amrex::MLLinOp::setDomainBC ( const Vector< Array< BCType, AMREX_SPACEDIM > > &  lobc,
const Vector< Array< BCType, AMREX_SPACEDIM > > &  hibc 
)
noexcept

◆ setDomainBCLoc()

void amrex::MLLinOp::setDomainBCLoc ( const Array< Real, AMREX_SPACEDIM > &  lo_bcloc,
const Array< Real, AMREX_SPACEDIM > &  hi_bcloc 
)
noexcept

◆ setEnforceSingularSolvable()

void amrex::MLLinOp::setEnforceSingularSolvable ( bool  o)
inlinenoexcept

◆ setLevelBC()

virtual void amrex::MLLinOp::setLevelBC ( int  ,
const MultiFab ,
const MultiFab = nullptr,
const MultiFab = nullptr,
const MultiFab = nullptr 
)
inlinevirtual

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.

Reimplemented in amrex::MLCellLinOp.

◆ setMaxOrder()

void amrex::MLLinOp::setMaxOrder ( int  o)
inlinenoexcept

◆ setVerbose()

void amrex::MLLinOp::setVerbose ( int  v)
inlinenoexcept

◆ smooth()

virtual void amrex::MLLinOp::smooth ( int  ,
int  ,
MultiFab ,
const MultiFab ,
bool  = false 
) const
inlinevirtual

◆ solutionResidual()

virtual void amrex::MLLinOp::solutionResidual ( int  ,
MultiFab ,
MultiFab ,
const MultiFab ,
const MultiFab = nullptr 
)
inlinevirtual

Reimplemented in amrex::MLNodeLinOp, and amrex::MLCellLinOp.

◆ supportInhomogNeumannBC()

virtual bool amrex::MLLinOp::supportInhomogNeumannBC ( ) const
inlineprotectedvirtualnoexcept

Reimplemented in amrex::MLCellABecLap.

◆ supportNSolve()

virtual bool amrex::MLLinOp::supportNSolve ( ) const
inlinevirtual

Reimplemented in amrex::MLPoisson, and amrex::MLABecLaplacian.

◆ supportRobinBC()

virtual bool amrex::MLLinOp::supportRobinBC ( ) const
inlineprotectedvirtualnoexcept

Reimplemented in amrex::MLABecLaplacian.

◆ unapplyMetricTerm()

virtual void amrex::MLLinOp::unapplyMetricTerm ( int  ,
int  ,
MultiFab  
) const
inlinevirtual

Reimplemented in amrex::MLNodeLinOp, and amrex::MLCellLinOp.

◆ unimposeNeumannBC()

virtual void amrex::MLLinOp::unimposeNeumannBC ( int  ,
Any  
) const
inlinevirtual

Reimplemented in amrex::MLNodeLaplacian.

◆ update()

virtual void amrex::MLLinOp::update ( )
inlinevirtual

◆ xdoty()

virtual Real amrex::MLLinOp::xdoty ( int  ,
int  ,
const MultiFab ,
const MultiFab ,
bool   
) const
inlinevirtual

Reimplemented in amrex::MLNodeLinOp, and amrex::MLCellLinOp.

Friends And Related Function Documentation

◆ MLABecLaplacian

friend class MLABecLaplacian
friend

◆ MLCGSolver

friend class MLCGSolver
friend

◆ MLMG

friend class MLMG
friend

◆ MLPoisson

friend class MLPoisson
friend

Member Data Documentation

◆ enforceSingularSolvable

bool amrex::MLLinOp::enforceSingularSolvable = true
protected

◆ info

LPInfo amrex::MLLinOp::info
protected

◆ m_amr_ref_ratio

Vector<int> amrex::MLLinOp::m_amr_ref_ratio
protected

◆ m_bottom_comm

MPI_Comm amrex::MLLinOp::m_bottom_comm = MPI_COMM_NULL
protected

◆ m_coarse_bc_loc

RealVect amrex::MLLinOp::m_coarse_bc_loc
protected

◆ m_coarse_data_crse_ratio

int amrex::MLLinOp::m_coarse_data_crse_ratio = -1
protected

◆ m_coarse_data_for_bc

const MultiFab* amrex::MLLinOp::m_coarse_data_for_bc = nullptr
protected

◆ m_default_comm

MPI_Comm amrex::MLLinOp::m_default_comm = MPI_COMM_NULL
protected

◆ m_dmap

Vector<Vector<DistributionMapping> > amrex::MLLinOp::m_dmap
protected

◆ m_do_agglomeration

bool amrex::MLLinOp::m_do_agglomeration = false
protected

◆ m_do_consolidation

bool amrex::MLLinOp::m_do_consolidation = false
protected

◆ m_do_semicoarsening

bool amrex::MLLinOp::m_do_semicoarsening = false
protected

◆ m_domain_bloc_hi

Array<Real, AMREX_SPACEDIM> amrex::MLLinOp::m_domain_bloc_hi {{AMREX_D_DECL(0.,0.,0.)}}
protected

◆ m_domain_bloc_lo

Array<Real, AMREX_SPACEDIM> amrex::MLLinOp::m_domain_bloc_lo {{AMREX_D_DECL(0.,0.,0.)}}
protected

◆ m_domain_covered

Vector<int> amrex::MLLinOp::m_domain_covered
protected

◆ m_factory

Vector<Vector<std::unique_ptr<FabFactory<FArrayBox> > > > amrex::MLLinOp::m_factory
protected

◆ m_geom

Vector<Vector<Geometry> > amrex::MLLinOp::m_geom
protected

first Vector is for amr level and second is mg level

◆ m_grids

Vector<Vector<BoxArray> > amrex::MLLinOp::m_grids
protected

◆ m_hibc

Vector<Array<BCType, AMREX_SPACEDIM> > amrex::MLLinOp::m_hibc
protected

◆ m_hibc_orig

Vector<Array<BCType, AMREX_SPACEDIM> > amrex::MLLinOp::m_hibc_orig
protected

◆ m_ixtype

IntVect amrex::MLLinOp::m_ixtype
protected

◆ m_lobc

Vector<Array<BCType, AMREX_SPACEDIM> > amrex::MLLinOp::m_lobc
protected

◆ m_lobc_orig

Vector<Array<BCType, AMREX_SPACEDIM> > amrex::MLLinOp::m_lobc_orig
protected

◆ m_needs_coarse_data_for_bc

bool amrex::MLLinOp::m_needs_coarse_data_for_bc
protected

◆ m_num_amr_levels

int amrex::MLLinOp::m_num_amr_levels
protected

◆ m_num_mg_levels

Vector<int> amrex::MLLinOp::m_num_mg_levels
protected

◆ m_parent

const MLLinOp* amrex::MLLinOp::m_parent = nullptr
protected

◆ m_raii_comm

std::unique_ptr<CommContainer> amrex::MLLinOp::m_raii_comm
protected

◆ maxorder

int amrex::MLLinOp::maxorder = 3
protected

◆ mg_box_min_width

constexpr int amrex::MLLinOp::mg_box_min_width = 2
staticconstexprprotected

◆ mg_coarsen_ratio

constexpr int amrex::MLLinOp::mg_coarsen_ratio = 2
staticconstexprprotected

◆ mg_coarsen_ratio_vec

Vector<IntVect> amrex::MLLinOp::mg_coarsen_ratio_vec
protected

◆ mg_domain_min_width

constexpr int amrex::MLLinOp::mg_domain_min_width = 2
staticconstexprprotected

◆ verbose

int amrex::MLLinOp::verbose = 0
protected

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