MLCurlCurl ()=default
MLCurlCurl (const Vector < Geometry > &a_geom, const Vector < BoxArray > &a_grids, const Vector < DistributionMapping > &a_dmap, const LPInfo &a_info=LPInfo (), int a_coord=0)
void define (const Vector < Geometry > &a_geom, const Vector < BoxArray > &a_grids, const Vector < DistributionMapping > &a_dmap, const LPInfo &a_info=LPInfo (), int a_coord=0)
void setScalars (RT a_alpha, RT a_beta) noexcept
void setBeta (const Vector < Array < MultiFab const *, 3 > > &a_bcoefs)
void setAlpha (const Vector < MultiFab const * > &a_acoeffs)
void prepareRHS (Vector < MF * > const &rhs) const
void setDirichletNodesToZero (int amrlev, int mglev, MF &a_mf) const override
std::string name () const override
bool setUsePCG (bool flag)
void setLevelBC (int amrlev, const MF *levelbcdata, const MF *robinbc_a=nullptr, const MF *robinbc_b=nullptr, const MF *robinbc_f=nullptr) override
void restriction (int amrlev, int cmglev, MF &crse , MF &fine ) const override
void interpolation (int amrlev, int fmglev, MF &fine , const MF &crse ) const override
void apply (int amrlev, int mglev, MF &out, MF &in, BCMode bc_mode, StateMode s_mode, const MLMGBndryT < MF > *bndry=nullptr) const override
void smooth (int amrlev, int mglev, MF &sol, const MF &rhs, bool skip_fillboundary, int niter) const override
void solutionResidual (int amrlev, MF &resid, MF &x , const MF &b, 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) override
bool needsUpdate () const override
Does it need update if it's reused?
void update () override
Update for reuse.
void prepareForSolve () override
void preparePrecond () override
bool isSingular (int ) const override
Is it singular on given AMR level?
bool isBottomSingular () const override
Is the bottom of MG singular?
RT xdoty (int amrlev, int mglev, const MF &x , const MF &y , bool local) const override
RT normInf (int amrlev, MF const &mf, bool local) const override
void averageDownAndSync (Vector < MF > &sol) const override
IntVect getNGrowVectRestriction () const override
void make (Vector < Vector < MF > > &mf, IntVect const &ng) const override
MF make (int amrlev, int mglev, IntVect const &ng) const override
MF makeAlias (MF const &mf) const override
MF makeCoarseMG (int amrlev, int mglev, IntVect const &ng) const override
MF makeCoarseAmr (int famrlev, IntVect const &ng) const override
void smooth4 (int amrlev, int mglev, MF &sol, MF const &rhs, int color) const
void compresid (int amrlev, int mglev, MF &resid, MF const &b) const
void applyPhysBC (int amrlev, int mglev, MultiFab &mf, CurlCurlStateType type) const
MLLinOpT ()=default
MLLinOpT (const MLLinOpT < Array < MultiFab , 3 > > &)=delete
MLLinOpT (MLLinOpT < Array < MultiFab , 3 > > &&)=delete
virtual ~MLLinOpT ()=default
MLLinOpT < Array < MultiFab , 3 > > & operator= (const MLLinOpT < Array < MultiFab , 3 > > &)=delete
MLLinOpT < Array < MultiFab , 3 > > & operator= (MLLinOpT < Array < MultiFab , 3 > > &&)=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)
void setDomainBC (const Array < BCType , 3 > &lobc, const Array < BCType , 3 > &hibc) noexcept
Boundary of the whole domain.
void setDomainBC (const Vector < Array < BCType , 3 > > &lobc, const Vector < Array < BCType , 3 > > &hibc)
Boundary of the whole domain.
void setDomainBCLoc (const Array < Real , 3 > &lo_bcloc, const Array < Real , 3 > &hi_bcloc) noexcept
Set location of domain boundaries.
bool needsCoarseDataForBC () const noexcept
Needs coarse data for bc?
void setCoarseFineBC (const Array < MultiFab , 3 > *crse , int crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
Set coarse/fine boundary conditions. For cell-centered solves only.
void setCoarseFineBC (const Array < MultiFab , 3 > *crse , IntVect const &crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
void setCoarseFineBC (const AMF *crse , int crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
void setCoarseFineBC (const AMF *crse , IntVect const &crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
virtual void setLevelBC (int, const Array < MultiFab , 3 > *, const Array < MultiFab , 3 > *=nullptr, const Array < MultiFab , 3 > *=nullptr, const Array < MultiFab , 3 > *=nullptr)=0
Set boundary conditions for given level. For cell-centered solves only.
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.
void setMaxOrder (int o) noexcept
Set order of interpolation at coarse/fine boundary.
int getMaxOrder () const noexcept
Get order of interpolation at coarse/fine boundary.
void setEnforceSingularSolvable (bool o) noexcept
bool getEnforceSingularSolvable () const noexcept
virtual BottomSolver getDefaultBottomSolver () const
virtual int getNComp () const
Return number of components.
virtual int getNGrow (int=0, int=0) const
virtual void restriction (int amrlev, int cmglev, Array < MultiFab , 3 > &crse , Array < MultiFab , 3 > &fine ) const=0
Restriction onto coarse MG level.
virtual void interpolation (int amrlev, int fmglev, Array < MultiFab , 3 > &fine , const Array < MultiFab , 3 > &crse ) const=0
Add interpolated coarse MG level data to fine MG level data.
virtual void interpAssign (int amrlev, int fmglev, Array < MultiFab , 3 > &fine , Array < MultiFab , 3 > &crse ) const
Overwrite fine MG level data with interpolated coarse data.
virtual void interpolationAmr (int famrlev, Array < MultiFab , 3 > &fine , const Array < MultiFab , 3 > &crse , IntVect const &nghost) const
Interpolation between AMR levels.
virtual void averageDownSolutionRHS (int camrlev, Array < MultiFab , 3 > &crse_sol, Array < MultiFab , 3 > &crse_rhs, const Array < MultiFab , 3 > &fine_sol, const Array < MultiFab , 3 > &fine_rhs)
Average-down data from fine AMR level to coarse AMR level.
virtual void apply (int amrlev, int mglev, Array < MultiFab , 3 > &out, Array < MultiFab , 3 > &in, BCMode bc_mode, StateMode s_mode, const MLMGBndryT < Array < MultiFab , 3 > > *bndry=nullptr) const=0
Apply the linear operator, out = L(in)
virtual void smooth (int amrlev, int mglev, Array < MultiFab , 3 > &sol, const Array < MultiFab , 3 > &rhs, bool skip_fillboundary, int niter) const=0
Smooth.
virtual void normalize (int amrlev, int mglev, Array < MultiFab , 3 > &mf) const
Divide mf by the diagonal component of the operator. Used by bicgstab.
virtual void solutionResidual (int amrlev, Array < MultiFab , 3 > &resid, Array < MultiFab , 3 > &x, const Array < MultiFab , 3 > &b, const Array < MultiFab , 3 > *crse_bcdata=nullptr)=0
Compute residual for solution.
virtual void prepareForFluxes (int, const Array < MultiFab , 3 > *=nullptr)
virtual void correctionResidual (int amrlev, int mglev, Array < MultiFab , 3 > &resid, Array < MultiFab , 3 > &x, const Array < MultiFab , 3 > &b, BCMode bc_mode, const Array < MultiFab , 3 > *crse_bcdata=nullptr)=0
Compute residual for the residual-correction form, resid = b - L(x)
virtual void reflux (int crse_amrlev, Array < MultiFab , 3 > &res, const Array < MultiFab , 3 > &crse_sol, const Array < MultiFab , 3 > &crse_rhs, Array < MultiFab , 3 > &fine_res, Array < MultiFab , 3 > &fine_sol, const Array < MultiFab , 3 > &fine_rhs) const
Reflux at AMR coarse/fine boundary.
virtual void compFlux (int amrlev, const Array < Array < MultiFab , 3 > *, 3 > &fluxes, Array < MultiFab , 3 > &sol, Location loc) const
Compute fluxes.
virtual void compGrad (int amrlev, const Array < Array < MultiFab , 3 > *, 3 > &grad, Array < MultiFab , 3 > &sol, Location loc) const
Compute gradients of the solution.
virtual void applyMetricTerm (int, int, Array < MultiFab , 3 > &) const
apply metric terms if there are any
virtual void unapplyMetricTerm (int, int, Array < MultiFab , 3 > &) const
unapply metric terms if there are any
virtual void unimposeNeumannBC (int, Array < MultiFab , 3 > &) const
This is needed for our nodal projection solver.
virtual void applyInhomogNeumannTerm (int, Array < MultiFab , 3 > &) const
Extra terms introduced when we treat inhomogeneous Nuemann BC as homogeneous.
virtual void applyOverset (int, Array < MultiFab , 3 > &) const
for overset solver only
virtual bool scaleRHS (int, Array < MultiFab , 3 > *) const
scale RHS to fix solvability
virtual Vector < RT > getSolvabilityOffset (int, int, Array < MultiFab , 3 > const &) const
get offset for fixing solvability
virtual void fixSolvabilityByOffset (int, int, Array < MultiFab , 3 > &, Vector < RT > const &) const
fix solvability by subtracting offset from RHS
virtual void setDirichletNodesToZero (int, int, Array < MultiFab , 3 > &) const
virtual RT xdoty (int amrlev, int mglev, const Array < MultiFab , 3 > &x, const Array < MultiFab , 3 > &y, bool local) const=0
x dot y, used by the bottom solver
virtual RT dotProductPrecond (Vector < Array < MultiFab , 3 > const * > const &x, Vector < Array < MultiFab , 3 > const * > const &y) const
virtual RT norm2Precond (Vector < Array < MultiFab , 3 > const * > const &x) const
virtual std::unique_ptr< MLLinOpT < Array < MultiFab , 3 > > > makeNLinOp (int) const
virtual void getFluxes (const Vector < Array < Array < MultiFab , 3 > *, 3 > > &, const Vector < Array < MultiFab , 3 > * > &, Location ) const
virtual void getFluxes (const Vector < Array < MultiFab , 3 > * > &, const Vector < Array < MultiFab , 3 > * > &) const
virtual void getEBFluxes (const Vector < Array < MultiFab , 3 > * > &, const Vector < Array < MultiFab , 3 > * > &) const
virtual bool supportNSolve () const
virtual void copyNSolveSolution (Array < MultiFab , 3 > &, Array < MultiFab , 3 > const &) const
virtual void postSolve (Vector < Array < MultiFab , 3 > * > const &) const
virtual RT normInf (int amrlev, Array < MultiFab , 3 > const &mf, bool local) const=0
virtual void averageDownAndSync (Vector < Array < MultiFab , 3 > > &sol) const=0
virtual void avgDownResAmr (int clev, Array < MultiFab , 3 > &cres, Array < MultiFab , 3 > const &fres) const
virtual void avgDownResMG (int clev, Array < MultiFab , 3 > &cres, Array < MultiFab , 3 > const &fres) const
virtual void beginPrecondBC ()
virtual void endPrecondBC ()
bool isMFIterSafe (int amrlev, int mglev1, int mglev2) const
int NAMRLevels () const noexcept
Return the number of AMR levels.
int NMGLevels (int amrlev) const noexcept
Return the number of MG levels at given AMR level.
const Geometry & Geom (int amr_lev, int mglev=0) const noexcept
curl (alpha curl E) + beta E = rhs
Here E is an Array of 3 MultiFabs on the staggered grid. The coefficient can be supplied either as a positive scalar (via setScalars) or as a nodal MultiFab (via setAlpha), and can be either a non-negative scalar or an edge-centered MultiFab (via setBeta).
It's the caller's responsibility to make sure rhs has consistent nodal data. If needed, one could call prepareRHS for this.
The smoother is based on the 4-color Gauss-Seidel smoother of Li et. al. 2020. "An Efficient Preconditioner for 3-D Finite Difference
Modeling of the Electromagnetic Diffusion Process in the Frequency
Domain", IEEE Transactions on Geoscience and Remote Sensing, 58, 500-509.