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

#include <AMReX_MLCellLinOp.H>

Inheritance diagram for amrex::MLCellLinOpT< MF >:
amrex::MLLinOpT< MF > amrex::MLCellABecLapT< MF > amrex::MLABecLaplacianT< MF > amrex::MLALaplacianT< MF > amrex::MLEBABecLap amrex::MLPoissonT< MF > amrex::MLTensorOp amrex::MLEBTensorOp

Classes

struct  BCTL
 
class  BndryCondLoc
 

Public Types

using FAB = typename MF::fab_type
 
using RT = typename MF::value_type
 
using BCType = LinOpBCType
 
using BCMode = typename MLLinOpT< MF >::BCMode
 
using StateMode = typename MLLinOpT< MF >::StateMode
 
using Location = typename MLLinOpT< MF >::Location
 
- Public Types inherited from amrex::MLLinOpT< MF >
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

 MLCellLinOpT ()
 
 ~MLCellLinOpT () override=default
 
 MLCellLinOpT (const MLCellLinOpT< MF > &)=delete
 
 MLCellLinOpT (MLCellLinOpT< MF > &&)=delete
 
MLCellLinOpT< MF > & operator= (const MLCellLinOpT< MF > &)=delete
 
MLCellLinOpT< MF > & operator= (MLCellLinOpT< MF > &&)=delete
 
void define (const Vector< Geometry > &a_geom, const Vector< BoxArray > &a_grids, const Vector< DistributionMapping > &a_dmap, const LPInfo &a_info=LPInfo(), const Vector< FabFactory< FAB > const * > &a_factory={})
 
void setLevelBC (int amrlev, const MF *levelbcdata, const MF *robinbc_a=nullptr, const MF *robinbc_b=nullptr, const MF *robinbc_f=nullptr) final
 Set boundary conditions for given level. For cell-centered solves only. More...
 
bool needsUpdate () const override
 Does it need update if it's reused? More...
 
void update () override
 Update for reuse. More...
 
void setGaussSeidel (bool flag) noexcept
 
virtual bool isCrossStencil () const
 
virtual bool isTensorOp () const
 
void updateSolBC (int amrlev, const MF &crse_bcdata) const
 
void updateCorBC (int amrlev, const MF &crse_bcdata) const
 
virtual void applyBC (int amrlev, int mglev, MF &in, BCMode bc_mode, StateMode s_mode, const MLMGBndryT< MF > *bndry=nullptr, bool skip_fillboundary=false) const
 
BoxArray makeNGrids (int grid_size) const
 
void restriction (int, int, MF &crse, MF &fine) const override
 Restriction onto coarse MG level. More...
 
void interpolation (int amrlev, int fmglev, MF &fine, const MF &crse) const override
 Add interpolated coarse MG level data to fine MG level data. More...
 
void interpAssign (int amrlev, int fmglev, MF &fine, MF &crse) const override
 Overwrite fine MG level data with interpolated coarse data. More...
 
void interpolationAmr (int famrlev, MF &fine, const MF &crse, IntVect const &nghost) const override
 Interpolation between AMR levels. More...
 
void averageDownSolutionRHS (int camrlev, MF &crse_sol, MF &crse_rhs, const MF &fine_sol, const MF &fine_rhs) override
 Average-down data from fine AMR level to coarse AMR level. More...
 
void apply (int amrlev, int mglev, MF &out, MF &in, BCMode bc_mode, StateMode s_mode, const MLMGBndryT< MF > *bndry=nullptr) const override
 Apply the linear operator, out = L(in) More...
 
void smooth (int amrlev, int mglev, MF &sol, const MF &rhs, bool skip_fillboundary=false) const final
 Smooth. More...
 
void solutionResidual (int amrlev, MF &resid, MF &x, const MF &b, const MF *crse_bcdata=nullptr) override
 Compute residual for solution. More...
 
void prepareForFluxes (int amrlev, const MF *crse_bcdata=nullptr) override
 
void correctionResidual (int amrlev, int mglev, MF &resid, MF &x, const MF &b, BCMode bc_mode, const MF *crse_bcdata=nullptr) final
 Compute residual for the residual-correction form, resid = b - L(x) More...
 
