Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_MLCurlCurl.H
Go to the documentation of this file.
1#ifndef AMREX_ML_CURL_CURL_H_
2#define AMREX_ML_CURL_CURL_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_MLLinOp.H>
6#include <AMReX_MLCurlCurl_K.H>
7
8namespace amrex {
9
27 : public MLLinOpT<Array<MultiFab,3> >
28{
29public:
31 using RT = typename MLLinOpT<MF>::RT;
32 using BCType = typename MLLinOpT<MF>::BCType;
33 using BCMode = typename MLLinOpT<MF>::BCMode;
36
37 MLCurlCurl () = default;
38 MLCurlCurl (const Vector<Geometry>& a_geom,
39 const Vector<BoxArray>& a_grids,
40 const Vector<DistributionMapping>& a_dmap,
41 const LPInfo& a_info = LPInfo(),
42 int a_coord = 0);
43
44 void define (const Vector<Geometry>& a_geom,
45 const Vector<BoxArray>& a_grids,
46 const Vector<DistributionMapping>& a_dmap,
47 const LPInfo& a_info = LPInfo(),
48 int a_coord = 0);
49
53 void setScalars (RT a_alpha, RT a_beta) noexcept;
54
58 void setBeta (const Vector<Array<MultiFab const*,3>>& a_bcoefs);
59
63 void setAlpha (const Vector<MultiFab const*>& a_acoeffs);
64
67 void prepareRHS (Vector<MF*> const& rhs) const;
68
69 void setDirichletNodesToZero (int amrlev, int mglev, MF& a_mf) const override;
70
71 [[nodiscard]] std::string name () const override {
72 return std::string("curl of curl");
73 }
74
75 bool setUsePCG (bool flag) { return std::exchange(m_use_pcg, flag); }
76
77 void setLevelBC (int amrlev, const MF* levelbcdata,
78 const MF* robinbc_a = nullptr,
79 const MF* robinbc_b = nullptr,
80 const MF* robinbc_f = nullptr) override;
81
82 void restriction (int amrlev, int cmglev, MF& crse, MF& fine) const override;
83
84 void interpolation (int amrlev, int fmglev, MF& fine, const MF& crse) const override;
85
86 void apply (int amrlev, int mglev, MF& out, MF& in, BCMode bc_mode,
87 StateMode s_mode, const MLMGBndryT<MF>* bndry=nullptr) const override;
88
89 void smooth (int amrlev, int mglev, MF& sol, const MF& rhs,
90 bool skip_fillboundary, int niter) const override;
91
92 void solutionResidual (int amrlev, MF& resid, MF& x, const MF& b,
93 const MF* crse_bcdata=nullptr) override;
94
95 void correctionResidual (int amrlev, int mglev, MF& resid, MF& x,
96 const MF& b, BCMode bc_mode,
97 const MF* crse_bcdata=nullptr) override;
98
99 [[nodiscard]] bool needsUpdate () const override {
100 return (m_needs_update || MLLinOpT<Array<MultiFab,3>>::needsUpdate());
101 }
102 void update () override;
103
104 void prepareForSolve () override;
105
106 void preparePrecond () override;
107
108 [[nodiscard]] bool isSingular (int /*amrlev*/) const override { return false; }
109 [[nodiscard]] bool isBottomSingular () const override { return false; }
110
111 RT xdoty (int amrlev, int mglev, const MF& x, const MF& y, bool local) const override;
112
113 [[nodiscard]] RT normInf (int amrlev, MF const& mf, bool local) const override;
114
115 void averageDownAndSync (Vector<MF>& sol) const override;
116
117 [[nodiscard]] IntVect getNGrowVectRestriction () const override {
118 return IntVect(1);
119 }
120
121 void make (Vector<Vector<MF> >& mf, IntVect const& ng) const override;
122
123 [[nodiscard]] MF make (int amrlev, int mglev, IntVect const& ng) const override;
124
125 [[nodiscard]] MF makeAlias (MF const& mf) const override;
126
127 [[nodiscard]] MF makeCoarseMG (int amrlev, int mglev, IntVect const& ng) const override;
128
129 [[nodiscard]] MF makeCoarseAmr (int famrlev, IntVect const& ng) const override;
130
131// public for cuda
132
133#if (AMREX_SPACEDIM > 1)
134 void smooth4 (int amrlev, int mglev, MF& sol, MF const& rhs, int color) const;
135#else
136 void smooth1D (int amrlev, int mglev, MF& sol, MF const& rhs, int color) const;
137#endif
138
139 void compresid (int amrlev, int mglev, MF& resid, MF const& b) const;
140
141 void applyPhysBC (int amrlev, int mglev, MultiFab& mf, CurlCurlStateType type) const;
142
143private:
144
145 void applyBC (int amrlev, int mglev, MF& in, CurlCurlStateType type) const;
146
147 [[nodiscard]] iMultiFab const& getDotMask (int amrlev, int mglev, int idim) const;
148
149 [[nodiscard]] CurlCurlDirichletInfo getDirichletInfo (int amrlev, int mglev) const;
150 [[nodiscard]] CurlCurlSymmetryInfo getSymmetryInfo (int amrlev, int mglev) const;
151
152 void update_lusolver ();
153
154 void set_curvilinear_domain_bc ();
155 void clearAlphaMultiFab ();
156 void clearBetaMultiFab ();
157
158 int m_coord = 0;
159
160 RT m_alpha = std::numeric_limits<RT>::lowest();
161 RT m_beta = std::numeric_limits<RT>::lowest();
162
163 Array<IntVect,3> m_etype
164#if (AMREX_SPACEDIM == 3)
165 {IntVect(0,1,1), IntVect(1,0,1), IntVect(1,1,0)};
166#elif (AMREX_SPACEDIM == 2)
167 {IntVect(0,1), IntVect(1,0), IntVect(1,1)};
168#else
169 {IntVect(0), IntVect(1), IntVect(1)};
170#endif
171
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;
180};
181
182}
183
184#endif
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
MLCurlCurl()=default
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