1 #ifndef SOLVER_NONLOCAL_LINEAR
2 #define SOLVER_NONLOCAL_LINEAR
3 #include "Operator/Operator.H"
4 #include <AMReX_MLMG.H>
41 this->
mlmg =
new amrex::MLMG(a_lp);
48 this->
linop =
nullptr;
56 Real a_tol_rel, Real a_tol_abs,
bool copyrhs =
false,
57 const char* checkpoint_file =
nullptr)
60 amrex::Vector<amrex::MultiFab*> rhs_tmp(a_rhs.size());
61 amrex::Vector<amrex::MultiFab*> zero_tmp(a_rhs.size());
62 for (
int i = 0; i < rhs_tmp.size(); i++)
64 rhs_tmp[i] =
new amrex::MultiFab(a_rhs[i]->boxArray(), a_rhs[i]->DistributionMap(), a_rhs[i]->nComp(), a_rhs[i]->nGrow());
65 zero_tmp[i] =
new amrex::MultiFab(a_rhs[i]->boxArray(), a_rhs[i]->DistributionMap(), a_rhs[i]->nComp(), a_rhs[i]->nGrow());
66 rhs_tmp[i]->setVal(0.0);
67 zero_tmp[i]->setVal(0.0);
72 mlmg->apply(rhs_tmp, zero_tmp);
74 for (
int lev = 0; lev < rhs_tmp.size(); lev++)
76 amrex::Box domain =
linop->
Geom(lev).Domain();
77 domain.convert(amrex::IntVect::TheNodeVector());
78 const Dim3 lo = amrex::lbound(domain), hi = amrex::ubound(domain);
79 for (MFIter mfi(*rhs_tmp[lev], amrex::TilingIfNotGPU());mfi.isValid();++mfi)
81 amrex::Box bx = mfi.growntilebox(rhs_tmp[lev]->nGrow());
83 amrex::Array4<amrex::Real>
const& rhstmp = rhs_tmp[lev]->array(mfi);
84 for (
int n = 0; n < rhs_tmp[lev]->nComp(); n++)
86 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
88 bool AMREX_D_DECL(xmin = (i == lo.x), ymin = (j == lo.y), zmin = (k == lo.z)),
89 AMREX_D_DECL(xmax = (i == hi.x), ymax = (j == hi.y), zmax = (k == hi.z));
90 if (AMREX_D_TERM(xmax || xmin, || ymax || ymin, || zmax || zmin))
91 rhstmp(i, j, k, n) = 0.0;
93 rhstmp(i, j, k, n) *= -1.0;
103 for (
int lev = 0; lev < rhs_tmp.size(); lev++)
106 amrex::Add(*rhs_tmp[lev], *a_rhs[lev], 0, 0, rhs_tmp[lev]->nComp(), rhs_tmp[lev]->nGrow());
108 amrex::Copy(*a_rhs[lev], *rhs_tmp[lev], 0, 0, rhs_tmp[lev]->nComp(), rhs_tmp[lev]->nGrow());
117 retval =
mlmg->solve(GetVecOfPtrs(a_sol), GetVecOfConstPtrs(rhs_tmp), a_tol_rel, a_tol_abs, checkpoint_file);
119 catch (
const std::exception& e)
124 if (a_sol[0]->contains_nan())
135 Real a_tol_rel, Real a_tol_abs,
const char* checkpoint_file =
nullptr)
141 retval =
mlmg->solve(GetVecOfPtrs(a_sol), GetVecOfConstPtrs(a_rhs), a_tol_rel, a_tol_abs, checkpoint_file);
143 catch (
const std::exception& e)
157 retval =
mlmg->solve(GetVecOfPtrs(a_sol), GetVecOfConstPtrs(a_rhs),
tol_rel,
tol_abs);
159 catch (
const std::exception& e)
164 if (a_sol[0]->contains_nan())
175 mlmg->apply(GetVecOfPtrs(a_rhs), GetVecOfPtrs(a_sol));
189 const amrex::Vector<amrex::MultiFab const*>& a_rhs_mf)
191 int nlevs = a_sol_mf.size();
192 int ncomps = a_sol_mf[0]->nComp();
194 amrex::Vector<amrex::Geometry> geom;
195 amrex::Vector<int> iter;
196 amrex::Vector<amrex::IntVect> refratio;
197 amrex::Vector<std::string> names;
198 for (
int i = 0; i < nlevs; i++)
202 if (i > 0) refratio.push_back(amrex::IntVect(2));
204 for (
int n = 0; n < ncomps; n++)
206 names.push_back(
"var" + std::to_string(n));
210 WriteMultiLevelPlotfile(outputdir +
"/mlmg_sol", nlevs,
211 amrex::GetVecOfConstPtrs(a_sol_mf),
212 names, geom, 0, iter, refratio);
213 WriteMultiLevelPlotfile(outputdir +
"/mlmg_rhs", nlevs,
214 amrex::GetVecOfConstPtrs(a_rhs_mf),
215 names, geom, 0, iter, refratio);
218 for (
int lev = 0; lev < nlevs; lev++)
220 res_mf.
Define(lev, a_sol_mf[lev]->boxArray(), a_sol_mf[lev]->DistributionMap(),
221 ncomps, a_sol_mf[lev]->nGrow());
224 mlmg->compResidual(amrex::GetVecOfPtrs(res_mf), a_sol_mf, a_rhs_mf);
226 WriteMultiLevelPlotfile(outputdir +
"/mlmg_res", nlevs,
227 amrex::GetVecOfConstPtrs(res_mf),
228 names, geom, 0, iter, refratio);
260 mlmg.setBottomSolver(MLMG::BottomSolver::bicgstab);
261 mlmg.setCFStrategy(MLMG::CFStrategy::ghostnodes);
262 mlmg.setFinalFillBC(
false);
263 mlmg.setMaxFmgIter(100000000);
274 else mlmg.setBottomVerbose(0);
283 else if (
bottom_solver ==
"bicgstab")
mlmg.setBottomSolver(MLMG::BottomSolver::bicgstab);
284 else if (
bottom_solver ==
"smoother")
mlmg.setBottomSolver(MLMG::BottomSolver::smoother);
310 if (pp.contains(
"max_fixed_iter"))
311 Util::Abort(
INFO,
"max_fixed_iter is depricated. Use fixed_iter instead.");
334 if (pp.contains(
"cg_tol_rel"))
335 Util::Abort(
INFO,
"cg_tol_rel is depricated. Use bottom_tol_rel instead.");
336 if (pp.contains(
"cg_tol_abs"))
337 Util::Abort(
INFO,
"cg_tol_abs is depricated. Use bottom_tol_abs instead.");
380 pp.add(
"signal_handling", 0);
381 pp.add(
"throw_exception", 1);