void reflux (int crse_amrlev, MF &res, const MF &crse_sol, const MF &, MF &, MF &fine_sol, const MF &) const final
 Reflux at AMR coarse/fine boundary. More...
 
void compFlux (int amrlev, const Array< MF *, AMREX_SPACEDIM > &fluxes, MF &sol, Location loc) const override
 Compute fluxes. More...
 
void compGrad (int amrlev, const Array< MF *, AMREX_SPACEDIM > &grad, MF &sol, Location loc) const override
 Compute gradients of the solution. More...
 
void applyMetricTerm (int amrlev, int mglev, MF &rhs) const final
 apply metric terms if there are any More...
 
void unapplyMetricTerm (int amrlev, int mglev, MF &rhs) const final
 unapply metric terms if there are any More...
 
Vector< RTgetSolvabilityOffset (int amrlev, int mglev, MF const &rhs) const override
 get offset for fixing solvability More...
 
void fixSolvabilityByOffset (int amrlev, int mglev, MF &rhs, Vector< RT > const &offset) const override
 
void prepareForSolve () override
 
RT xdoty (int amrlev, int mglev, const MF &x, const MF &y, bool local) const final
 x dot y, used by the bottom solver More...
 
virtual void Fapply (int amrlev, int mglev, MF &out, const MF &in) const =0
 
virtual void Fsmooth (int amrlev, int mglev, MF &sol, const MF &rhs, int redblack) const =0
 
virtual void FFlux (int amrlev, const MFIter &mfi, const Array< FAB *, AMREX_SPACEDIM > &flux, const FAB &sol, Location loc, int face_only=0) const =0
 
virtual void addInhomogNeumannFlux (int, const Array< MF *, AMREX_SPACEDIM > &, MF const &, bool) const
 
RT normInf (int amrlev, MF const &mf, bool local) const override
 
void averageDownAndSync (Vector< MF > &sol) const override
 
void avgDownResAmr (int clev, MF &cres, MF const &fres) const override
 
- Public Member Functions inherited from amrex::MLLinOpT< MF >
 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
 
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 void normalize (int amrlev, int mglev, MF &mf) const
 Divide mf by the diagonal component of the operator. Used by bicgstab. 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 void fixSolvabilityByOffset (int, int, MF &, Vector< RT > const &) const
 fix solvability by subtracting offset from RHS More...
 
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 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 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< std::unique_ptr< MF > > m_robin_bcval
 
- Public Attributes inherited from amrex::MLLinOpT< MF >
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 Types

using RealTuple = Array< RT, 2 *BL_SPACEDIM >
 
using BCTuple = Array< BoundCond, 2 *BL_SPACEDIM >
 

Protected Attributes

bool m_has_metric_term = false
 
Vector< std::unique_ptr< MLMGBndryT< MF > > > m_bndry_sol
 
Vector< std::unique_ptr< BndryRegisterT< MF > > > m_crse_sol_br
 
Vector< std::unique_ptr< MLMGBndryT< MF > > > m_bndry_cor
 
Vector< std::unique_ptr< BndryRegisterT< MF > > > m_crse_cor_br
 
Vector< Vector< std::unique_ptr< BndryCondLoc > > > m_bcondloc
 
Vector< Vector< BndryRegisterT< MF > > > m_undrrelxr
 
Vector< Vector< Array< MultiMask, 2 *AMREX_SPACEDIM > > > m_maskvals
 
Vector< std::unique_ptr< iMultiFab > > m_norm_fine_mask
 
Vector< YAFluxRegisterT< MF > > m_fluxreg
 
bool m_use_gauss_seidel = true
 
- Protected Attributes inherited from amrex::MLLinOpT< MF >
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
 

Private Member Functions

void defineAuxData ()
 
void defineBC ()
 
void computeVolInv () const
 

Private Attributes

Vector< Vector< RT > > m_volinv
 
int m_interpbndry_halfwidth = 2
 

Additional Inherited Members

- Protected Member Functions inherited from amrex::MLLinOpT< MF >
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
 
- Static Protected Attributes inherited from amrex::MLLinOpT< MF >
static constexpr int mg_coarsen_ratio = 2
 
