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
26template <typename MF>
28{
29public:
30 using VEC = Vector<MF>;
31 using MG = MLMGT<MF>;
32 using RT = typename MG::RT; // double or float
34
40 explicit GMRESMLMGT (MG& mlmg);
41
50 void solve (MF& a_sol, MF const& a_rhs, RT a_tol_rel, RT a_tol_abs);
51
60 void solve (Vector<MF*> const& a_sol, Vector<MF const*> const& a_rhs, RT a_tol_rel, RT a_tol_abs);
61
63 void setVerbose (int v) {
64 m_verbose = v;
65 m_gmres.setVerbose(v);
66 }
67
69 void setMaxIters (int niters) { m_gmres.setMaxIters(niters); }
70
72 [[nodiscard]] int getNumIters () const { return m_gmres.getNumIters(); }
73
75 [[nodiscard]] RT getResidualNorm () const { return m_gmres.getResidualNorm(); }
76
78 GM& getGMRES () { return m_gmres; }
79
81 VEC makeVecRHS () const;
82
85 VEC makeVecLHS () const;
86
88 RT norm2 (VEC const& mf) const;
89
91 static void scale (VEC& mf, RT scale_factor);
92
94 RT dotProduct (VEC const& mf1, VEC const& mf2) const;
95
101 static void setToZero (VEC& lhs);
102
109 static void assign (VEC& lhs, VEC const& rhs);
110
118 static void increment (VEC& lhs, VEC const& rhs, RT a);
119
129 static void linComb (VEC& lhs, RT a, VEC const& rhs_a, RT b, VEC const& rhs_b);
130
137 void apply (VEC& lhs, VEC const& rhs) const;
138
145 void precond (VEC& lhs, VEC const& rhs) const;
146
153 bool usePrecond (bool new_flag) { return std::exchange(m_use_precond, new_flag); }
154
160 void setPrecondNumIters (int precond_niters) { m_precond_niters = precond_niters; }
161
162private:
163 GM m_gmres;
164 MG* m_mlmg;
165 MLLinOpT<MF>* m_linop;
166 int m_verbose = 0;
167 int m_nlevels = 0;
168 bool m_use_precond = true;
169 int m_precond_niters = 1;
170};
171
172template <typename MF>
174 : m_mlmg(&mlmg), m_linop(&mlmg.getLinOp()), m_nlevels(m_linop->NAMRLevels())
175{
176 m_mlmg->preparePrecond();
177 m_gmres.define(*this);
178}
179
180template <typename MF>
182{
183 VEC vmf(m_nlevels);
184 for (int ilev = 0; ilev < m_nlevels; ++ilev) {
185 vmf[ilev] = m_linop->make(ilev, 0, IntVect(0));
186 }
187 return vmf;
188}
189
190template <typename MF>
192{
193 VEC vmf(m_nlevels);
194 for (int ilev = 0; ilev < m_nlevels; ++ilev) {
195 vmf[ilev] = m_linop->make(ilev, 0, IntVect(1));
196 setBndry(vmf[ilev], RT(0), 0, nComp(vmf[ilev]));
197 }
198 return vmf;
199}
200
201template <typename MF>
202auto GMRESMLMGT<MF>::norm2 (VEC const& mf) const -> RT
203{
204 return m_linop->norm2Precond(GetVecOfConstPtrs(mf));
205}
206
207template <typename MF>
208void GMRESMLMGT<MF>::scale (VEC& mf, RT scale_factor)
209{
210 for (auto& xmf : mf) {
211 Scale(xmf, scale_factor, 0, nComp(xmf), 0);
212 }
213}
214
215template <typename MF>
216auto GMRESMLMGT<MF>::dotProduct (VEC const& mf1, VEC const& mf2) const -> RT
217{
218 return m_linop->dotProductPrecond(GetVecOfConstPtrs(mf1), GetVecOfConstPtrs(mf2));
219}
220
221template <typename MF>
223{
224 for (auto& xmf : lhs) {
225 setVal(xmf, RT(0.0));
226 }
227}
228
229template <typename MF>
230void GMRESMLMGT<MF>::assign (VEC& lhs, VEC const& rhs)
231{
232 auto nlevels = int(lhs.size());
233 for (int ilev = 0; ilev < nlevels; ++ilev) {
234 LocalCopy(lhs[ilev], rhs[ilev], 0, 0, nComp(lhs[ilev]), IntVect(0));
235 }
236}
237
238template <typename MF>
239void GMRESMLMGT<MF>::increment (VEC& lhs, VEC const& rhs, RT a)
240{
241 auto nlevels = int(lhs.size());
242 for (int ilev = 0; ilev < nlevels; ++ilev) {
243 Saxpy(lhs[ilev], a, rhs[ilev], 0, 0, nComp(lhs[ilev]), IntVect(0));
244 }
245}
246
247template <typename MF>
248void GMRESMLMGT<MF>::linComb (VEC& lhs, RT a, VEC const& rhs_a, RT b, VEC const& rhs_b)
249{
250 auto nlevels = int(lhs.size());
251 for (int ilev = 0; ilev < nlevels; ++ilev) {
252 LinComb(lhs[ilev], a, rhs_a[ilev], 0, b, rhs_b[ilev], 0, 0, nComp(lhs[ilev]), IntVect(0));
253 }
254}
255
256template <typename MF>
257void GMRESMLMGT<MF>::apply (VEC& lhs, VEC const& rhs) const
258{
259 m_mlmg->applyPrecond(GetVecOfPtrs(lhs), GetVecOfPtrs(const_cast<VEC&>(rhs)));
260}
261
262template <typename MF>
263void GMRESMLMGT<MF>::precond (VEC& lhs, VEC const& rhs) const
264{
265 if (m_use_precond) {
266 m_mlmg->setPrecondIter(m_precond_niters);
267 setToZero(lhs);
268 m_mlmg->precond(GetVecOfPtrs(lhs), GetVecOfConstPtrs(rhs), 0, 0);
269 } else {
270 for (int ilev = 0; ilev < m_nlevels; ++ilev) {
271 LocalCopy(lhs[ilev], rhs[ilev], 0, 0, nComp(lhs[ilev]), IntVect(0));
272 }
273 }
274}
275
276template <typename MF>
277void GMRESMLMGT<MF>::solve (MF& a_sol, MF const& a_rhs, RT a_tol_rel, RT a_tol_abs)
278{
279 AMREX_ALWAYS_ASSERT(m_nlevels == 1);
280 this->solve({&a_sol}, {&a_rhs}, a_tol_rel, a_tol_abs);
281}
282
283template <typename MF>
284void GMRESMLMGT<MF>::solve (Vector<MF*> const& a_sol, Vector<MF const*> const& a_rhs, RT a_tol_rel, RT a_tol_abs)
285{
286 m_mlmg->incPrintIdentation();
287 auto mlmg_verbose = m_mlmg->getVerbose();
288 auto mlmg_bottom_verbose = m_mlmg->getBottomVerbose();
289 m_mlmg->setVerbose(m_verbose);
290 auto mlmg_bottom_solver = m_mlmg->getBottomSolver();
291
292 if (mlmg_bottom_solver != BottomSolver::smoother &&
293 mlmg_bottom_solver != BottomSolver::hypre &&
294 mlmg_bottom_solver != BottomSolver::petsc)
295 {
296 m_mlmg->setBottomSolver(BottomSolver::smoother);
297 }
298
299 auto res = makeVecLHS();
300 auto cor = makeVecLHS();
301
302 m_mlmg->apply(GetVecOfPtrs(res), a_sol); // res = L(sol)
303 // res = L(sol) - rhs
304 bool need_to_scale_rhs = m_linop->scaleRHS(0,nullptr);
305 for (int ilev = 0; ilev < m_nlevels; ++ilev) {
306 MF const* prhs;
307 if (need_to_scale_rhs) {
308 LocalCopy(cor[ilev], *a_rhs[ilev], 0, 0, nComp(cor[ilev]), IntVect(0));
309 auto r = m_linop->scaleRHS(ilev, &(cor[ilev]));
311 prhs = &(cor[ilev]);
312 } else {
313 prhs = a_rhs[ilev];
314 }
315 Saxpy(res[ilev], RT(-1), *prhs, 0, 0, nComp(res[ilev]), IntVect(0));
316 }
317 for (int ilev = 0; ilev < m_nlevels; ++ilev) {
318 m_linop->setDirichletNodesToZero(ilev,0,res[ilev]);
319 }
320 m_gmres.solve(cor, res, a_tol_rel, a_tol_abs); // L(cor) = res
321 // sol = sol - cor
322 for (int ilev = 0; ilev < m_nlevels; ++ilev) {
323 Saxpy(*a_sol[ilev], RT(-1), cor[ilev], 0, 0, nComp(*a_sol[ilev]), IntVect(0));
324 }
325
326 m_mlmg->setBottomSolver(mlmg_bottom_solver);
327 m_mlmg->setVerbose(mlmg_verbose);
328 m_mlmg->setBottomVerbose(mlmg_bottom_verbose);
329 m_mlmg->decPrintIdentation();
330}
331
334
335}
336
337#endif
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
Solve using GMRES with multigrid as preconditioner.
Definition AMReX_GMRES_MLMG.H:28
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:277
static void assign(VEC &lhs, VEC const &rhs)
Copy rhs into lhs level by level.
Definition AMReX_GMRES_MLMG.H:230
VEC makeVecRHS() const
Return a MultiFab vector without ghost cells for RHS storage.
Definition AMReX_GMRES_MLMG.H:181
GMRESMLMGT(MG &mlmg)
Wrap an existing MLMG hierarchy as a GMRES operator/preconditioner.
Definition AMReX_GMRES_MLMG.H:173
int getNumIters() const
Number of iterations executed by the last solve().
Definition AMReX_GMRES_MLMG.H:72
void setPrecondNumIters(int precond_niters)
Adjust how many MLMG smoothing iterations run inside each GMRES iteration.
Definition AMReX_GMRES_MLMG.H:160
Vector< MF > VEC
Definition AMReX_GMRES_MLMG.H:30
static void setToZero(VEC &lhs)
Reset lhs to zero on every AMR level.
Definition AMReX_GMRES_MLMG.H:222
GM & getGMRES()
Direct access to the underlying GMRES driver for fine-grained tuning.
Definition AMReX_GMRES_MLMG.H:78
void setVerbose(int v)
Set verbosity level v for both GMRES and MLMG components.
Definition AMReX_GMRES_MLMG.H:63
typename MG::RT RT
Definition AMReX_GMRES_MLMG.H:32
RT dotProduct(VEC const &mf1, VEC const &mf2) const
Return the preconditioner-aware dot product between mf1 and mf2.
Definition AMReX_GMRES_MLMG.H:216
GMRES< Vector< MF >, GMRESMLMGT< MF > > GM
Definition AMReX_GMRES_MLMG.H:33
VEC makeVecLHS() const
Definition AMReX_GMRES_MLMG.H:191
RT norm2(VEC const &mf) const
Return the 2-norm of mf using the multigrid operator's weighting.
Definition AMReX_GMRES_MLMG.H:202
static void scale(VEC &mf, RT scale_factor)
Scale each component of mf by scale_factor.
Definition AMReX_GMRES_MLMG.H:208
MLMGT< MF > MG
Definition AMReX_GMRES_MLMG.H:31
static void linComb(VEC &lhs, RT a, VEC const &rhs_a, RT b, VEC const &rhs_b)
Form a linear combination: .
Definition AMReX_GMRES_MLMG.H:248
static void increment(VEC &lhs, VEC const &rhs, RT a)
Add a scaled vector: .
Definition AMReX_GMRES_MLMG.H:239
RT getResidualNorm() const
Final residual 2-norm produced by the last solve().
Definition AMReX_GMRES_MLMG.H:75
bool usePrecond(bool new_flag)
Enable or disable MLMG as a preconditioner.
Definition AMReX_GMRES_MLMG.H:153
void apply(VEC &lhs, VEC const &rhs) const
Apply the operator supplied to the constructor.
Definition AMReX_GMRES_MLMG.H:257
void setMaxIters(int niters)
Cap the number of GMRES iterations executed per solve to niters.
Definition AMReX_GMRES_MLMG.H:69
void precond(VEC &lhs, VEC const &rhs) const
Apply the MLMG-based preconditioner: .
Definition AMReX_GMRES_MLMG.H:263
int getNumIters() const
Number of iterations executed by the last solve().
Definition AMReX_GMRES.H:134
void define(M &linop)
Bind the solver to a linear operator.
Definition AMReX_GMRES.H:208
void setVerbose(int v)
Set verbosity level v (0 = silent).
Definition AMReX_GMRES.H:117
void setMaxIters(int niters)
Cap the number of iterations performed by solve().
Definition AMReX_GMRES.H:131
RT getResidualNorm() const
Final residual 2-norm from the last solve().
Definition AMReX_GMRES.H:140
Definition AMReX_MLLinOp.H:102
Definition AMReX_MLMG.H:17
void preparePrecond()
Definition AMReX_MLMG.H:1229
typename MLLinOpT< MF >::RT RT
Definition AMReX_MLMG.H:32
Long size() const noexcept
Definition AMReX_Vector.H:53
Definition AMReX_Amr.cpp:49
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
int nComp(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2851
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:1966
void LocalCopy(DMF &dst, SMF const &src, int scomp, int dcomp, int ncomp, IntVect const &nghost)
dst = src
Definition AMReX_FabArrayUtility.H:1950
Vector< const T * > GetVecOfConstPtrs(const Vector< T > &a)
Definition AMReX_Vector.H:93
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:2009
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
void setBndry(MF &dst, typename MF::value_type val, int scomp, int ncomp)
dst = val in ghost cells.
Definition AMReX_FabArrayUtility.H:1934
Vector< T * > GetVecOfPtrs(Vector< T > &a)
Definition AMReX_Vector.H:64
void Scale(MF &dst, typename MF::value_type val, int scomp, int ncomp, int nghost)
dst *= val
Definition AMReX_FabArrayUtility.H:1941
const int[]
Definition AMReX_BLProfiler.cpp:1664
void setVal(MF &dst, typename MF::value_type val)
dst = val
Definition AMReX_FabArrayUtility.H:1927