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

#include <AMReX_MLCellLinOp.H>

Inheritance diagram for amrex::MLCellLinOp:
amrex::MLLinOp amrex::MLCellABecLap amrex::MLABecLaplacian amrex::MLALaplacian amrex::MLEBABecLap amrex::MLPoisson amrex::MLTensorOp amrex::MLEBTensorOp

Classes

struct  BCTL
 
class  BndryCondLoc
 

Public Member Functions

 MLCellLinOp ()
 
virtual ~MLCellLinOp ()
 
 MLCellLinOp (const MLCellLinOp &)=delete
 
 MLCellLinOp (MLCellLinOp &&)=delete
 
MLCellLinOpoperator= (const MLCellLinOp &)=delete
 
MLCellLinOpoperator= (MLCellLinOp &&)=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< FArrayBox > const * > &a_factory={})
 
virtual void setLevelBC (int amrlev, const MultiFab *levelbcdata, const MultiFab *robinbc_a=nullptr, const MultiFab *robinbc_b=nullptr, const MultiFab *robinbc_f=nullptr) final override
 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...
 
virtual bool needsUpdate () const override
 
virtual void update () override
 
virtual bool isCrossStencil () const
 
virtual bool isTensorOp () const
 
void updateSolBC (int amrlev, const MultiFab &crse_bcdata) const
 
void updateCorBC (int amrlev, const MultiFab &crse_bcdata) const
 
virtual void applyBC (int amrlev, int mglev, MultiFab &in, BCMode bc_mode, StateMode s_mode, const MLMGBndry *bndry=nullptr, bool skip_fillboundary=false) const
 
BoxArray makeNGrids (int grid_size) const
 
virtual void restriction (int, int, MultiFab &crse, MultiFab &fine) const override
 
virtual void interpolation (int amrlev, int fmglev, MultiFab &fine, const MultiFab &crse) const override
 
virtual void interpAssign (int amrlev, int fmglev, MultiFab &fine, MultiFab &crse) const override
 
virtual void averageDownSolutionRHS (int camrlev, MultiFab &crse_sol, MultiFab &crse_rhs, const MultiFab &fine_sol, const MultiFab &fine_rhs) override
 
virtual void apply (int amrlev, int mglev, MultiFab &out, MultiFab &in, BCMode bc_mode, StateMode s_mode, const MLMGBndry *bndry=nullptr) const override
 
virtual void smooth (int amrlev, int mglev, MultiFab &sol, const MultiFab &rhs, bool skip_fillboundary=false) const final override
 
virtual void solutionResidual (int amrlev, MultiFab &resid, MultiFab &x, const MultiFab &b, const MultiFab *crse_bcdata=nullptr) override
 
virtual void correctionResidual (int amrlev, int mglev, MultiFab &resid, MultiFab &x, const MultiFab &b, BCMode bc_mode, const MultiFab *crse_bcdata=nullptr) final override
 
virtual void reflux (int crse_amrlev, MultiFab &res, const MultiFab &crse_sol, const MultiFab &, MultiFab &, MultiFab &fine_sol, const MultiFab &) const final override
 
virtual void compFlux (int amrlev, const Array< MultiFab *, AMREX_SPACEDIM > &fluxes, MultiFab &sol, Location loc) const override
 
virtual void compGrad (int amrlev, const Array< MultiFab *, AMREX_SPACEDIM > &grad, MultiFab &sol, Location loc) const override
 
virtual void applyMetricTerm (int amrlev, int mglev, Any &rhs) const final override
 
virtual void unapplyMetricTerm (int amrlev, int mglev, MultiFab &rhs) const final override
 
virtual Vector< Real > getSolvabilityOffset (int amrlev, int mglev, Any const &rhs) const override
 
virtual void fixSolvabilityByOffset (int amrlev, int mglev, Any &rhs, Vector< Real > const &offset) const override
 
virtual void prepareForSolve () override
 
virtual Real xdoty (int amrlev, int mglev, const MultiFab &x, const MultiFab &y, bool local) const final override
 
virtual void Fapply (int amrlev, int mglev, MultiFab &out, const MultiFab &in) const =0
 
virtual void Fsmooth (int amrlev, int mglev, MultiFab &sol, const MultiFab &rsh, int redblack) const =0
 
virtual void FFlux (int amrlev, const MFIter &mfi, const Array< FArrayBox *, AMREX_SPACEDIM > &flux, const FArrayBox &sol, Location loc, const int face_only=0) const =0
 
void applyMetricTermToMF (int amrlev, int mglev, MultiFab &rhs) const
 
virtual Real AnyNormInfMask (int amrlev, Any const &a, bool local) const override
 
virtual void AnyAvgDownResAmr (int clev, Any &cres, Any const &fres) const override
 