static constexpr int mg_box_min_width = 2
 

Member Typedef Documentation

◆ BCMode

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

◆ BCTuple

template<typename MF >
using amrex::MLCellLinOpT< MF >::BCTuple = Array<BoundCond,2*BL_SPACEDIM>
protected

◆ BCType

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

◆ FAB

template<typename MF >
using amrex::MLCellLinOpT< MF >::FAB = typename MF::fab_type

◆ Location

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

◆ RealTuple

template<typename MF >
using amrex::MLCellLinOpT< MF >::RealTuple = Array<RT,2*BL_SPACEDIM>
protected

◆ RT

template<typename MF >
using amrex::MLCellLinOpT< MF >::RT = typename MF::value_type

◆ StateMode

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

Constructor & Destructor Documentation

◆ MLCellLinOpT() [1/3]

template<typename MF >
amrex::MLCellLinOpT< MF >::MLCellLinOpT

◆ ~MLCellLinOpT()

template<typename MF >
amrex::MLCellLinOpT< MF >::~MLCellLinOpT ( )
overridedefault

◆ MLCellLinOpT() [2/3]

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

◆ MLCellLinOpT() [3/3]

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

Member Function Documentation

◆ addInhomogNeumannFlux()

template<typename MF >
virtual void amrex::MLCellLinOpT< MF >::addInhomogNeumannFlux ( int  ,
const Array< MF *, AMREX_SPACEDIM > &  ,
MF const &  ,
bool   
) const
inlinevirtual

Reimplemented in amrex::MLCellABecLapT< MF >.

◆ apply()

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

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

Implements amrex::MLLinOpT< MF >.

◆ applyBC()

template<typename MF >
void amrex::MLCellLinOpT< MF >::applyBC ( int  amrlev,
int  mglev,
MF &  in,
BCMode  bc_mode,
StateMode  s_mode,
const MLMGBndryT< MF > *  bndry = nullptr,
bool  skip_fillboundary = false 
) const
virtual

◆ applyMetricTerm()

template<typename MF >
void amrex::MLCellLinOpT< MF >::applyMetricTerm ( int  ,
int  ,
MF &   
) const
finalvirtual

apply metric terms if there are any

Reimplemented from amrex::MLLinOpT< MF >.

◆ averageDownAndSync()

template<typename MF >
void amrex::MLCellLinOpT< MF >::averageDownAndSync ( Vector< MF > &  sol) const
overridevirtual

Implements amrex::MLLinOpT< MF >.

◆ averageDownSolutionRHS()

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

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 from amrex::MLLinOpT< MF >.

◆ avgDownResAmr()

template<typename MF >
void amrex::MLCellLinOpT< MF >::avgDownResAmr ( int  clev,
MF &  cres,
MF const &  fres 
) const
overridevirtual

Reimplemented from amrex::MLLinOpT< MF >.

◆ compFlux()

template<typename MF >
void amrex::MLCellLinOpT< MF >::compFlux ( int  ,
const Array< MF *, AMREX_SPACEDIM > &  ,
MF &  ,
Location   
) const
overridevirtual

Compute fluxes.

Parameters
amrlevAMR level
fluxesfluxes
solsolution
loclocation of the fluxes

Reimplemented from amrex::MLLinOpT< MF >.

◆ compGrad()

template<typename MF >
void amrex::MLCellLinOpT< MF >::compGrad ( int  ,
const Array< MF *, AMREX_SPACEDIM > &  ,
MF &  ,
Location   
) const
overridevirtual

Compute gradients of the solution.

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

Reimplemented from amrex::MLLinOpT< MF >.

◆ computeVolInv()

template<typename MF >
void amrex::MLCellLinOpT< MF >::computeVolInv
private

◆ correctionResidual()

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

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.

Implements amrex::MLLinOpT< MF >.

◆ define()

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

◆ defineAuxData()

template<typename MF >
void amrex::MLCellLinOpT< MF >::defineAuxData
private

◆ defineBC()

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

◆ Fapply()

template<typename MF >
virtual void amrex::MLCellLinOpT< MF >::Fapply ( int  amrlev,
int  mglev,
MF &  out,
const MF &  in 
) const
pure virtual

