1#ifndef AMREX_ML_CURL_CURL_H_
2#define AMREX_ML_CURL_CURL_H_
3#include <AMReX_Config.H>
6#include <AMReX_MLCurlCurl_K.H>
27 :
public MLLinOpT<Array<MultiFab,3> >
71 [[nodiscard]] std::string
name ()
const override {
72 return std::string(
"curl of curl");
75 bool setUsePCG (
bool flag) {
return std::exchange(m_use_pcg, flag); }
78 const MF* robinbc_a =
nullptr,
79 const MF* robinbc_b =
nullptr,
80 const MF* robinbc_f =
nullptr)
override;
89 void smooth (
int amrlev,
int mglev,
MF& sol,
const MF& rhs,
90 bool skip_fillboundary,
int niter)
const override;
93 const MF* crse_bcdata=
nullptr)
override;
97 const MF* crse_bcdata=
nullptr)
override;
108 [[nodiscard]]
bool isSingular (
int )
const override {
return false; }
111 RT xdoty (
int amrlev,
int mglev,
const MF&
x,
const MF&
y,
bool local)
const override;
113 [[nodiscard]]
RT normInf (
int amrlev,
MF const& mf,
bool local)
const override;
123 [[nodiscard]]
MF make (
int amrlev,
int mglev,
IntVect const& ng)
const override;
133#if (AMREX_SPACEDIM > 1)
134 void smooth4 (
int amrlev,
int mglev,
MF& sol,
MF const& rhs,
int color)
const;
136 void smooth1D (
int amrlev,
int mglev,
MF& sol,
MF const& rhs,
int color)
const;
139 void compresid (
int amrlev,
int mglev,
MF& resid,
MF const& b)
const;
145 void applyBC (
int amrlev,
int mglev,
MF& in, CurlCurlStateType type)
const;
147 [[nodiscard]]
iMultiFab const& getDotMask (
int amrlev,
int mglev,
int idim)
const;
149 [[nodiscard]] CurlCurlDirichletInfo getDirichletInfo (
int amrlev,
int mglev)
const;
150 [[nodiscard]] CurlCurlSymmetryInfo getSymmetryInfo (
int amrlev,
int mglev)
const;
152 void update_lusolver ();
154 void set_curvilinear_domain_bc ();
155 void clearAlphaMultiFab ();
156 void clearBetaMultiFab ();
160 RT m_alpha = std::numeric_limits<RT>::lowest();
161 RT m_beta = std::numeric_limits<RT>::lowest();
164#if (AMREX_SPACEDIM == 3)
166#elif (AMREX_SPACEDIM == 2)
172 mutable Vector<Vector<Array<std::unique_ptr<iMultiFab>,3>>> m_dotmask;
173 static constexpr int m_ncomp = 1;
174 Vector<Vector<std::unique_ptr<Gpu::DeviceScalar
175 <LUSolver<AMREX_SPACEDIM*2,RT>>>>> m_lusolver;
176 Vector<Vector<Array<std::unique_ptr<MultiFab>,3>>> m_bcoefs;
177 Vector<Vector<Array<std::unique_ptr<MultiFab>,3>>> m_acoefs;
178 bool m_use_pcg =
false;
179 bool m_needs_update =
true;
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
curl (alpha curl E) + beta E = rhs
Definition AMReX_MLCurlCurl.H:28
void setLevelBC(int amrlev, const MF *levelbcdata, const MF *robinbc_a=nullptr, const MF *robinbc_b=nullptr, const MF *robinbc_f=nullptr) override
Definition AMReX_MLCurlCurl.cpp:297
void setAlpha(const Vector< MultiFab const * > &a_acoeffs)
Definition AMReX_MLCurlCurl.cpp:134
void preparePrecond() override
Definition AMReX_MLCurlCurl.cpp:906
void restriction(int amrlev, int cmglev, MF &crse, MF &fine) const override
Definition AMReX_MLCurlCurl.cpp:304
void setDirichletNodesToZero(int amrlev, int mglev, MF &a_mf) const override
Definition AMReX_MLCurlCurl.cpp:239
void interpolation(int amrlev, int fmglev, MF &fine, const MF &crse) const override
Definition AMReX_MLCurlCurl.cpp:340
bool isSingular(int) const override
Is it singular on given AMR level?
Definition AMReX_MLCurlCurl.H:108
void update() override
Update for reuse.
Definition AMReX_MLCurlCurl.cpp:1279
typename MLLinOpT< MF >::BCType BCType
Definition AMReX_MLCurlCurl.H:32
typename MLLinOpT< MF >::RT RT
Definition AMReX_MLCurlCurl.H:31
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)
Definition AMReX_MLCurlCurl.cpp:16
void smooth(int amrlev, int mglev, MF &sol, const MF &rhs, bool skip_fillboundary, int niter) const override
Definition AMReX_MLCurlCurl.cpp:531
RT normInf(int amrlev, MF const &mf, bool local) const override
Definition AMReX_MLCurlCurl.cpp:959
void prepareForSolve() override
Definition AMReX_MLCurlCurl.cpp:900
void averageDownAndSync(Vector< MF > &sol) const override
Definition AMReX_MLCurlCurl.cpp:964
void solutionResidual(int amrlev, MF &resid, MF &x, const MF &b, const MF *crse_bcdata=nullptr) override
Definition AMReX_MLCurlCurl.cpp:746
void setScalars(RT a_alpha, RT a_beta) noexcept
Definition AMReX_MLCurlCurl.cpp:56
void make(Vector< Vector< MF > > &mf, IntVect const &ng) const override
Definition AMReX_MLCurlCurl.cpp:976
void smooth4(int amrlev, int mglev, MF &sol, MF const &rhs, int color) const
Definition AMReX_MLCurlCurl.cpp:643
IntVect getNGrowVectRestriction() const override
Definition AMReX_MLCurlCurl.H:117
void setBeta(const Vector< Array< MultiFab const *, 3 > > &a_bcoefs)
Definition AMReX_MLCurlCurl.cpp:66
std::string name() const override
Definition AMReX_MLCurlCurl.H:71
typename MLLinOpT< MF >::Location Location
Definition AMReX_MLCurlCurl.H:35
void correctionResidual(int amrlev, int mglev, MF &resid, MF &x, const MF &b, BCMode bc_mode, const MF *crse_bcdata=nullptr) override
Definition AMReX_MLCurlCurl.cpp:755
typename MLLinOpT< MF >::StateMode StateMode
Definition AMReX_MLCurlCurl.H:34
void compresid(int amrlev, int mglev, MF &resid, MF const &b) const
Definition AMReX_MLCurlCurl.cpp:764
void apply(int amrlev, int mglev, MF &out, MF &in, BCMode bc_mode, StateMode s_mode, const MLMGBndryT< MF > *bndry=nullptr) const override
Definition AMReX_MLCurlCurl.cpp:374
MF makeAlias(MF const &mf) const override
Definition AMReX_MLCurlCurl.cpp:994
RT xdoty(int amrlev, int mglev, const MF &x, const MF &y, bool local) const override
Definition AMReX_MLCurlCurl.cpp:944
bool setUsePCG(bool flag)
Definition AMReX_MLCurlCurl.H:75
MF makeCoarseMG(int amrlev, int mglev, IntVect const &ng) const override
Definition AMReX_MLCurlCurl.cpp:1004
bool needsUpdate() const override
Does it need update if it's reused?
Definition AMReX_MLCurlCurl.H:99
typename MLLinOpT< MF >::BCMode BCMode
Definition AMReX_MLCurlCurl.H:33
void prepareRHS(Vector< MF * > const &rhs) const
Definition AMReX_MLCurlCurl.cpp:230
Array< MultiFab, 3 > MF
Definition AMReX_MLCurlCurl.H:30
void applyPhysBC(int amrlev, int mglev, MultiFab &mf, CurlCurlStateType type) const
Definition AMReX_MLCurlCurl.cpp:1055
MF makeCoarseAmr(int famrlev, IntVect const &ng) const override
Definition AMReX_MLCurlCurl.cpp:1019
bool isBottomSingular() const override
Is the bottom of MG singular?
Definition AMReX_MLCurlCurl.H:109
Definition AMReX_MLLinOp.H:102
typename FabDataType< MF >::value_type RT
Definition AMReX_MLLinOp.H:113
Definition AMReX_MLMGBndry.H:12
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
A Collection of IArrayBoxes.
Definition AMReX_iMultiFab.H:34
std::array< T, N > Array
Definition AMReX_Array.H:26
Definition AMReX_Amr.cpp:49
LinOpBCType
Definition AMReX_LO_BCTYPES.H:27
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
Definition AMReX_MLLinOp.H:36
StateMode
Definition AMReX_MLLinOp.H:89
BCMode
Definition AMReX_MLLinOp.H:88
Location
Definition AMReX_MLLinOp.H:90