1#ifndef AMREX_GMRES_MLMG_H_
2#define AMREX_GMRES_MLMG_H_
3#include <AMReX_Config.H>
39 void solve (MF& a_sol, MF
const& a_rhs,
RT a_tol_rel,
RT a_tol_abs);
106template <
typename MF>
108 : m_mlmg(&mlmg), m_linop(&mlmg.getLinOp()), m_nlevels(m_linop->NAMRLevels())
114template <
typename MF>
118 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
119 vmf[ilev] = m_linop->make(ilev, 0,
IntVect(0));
124template <
typename MF>
128 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
129 vmf[ilev] = m_linop->make(ilev, 0,
IntVect(1));
135template <
typename MF>
141template <
typename MF>
144 for (
auto& xmf : mf) {
149template <
typename MF>
155template <
typename MF>
158 for (
auto& xmf : lhs) {
163template <
typename MF>
166 auto nlevels = int(lhs.
size());
167 for (
int ilev = 0; ilev < nlevels; ++ilev) {
172template <
typename MF>
175 auto nlevels = int(lhs.
size());
176 for (
int ilev = 0; ilev < nlevels; ++ilev) {
181template <
typename MF>
184 auto nlevels = int(lhs.
size());
185 for (
int ilev = 0; ilev < nlevels; ++ilev) {
190template <
typename MF>
196template <
typename MF>
200 m_mlmg->setPrecondIter(m_precond_niters);
204 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
210template <
typename MF>
214 this->solve({&a_sol}, {&a_rhs}, a_tol_rel, a_tol_abs);
217template <
typename MF>
220 m_mlmg->incPrintIdentation();
221 auto mlmg_verbose = m_mlmg->getVerbose();
222 auto mlmg_bottom_verbose = m_mlmg->getBottomVerbose();
223 m_mlmg->setVerbose(m_verbose);
224 auto mlmg_bottom_solver = m_mlmg->getBottomSolver();
233 auto res = makeVecLHS();
234 auto cor = makeVecLHS();
238 bool need_to_scale_rhs = m_linop->scaleRHS(0,
nullptr);
239 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
241 if (need_to_scale_rhs) {
243 auto r = m_linop->scaleRHS(ilev, &(cor[ilev]));
251 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
252 m_linop->setDirichletNodesToZero(ilev,0,res[ilev]);
254 m_gmres.solve(cor, res, a_tol_rel, a_tol_abs);
256 for (
int ilev = 0; ilev < m_nlevels; ++ilev) {
260 m_mlmg->setBottomSolver(mlmg_bottom_solver);
261 m_mlmg->setVerbose(mlmg_verbose);
262 m_mlmg->setBottomVerbose(mlmg_bottom_verbose);
263 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:22
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:211
static void assign(VEC &lhs, VEC const &rhs)
lhs = rhs
Definition AMReX_GMRES_MLMG.H:164
VEC makeVecRHS() const
Make MultiFab without ghost cells.
Definition AMReX_GMRES_MLMG.H:115
int m_precond_niters
Definition AMReX_GMRES_MLMG.H:103
GMRESMLMGT(MG &mlmg)
Definition AMReX_GMRES_MLMG.H:107
GM m_gmres
Definition AMReX_GMRES_MLMG.H:97
int getNumIters() const
Gets the number of iterations.
Definition AMReX_GMRES_MLMG.H:53
void setPrecondNumIters(int precond_niters)
Set the number of MLMG preconditioner iterations per GMRES iteration.
Definition AMReX_GMRES_MLMG.H:94
Vector< MF > VEC
Definition AMReX_GMRES_MLMG.H:24
static void setToZero(VEC &lhs)
lhs = 0
Definition AMReX_GMRES_MLMG.H:156
GM & getGMRES()
Get the GMRES object.
Definition AMReX_GMRES_MLMG.H:59
void setVerbose(int v)
Sets verbosity.
Definition AMReX_GMRES_MLMG.H:44
typename MG::RT RT
Definition AMReX_GMRES_MLMG.H:26
RT dotProduct(VEC const &mf1, VEC const &mf2) const
Definition AMReX_GMRES_MLMG.H:150
VEC makeVecLHS() const
Make MultiFab with ghost cells and set ghost cells to zero.
Definition AMReX_GMRES_MLMG.H:125
RT norm2(VEC const &mf) const
Definition AMReX_GMRES_MLMG.H:136
static void scale(VEC &mf, RT scale_factor)
Definition AMReX_GMRES_MLMG.H:142
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:182
MLLinOpT< MF > * m_linop
Definition AMReX_GMRES_MLMG.H:99
static void increment(VEC &lhs, VEC const &rhs, RT a)
lhs += a*rhs
Definition AMReX_GMRES_MLMG.H:173
int m_verbose
Definition AMReX_GMRES_MLMG.H:100
bool m_use_precond
Definition AMReX_GMRES_MLMG.H:102
MG * m_mlmg
Definition AMReX_GMRES_MLMG.H:98
RT getResidualNorm() const
Gets the 2-norm of the residual.
Definition AMReX_GMRES_MLMG.H:56
int m_nlevels
Definition AMReX_GMRES_MLMG.H:101
bool usePrecond(bool new_flag)
Control whether or not to use MLMG as preconditioner.
Definition AMReX_GMRES_MLMG.H:91
void apply(VEC &lhs, VEC const &rhs) const
lhs = L(rhs)
Definition AMReX_GMRES_MLMG.H:191
void setMaxIters(int niters)
Sets the max number of iterations.
Definition AMReX_GMRES_MLMG.H:50
void precond(VEC &lhs, VEC const &rhs) const
Definition AMReX_GMRES_MLMG.H:197
int getNumIters() const
Gets the number of iterations.
Definition AMReX_GMRES.H:114
void define(M &linop)
Definition AMReX_GMRES.H:188
void setVerbose(int v)
Sets verbosity.
Definition AMReX_GMRES.H:105
void setMaxIters(int niters)
Sets the max number of iterations.
Definition AMReX_GMRES.H:111
RT getResidualNorm() const
Gets the 2-norm of the residual.
Definition AMReX_GMRES.H:120
Definition AMReX_MLLinOp.H:102
Definition AMReX_MLMG.H:17
void preparePrecond()
Definition AMReX_MLMG.H:1231
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:138
int nComp(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2854
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:30
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
void setVal(MF &dst, typename MF::value_type val)
dst = val
Definition AMReX_FabArrayUtility.H:1927