◆ FFlux()

template<typename MF >
virtual void amrex::MLCellLinOpT< MF >::FFlux ( int  amrlev,
const MFIter mfi,
const Array< FAB *, AMREX_SPACEDIM > &  flux,
const FAB sol,
Location  loc,
int  face_only = 0 
) const
pure virtual

◆ fixSolvabilityByOffset()

template<typename MF >
void amrex::MLCellLinOpT< MF >::fixSolvabilityByOffset ( int  amrlev,
int  mglev,
MF &  rhs,
Vector< RT > const &  offset 
) const
override

◆ Fsmooth()

template<typename MF >
virtual void amrex::MLCellLinOpT< MF >::Fsmooth ( int  amrlev,
int  mglev,
MF &  sol,
const MF &  rhs,
int  redblack 
) const
pure virtual

◆ getSolvabilityOffset()

template<typename MF >
auto amrex::MLCellLinOpT< MF >::getSolvabilityOffset ( int  ,
int  ,
MF const &   
) const
overridevirtual

get offset for fixing solvability

Reimplemented from amrex::MLLinOpT< MF >.

◆ interpAssign()

template<typename MF >
void amrex::MLCellLinOpT< MF >::interpAssign ( int  amrlev,
int  fmglev,
MF &  fine,
MF &  crse 
) const
overridevirtual

Overwrite fine MG level data with interpolated coarse data.

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

Reimplemented from amrex::MLLinOpT< MF >.

◆ interpolation()

template<typename MF >
void amrex::MLCellLinOpT< MF >::interpolation ( int  amrlev,
int  fmglev,
MF &  fine,
const MF &  crse 
) const
overridevirtual

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

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

Implements amrex::MLLinOpT< MF >.

◆ interpolationAmr()

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

Interpolation between AMR levels.

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

Reimplemented from amrex::MLLinOpT< MF >.

◆ isCrossStencil()

template<typename MF >
virtual bool amrex::MLCellLinOpT< MF >::isCrossStencil ( ) const
inlinevirtual

◆ isTensorOp()

template<typename MF >
virtual bool amrex::MLCellLinOpT< MF >::isTensorOp ( ) const
inlinevirtual

Reimplemented in amrex::MLTensorOp, and amrex::MLEBTensorOp.

◆ makeNGrids()

template<typename MF >
BoxArray amrex::MLCellLinOpT< MF >::makeNGrids ( int  grid_size) const

◆ needsUpdate()

template<typename MF >
bool amrex::MLCellLinOpT< MF >::needsUpdate ( ) const
inlineoverridevirtual

Does it need update if it's reused?

Reimplemented from amrex::MLLinOpT< MF >.

Reimplemented in amrex::MLEBABecLap, amrex::MLTensorOp, and amrex::MLEBTensorOp.

◆ normInf()

template<typename MF >
auto amrex::MLCellLinOpT< MF >::normInf ( int  amrlev,
MF const &  mf,
bool  local 
) const
overridevirtual

Implements amrex::MLLinOpT< MF >.

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ prepareForFluxes()

template<typename MF >
void amrex::MLCellLinOpT< MF >::prepareForFluxes ( int  amrlev,
const MF *  crse_bcdata = nullptr 
)
overridevirtual

Reimplemented from amrex::MLLinOpT< MF >.

◆ prepareForSolve()

template<typename MF >
void amrex::MLCellLinOpT< MF >::prepareForSolve
overridevirtual

◆ reflux()

template<typename MF >
void amrex::MLCellLinOpT< 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
finalvirtual

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 from amrex::MLLinOpT< MF >.

◆ restriction()

template<typename MF >
void amrex::MLCellLinOpT< MF >::restriction ( int  amrlev,
int  cmglev,
MF &  crse,
MF &  fine 
) const
overridevirtual

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.

Implements amrex::MLLinOpT< MF >.

◆ setGaussSeidel()

template<typename MF >
void amrex::MLCellLinOpT< MF >::setGaussSeidel ( bool  flag)
inlinenoexcept

◆ setLevelBC()

template<typename MF >
void amrex::MLCellLinOpT< MF >::setLevelBC ( int  ,
const MF *  ,
const MF *  = nullptr,
const MF *  = nullptr,
const MF *  = nullptr 
)
finalvirtual

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.

