1#ifndef AMREX_GMRES_MLMG_H_
2#define AMREX_GMRES_MLMG_H_
3#include <AMReX_Config.H>
50 void solve (MF& a_sol, MF
const& a_rhs,
RT a_tol_rel,
RT a_tol_abs);
153 bool usePrecond (
bool new_flag) {
return std::exchange(m_use_precond, new_flag); }
168 bool m_use_precond =
true;
169 int m_precond_niters = 1;
172template <
typename MF>
174 : m_mlmg(&mlmg), m_linop(&mlmg.getLinOp()), m_nlevels(m_linop->NAMRLevels())
180template <
typename MF>
184 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
185 vmf[ilev] = m_linop->make(ilev, 0,
IntVect(0));
190template <
typename MF>
194 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
195 vmf[ilev] = m_linop->make(ilev, 0,
IntVect(1));
201template <
typename MF>
207template <
typename MF>
210 for (
auto& xmf : mf) {
215template <
typename MF>
221template <
typename MF>
224 for (
auto& xmf : lhs) {
229template <
typename MF>
232 auto nlevels =
int(lhs.
size());
233 for (
int ilev = 0; ilev < nlevels; ++ilev) {
238template <
typename MF>
241 auto nlevels =
int(lhs.
size());
242 for (
int ilev = 0; ilev < nlevels; ++ilev) {
247template <
typename MF>
250 auto nlevels =
int(lhs.
size());
251 for (
int ilev = 0; ilev < nlevels; ++ilev) {
256template <
typename MF>
262template <
typename MF>
266 m_mlmg->setPrecondIter(m_precond_niters);
270 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
276template <
typename MF>
280 this->solve({&a_sol}, {&a_rhs}, a_tol_rel, a_tol_abs);
283template <
typename MF>
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();
299 auto res = makeVecLHS();
300 auto cor = makeVecLHS();
304 bool need_to_scale_rhs = m_linop->scaleRHS(0,
nullptr);
305 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
307 if (need_to_scale_rhs) {
309 auto r = m_linop->scaleRHS(ilev, &(cor[ilev]));
317 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
318 m_linop->setDirichletNodesToZero(ilev,0,res[ilev]);
320 m_gmres.solve(cor, res, a_tol_rel, a_tol_abs);
322 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
326 m_mlmg->setBottomSolver(mlmg_bottom_solver);
327 m_mlmg->setVerbose(mlmg_verbose);
328 m_mlmg->setBottomVerbose(mlmg_bottom_verbose);
329 m_mlmg->decPrintIdentation();
#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