1#ifndef AMREX_GMRES_MLMG_H_
2#define AMREX_GMRES_MLMG_H_
3#include <AMReX_Config.H>
38 void solve (MF& a_sol, MF
const& a_rhs,
RT a_tol_rel,
RT a_tol_abs);
105template <
typename MF>
107 : m_mlmg(&mlmg), m_linop(&mlmg.getLinOp()), m_nlevels(m_linop->NAMRLevels())
113template <
typename MF>
117 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
118 vmf[ilev] = m_linop->make(ilev, 0,
IntVect(0));
123template <
typename MF>
127 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
128 vmf[ilev] = m_linop->make(ilev, 0,
IntVect(1));
134template <
typename MF>
140template <
typename MF>
143 for (
auto& xmf : mf) {
148template <
typename MF>
154template <
typename MF>
157 for (
auto& xmf : lhs) {
162template <
typename MF>
165 auto nlevels =
int(lhs.
size());
166 for (
int ilev = 0; ilev < nlevels; ++ilev) {
171template <
typename MF>
174 auto nlevels =
int(lhs.
size());
175 for (
int ilev = 0; ilev < nlevels; ++ilev) {
180template <
typename MF>
183 auto nlevels =
int(lhs.
size());
184 for (
int ilev = 0; ilev < nlevels; ++ilev) {
189template <
typename MF>
195template <
typename MF>
199 m_mlmg->setPrecondIter(m_precond_niters);
203 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
209template <
typename MF>
213 this->solve({&a_sol}, {&a_rhs}, a_tol_rel, a_tol_abs);
216template <
typename MF>
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();
232 auto res = makeVecLHS();
233 auto cor = makeVecLHS();
237 bool need_to_scale_rhs = m_linop->scaleRHS(0,
nullptr);
238 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
240 if (need_to_scale_rhs) {
242 auto r = m_linop->scaleRHS(ilev, &(cor[ilev]));
250 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
251 m_linop->setDirichletNodesToZero(ilev,0,res[ilev]);
253 m_gmres.solve(cor, res, a_tol_rel, a_tol_abs);
255 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
259 m_mlmg->setBottomSolver(mlmg_bottom_solver);
260 m_mlmg->setVerbose(mlmg_verbose);
261 m_mlmg->setBottomVerbose(mlmg_bottom_verbose);
262 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: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