Implements amrex::MLLinOpT< MF >.

◆ smooth()

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

Smooth.

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

Implements amrex::MLLinOpT< MF >.

◆ solutionResidual()

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

Compute residual for solution.

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

Implements amrex::MLLinOpT< MF >.

◆ unapplyMetricTerm()

template<typename MF >
void amrex::MLCellLinOpT< MF >::unapplyMetricTerm ( int  ,
int  ,
MF &   
) const
finalvirtual

unapply metric terms if there are any

Reimplemented from amrex::MLLinOpT< MF >.

◆ update()

template<typename MF >
void amrex::MLCellLinOpT< MF >::update
overridevirtual

Update for reuse.

Reimplemented from amrex::MLLinOpT< MF >.

Reimplemented in amrex::MLEBABecLap, amrex::MLTensorOp, and amrex::MLEBTensorOp.

◆ updateCorBC()

template<typename MF >
void amrex::MLCellLinOpT< MF >::updateCorBC ( int  amrlev,
const MF &  crse_bcdata 
) const

◆ updateSolBC()

template<typename MF >
void amrex::MLCellLinOpT< MF >::updateSolBC ( int  amrlev,
const MF &  crse_bcdata 
) const

◆ xdoty()

template<typename MF >
auto amrex::MLCellLinOpT< MF >::xdoty ( int  amrlev,
int  mglev,
const MF &  x,
const MF &  y,
bool  local 
) const
finalvirtual

x dot y, used by the bottom solver

Implements amrex::MLLinOpT< MF >.

Member Data Documentation

◆ m_bcondloc

template<typename MF >
Vector<Vector<std::unique_ptr<BndryCondLoc> > > amrex::MLCellLinOpT< MF >::m_bcondloc
protected

◆ m_bndry_cor

template<typename MF >
Vector<std::unique_ptr<MLMGBndryT<MF> > > amrex::MLCellLinOpT< MF >::m_bndry_cor
protected

◆ m_bndry_sol

template<typename MF >
Vector<std::unique_ptr<MLMGBndryT<MF> > > amrex::MLCellLinOpT< MF >::m_bndry_sol
protected

◆ m_crse_cor_br

template<typename MF >
Vector<std::unique_ptr<BndryRegisterT<MF> > > amrex::MLCellLinOpT< MF >::m_crse_cor_br
protected

◆ m_crse_sol_br

template<typename MF >
Vector<std::unique_ptr<BndryRegisterT<MF> > > amrex::MLCellLinOpT< MF >::m_crse_sol_br
protected

◆ m_fluxreg

template<typename MF >
Vector<YAFluxRegisterT<MF> > amrex::MLCellLinOpT< MF >::m_fluxreg
mutableprotected

◆ m_has_metric_term

template<typename MF >
bool amrex::MLCellLinOpT< MF >::m_has_metric_term = false
protected

◆ m_interpbndry_halfwidth

template<typename MF >
int amrex::MLCellLinOpT< MF >::m_interpbndry_halfwidth = 2
private

◆ m_maskvals

template<typename MF >
Vector<Vector<Array<MultiMask,2*AMREX_SPACEDIM> > > amrex::MLCellLinOpT< MF >::m_maskvals
protected

◆ m_norm_fine_mask

template<typename MF >
Vector<std::unique_ptr<iMultiFab> > amrex::MLCellLinOpT< MF >::m_norm_fine_mask
protected

◆ m_robin_bcval

template<typename MF >
Vector<std::unique_ptr<MF> > amrex::MLCellLinOpT< MF >::m_robin_bcval

◆ m_undrrelxr

template<typename MF >
Vector<Vector<BndryRegisterT<MF> > > amrex::MLCellLinOpT< MF >::m_undrrelxr
mutableprotected

◆ m_use_gauss_seidel

template<typename MF >
bool amrex::MLCellLinOpT< MF >::m_use_gauss_seidel = true
protected

◆ m_volinv

template<typename MF >
Vector<Vector<RT> > amrex::MLCellLinOpT< MF >::m_volinv
mutableprivate

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