virtual void AnyInterpolationAmr (int famrlev, Any &fine, const Any &crse, IntVect const &) const override
 
virtual void AnyAverageDownAndSync (Vector< Any > &sol) const override
 
virtual void addInhomogNeumannFlux (int, const Array< MultiFab *, AMREX_SPACEDIM > &, MultiFab const &, bool) const
 
- Public Member Functions inherited from amrex::MLLinOp
 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...
 
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 void normalize (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 bool isSingular (int) const
 
virtual bool isBottomSingular () 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 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 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 AnyAverageDownSolutionRHS (int camrlev, Any &crse_sol, Any &crse_rhs, const Any &fine_sol, const Any &fine_rhs)
 
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
 

Protected Types

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

Protected Attributes

bool m_has_metric_term = false
 
Vector< std::unique_ptr< MLMGBndry > > m_bndry_sol
 
Vector< std::unique_ptr< BndryRegister > > m_crse_sol_br
 
Vector< std::unique_ptr< MLMGBndry > > m_bndry_cor
 
Vector< std::unique_ptr< BndryRegister > > m_crse_cor_br
 
Vector< Vector< std::unique_ptr< BndryCondLoc > > > m_bcondloc
 
Vector< std::unique_ptr< MultiFab > > m_robin_bcval
 
Vector< Vector< BndryRegister > > m_undrrelxr
 
Vector< Vector< Array< MultiMask, 2 *AMREX_SPACEDIM > > > m_maskvals
 
Vector< std::unique_ptr< iMultiFab > > m_norm_fine_mask
 
Vector< YAFluxRegisterm_fluxreg
 
- Protected Attributes inherited from amrex::MLLinOp
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
 

Private Member Functions

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

Private Attributes

Vector< Vector< Real > > m_volinv
 

Friends

class MLMG
 
class MLCGSolver
 

Additional Inherited Members

- Public Types inherited from amrex::MLLinOp
enum class  BCMode { Homogeneous , Inhomogeneous }
 
enum class  StateMode { Solution , Correction }
 
enum class  Location { FaceCenter , FaceCentroid , CellCenter , CellCentroid }
 
using BCType = LinOpBCType
 
- Static Public Member Functions inherited from amrex::MLLinOp
static void Initialize ()
 
static void Finalize ()
 
- Protected Member Functions inherited from amrex::MLLinOp
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
 
- Static Protected Attributes inherited from amrex::MLLinOp
static constexpr int mg_coarsen_ratio = 2
 
static constexpr int mg_box_min_width = 2
 
static constexpr int mg_domain_min_width = 2
 

Member Typedef Documentation

◆ BCTuple

◆ RealTuple

using amrex::MLCellLinOp::RealTuple = Array<Real,2*BL_SPACEDIM>
protected

Constructor & Destructor Documentation

◆ MLCellLinOp() [1/3]

amrex::MLCellLinOp::MLCellLinOp ( )

◆ ~MLCellLinOp()

amrex::MLCellLinOp::~MLCellLinOp ( )
virtual

◆ MLCellLinOp() [2/3]

amrex::MLCellLinOp::MLCellLinOp ( const MLCellLinOp )
delete

◆ MLCellLinOp() [3/3]

amrex::MLCellLinOp::MLCellLinOp ( MLCellLinOp &&  )
delete

Member Function Documentation

◆ addInhomogNeumannFlux()

virtual void amrex::MLCellLinOp::addInhomogNeumannFlux ( int  ,
const Array< MultiFab *, AMREX_SPACEDIM > &  ,
MultiFab const &  ,
bool   
) const
inlinevirtual

Reimplemented in amrex::MLCellABecLap.

◆ AnyAverageDownAndSync()

void amrex::MLCellLinOp::AnyAverageDownAndSync ( Vector< Any > &  sol) const
overridevirtual

Implements amrex::MLLinOp.

◆ AnyAvgDownResAmr()

void amrex::MLCellLinOp::AnyAvgDownResAmr ( int  clev,
Any cres,
Any const &  fres 
) const
overridevirtual

Implements amrex::MLLinOp.

◆ AnyInterpolationAmr()

void amrex::MLCellLinOp::AnyInterpolationAmr ( int  famrlev,
Any fine,
const Any crse,
IntVect const &   
) const
overridevirtual

Implements amrex::MLLinOp.

◆ AnyNormInfMask()

Real amrex::MLCellLinOp::AnyNormInfMask ( int  amrlev,
Any const &  a,
bool  local 
) const
overridevirtual

Implements amrex::MLLinOp.

◆ apply()

void amrex::MLCellLinOp::apply ( int  amrlev,
int  mglev,
MultiFab out,
MultiFab in,
BCMode  bc_mode,
StateMode  s_mode,
const MLMGBndry bndry = nullptr 
) const
overridevirtual

Reimplemented from amrex::MLLinOp.

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

◆ applyBC()

void amrex::MLCellLinOp::applyBC ( int  amrlev,
int  mglev,
MultiFab in,
BCMode  bc_mode,
StateMode  s_mode,
const MLMGBndry bndry = nullptr,
bool  skip_fillboundary = false 
) const
virtual

Reimplemented in amrex::MLEBABecLap.

◆ applyMetricTerm()

void amrex::MLCellLinOp::applyMetricTerm ( int  amrlev,
int  mglev,
Any rhs 
) const
finaloverridevirtual

Reimplemented from amrex::MLLinOp.

◆ applyMetricTermToMF()

void amrex::MLCellLinOp::applyMetricTermToMF ( int  amrlev,
int  mglev,
MultiFab rhs 
) const

◆ averageDownSolutionRHS()

void amrex::MLCellLinOp::averageDownSolutionRHS ( int  camrlev,
MultiFab crse_sol,
MultiFab crse_rhs,
const MultiFab fine_sol,
const MultiFab fine_rhs 
)
overridevirtual

Reimplemented from amrex::MLLinOp.

Reimplemented in amrex::MLEBABecLap.

◆ compFlux()

void amrex::MLCellLinOp::compFlux ( int  amrlev,
const Array< MultiFab *, AMREX_SPACEDIM > &  fluxes,
MultiFab sol,
Location  loc 
) const
overridevirtual

Reimplemented from amrex::MLLinOp.

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

◆ compGrad()

void amrex::MLCellLinOp::compGrad ( int  amrlev,
const Array< MultiFab *, AMREX_SPACEDIM > &  grad,
MultiFab sol,
Location  loc 
) const
overridevirtual

Reimplemented from amrex::MLLinOp.

Reimplemented in amrex::MLEBABecLap.

◆ computeVolInv()

void amrex::MLCellLinOp::computeVolInv ( ) const
private

◆ correctionResidual()

void amrex::MLCellLinOp::correctionResidual ( int  amrlev,
int  mglev,
MultiFab resid,
MultiFab x,
const MultiFab b,
BCMode  bc_mode,
const MultiFab crse_bcdata = nullptr 
)
finaloverridevirtual

Reimplemented from amrex::MLLinOp.

◆ define()

void amrex::MLCellLinOp::define ( const Vector< Geometry > &  a_geom,
const Vector< BoxArray > &  a_grids,
const Vector< DistributionMapping > &  a_dmap,
const LPInfo a_info = LPInfo(),
const Vector< FabFactory< FArrayBox > const * > &  a_factory = {} 
)

◆ defineAuxData()

void amrex::MLCellLinOp::defineAuxData ( )
private

◆ defineBC()

void amrex::MLCellLinOp::defineBC ( )
private

◆ Fapply()

virtual void amrex::MLCellLinOp::Fapply ( int  amrlev,
int  mglev,
MultiFab out,
const MultiFab in 
) const
pure virtual

◆ FFlux()

virtual void amrex::MLCellLinOp::FFlux ( int  amrlev,
const MFIter mfi,
const Array< FArrayBox *, AMREX_SPACEDIM > &  flux,
const FArrayBox sol,
Location  loc,
const int  face_only = 0 
) const
pure virtual

◆ fixSolvabilityByOffset()

void amrex::MLCellLinOp::fixSolvabilityByOffset ( int  amrlev,
int  mglev,
Any rhs,
Vector< Real > const &  offset 
) const
overridevirtual

Reimplemented from amrex::MLLinOp.

◆ Fsmooth()

virtual void amrex::MLCellLinOp::Fsmooth ( int  amrlev,
int  mglev,
MultiFab sol,
const MultiFab rsh,
int  redblack 
) const
pure virtual

◆ getSolvabilityOffset()

Vector< Real > amrex::MLCellLinOp::getSolvabilityOffset ( int  amrlev,
int  mglev,
Any const &  rhs 
) const
overridevirtual

Reimplemented from amrex::MLLinOp.

◆ interpAssign()

void amrex::MLCellLinOp::interpAssign ( int  amrlev,
int  fmglev,
MultiFab fine,
MultiFab crse 
) const
overridevirtual

Reimplemented from amrex::MLLinOp.

◆ interpolation()

void amrex::MLCellLinOp::interpolation ( int  amrlev,
int  fmglev,
MultiFab fine,
const MultiFab crse 
) const
overridevirtual

Reimplemented from amrex::MLLinOp.

Reimplemented in amrex::MLEBABecLap.

◆ isCrossStencil()

virtual bool amrex::MLCellLinOp::isCrossStencil ( ) const
inlinevirtual

◆ isTensorOp()

virtual bool amrex::MLCellLinOp::isTensorOp ( ) const
inlinevirtual

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

◆ makeNGrids()

BoxArray amrex::MLCellLinOp::makeNGrids ( int  grid_size) const

◆ needsUpdate()

virtual bool amrex::MLCellLinOp::needsUpdate ( ) const
inlineoverridevirtual

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ prepareForSolve()

void amrex::MLCellLinOp::prepareForSolve ( )
overridevirtual

◆ reflux()

void amrex::MLCellLinOp::reflux ( int  crse_amrlev,
MultiFab res,
const MultiFab crse_sol,
const MultiFab ,
MultiFab ,
MultiFab fine_sol,
const MultiFab  
) const
finaloverridevirtual

Reimplemented from amrex::MLLinOp.

◆ restriction()

void amrex::MLCellLinOp::restriction ( int  amrlev,
int  cmglev,
MultiFab crse,
MultiFab fine 
) const
overridevirtual

Reimplemented from amrex::MLLinOp.

Reimplemented in amrex::MLEBABecLap.

◆ setLevelBC()

void amrex::MLCellLinOp::setLevelBC ( int  ,
const MultiFab ,
const MultiFab = nullptr,
const MultiFab = nullptr,
const MultiFab = nullptr 
)
finaloverridevirtual

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

◆ smooth()

void amrex::MLCellLinOp::smooth ( int  amrlev,
int  mglev,
MultiFab sol,
const MultiFab rhs,
bool  skip_fillboundary = false 
) const
finaloverridevirtual

Reimplemented from amrex::MLLinOp.

◆ solutionResidual()

void amrex::MLCellLinOp::solutionResidual ( int  amrlev,
MultiFab resid,
MultiFab x,
const MultiFab b,
const MultiFab crse_bcdata = nullptr 
)
overridevirtual

Reimplemented from amrex::MLLinOp.

◆ unapplyMetricTerm()

void amrex::MLCellLinOp::unapplyMetricTerm ( int  amrlev,
int  mglev,
MultiFab rhs 
) const
finaloverridevirtual

Reimplemented from amrex::MLLinOp.

◆ update()

void amrex::MLCellLinOp::update ( )
overridevirtual

◆ updateCorBC()

void amrex::MLCellLinOp::updateCorBC ( int  amrlev,
const MultiFab crse_bcdata 
) const

◆ updateSolBC()

void amrex::MLCellLinOp::updateSolBC ( int  amrlev,
const MultiFab crse_bcdata 
) const

◆ xdoty()

Real amrex::MLCellLinOp::xdoty ( int  amrlev,
int  mglev,
const MultiFab x,
const MultiFab y,
bool  local 
) const
finaloverridevirtual

Reimplemented from amrex::MLLinOp.

Friends And Related Function Documentation

◆ MLCGSolver

friend class MLCGSolver
friend

◆ MLMG

friend class MLMG
friend

Member Data Documentation

◆ m_bcondloc

Vector<Vector<std::unique_ptr<BndryCondLoc> > > amrex::MLCellLinOp::m_bcondloc
protected

◆ m_bndry_cor

Vector<std::unique_ptr<MLMGBndry> > amrex::MLCellLinOp::m_bndry_cor
protected

◆ m_bndry_sol

Vector<std::unique_ptr<MLMGBndry> > amrex::MLCellLinOp::m_bndry_sol
protected

◆ m_crse_cor_br

Vector<std::unique_ptr<BndryRegister> > amrex::MLCellLinOp::m_crse_cor_br
protected

◆ m_crse_sol_br

Vector<std::unique_ptr<BndryRegister> > amrex::MLCellLinOp::m_crse_sol_br
protected

◆ m_fluxreg

Vector<YAFluxRegister> amrex::MLCellLinOp::m_fluxreg
mutableprotected

◆ m_has_metric_term

bool amrex::MLCellLinOp::m_has_metric_term = false
protected

◆ m_maskvals

Vector<Vector<Array<MultiMask,2*AMREX_SPACEDIM> > > amrex::MLCellLinOp::m_maskvals
protected

◆ m_norm_fine_mask

Vector<std::unique_ptr<iMultiFab> > amrex::MLCellLinOp::m_norm_fine_mask
protected

◆ m_robin_bcval

Vector<std::unique_ptr<MultiFab> > amrex::MLCellLinOp::m_robin_bcval
protected

◆ m_undrrelxr

Vector<Vector<BndryRegister> > amrex::MLCellLinOp::m_undrrelxr
mutableprotected

◆ m_volinv

Vector<Vector<Real> > amrex::MLCellLinOp::m_volinv
mutableprivate

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