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>
7
8namespace amrex {
9
25 : public MLLinOpT<Array<MultiFab,3> >
26{
27public:
29 using RT = typename MLLinOpT<MF>::RT;
30 using BCType = typename MLLinOpT<MF>::BCType;
31 using BCMode = typename MLLinOpT<MF>::BCMode;
34
35 MLCurlCurl () = default;
36 MLCurlCurl (const Vector<Geometry>& a_geom,
37 const Vector<BoxArray>& a_grids,
38 const Vector<DistributionMapping>& a_dmap,
39 const LPInfo& a_info = LPInfo(),
40 int a_coord = 0);
41
42 void define (const Vector<Geometry>& a_geom,
43 const Vector<BoxArray>& a_grids,
44 const Vector<DistributionMapping>& a_dmap,
45 const LPInfo& a_info = LPInfo(),
46 int a_coord = 0);
47
48 void setScalars (RT a_alpha, RT a_beta) noexcept;
49
51 void setBeta (const Vector<Array<MultiFab const*,3>>& a_bcoefs);
52
55 void prepareRHS (Vector<MF*> const& rhs) const;
56
57 void setDirichletNodesToZero (int amrlev, int mglev, MF& a_mf) const override;
58
59 [[nodiscard]] std::string name () const override {
60 return std::string("curl of curl");
61 }
62
63 bool setUsePCG (bool flag) { return std::exchange(m_use_pcg, flag); }
64
65 void setLevelBC (int amrlev, const MF* levelbcdata,
66 const MF* robinbc_a = nullptr,
67 const MF* robinbc_b = nullptr,
68 const MF* robinbc_f = nullptr) override;
69
70 void restriction (int amrlev, int cmglev, MF& crse, MF& fine) const override;
71
72 void interpolation (int amrlev, int fmglev, MF& fine, const MF& crse) const override;
73
74 void apply (int amrlev, int mglev, MF& out, MF& in, BCMode bc_mode,
75 StateMode s_mode, const MLMGBndryT<MF>* bndry=nullptr) const override;
76
77 void smooth (int amrlev, int mglev, MF& sol, const MF& rhs,
78 bool skip_fillboundary, int niter) const override;
79
80 void solutionResidual (int amrlev, MF& resid, MF& x, const MF& b,
81 const MF* crse_bcdata=nullptr) override;
82
83 void correctionResidual (int amrlev, int mglev, MF& resid, MF& x,
84 const MF& b, BCMode bc_mode,
85 const MF* crse_bcdata=nullptr) override;
86
87 [[nodiscard]] bool needsUpdate () const override {
89 }
90 void update () override;
91
92 void prepareForSolve () override;
93
94 void preparePrecond () override;
95
96 [[nodiscard]] bool isSingular (int /*amrlev*/) const override { return false; }
97 [[nodiscard]] bool isBottomSingular () const override { return false; }
98
99 RT xdoty (int amrlev, int mglev, const MF& x, const MF& y, bool local) const override;
100
101 [[nodiscard]] RT normInf (int amrlev, MF const& mf, bool local) const override;
102
103 void averageDownAndSync (Vector<MF>& sol) const override;
104
105 [[nodiscard]] IntVect getNGrowVectRestriction () const override {
106 return IntVect(1);
107 }
108
109 void make (Vector<Vector<MF> >& mf, IntVect const& ng) const override;
110
111 [[nodiscard]] MF make (int amrlev, int mglev, IntVect const& ng) const override;
112
113 [[nodiscard]] MF makeAlias (MF const& mf) const override;
114
115 [[nodiscard]] MF makeCoarseMG (int amrlev, int mglev, IntVect const& ng) const override;
116
117 [[nodiscard]] MF makeCoarseAmr (int famrlev, IntVect const& ng) const override;
118
119// public for cuda
120
121#if (AMREX_SPACEDIM > 1)
122 void smooth4 (int amrlev, int mglev, MF& sol, MF const& rhs, int color) const;
123#else
124 void smooth1D (int amrlev, int mglev, MF& sol, MF const& rhs, int color) const;
125#endif
126
127 void compresid (int amrlev, int mglev, MF& resid, MF const& b) const;
128
129 void applyPhysBC (int amrlev, int mglev, MultiFab& mf, CurlCurlStateType type) const;
130
131private:
132
133 void applyBC (int amrlev, int mglev, MF& in, CurlCurlStateType type) const;
134
135 [[nodiscard]] iMultiFab const& getDotMask (int amrlev, int mglev, int idim) const;
136
137 [[nodiscard]] CurlCurlDirichletInfo getDirichletInfo (int amrlev, int mglev) const;
138 [[nodiscard]] CurlCurlSymmetryInfo getSymmetryInfo (int amrlev, int mglev) const;
139
140 void update_lusolver ();
141
143
144 int m_coord = 0;
145
146 RT m_alpha = std::numeric_limits<RT>::lowest();
147 RT m_beta = std::numeric_limits<RT>::lowest();
148
150#if (AMREX_SPACEDIM == 3)
151 {IntVect(0,1,1), IntVect(1,0,1), IntVect(1,1,0)};
152#elif (AMREX_SPACEDIM == 2)
153 {IntVect(0,1), IntVect(1,0), IntVect(1,1)};
154#else
155 {IntVect(0), IntVect(1), IntVect(1)};
156#endif
157
159 static constexpr int m_ncomp = 1;
160 Vector<Vector<std::unique_ptr<Gpu::DeviceScalar
163 bool m_use_pcg = false;
164 bool m_needs_update = true;
165};
166
167}
168
169#endif
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
Definition AMReX_LUSolver.H:16
curl (alpha curl E) + beta E = rhs
Definition AMReX_MLCurlCurl.H:26
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:185
void preparePrecond() override
Definition AMReX_MLCurlCurl.cpp:661
void applyBC(int amrlev, int mglev, MF &in, CurlCurlStateType type) const
Definition AMReX_MLCurlCurl.cpp:766
void restriction(int amrlev, int cmglev, MF &crse, MF &fine) const override
Definition AMReX_MLCurlCurl.cpp:192
void setDirichletNodesToZero(int amrlev, int mglev, MF &a_mf) const override
Definition AMReX_MLCurlCurl.cpp:127
void interpolation(int amrlev, int fmglev, MF &fine, const MF &crse) const override
Definition AMReX_MLCurlCurl.cpp:225
bool isSingular(int) const override
Is it singular on given AMR level?
Definition AMReX_MLCurlCurl.H:96
CurlCurlSymmetryInfo getSymmetryInfo(int amrlev, int mglev) const
Definition AMReX_MLCurlCurl.cpp:974
void update() override
Update for reuse.
Definition AMReX_MLCurlCurl.cpp:1012
CurlCurlDirichletInfo getDirichletInfo(int amrlev, int mglev) const
Definition AMReX_MLCurlCurl.cpp:932
typename MLLinOpT< MF >::BCType BCType
Definition AMReX_MLCurlCurl.H:30
typename MLLinOpT< MF >::RT RT
Definition AMReX_MLCurlCurl.H:29
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:13
void smooth(int amrlev, int mglev, MF &sol, const MF &rhs, bool skip_fillboundary, int niter) const override
Definition AMReX_MLCurlCurl.cpp:360
RT normInf(int amrlev, MF const &mf, bool local) const override
Definition AMReX_MLCurlCurl.cpp:692
void prepareForSolve() override
Definition AMReX_MLCurlCurl.cpp:655
Vector< Vector< Array< std::unique_ptr< iMultiFab >, 3 > > > m_dotmask
Definition AMReX_MLCurlCurl.H:158
Vector< Vector< Array< std::unique_ptr< MultiFab >, 3 > > > m_bcoefs
Definition AMReX_MLCurlCurl.H:162
void averageDownAndSync(Vector< MF > &sol) const override
Definition AMReX_MLCurlCurl.cpp:697
void solutionResidual(int amrlev, MF &resid, MF &x, const MF &b, const MF *crse_bcdata=nullptr) override
Definition AMReX_MLCurlCurl.cpp:502
void set_curvilinear_domain_bc()
Definition AMReX_MLCurlCurl.cpp:666
void setScalars(RT a_alpha, RT a_beta) noexcept
Definition AMReX_MLCurlCurl.cpp:48
void make(Vector< Vector< MF > > &mf, IntVect const &ng) const override
Definition AMReX_MLCurlCurl.cpp:709
RT m_beta
Definition AMReX_MLCurlCurl.H:147
void smooth4(int amrlev, int mglev, MF &sol, MF const &rhs, int color) const
Definition AMReX_MLCurlCurl.cpp:439
IntVect getNGrowVectRestriction() const override
Definition AMReX_MLCurlCurl.H:105
void setBeta(const Vector< Array< MultiFab const *, 3 > > &a_bcoefs)
This is needed only if there is variable beta coefficient.
Definition AMReX_MLCurlCurl.cpp:56
static constexpr int m_ncomp
Definition AMReX_MLCurlCurl.H:159
MLCurlCurl()=default
RT m_alpha
Definition AMReX_MLCurlCurl.H:146
bool m_use_pcg
Definition AMReX_MLCurlCurl.H:163
std::string name() const override
Definition AMReX_MLCurlCurl.H:59
iMultiFab const & getDotMask(int amrlev, int mglev, int idim) const
Definition AMReX_MLCurlCurl.cpp:921
bool m_needs_update
Definition AMReX_MLCurlCurl.H:164
typename MLLinOpT< MF >::Location Location
Definition AMReX_MLCurlCurl.H:33
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:511
typename MLLinOpT< MF >::StateMode StateMode
Definition AMReX_MLCurlCurl.H:32
void compresid(int amrlev, int mglev, MF &resid, MF const &b) const
Definition AMReX_MLCurlCurl.cpp:520
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:256
MF makeAlias(MF const &mf) const override
Definition AMReX_MLCurlCurl.cpp:727
RT xdoty(int amrlev, int mglev, const MF &x, const MF &y, bool local) const override
Definition AMReX_MLCurlCurl.cpp:677
void update_lusolver()
Definition AMReX_MLCurlCurl.cpp:566
bool setUsePCG(bool flag)
Definition AMReX_MLCurlCurl.H:63
MF makeCoarseMG(int amrlev, int mglev, IntVect const &ng) const override
Definition AMReX_MLCurlCurl.cpp:737
int m_coord
Definition AMReX_MLCurlCurl.H:144
bool needsUpdate() const override
Does it need update if it's reused?
Definition AMReX_MLCurlCurl.H:87
typename MLLinOpT< MF >::BCMode BCMode
Definition AMReX_MLCurlCurl.H:31
void prepareRHS(Vector< MF * > const &rhs) const
Definition AMReX_MLCurlCurl.cpp:118
Array< MultiFab, 3 > MF
Definition AMReX_MLCurlCurl.H:28
void applyPhysBC(int amrlev, int mglev, MultiFab &mf, CurlCurlStateType type) const
Definition AMReX_MLCurlCurl.cpp:788
Vector< Vector< std::unique_ptr< Gpu::DeviceScalar< LUSolver< 3 *2, RT > > > > > m_lusolver
Definition AMReX_MLCurlCurl.H:161
Array< IntVect, 3 > m_etype
Definition AMReX_MLCurlCurl.H:151
MF makeCoarseAmr(int famrlev, IntVect const &ng) const override
Definition AMReX_MLCurlCurl.cpp:752
bool isBottomSingular() const override
Is the bottom of MG singular?
Definition AMReX_MLCurlCurl.H:97
Definition AMReX_MLLinOp.H:98
LinOpBCType BCType
Definition AMReX_MLLinOp.H:111
typename FabDataType< MF >::value_type RT
Definition AMReX_MLLinOp.H:109
Definition AMReX_MLMGBndry.H:12
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:38
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
Definition AMReX_iMultiFab.H:32
Definition AMReX_Amr.cpp:49
CurlCurlStateType
Definition AMReX_MLCurlCurl_K.H:333
IntVectND< 3 > IntVect
Definition AMReX_BaseFwd.H:30
std::array< T, N > Array
Definition AMReX_Array.H:24
Definition AMReX_MLCurlCurl_K.H:185
Definition AMReX_MLCurlCurl_K.H:266
Definition AMReX_GpuMemory.H:56
Definition AMReX_MLLinOp.H:35
StateMode
Definition AMReX_MLLinOp.H:86
BCMode
Definition AMReX_MLLinOp.H:85
Location
Definition AMReX_MLLinOp.H:87