Block-Structured AMR Software Framework
 
Loading...
Searching...
No Matches
AMReX_GMRES_MLMG.H
Go to the documentation of this file.
1#ifndef AMREX_GMRES_MLMG_H_
2#define AMREX_GMRES_MLMG_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_GMRES.H>
6#include <AMReX_MLMG.H>
7#include <AMReX_TypeTraits.H>
8#include <utility>
9
10namespace amrex {
11
19template <typename MF>
21{
22public:
23 using VEC = Vector<MF>;
24 using MG = MLMGT<MF>;
25 using RT = typename MG::RT; // double or float
27
28 explicit GMRESMLMGT (MG& mlmg);
29
38 void solve (MF& a_sol, MF const& a_rhs, RT a_tol_rel, RT a_tol_abs);
39
40 void solve (Vector<MF*> const& a_sol, Vector<MF const*> const& a_rhs, RT a_tol_rel, RT a_tol_abs);
41
43 void setVerbose (int v) {
44 m_verbose = v;
46 }
47
49 void setMaxIters (int niters) { m_gmres.setMaxIters(niters); }
50
52 [[nodiscard]] int getNumIters () const { return m_gmres.getNumIters(); }
53
55 [[nodiscard]] RT getResidualNorm () const { return m_gmres.getResidualNorm(); }
56
58 GM& getGMRES () { return m_gmres; }
59
61 VEC makeVecRHS () const;
62
64 VEC makeVecLHS () const;
65
66 RT norm2 (VEC const& mf) const;
67
68 static void scale (VEC& mf, RT scale_factor);
69
70 RT dotProduct (VEC const& mf1, VEC const& mf2) const;
71
73 static void setToZero (VEC& lhs);
74
76 static void assign (VEC& lhs, VEC const& rhs);
77
79 static void increment (VEC& lhs, VEC const& rhs, RT a);
80
82 static void linComb (VEC& lhs, RT a, VEC const& rhs_a, RT b, VEC const& rhs_b);
83
85 void apply (VEC& lhs, VEC const& rhs) const;
86
87 void precond (VEC& lhs, VEC const& rhs) const;
88
90 bool usePrecond (bool new_flag) { return std::exchange(m_use_precond, new_flag); }
91
93 void setPrecondNumIters (int precond_niters) { m_precond_niters = precond_niters; }
94
95private:
99 int m_verbose = 0;
100 int m_nlevels = 0;
101 bool m_use_precond = true;
103};
104
105template <typename MF>
107 : m_mlmg(&mlmg), m_linop(&mlmg.getLinOp()), m_nlevels(m_linop->NAMRLevels())
108{
110 m_gmres.define(*this);
111}
112
113template <typename MF>
115{
116 VEC vmf(m_nlevels);
117 for (int ilev = 0; ilev < m_nlevels; ++ilev) {
118 vmf[ilev] = m_linop->make(ilev, 0, IntVect(0));
119 }
120 return vmf;
121}
122
123template <typename MF>
125{
126 VEC vmf(m_nlevels);
127 for (int ilev = 0; ilev < m_nlevels; ++ilev) {
128 vmf[ilev] = m_linop->make(ilev, 0, IntVect(1));
129 setBndry(vmf[ilev], RT(0), 0, nComp(vmf[ilev]));
130 }
131 return vmf;
132}
133
134template <typename MF>
135auto GMRESMLMGT<MF>::norm2 (VEC const& mf) const -> RT
136{
137 return m_linop->norm2Precond(GetVecOfConstPtrs(mf));
138}
139
140template <typename MF>
141void GMRESMLMGT<MF>::scale (VEC& mf, RT scale_factor)
142{
143 for (auto& xmf : mf) {
144 Scale(xmf, scale_factor, 0, nComp(xmf), 0);
145 }
146}
147
148template <typename MF>
149auto GMRESMLMGT<MF>::dotProduct (VEC const& mf1, VEC const& mf2) const -> RT
150{
151 return m_linop->dotProductPrecond(GetVecOfConstPtrs(mf1), GetVecOfConstPtrs(mf2));
152}
153
154template <typename MF>
156{
157 for (auto& xmf : lhs) {
158 setVal(xmf, RT(0.0));
159 }
160}
161
162template <typename MF>
163void GMRESMLMGT<MF>::assign (VEC& lhs, VEC const& rhs)
164{
165 auto nlevels = int(lhs.size());
166 for (int ilev = 0; ilev < nlevels; ++ilev) {
167 LocalCopy(lhs[ilev], rhs[ilev], 0, 0, nComp(lhs[ilev]), IntVect(0));
168 }
169}
170
171template <typename MF>
172void GMRESMLMGT<MF>::increment (VEC& lhs, VEC const& rhs, RT a)
173{
174 auto nlevels = int(lhs.size());
175 for (int ilev = 0; ilev < nlevels; ++ilev) {
176 Saxpy(lhs[ilev], a, rhs[ilev], 0, 0, nComp(lhs[ilev]), IntVect(0));
177 }
178}
179
180template <typename MF>
181void GMRESMLMGT<MF>::linComb (VEC& lhs, RT a, VEC const& rhs_a, RT b, VEC const& rhs_b)
182{
183 auto nlevels = int(lhs.size());
184 for (int ilev = 0; ilev < nlevels; ++ilev) {
185 LinComb(lhs[ilev], a, rhs_a[ilev], 0, b, rhs_b[ilev], 0, 0, nComp(lhs[ilev]), IntVect(0));
186 }
187}
188
189template <typename MF>
190void GMRESMLMGT<MF>::apply (VEC& lhs, VEC const& rhs) const
191{
192 m_mlmg->applyPrecond(GetVecOfPtrs(lhs), GetVecOfPtrs(const_cast<VEC&>(rhs)));
193}
194
195template <typename MF>
196void GMRESMLMGT<MF>::precond (VEC& lhs, VEC const& rhs) const
197{
198 if (m_use_precond) {
199 m_mlmg->setPrecondIter(m_precond_niters);
200 setToZero(lhs);
201 m_mlmg->precond(GetVecOfPtrs(lhs), GetVecOfConstPtrs(rhs), 0, 0);
202 } else {
203 for (int ilev = 0; ilev < m_nlevels; ++ilev) {
204 LocalCopy(lhs[ilev], rhs[ilev], 0, 0, nComp(lhs[ilev]), IntVect(0));
205 }
206 }
207}
208
209template <typename MF>
210void GMRESMLMGT<MF>::solve (MF& a_sol, MF const& a_rhs, RT a_tol_rel, RT a_tol_abs)
211{
212 AMREX_ALWAYS_ASSERT(m_nlevels == 1);
213 this->solve({&a_sol}, {&a_rhs}, a_tol_rel, a_tol_abs);
214}
215
216template <typename MF>
217void GMRESMLMGT<MF>::solve (Vector<MF*> const& a_sol, Vector<MF const*> const& a_rhs, RT a_tol_rel, RT a_tol_abs)
218{
219 m_mlmg->incPrintIdentation();
220 auto mlmg_verbose = m_mlmg->getVerbose();
221 auto mlmg_bottom_verbose = m_mlmg->getBottomVerbose();
222 m_mlmg->setVerbose(m_verbose);
223 auto mlmg_bottom_solver = m_mlmg->getBottomSolver();
224
225 if (mlmg_bottom_solver != BottomSolver::smoother &&
226 mlmg_bottom_solver != BottomSolver::hypre &&
227 mlmg_bottom_solver != BottomSolver::petsc)
228 {
229 m_mlmg->setBottomSolver(BottomSolver::smoother);
230 }
231
232 auto res = makeVecLHS();
233 auto cor = makeVecLHS();
234
235 m_mlmg->apply(GetVecOfPtrs(res), a_sol); // res = L(sol)
236 // res = L(sol) - rhs
237 bool need_to_scale_rhs = m_linop->scaleRHS(0,nullptr);
238 for (int ilev = 0; ilev < m_nlevels; ++ilev) {
239 MF const* prhs;
240 if (need_to_scale_rhs) {
241 LocalCopy(cor[ilev], *a_rhs[ilev], 0, 0, nComp(cor[ilev]), IntVect(0));
242 auto r = m_linop->scaleRHS(ilev, &(cor[ilev]));
244 prhs = &(cor[ilev]);
245 } else {
246 prhs = a_rhs[ilev];
247 }
248 Saxpy(res[ilev], RT(-1), *prhs, 0, 0, nComp(res[ilev]), IntVect(0));
249 }
250 for (int ilev = 0; ilev < m_nlevels; ++ilev) {
251 m_linop->setDirichletNodesToZero(ilev,0,res[ilev]);
252 }
253 m_gmres.solve(cor, res, a_tol_rel, a_tol_abs); // L(cor) = res
254 // sol = sol - cor
255 for (int ilev = 0; ilev < m_nlevels; ++ilev) {
256 Saxpy(*a_sol[ilev], RT(-1), cor[ilev], 0, 0, nComp(*a_sol[ilev]), IntVect(0));
257 }
258
259 m_mlmg->setBottomSolver(mlmg_bottom_solver);
260 m_mlmg->setVerbose(mlmg_verbose);
261 m_mlmg->setBottomVerbose(mlmg_bottom_verbose);
262 m_mlmg->decPrintIdentation();
263}
264
266
267}
268
269#endif
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
Solve using GMRES with multigrid as preconditioner.
Definition AMReX_GMRES_MLMG.H:21
void solve(MF &a_sol, MF const &a_rhs, RT a_tol_rel, RT a_tol_abs)
Solve the linear system.
Definition AMReX_GMRES_MLMG.H:210
static void assign(VEC &lhs, VEC const &rhs)
lhs = rhs
Definition AMReX_GMRES_MLMG.H:163
VEC makeVecRHS() const
Make MultiFab without ghost cells.
Definition AMReX_GMRES_MLMG.H:114
int m_precond_niters
Definition AMReX_GMRES_MLMG.H:102
GMRESMLMGT(MG &mlmg)
Definition AMReX_GMRES_MLMG.H:106
GM m_gmres
Definition AMReX_GMRES_MLMG.H:96
int getNumIters() const
Gets the number of iterations.
Definition AMReX_GMRES_MLMG.H:52
void setPrecondNumIters(int precond_niters)
Set the number of MLMG preconditioner iterations per GMRES iteration.
Definition AMReX_GMRES_MLMG.H:93
Vector< MF > VEC
Definition AMReX_GMRES_MLMG.H:23
static void setToZero(VEC &lhs)
lhs = 0
Definition AMReX_GMRES_MLMG.H:155
GM & getGMRES()
Get the GMRES object.
Definition AMReX_GMRES_MLMG.H:58
void setVerbose(int v)
Sets verbosity.
Definition AMReX_GMRES_MLMG.H:43
typename MG::RT RT
Definition AMReX_GMRES_MLMG.H:25
RT dotProduct(VEC const &mf1, VEC const &mf2) const
Definition AMReX_GMRES_MLMG.H:149
VEC makeVecLHS() const
Make MultiFab with ghost cells and set ghost cells to zero.
Definition AMReX_GMRES_MLMG.H:124
RT norm2(VEC const &mf) const
Definition AMReX_GMRES_MLMG.H:135
static void scale(VEC &mf, RT scale_factor)
Definition AMReX_GMRES_MLMG.H:141
static void linComb(VEC &lhs, RT a, VEC const &rhs_a, RT b, VEC const &rhs_b)
lhs = a*rhs_a + b*rhs_b
Definition AMReX_GMRES_MLMG.H:181
MLLinOpT< MF > * m_linop
Definition AMReX_GMRES_MLMG.H:98
static void increment(VEC &lhs, VEC const &rhs, RT a)
lhs += a*rhs
Definition AMReX_GMRES_MLMG.H:172
int m_verbose
Definition AMReX_GMRES_MLMG.H:99
bool m_use_precond
Definition AMReX_GMRES_MLMG.H:101
MG * m_mlmg
Definition AMReX_GMRES_MLMG.H:97
RT getResidualNorm() const
Gets the 2-norm of the residual.
Definition AMReX_GMRES_MLMG.H:55
int m_nlevels
Definition AMReX_GMRES_MLMG.H:100
bool usePrecond(bool new_flag)
Control whether or not to use MLMG as preconditioner.
Definition AMReX_GMRES_MLMG.H:90
void apply(VEC &lhs, VEC const &rhs) const
lhs = L(rhs)
Definition AMReX_GMRES_MLMG.H:190
void setMaxIters(int niters)
Sets the max number of iterations.
Definition AMReX_GMRES_MLMG.H:49
void precond(VEC &lhs, VEC const &rhs) const
Definition AMReX_GMRES_MLMG.H:196
int getNumIters() const
Gets the number of iterations.
Definition AMReX_GMRES.H:113
void define(M &linop)
Definition AMReX_GMRES.H:187
void setVerbose(int v)
Sets verbosity.
Definition AMReX_GMRES.H:104
void setMaxIters(int niters)
Sets the max number of iterations.
Definition AMReX_GMRES.H:110
RT getResidualNorm() const
Gets the 2-norm of the residual.
Definition AMReX_GMRES.H:119
Definition AMReX_MLLinOp.H:98
Definition AMReX_MLMG.H:12
void preparePrecond()
Definition AMReX_MLMG.H:1180
typename MLLinOpT< MF >::RT RT
Definition AMReX_MLMG.H:27
Long size() const noexcept
Definition AMReX_Vector.H:50
Definition AMReX_Amr.cpp:49
int nComp(FabArrayBase const &fa)
void Saxpy(MF &dst, typename MF::value_type a, MF const &src, int scomp, int dcomp, int ncomp, IntVect const &nghost)
dst += a * src
Definition AMReX_FabArrayUtility.H:1847
void LocalCopy(DMF &dst, SMF const &src, int scomp, int dcomp, int ncomp, IntVect const &nghost)
dst = src
Definition AMReX_FabArrayUtility.H:1831
Vector< const T * > GetVecOfConstPtrs(const Vector< T > &a)
Definition AMReX_Vector.H:90
IntVectND< AMREX_SPACEDIM > IntVect
Definition AMReX_BaseFwd.H:30
void LinComb(MF &dst, typename MF::value_type a, MF const &src_a, int acomp, typename MF::value_type b, MF const &src_b, int bcomp, int dcomp, int ncomp, IntVect const &nghost)
dst = a*src_a + b*src_b
Definition AMReX_FabArrayUtility.H:1863
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:111
void setBndry(MF &dst, typename MF::value_type val, int scomp, int ncomp)
dst = val in ghost cells.
Definition AMReX_FabArrayUtility.H:1815
Vector< T * > GetVecOfPtrs(Vector< T > &a)
Definition AMReX_Vector.H:61
void Scale(MF &dst, typename MF::value_type val, int scomp, int ncomp, int nghost)
dst *= val
Definition AMReX_FabArrayUtility.H:1822
const int[]
Definition AMReX_BLProfiler.cpp:1664
void setVal(MF &dst, typename MF::value_type val)
dst = val
Definition AMReX_FabArrayUtility.H:1808