1#ifndef AMREX_ML_LINOP_H_
2#define AMREX_ML_LINOP_H_
3#include <AMReX_Config.H>
5#if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
10#if defined(AMREX_USE_PETSC) && (AMREX_SPACEDIM > 1)
90template <
typename T>
class MLMGT;
91template <
typename T>
class MLCGSolverT;
92template <
typename T>
class MLPoissonT;
93template <
typename T>
class MLABecLaplacianT;
94template <
typename T>
class GMRESMLMGT;
101 template <
typename T>
friend class MLMGT;
129 bool eb_limit_coarsening =
true);
131 [[nodiscard]]
virtual std::string
name ()
const {
return std::string(
"Unspecified"); }
201 LinOpBCType bc_type = LinOpBCType::Dirichlet)
noexcept;
204 LinOpBCType bc_type = LinOpBCType::Dirichlet)
noexcept;
206 template <
typename AMF, std::enable_if_t<!std::is_same_v<MF,AMF>,
int> = 0>
208 LinOpBCType bc_type = LinOpBCType::Dirichlet)
noexcept;
210 template <
typename AMF, std::enable_if_t<!std::is_same_v<MF,AMF>,
int> = 0>
212 LinOpBCType bc_type = LinOpBCType::Dirichlet)
noexcept;
235 const MF* =
nullptr) = 0;
237 template <
typename AMF, std::enable_if_t<!std::is_same_v<MF,AMF>,
int> = 0>
239 const AMF* robinbc_a =
nullptr,
240 const AMF* robinbc_b =
nullptr,
241 const AMF* robinbc_f =
nullptr);
261 [[nodiscard]]
virtual int getNComp ()
const {
return 1; }
263 [[nodiscard]]
virtual int getNGrow (
int = 0,
int = 0)
const {
return 0; }
302 amrex::Abort(
"MLLinOpT::interpAssign: Must be implemented for FMG cycle");
317 amrex::Abort(
"MLLinOpT::interpolationAmr: Must be implemented for composite solves across multiple AMR levels");
330 const MF& fine_sol,
const MF& fine_rhs)
333 amrex::Abort(
"MLLinOpT::averageDownSolutionRHS: Must be implemented for composite solves across multiple AMR levels");
347 virtual void apply (
int amrlev,
int mglev, MF& out, MF& in,
BCMode bc_mode,
360 virtual void smooth (
int amrlev,
int mglev, MF& sol,
const MF& rhs,
361 bool skip_fillboundary,
int niter)
const = 0;
364 virtual void normalize (
int amrlev,
int mglev, MF& mf)
const {
378 const MF* crse_bcdata=
nullptr) = 0;
394 BCMode bc_mode,
const MF* crse_bcdata=
nullptr) = 0;
408 MF& res,
const MF& crse_sol,
const MF& crse_rhs,
409 MF& fine_res, MF& fine_sol,
const MF& fine_rhs)
const
413 amrex::Abort(
"MLLinOpT::reflux: Must be implemented for composite solves across multiple AMR levels");
427 amrex::Abort(
"AMReX_MLLinOp::compFlux::How did we get here?");
441 amrex::Abort(
"AMReX_MLLinOp::compGrad::How did we get here?");
459 [[nodiscard]]
virtual bool scaleRHS (
int , MF* )
const {
465 MF
const& )
const {
return {}; }
480 amrex::Warning(
"This function might need to be implemented for GMRES to work with this LinOp.");
484 [[nodiscard]]
virtual bool isSingular (
int amrlev)
const = 0;
489 virtual RT xdoty (
int amrlev,
int mglev,
const MF&
x,
const MF& y,
bool local)
const = 0;
495 virtual std::unique_ptr<MLLinOpT<MF>>
makeNLinOp (
int )
const
497 amrex::Abort(
"MLLinOp::makeNLinOp: N-Solve not supported");
504 amrex::Abort(
"MLLinOp::getFluxes: How did we get here?");
508 amrex::Abort(
"MLLinOp::getFluxes: How did we get here?");
514 amrex::Abort(
"MLLinOp::getEBFluxes: How did we get here?");
518#if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
520 amrex::Abort(
"MLLinOp::makeHypre: How did we get here?");
523 [[nodiscard]]
virtual std::unique_ptr<HypreNodeLap> makeHypreNodeLap(
525 const std::string& )
const
527 amrex::Abort(
"MLLinOp::makeHypreNodeLap: How did we get here?");
532#if defined(AMREX_USE_PETSC) && (AMREX_SPACEDIM > 1)
533 [[nodiscard]]
virtual std::unique_ptr<PETScABecLap> makePETSc ()
const {
534 amrex::Abort(
"MLLinOp::makePETSc: How did we get here?");
545 [[nodiscard]]
virtual RT normInf (
int amrlev, MF
const& mf,
bool local)
const = 0;
552 amrex::Abort(
"MLLinOpT::avgDownResAmr: Must be implemented for composite solves across multiple AMR levels");
561 [[nodiscard]]
bool isMFIterSafe (
int amrlev,
int mglev1,
int mglev2)
const;
569 [[nodiscard]]
const Geometry&
Geom (
int amr_lev,
int mglev=0) const noexcept {
return m_geom[amr_lev][mglev]; }
693 [[nodiscard]]
virtual MF
make (
int amrlev,
int mglev,
IntVect const& ng)
const;
695 [[nodiscard]]
virtual MF
makeAlias (MF
const& mf)
const;
701 [[nodiscard]]
virtual std::unique_ptr<FabFactory<FAB> >
makeFactory (
int ,
int )
const {
702 return std::make_unique<DefaultFabFactory<FAB>>();
711 template <
typename T>
715 return Array4<T>(a.dataPtr(), {a.begin.y,a.begin.z,0}, {a.end.y,a.end.z,1}, a.nComp());
717 return Array4<T>(a.dataPtr(), {a.begin.x,a.begin.z,0}, {a.end.x,a.end.z,1}, a.nComp());
719 return Array4<T>(a.dataPtr(), {a.begin.x,a.begin.y,0}, {a.end.x,a.end.y,1}, a.nComp());
725 template <
typename T>
726 [[nodiscard]] T
get_d0 (T
const& dx, T
const& dy, T
const&)
const noexcept
735 template <
typename T>
736 [[nodiscard]] T
get_d1 (T
const&, T
const& dy, T
const& dz)
const noexcept
754 int ratio,
int strategy);
767template <
typename MF>
774 [[maybe_unused]]
bool eb_limit_coarsening)
783 if (info.con_grid_size <= 0) { info.con_grid_size =
AMREX_D_PICK(32, 16, 8); }
793 if (!a_factory.empty() && eb_limit_coarsening) {
796 info.max_coarsening_level = std::min(info.max_coarsening_level,
797 f->maxCoarseningLevel());
801 defineGrids(a_geom, a_grids, a_dmap, a_factory);
805template <
typename MF>
815 if ( ! a_factory.empty() ) {
817 if (ebf && !(ebf->isAllRegular())) {
818 mg_domain_min_width = 4;
823 m_num_amr_levels = 0;
824 for (
int amrlev = 0; amrlev < a_geom.
size(); amrlev++) {
825 if (!a_grids[amrlev].empty()) {
830 m_amr_ref_ratio.resize(m_num_amr_levels);
831 m_num_mg_levels.resize(m_num_amr_levels);
833 m_geom.resize(m_num_amr_levels);
834 m_grids.resize(m_num_amr_levels);
835 m_dmap.resize(m_num_amr_levels);
836 m_factory.resize(m_num_amr_levels);
840 const RealBox& rb = a_geom[0].ProbDomain();
841 const int coord = a_geom[0].Coord();
844 IntVect mg_coarsen_ratio_v(mg_coarsen_ratio);
845 IntVect mg_box_min_width_v(mg_box_min_width);
846 IntVect mg_domain_min_width_v(mg_domain_min_width);
847 if (hasHiddenDimension()) {
849 "Hidden direction only supported for 3d level solve");
850 mg_coarsen_ratio_v[info.hidden_direction] = 1;
851 mg_box_min_width_v[info.hidden_direction] = 0;
852 mg_domain_min_width_v[info.hidden_direction] = 0;
856 for (
int amrlev = m_num_amr_levels-1; amrlev > 0; --amrlev)
858 m_num_mg_levels[amrlev] = 1;
859 m_geom[amrlev].push_back(a_geom[amrlev]);
860 m_grids[amrlev].push_back(a_grids[amrlev]);
861 m_dmap[amrlev].push_back(a_dmap[amrlev]);
862 if (amrlev < a_factory.size()) {
863 m_factory[amrlev].emplace_back(a_factory[amrlev]->clone());
868 IntVect rr = mg_coarsen_ratio_v;
869 const Box& dom = a_geom[amrlev].Domain();
870 for (
int i = 0; i < 2; ++i)
875 if (cdom == a_geom[amrlev-1].Domain()) {
break; }
877 ++(m_num_mg_levels[amrlev]);
879 m_geom[amrlev].emplace_back(cdom, rb, coord, is_per);
881 m_grids[amrlev].push_back(a_grids[amrlev]);
883 m_grids[amrlev].back().coarsen(rr);
885 m_dmap[amrlev].push_back(a_dmap[amrlev]);
887 rr *= mg_coarsen_ratio_v;
890#if (AMREX_SPACEDIM > 1)
891 if (hasHiddenDimension()) {
892 m_amr_ref_ratio[amrlev-1] = rr[AMREX_SPACEDIM-info.hidden_direction];
896 m_amr_ref_ratio[amrlev-1] = rr[0];
901 m_num_mg_levels[0] = 1;
902 m_geom[0].push_back(a_geom[0]);
903 m_grids[0].push_back(a_grids[0]);
904 m_dmap[0].push_back(a_dmap[0]);
905 if (!a_factory.empty()) {
906 m_factory[0].emplace_back(a_factory[0]->clone());
911 m_domain_covered.
resize(m_num_amr_levels,
false);
912 auto npts0 = m_grids[0][0].numPts();
913 m_domain_covered[0] = (npts0 == compactify(m_geom[0][0].Domain()).numPts());
914 for (
int amrlev = 1; amrlev < m_num_amr_levels; ++amrlev)
916 if (!m_domain_covered[amrlev-1]) {
break; }
917 m_domain_covered[amrlev] = (m_grids[amrlev][0].numPts() ==
918 compactify(m_geom[amrlev][0].Domain()).numPts());
922 bool aggable =
false;
924 if (m_grids[0][0].size() > 1 && info.do_agglomeration)
926 if (m_domain_covered[0])
928 aggbox = m_geom[0][0].Domain();
929 if (hasHiddenDimension()) {
930 aggbox.
makeSlab(hiddenDirection(), m_grids[0][0][0].smallEnd(hiddenDirection()));
936 aggbox = m_grids[0][0].minimalBox();
937 aggable = (aggbox.
numPts() == npts0);
943 int agg_lev = 0, con_lev = 0;
946 && info.semicoarsening_direction >= -1
947 && info.semicoarsening_direction < AMREX_SPACEDIM );
949 if (info.do_agglomeration && aggable)
951 Box dbx = m_geom[0][0].Domain();
953 Real
const nbxs =
static_cast<Real
>(m_grids[0][0].size());
954 Real
const threshold_npts =
static_cast<Real
>(
AMREX_D_TERM(info.agg_grid_size,
956 *info.agg_grid_size));
963 for (
int lev = 0; lev < info.max_coarsening_level; ++lev)
965 IntVect rr_level = mg_coarsen_ratio_v;
966 bool const do_semicoarsening_level = info.do_semicoarsening
967 && numsclevs < info.max_semicoarsening_level;
968 if (do_semicoarsening_level
969 && info.semicoarsening_direction != -1)
971 rr_level[info.semicoarsening_direction] = 1;
974 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
976 rr_dir[idim] = rr_level[idim];
977 is_coarsenable[idim] = dbx.
coarsenable(rr_dir, mg_domain_min_width_v)
979 if (!is_coarsenable[idim] && do_semicoarsening_level
980 && info.semicoarsening_direction == -1)
982 is_coarsenable[idim] =
true;
989 if (do_semicoarsening_level && info.semicoarsening_direction == -1) {
991 int n_ones =
AMREX_D_TERM(
static_cast<int>(rr_level[0] == 1),
992 +
static_cast<int>(rr_level[1] == 1),
993 +
static_cast<int>(rr_level[2] == 1));
994 if (n_ones > 1) {
break; }
996 if (rr_level != mg_coarsen_ratio_v) {
1000 accum_coarsen_ratio.push_back(accum_coarsen_ratio.back()*rr_level);
1001 domainboxes.push_back(dbx.
coarsen(rr_level));
1002 boundboxes.push_back(bbx.
coarsen(rr_level));
1003 bool to_agg = (bbx.
d_numPts() / nbxs) < 0.999*threshold_npts;
1004 agg_flag.push_back(to_agg);
1007 for (
int lev = 1, nlevs =
static_cast<int>(domainboxes.size()); lev < nlevs; ++lev) {
1008 if (!agged && !agg_flag[lev] &&
1009 a_grids[0].coarsenable(accum_coarsen_ratio[lev], mg_box_min_width_v))
1011 m_grids[0].push_back(
amrex::coarsen(a_grids[0], accum_coarsen_ratio[lev]));
1012 m_dmap[0].push_back(a_dmap[0]);
1014 IntVect cr = domainboxes[lev-1].length() / domainboxes[lev].length();
1015 if (!m_grids[0].back().coarsenable(cr)) {
1018 m_grids[0].emplace_back(boundboxes[lev]);
1019 IntVect max_grid_size(info.agg_grid_size);
1020 if (info.do_semicoarsening && info.max_semicoarsening_level >= lev
1021 && info.semicoarsening_direction != -1)
1024 AMREX_D_TERM(
int mgs_0 = (max_grid_size[0]+blen[0]-1) / blen[0];,
1025 int mgs_1 = (max_grid_size[1]+blen[1]-1) / blen[1];,
1026 int mgs_2 = (max_grid_size[2]+blen[2]-1) / blen[2]);
1027 max_grid_size[info.semicoarsening_direction]
1030 m_grids[0].back().maxSize(max_grid_size);
1037 m_geom[0].emplace_back(domainboxes[lev],rb,coord,is_per);
1042 Long consolidation_threshold = 0;
1043 Real avg_npts = 0.0;
1044 if (info.do_consolidation) {
1046 consolidation_threshold =
AMREX_D_TERM(Long(info.con_grid_size),
1047 *info.con_grid_size,
1048 *info.con_grid_size);
1051 Box const& dom0 = a_geom[0].Domain();
1054 for (
int lev = 0; lev < info.max_coarsening_level; ++lev)
1056 IntVect rr_level = mg_coarsen_ratio_v;
1057 bool do_semicoarsening_level = info.do_semicoarsening
1058 && numsclevs < info.max_semicoarsening_level;
1059 if (do_semicoarsening_level
1060 && info.semicoarsening_direction != -1)
1062 rr_level[info.semicoarsening_direction] = 1;
1065 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
1067 rr_dir[idim] = rr_vec[idim] * rr_level[idim];
1068 is_coarsenable[idim] = dom0.
coarsenable(rr_dir, mg_domain_min_width_v)
1069 && a_grids[0].coarsenable(rr_dir, mg_box_min_width_v);
1070 if (!is_coarsenable[idim] && do_semicoarsening_level
1071 && info.semicoarsening_direction == -1)
1073 is_coarsenable[idim] =
true;
1080 if (do_semicoarsening_level && info.semicoarsening_direction == -1) {
1082 int n_ones =
AMREX_D_TERM(
static_cast<int>(rr_level[0] == 1),
1083 +
static_cast<int>(rr_level[1] == 1),
1084 +
static_cast<int>(rr_level[2] == 1));
1085 if (n_ones > 1) {
break; }
1087 if (rr_level != mg_coarsen_ratio_v) {
1092 m_geom[0].emplace_back(
amrex::coarsen(dom0, rr_vec), rb, coord, is_per);
1095 if (info.do_consolidation)
1097 if (avg_npts/
static_cast<Real
>(
AMREX_D_TERM(rr_vec[0], *rr_vec[1], *rr_vec[2]))
1098 < Real(0.999)*
static_cast<Real
>(consolidation_threshold))
1101 con_lev = m_dmap[0].
size();
1106 m_dmap[0].push_back(m_dmap[0].back());
1111 m_dmap[0].push_back(a_dmap[0]);
1116 m_num_mg_levels[0] = m_grids[0].size();
1118 for (
int mglev = 0; mglev < m_num_mg_levels[0] - 1; mglev++){
1119 const Box& fine_domain = m_geom[0][mglev].Domain();
1120 const Box& crse_domain = m_geom[0][mglev+1].Domain();
1121 mg_coarsen_ratio_vec.push_back(fine_domain.
length()/crse_domain.
length());
1124 for (
int amrlev = 0; amrlev < m_num_amr_levels; ++amrlev) {
1125 if (AMRRefRatio(amrlev) == 4 && mg_coarsen_ratio_vec.empty()) {
1126 mg_coarsen_ratio_vec.push_back(
IntVect(2));
1132 makeAgglomeratedDMap(m_grids[0], m_dmap[0]);
1136 makeConsolidatedDMap(m_grids[0], m_dmap[0], info.con_ratio, info.con_strategy);
1141 m_bottom_comm = makeSubCommunicator(m_dmap[0].back());
1145 m_bottom_comm = m_default_comm;
1148 m_do_agglomeration = agged;
1149 m_do_consolidation = coned;
1153 Print() <<
"MLLinOp::defineGrids(): agglomerated AMR level 0 starting at MG level "
1154 << agg_lev <<
" of " << m_num_mg_levels[0] <<
"\n";
1156 Print() <<
"MLLinOp::defineGrids(): consolidated AMR level 0 starting at MG level "
1157 << con_lev <<
" of " << m_num_mg_levels[0]
1158 <<
" (ratio = " << info.con_ratio <<
")" <<
"\n";
1160 Print() <<
"MLLinOp::defineGrids(): no agglomeration or consolidation of AMR level 0\n";
1164 for (
int amrlev = 0; amrlev < m_num_amr_levels; ++amrlev)
1166 for (
int mglev = 1; mglev < m_num_mg_levels[amrlev]; ++mglev)
1168 m_factory[amrlev].emplace_back(makeFactory(amrlev,mglev));
1172 for (
int amrlev = 1; amrlev < m_num_amr_levels; ++amrlev)
1175 "MLLinOp: grids not coarsenable between AMR levels");
1179template <
typename MF>
1183 m_needs_coarse_data_for_bc = !m_domain_covered[0];
1185 levelbc_raii.resize(m_num_amr_levels);
1186 robin_a_raii.resize(m_num_amr_levels);
1187 robin_b_raii.resize(m_num_amr_levels);
1188 robin_f_raii.resize(m_num_amr_levels);
1191template <
typename MF>
1196 const int ncomp = getNComp();
1201template <
typename MF>
1206 const int ncomp = getNComp();
1208 "MLLinOp::setDomainBC: wrong size");
1211 m_lobc_orig = m_lobc;
1212 m_hibc_orig = m_hibc;
1213 for (
int icomp = 0; icomp < ncomp; ++icomp) {
1214 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
1215 if (m_geom[0][0].isPeriodic(idim)) {
1217 m_hibc[icomp][idim] == BCType::Periodic);
1220 m_hibc[icomp][idim] != BCType::Periodic);
1223 if (m_lobc[icomp][idim] == LinOpBCType::inhomogNeumann ||
1224 m_lobc[icomp][idim] == LinOpBCType::Robin)
1226 m_lobc[icomp][idim] = LinOpBCType::Neumann;
1229 if (m_hibc[icomp][idim] == LinOpBCType::inhomogNeumann ||
1230 m_hibc[icomp][idim] == LinOpBCType::Robin)
1232 m_hibc[icomp][idim] = LinOpBCType::Neumann;
1237 if (hasHiddenDimension()) {
1238 const int hd = hiddenDirection();
1239 for (
int n = 0; n < ncomp; ++n) {
1240 m_lobc[n][hd] = LinOpBCType::Neumann;
1241 m_hibc[n][hd] = LinOpBCType::Neumann;
1245 if (hasInhomogNeumannBC() && !supportInhomogNeumannBC()) {
1246 amrex::Abort(
"Inhomogeneous Neumann BC not supported");
1248 if (hasRobinBC() && !supportRobinBC()) {
1253template <
typename MF>
1257 int ncomp = m_lobc_orig.size();
1258 for (
int n = 0; n < ncomp; ++n) {
1259 for (
int idim = 0; idim <AMREX_SPACEDIM; ++idim) {
1260 if (m_lobc_orig[n][idim] == bct || m_hibc_orig[n][idim] == bct) {
1268template <
typename MF>
1272 return hasBC(BCType::inhomogNeumann);
1275template <
typename MF>
1279 return hasBC(BCType::Robin);
1282template <
typename MF>
1286#if (AMREX_SPACEDIM == 3)
1287 if (info.hasHiddenDimension()) {
1289 const auto& hi =
b.bigEnd();
1290 if (info.hidden_direction == 0) {
1292 }
else if (info.hidden_direction == 1) {
1304template <
typename MF>
1312 for (
int i = 1, N=
static_cast<int>(ba.
size()); i < N; ++i)
1322 for (
int iproc = 0; iproc < nprocs; ++iproc) {
1324 for (
int ibox : sfc[iproc]) {
1328 dm[i].define(std::move(pmap));
1333template <
typename MF>
1337 int ratio,
int strategy)
1339 BL_PROFILE(
"MLLinOp::makeConsolidatedDMap()");
1343 for (
int i = 1, N=
static_cast<int>(ba.
size()); i < N; ++i)
1350 const auto& pmap_fine = dm[i-1].ProcessorMap();
1353 if (strategy == 1) {
1354 for (
auto&
x: pmap) {
1357 }
else if (strategy == 2) {
1358 int nprocs_con =
static_cast<int>(std::ceil(
static_cast<Real
>(nprocs)
1359 /
static_cast<Real
>(factor)));
1360 for (
auto&
x: pmap) {
1361 auto d = std::div(
x,nprocs_con);
1364 }
else if (strategy == 3) {
1365 if (factor == ratio) {
1367 for (
int iproc = 0; iproc < nprocs; ++iproc) {
1368 for (
int ibox : sfc[iproc]) {
1373 for (
auto&
x: pmap) {
1379 dm[i].define(std::move(pmap));
1383 dm[i].define(std::move(pmap_g));
1389template <
typename MF>
1393 BL_PROFILE(
"MLLinOp::makeSubCommunicator()");
1398 std::sort(newgrp_ranks.begin(), newgrp_ranks.end());
1399 auto last = std::unique(newgrp_ranks.begin(), newgrp_ranks.end());
1400 newgrp_ranks.erase(last, newgrp_ranks.end());
1404 MPI_Comm_group(m_default_comm, &defgrp);
1406 MPI_Group_incl(defgrp,
static_cast<int>(newgrp_ranks.
size()), newgrp_ranks.data(), &newgrp);
1410 newgrp_ranks.data(),
static_cast<int>(newgrp_ranks.
size()));
1411 MPI_Group_incl(defgrp,
static_cast<int>(local_newgrp_ranks.
size()), local_newgrp_ranks.data(), &newgrp);
1414 MPI_Comm_create(m_default_comm, newgrp, &newcomm);
1416 m_raii_comm = std::make_unique<CommContainer>(newcomm);
1418 MPI_Group_free(&defgrp);
1419 MPI_Group_free(&newgrp);
1424 return m_default_comm;
1428template <
typename MF>
1433 m_domain_bloc_lo = lo_bcloc;
1434 m_domain_bloc_hi = hi_bcloc;
1437template <
typename MF>
1440 LinOpBCType bc_type)
noexcept
1442 setCoarseFineBC(
crse,
IntVect(crse_ratio), bc_type);
1445template <
typename MF>
1448 LinOpBCType bc_type)
noexcept
1450 m_coarse_data_for_bc =
crse;
1451 m_coarse_data_crse_ratio = crse_ratio;
1452 m_coarse_fine_bc_type = bc_type;
1455template <
typename MF>
1456template <
typename AMF, std::enable_if_t<!std::is_same_v<MF,AMF>,
int>>
1459 LinOpBCType bc_type)
noexcept
1461 setCoarseFineBC(
crse,
IntVect(crse_ratio), bc_type);
1464template <
typename MF>
1465template <
typename AMF, std::enable_if_t<!std::is_same_v<MF,AMF>,
int>>
1468 LinOpBCType bc_type)
noexcept
1470 m_coarse_data_for_bc_raii = MF(
crse->boxArray(),
crse->DistributionMap(),
1472 m_coarse_data_for_bc_raii.LocalCopy(*
crse, 0, 0,
crse->nComp(),
1474 m_coarse_data_for_bc = &m_coarse_data_for_bc_raii;
1475 m_coarse_data_crse_ratio = crse_ratio;
1476 m_coarse_fine_bc_type = bc_type;
1479template <
typename MF>
1484 mf.resize(m_num_amr_levels);
1485 for (
int alev = 0; alev < m_num_amr_levels; ++alev) {
1486 mf[alev].resize(m_num_mg_levels[alev]);
1487 for (
int mlev = 0; mlev < m_num_mg_levels[alev]; ++mlev) {
1488 mf[alev][mlev] = make(alev, mlev, ng);
1493template <
typename MF>
1497 if constexpr (IsMultiFabLike_v<MF>) {
1499 m_dmap[amrlev][mglev], getNComp(), ng,
MFInfo(),
1500 *m_factory[amrlev][mglev]);
1508template <
typename MF>
1512 if constexpr (IsMultiFabLike_v<MF>) {
1516 amrex::Abort(
"MLLinOpT::makeAlias: how did we get here?");
1521template <
typename MF>
1525 if constexpr (IsMultiFabLike_v<MF>) {
1526 BoxArray cba = m_grids[amrlev][mglev];
1527 IntVect ratio = (amrlev > 0) ?
IntVect(2) : mg_coarsen_ratio_vec[mglev];
1530 return MF(cba, m_dmap[amrlev][mglev], getNComp(), ng);
1533 amrex::Abort(
"MLLinOpT::makeCoarseMG: how did we get here?");
1538template <
typename MF>
1542 if constexpr (IsMultiFabLike_v<MF>) {
1543 BoxArray cba = m_grids[famrlev][0];
1544 IntVect ratio(AMRRefRatio(famrlev-1));
1547 return MF(cba, m_dmap[famrlev][0], getNComp(), ng);
1550 amrex::Abort(
"MLLinOpT::makeCoarseAmr: how did we get here?");
1555template <
typename MF>
1559 if (new_size <= 0 || new_size >= m_num_mg_levels[0]) {
return; }
1561 m_num_mg_levels[0] = new_size;
1563 m_geom[0].resize(new_size);
1564 m_grids[0].resize(new_size);
1565 m_dmap[0].resize(new_size);
1566 m_factory[0].resize(new_size);
1568 if (m_bottom_comm != m_default_comm) {
1569 m_bottom_comm = makeSubCommunicator(m_dmap[0].back());
1573template <
typename MF>
1579 const int ncomp = this->getNComp();
1581 if (!fres.isAllRegular()) {
1582 if constexpr (std::is_same<MF,MultiFab>()) {
1584 mg_coarsen_ratio_vec[clev-1]);
1586 amrex::Abort(
"EB_average_down only works with MultiFab");
1594 amrex::Abort(
"For non-FabArray, MLLinOpT<MF>::avgDownResMG should be overridden.");
1598template <
typename MF>
1602 return m_dmap[amrlev][mglev1] == m_dmap[amrlev][mglev2]
1606template <
typename MF>
1607template <
typename AMF, std::enable_if_t<!std::is_same_v<MF,AMF>,
int>>
1610 const AMF* robinbc_a,
const AMF* robinbc_b,
1611 const AMF* robinbc_f)
1613 const int ncomp = this->getNComp();
1615 levelbc_raii[amrlev] = std::make_unique<MF>(levelbcdata->boxArray(),
1616 levelbcdata->DistributionMap(),
1617 ncomp, levelbcdata->nGrowVect());
1618 levelbc_raii[amrlev]->LocalCopy(*levelbcdata, 0, 0, ncomp,
1619 levelbcdata->nGrowVect());
1621 levelbc_raii[amrlev].reset();
1625 robin_a_raii[amrlev] = std::make_unique<MF>(robinbc_a->boxArray(),
1626 robinbc_a->DistributionMap(),
1627 ncomp, robinbc_a->nGrowVect());
1628 robin_a_raii[amrlev]->LocalCopy(*robinbc_a, 0, 0, ncomp,
1629 robinbc_a->nGrowVect());
1631 robin_a_raii[amrlev].reset();
1635 robin_b_raii[amrlev] = std::make_unique<MF>(robinbc_b->boxArray(),
1636 robinbc_b->DistributionMap(),
1637 ncomp, robinbc_b->nGrowVect());
1638 robin_b_raii[amrlev]->LocalCopy(*robinbc_b, 0, 0, ncomp,
1639 robinbc_b->nGrowVect());
1641 robin_b_raii[amrlev].reset();
1645 robin_f_raii[amrlev] = std::make_unique<MF>(robinbc_f->boxArray(),
1646 robinbc_f->DistributionMap(),
1647 ncomp, robinbc_f->nGrowVect());
1648 robin_f_raii[amrlev]->LocalCopy(*robinbc_f, 0, 0, ncomp,
1649 robinbc_f->nGrowVect());
1651 robin_f_raii[amrlev].reset();
1654 this->setLevelBC(amrlev, levelbc_raii[amrlev].
get(), robin_a_raii[amrlev].
get(),
1655 robin_b_raii[amrlev].
get(), robin_f_raii[amrlev].
get());
1658template <
typename MF>
1663 return xdoty(0,0,*
x[0],*y[0],
false);
1666template <
typename MF>
1671 auto r = xdoty(0,0,*
x[0],*
x[0],
false);
1672 return std::sqrt(
r);
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:37
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:129
#define AMREX_D_PICK(a, b, c)
Definition AMReX_SPACE.H:151
int MPI_Comm
Definition AMReX_ccse-mpi.H:47
int MPI_Group
Definition AMReX_ccse-mpi.H:48
static constexpr int MPI_COMM_NULL
Definition AMReX_ccse-mpi.H:55
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:550
static bool SameRefs(const BoxArray &lhs, const BoxArray &rhs)
whether two BoxArrays share the same data
Definition AMReX_BoxArray.H:823
BoxArray & coarsen(int refinement_ratio)
Coarsen each Box in the BoxArray to the specified ratio.
BoxArray & convert(IndexType typ)
Apply Box::convert(IndexType) to each Box in the BoxArray.
AMREX_GPU_HOST_DEVICE BoxND & coarsen(int ref_ratio) noexcept
Coarsen BoxND by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/rati...
Definition AMReX_Box.H:708
AMREX_GPU_HOST_DEVICE BoxND & makeSlab(int direction, int slab_index) noexcept
Definition AMReX_Box.H:779
AMREX_GPU_HOST_DEVICE double d_numPts() const noexcept
Returns the number of points contained in the BoxND. This is intended for use only in diagnostic mess...
Definition AMReX_Box.H:366
AMREX_GPU_HOST_DEVICE IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:146
AMREX_GPU_HOST_DEVICE bool coarsenable(const IntVectND< dim > &refrat, const IntVectND< dim > &min_width) const noexcept
Definition AMReX_Box.H:745
AMREX_GPU_HOST_DEVICE const IntVectND< dim > & smallEnd() const &noexcept
Get the smallend of the BoxND.
Definition AMReX_Box.H:105
AMREX_GPU_HOST_DEVICE Long numPts() const noexcept
Returns the number of points contained in the BoxND.
Definition AMReX_Box.H:346
Definition AMReX_FabFactory.H:76
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:41
static DistributionMapping makeSFC(const MultiFab &weight, bool sort=true)
Definition AMReX_DistributionMapping.cpp:1764
const Vector< int > & ProcessorMap() const noexcept
Returns a constant reference to the mapping of boxes in the underlying BoxArray to the CPU that holds...
Definition AMReX_DistributionMapping.cpp:47
Definition AMReX_EBFabFactory.H:24
Definition AMReX_FabFactory.H:50
Solve using GMRES with multigrid as preconditioner.
Definition AMReX_GMRES_MLMG.H:21
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:73
Interface
Definition AMReX_Hypre.H:21
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE constexpr std::size_t size() noexcept
Definition AMReX_IntVect.H:723
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE IntVectND< new_dim > resize(int fill_extra=0) const noexcept
Returns a new IntVectND of size new_dim by either shrinking or expanding this IntVectND.
Definition AMReX_IntVect.H:766
Definition AMReX_MLABecLaplacian.H:15
Definition AMReX_MLCGSolver.H:12
Definition AMReX_MLLinOp.H:98
const MF * m_coarse_data_for_bc
Definition AMReX_MLLinOp.H:638
void setDomainBC(const Array< BCType, AMREX_SPACEDIM > &lobc, const Array< BCType, AMREX_SPACEDIM > &hibc) noexcept
Boundary of the whole domain.
Definition AMReX_MLLinOp.H:1193
Vector< Vector< std::unique_ptr< FabFactory< FAB > > > > m_factory
Definition AMReX_MLLinOp.H:611
virtual void avgDownResMG(int clev, MF &cres, MF const &fres) const
Definition AMReX_MLLinOp.H:1575
int NAMRLevels() const noexcept
Return the number of AMR levels.
Definition AMReX_MLLinOp.H:564
bool m_do_consolidation
Definition AMReX_MLLinOp.H:602
virtual void compGrad(int, const Array< MF *, AMREX_SPACEDIM > &, MF &, Location) const
Compute gradients of the solution.
Definition AMReX_MLLinOp.H:438
bool isCellCentered() const noexcept
Definition AMReX_MLLinOp.H:685
IntVect m_ixtype
Definition AMReX_MLLinOp.H:599
void setVerbose(int v) noexcept
Set verbosity.
Definition AMReX_MLLinOp.H:244
Array< Real, AMREX_SPACEDIM > m_domain_bloc_hi
Definition AMReX_MLLinOp.H:632
bool isMFIterSafe(int amrlev, int mglev1, int mglev2) const
Definition AMReX_MLLinOp.H:1600
RealVect m_coarse_bc_loc
Definition AMReX_MLLinOp.H:637
virtual bool needsUpdate() const
Does it need update if it's reused?
Definition AMReX_MLLinOp.H:266
virtual void interpolation(int amrlev, int fmglev, MF &fine, const MF &crse) const =0
Add interpolated coarse MG level data to fine MG level data.
virtual void setLevelBC(int, const MF *, const MF *=nullptr, const MF *=nullptr, const MF *=nullptr)=0
Set boundary conditions for given level. For cell-centered solves only.
virtual MF make(int amrlev, int mglev, IntVect const &ng) const
Definition AMReX_MLLinOp.H:1495
FabFactory< FAB > const * Factory(int amr_lev, int mglev=0) const noexcept
Definition AMReX_MLLinOp.H:649
T get_d0(T const &dx, T const &dy, T const &) const noexcept
Definition AMReX_MLLinOp.H:726
MPI_Comm BottomCommunicator() const noexcept
Definition AMReX_MLLinOp.H:676
void setEnforceSingularSolvable(bool o) noexcept
Definition AMReX_MLLinOp.H:253
MPI_Comm Communicator() const noexcept
Definition AMReX_MLLinOp.H:677
virtual void getFluxes(const Vector< Array< MF *, AMREX_SPACEDIM > > &, const Vector< MF * > &, Location) const
Definition AMReX_MLLinOp.H:501
int mg_domain_min_width
Definition AMReX_MLLinOp.H:583
void setMaxOrder(int o) noexcept
Set order of interpolation at coarse/fine boundary.
Definition AMReX_MLLinOp.H:247
Vector< Array< BCType, AMREX_SPACEDIM > > m_hibc_orig
Definition AMReX_MLLinOp.H:577
virtual void interpAssign(int amrlev, int fmglev, MF &fine, MF &crse) const
Overwrite fine MG level data with interpolated coarse data.
Definition AMReX_MLLinOp.H:299
virtual std::string name() const
Definition AMReX_MLLinOp.H:131
void setCoarseFineBC(const AMF *crse, IntVect const &crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
Definition AMReX_MLLinOp.H:1467
bool doAgglomeration() const noexcept
Definition AMReX_MLLinOp.H:681
Vector< std::unique_ptr< MF > > levelbc_raii
Definition AMReX_MLLinOp.H:761
MF m_coarse_data_for_bc_raii
Definition AMReX_MLLinOp.H:639
virtual void setDirichletNodesToZero(int, int, MF &) const
Definition AMReX_MLLinOp.H:477
MLLinOpT< MF > & operator=(const MLLinOpT< MF > &)=delete
std::unique_ptr< CommContainer > m_raii_comm
Definition AMReX_MLLinOp.H:629
bool m_do_semicoarsening
Definition AMReX_MLLinOp.H:604
bool hasRobinBC() const noexcept
Definition AMReX_MLLinOp.H:1277
Array< Real, AMREX_SPACEDIM > m_domain_bloc_lo
Definition AMReX_MLLinOp.H:631
virtual void resizeMultiGrid(int new_size)
Definition AMReX_MLLinOp.H:1557
Vector< Vector< BoxArray > > m_grids
Definition AMReX_MLLinOp.H:609
virtual MF makeCoarseAmr(int famrlev, IntVect const &ng) const
Definition AMReX_MLLinOp.H:1540
Vector< std::unique_ptr< MF > > robin_f_raii
Definition AMReX_MLLinOp.H:764
bool m_do_agglomeration
Definition AMReX_MLLinOp.H:601
void setCoarseFineBC(const AMF *crse, int crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
Definition AMReX_MLLinOp.H:1458
virtual MF makeAlias(MF const &mf) const
Definition AMReX_MLLinOp.H:1510
void setDomainBC(const Vector< Array< BCType, AMREX_SPACEDIM > > &lobc, const Vector< Array< BCType, AMREX_SPACEDIM > > &hibc) noexcept
Boundary of the whole domain.
Definition AMReX_MLLinOp.H:1203
Array4< T > compactify(Array4< T > const &a) const noexcept
Definition AMReX_MLLinOp.H:712
static constexpr int mg_coarsen_ratio
Definition AMReX_MLLinOp.H:581
virtual void copyNSolveSolution(MF &, MF const &) const
Definition AMReX_MLLinOp.H:541
virtual void solutionResidual(int amrlev, MF &resid, MF &x, const MF &b, const MF *crse_bcdata=nullptr)=0
Compute residual for solution.
int getMaxOrder() const noexcept
Get order of interpolation at coarse/fine boundary.
Definition AMReX_MLLinOp.H:249
virtual int getNComp() const
Return number of components.
Definition AMReX_MLLinOp.H:261
virtual void getFluxes(const Vector< MF * > &, const Vector< MF * > &) const
Definition AMReX_MLLinOp.H:506
void defineBC()
Definition AMReX_MLLinOp.H:1181
void setCoarseFineBCLocation(const RealVect &cloc) noexcept
Definition AMReX_MLLinOp.H:679
Vector< int > m_amr_ref_ratio
Definition AMReX_MLLinOp.H:594
MPI_Comm m_default_comm
Definition AMReX_MLLinOp.H:614
static void makeAgglomeratedDMap(const Vector< BoxArray > &ba, Vector< DistributionMapping > &dm)
Definition AMReX_MLLinOp.H:1306
bool isBottomActive() const noexcept
Definition AMReX_MLLinOp.H:674
virtual Vector< RT > getSolvabilityOffset(int, int, MF const &) const
get offset for fixing solvability
Definition AMReX_MLLinOp.H:464
void defineGrids(const Vector< Geometry > &a_geom, const Vector< BoxArray > &a_grids, const Vector< DistributionMapping > &a_dmap, const Vector< FabFactory< FAB > const * > &a_factory)
Definition AMReX_MLLinOp.H:807
virtual BottomSolver getDefaultBottomSolver() const
Definition AMReX_MLLinOp.H:258
typename FabDataType< MF >::fab_type FAB
Definition AMReX_MLLinOp.H:108
virtual RT normInf(int amrlev, MF const &mf, bool local) const =0
Vector< int > m_num_mg_levels
Definition AMReX_MLLinOp.H:596
bool hasBC(BCType bct) const noexcept
Definition AMReX_MLLinOp.H:1255
static void makeConsolidatedDMap(const Vector< BoxArray > &ba, Vector< DistributionMapping > &dm, int ratio, int strategy)
Definition AMReX_MLLinOp.H:1335
Vector< Vector< DistributionMapping > > m_dmap
Definition AMReX_MLLinOp.H:610
int verbose
Definition AMReX_MLLinOp.H:587
virtual void checkPoint(std::string const &) const
Definition AMReX_MLLinOp.H:757
IntVect m_coarse_data_crse_ratio
Definition AMReX_MLLinOp.H:636
LinOpBCType BCType
Definition AMReX_MLLinOp.H:111
virtual void correctionResidual(int amrlev, int mglev, MF &resid, MF &x, const MF &b, BCMode bc_mode, const MF *crse_bcdata=nullptr)=0
Compute residual for the residual-correction form, resid = b - L(x)
virtual void compFlux(int, const Array< MF *, AMREX_SPACEDIM > &, MF &, Location) const
Compute fluxes.
Definition AMReX_MLLinOp.H:424
const Vector< int > & AMRRefRatio() const noexcept
Return AMR refinement ratios.
Definition AMReX_MLLinOp.H:644
MLLinOpT(MLLinOpT< MF > &&)=delete
Vector< Array< BCType, AMREX_SPACEDIM > > m_lobc_orig
Definition AMReX_MLLinOp.H:576
void setCoarseFineBC(const MF *crse, int crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
Set coarse/fine boundary conditions. For cell-centered solves only.
Definition AMReX_MLLinOp.H:1439
virtual void apply(int amrlev, int mglev, MF &out, MF &in, BCMode bc_mode, StateMode s_mode, const MLMGBndryT< MF > *bndry=nullptr) const =0
Apply the linear operator, out = L(in)
GpuArray< BCType, AMREX_SPACEDIM > HiBC(int icomp=0) const noexcept
Definition AMReX_MLLinOp.H:658
void setDomainBCLoc(const Array< Real, AMREX_SPACEDIM > &lo_bcloc, const Array< Real, AMREX_SPACEDIM > &hi_bcloc) noexcept
Set location of domain boundaries.
Definition AMReX_MLLinOp.H:1430
bool needsCoarseDataForBC() const noexcept
Needs coarse data for bc?
Definition AMReX_MLLinOp.H:177
typename FabDataType< MF >::value_type RT
Definition AMReX_MLLinOp.H:109
virtual void update()
Update for reuse.
Definition AMReX_MLLinOp.H:268
bool m_precond_mode
Definition AMReX_MLLinOp.H:641
virtual std::unique_ptr< FabFactory< FAB > > makeFactory(int, int) const
Definition AMReX_MLLinOp.H:701
virtual void applyMetricTerm(int, int, MF &) const
apply metric terms if there are any
Definition AMReX_MLLinOp.H:445
virtual void unapplyMetricTerm(int, int, MF &) const
unapply metric terms if there are any
Definition AMReX_MLLinOp.H:447
bool hasHiddenDimension() const noexcept
Definition AMReX_MLLinOp.H:707
virtual bool isBottomSingular() const =0
Is the bottom of MG singular?
virtual void reflux(int crse_amrlev, MF &res, const MF &crse_sol, const MF &crse_rhs, MF &fine_res, MF &fine_sol, const MF &fine_rhs) const
Reflux at AMR coarse/fine boundary.
Definition AMReX_MLLinOp.H:407
virtual IntVect getNGrowVectRestriction() const
Definition AMReX_MLLinOp.H:687
Vector< std::unique_ptr< MF > > robin_b_raii
Definition AMReX_MLLinOp.H:763
static constexpr int mg_box_min_width
Definition AMReX_MLLinOp.H:582
virtual RT dotProductPrecond(Vector< MF const * > const &x, Vector< MF const * > const &y) const
Definition AMReX_MLLinOp.H:1660
virtual MF makeCoarseMG(int amrlev, int mglev, IntVect const &ng) const
Definition AMReX_MLLinOp.H:1523
virtual void fixSolvabilityByOffset(int, int, MF &, Vector< RT > const &) const
fix solvability by subtracting offset from RHS
Definition AMReX_MLLinOp.H:468
const Geometry & Geom(int amr_lev, int mglev=0) const noexcept
Definition AMReX_MLLinOp.H:569
virtual void endPrecondBC()
Definition AMReX_MLLinOp.H:559
Vector< Array< BCType, AMREX_SPACEDIM > > m_lobc
Definition AMReX_MLLinOp.H:572
int hiddenDirection() const noexcept
Definition AMReX_MLLinOp.H:708
Vector< int > m_domain_covered
Definition AMReX_MLLinOp.H:612
const MLLinOpT< MF > * m_parent
Definition AMReX_MLLinOp.H:597
bool doSemicoarsening() const noexcept
Definition AMReX_MLLinOp.H:683
virtual bool supportNSolve() const
Definition AMReX_MLLinOp.H:539
virtual bool supportRobinBC() const noexcept
Definition AMReX_MLLinOp.H:668
virtual void normalize(int amrlev, int mglev, MF &mf) const
Divide mf by the diagonal component of the operator. Used by bicgstab.
Definition AMReX_MLLinOp.H:364
virtual void applyInhomogNeumannTerm(int, MF &) const
Extra terms introduced when we treat inhomogeneous Nuemann BC as homogeneous.
Definition AMReX_MLLinOp.H:453
virtual void avgDownResAmr(int clev, MF &cres, MF const &fres) const
Definition AMReX_MLLinOp.H:549
Vector< Vector< Geometry > > m_geom
first Vector is for amr level and second is mg level
Definition AMReX_MLLinOp.H:608
virtual RT norm2Precond(Vector< MF const * > const &x) const
Definition AMReX_MLLinOp.H:1668
MLLinOpT(const MLLinOpT< MF > &)=delete
virtual bool scaleRHS(int, MF *) const
scale RHS to fix solvability
Definition AMReX_MLLinOp.H:459
virtual void averageDownAndSync(Vector< MF > &sol) const =0
void setCoarseFineBC(const MF *crse, IntVect const &crse_ratio, LinOpBCType bc_type=LinOpBCType::Dirichlet) noexcept
Definition AMReX_MLLinOp.H:1447
GpuArray< BCType, AMREX_SPACEDIM > LoBC(int icomp=0) const noexcept
Definition AMReX_MLLinOp.H:653
Box compactify(Box const &b) const noexcept
Definition AMReX_MLLinOp.H:1284
bool m_needs_coarse_data_for_bc
Definition AMReX_MLLinOp.H:634
int maxorder
Definition AMReX_MLLinOp.H:589
void define(const Vector< Geometry > &a_geom, const Vector< BoxArray > &a_grids, const Vector< DistributionMapping > &a_dmap, const LPInfo &a_info, const Vector< FabFactory< FAB > const * > &a_factory, bool eb_limit_coarsening=true)
Definition AMReX_MLLinOp.H:769
Vector< IntVect > mg_coarsen_ratio_vec
Definition AMReX_MLLinOp.H:605
MF MFType
Definition AMReX_MLLinOp.H:107
virtual void preparePrecond()
Definition AMReX_MLLinOp.H:473
virtual ~MLLinOpT()=default
virtual void averageDownSolutionRHS(int camrlev, MF &crse_sol, MF &crse_rhs, const MF &fine_sol, const MF &fine_rhs)
Average-down data from fine AMR level to coarse AMR level.
Definition AMReX_MLLinOp.H:329
LPInfo info
Definition AMReX_MLLinOp.H:585
MPI_Comm makeSubCommunicator(const DistributionMapping &dm)
Definition AMReX_MLLinOp.H:1391
Vector< Array< BCType, AMREX_SPACEDIM > > m_hibc
Definition AMReX_MLLinOp.H:573
int NMGLevels(int amrlev) const noexcept
Return the number of MG levels at given AMR level.
Definition AMReX_MLLinOp.H:567
T get_d1(T const &, T const &dy, T const &dz) const noexcept
Definition AMReX_MLLinOp.H:736
bool enforceSingularSolvable
Definition AMReX_MLLinOp.H:591
void setLevelBC(int amrlev, const AMF *levelbcdata, const AMF *robinbc_a=nullptr, const AMF *robinbc_b=nullptr, const AMF *robinbc_f=nullptr)
Definition AMReX_MLLinOp.H:1609
virtual RT xdoty(int amrlev, int mglev, const MF &x, const MF &y, bool local) const =0
x dot y, used by the bottom solver
virtual void interpolationAmr(int famrlev, MF &fine, const MF &crse, IntVect const &nghost) const
Interpolation between AMR levels.
Definition AMReX_MLLinOp.H:313
bool doConsolidation() const noexcept
Definition AMReX_MLLinOp.H:682
virtual std::unique_ptr< MLLinOpT< MF > > makeNLinOp(int) const
Definition AMReX_MLLinOp.H:495
virtual bool supportInhomogNeumannBC() const noexcept
Definition AMReX_MLLinOp.H:669
virtual void prepareForFluxes(int, const MF *=nullptr)
Definition AMReX_MLLinOp.H:380
LinOpBCType m_coarse_fine_bc_type
Definition AMReX_MLLinOp.H:635
virtual bool isSingular(int amrlev) const =0
Is it singular on given AMR level?
virtual void postSolve(Vector< MF * > const &) const
Definition AMReX_MLLinOp.H:543
int AMRRefRatio(int amr_lev) const noexcept
Return AMR refinement ratio at given AMR level.
Definition AMReX_MLLinOp.H:647
MPI_Comm m_bottom_comm
Definition AMReX_MLLinOp.H:615
virtual void restriction(int amrlev, int cmglev, MF &crse, MF &fine) const =0
Restriction onto coarse MG level.
virtual void beginPrecondBC()
Definition AMReX_MLLinOp.H:558
virtual void applyOverset(int, MF &) const
for overset solver only
Definition AMReX_MLLinOp.H:456
virtual void unimposeNeumannBC(int, MF &) const
This is needed for our nodal projection solver.
Definition AMReX_MLLinOp.H:450
bool getEnforceSingularSolvable() const noexcept
Definition AMReX_MLLinOp.H:256
virtual void prepareForSolve()=0
virtual void smooth(int amrlev, int mglev, MF &sol, const MF &rhs, bool skip_fillboundary, int niter) const =0
Smooth.
virtual void make(Vector< Vector< MF > > &mf, IntVect const &ng) const
Definition AMReX_MLLinOp.H:1481
bool hasInhomogNeumannBC() const noexcept
Definition AMReX_MLLinOp.H:1270
Vector< std::unique_ptr< MF > > robin_a_raii
Definition AMReX_MLLinOp.H:762
virtual int getNGrow(int=0, int=0) const
Definition AMReX_MLLinOp.H:263
int m_num_amr_levels
Definition AMReX_MLLinOp.H:593
Definition AMReX_MLMGBndry.H:12
Definition AMReX_MLMG.H:12
Definition AMReX_MLPoisson.H:16
This class provides the user with a few print options.
Definition AMReX_Print.H:35
A Box with real dimensions. A RealBox is OK iff volume >= 0.
Definition AMReX_RealBox.H:21
A Real vector in SpaceDim-dimensional space.
Definition AMReX_RealVect.H:32
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:27
Long size() const noexcept
Definition AMReX_Vector.H:50
bool notInLaunchRegion() noexcept
Definition AMReX_GpuControl.H:87
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition AMReX_ParallelContext.H:70
int local_to_global_rank(int rank) noexcept
translate between local rank and global rank
Definition AMReX_ParallelContext.H:95
int global_to_local_rank(int rank) noexcept
Definition AMReX_ParallelContext.H:98
int NProcsSub() noexcept
number of ranks in current frame
Definition AMReX_ParallelContext.H:74
MPI_Comm Communicator() noexcept
Definition AMReX_ParallelDescriptor.H:210
Definition AMReX_Amr.cpp:49
@ make_alias
Definition AMReX_MakeType.H:7
void average_down(const MultiFab &S_fine, MultiFab &S_crse, const Geometry &fgeom, const Geometry &cgeom, int scomp, int ncomp, int rr)
Definition AMReX_MultiFabUtil.cpp:314
BoxND< AMREX_SPACEDIM > Box
Definition AMReX_BaseFwd.H:27
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > enclosedCells(const BoxND< dim > &b, int dir) noexcept
Returns a BoxND with CELL based coordinates in direction dir that is enclosed by b....
Definition AMReX_Box.H:1463
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Returns a BoxND with different type.
Definition AMReX_Box.H:1435
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/rati...
Definition AMReX_Box.H:1304
void EB_average_down(const MultiFab &S_fine, MultiFab &S_crse, const MultiFab &vol_fine, const MultiFab &vfrac_fine, int scomp, int ncomp, const IntVect &ratio)
Definition AMReX_EBMultiFabUtil.cpp:336
IntVectND< AMREX_SPACEDIM > IntVect
Definition AMReX_BaseFwd.H:30
AMREX_GPU_HOST_DEVICE constexpr GpuTupleElement< I, GpuTuple< Ts... > >::type & get(GpuTuple< Ts... > &tup) noexcept
Definition AMReX_Tuple.H:179
BottomSolver
Definition AMReX_MLLinOp.H:30
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:127
void Warning(const std::string &msg)
Print out warning message to cerr.
Definition AMReX.cpp:236
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230
int verbose
Definition AMReX_DistributionMapping.cpp:36
std::unique_ptr< Hypre > makeHypre(const BoxArray &grids, const DistributionMapping &dmap, const Geometry &geom, MPI_Comm comm_, Hypre::Interface interface, const iMultiFab *overset_mask)
Definition AMReX_Hypre.cpp:12
std::array< T, N > Array
Definition AMReX_Array.H:24
Definition AMReX_Array4.H:61
Definition AMReX_FabDataType.H:9
Definition AMReX_Array.H:34
Definition AMReX_TypeTraits.H:29
Definition AMReX_MLLinOp.H:35
LPInfo & setConsolidationRatio(int x) noexcept
Definition AMReX_MLLinOp.H:54
int con_strategy
Definition AMReX_MLLinOp.H:42
bool do_semicoarsening
Definition AMReX_MLLinOp.H:38
bool has_metric_term
Definition AMReX_MLLinOp.H:43
LPInfo & setConsolidationGridSize(int x) noexcept
Definition AMReX_MLLinOp.H:53
int max_semicoarsening_level
Definition AMReX_MLLinOp.H:45
bool hasHiddenDimension() const noexcept
Definition AMReX_MLLinOp.H:62
int con_ratio
Definition AMReX_MLLinOp.H:41
bool do_consolidation
Definition AMReX_MLLinOp.H:37
int con_grid_size
Definition AMReX_MLLinOp.H:40
LPInfo & setSemicoarsening(bool x) noexcept
Definition AMReX_MLLinOp.H:51
LPInfo & setHiddenDirection(int n) noexcept
Definition AMReX_MLLinOp.H:60
LPInfo & setConsolidation(bool x) noexcept
Definition AMReX_MLLinOp.H:50
LPInfo & setSemicoarseningDirection(int n) noexcept
Definition AMReX_MLLinOp.H:59
LPInfo & setMaxSemicoarseningLevel(int n) noexcept
Definition AMReX_MLLinOp.H:58
bool do_agglomeration
Definition AMReX_MLLinOp.H:36
LPInfo & setMetricTerm(bool x) noexcept
Definition AMReX_MLLinOp.H:56
static constexpr int getDefaultConsolidationGridSize()
Definition AMReX_MLLinOp.H:74
LPInfo & setConsolidationStrategy(int x) noexcept
Definition AMReX_MLLinOp.H:55
int max_coarsening_level
Definition AMReX_MLLinOp.H:44
int agg_grid_size
Definition AMReX_MLLinOp.H:39
static constexpr int getDefaultAgglomerationGridSize()
Definition AMReX_MLLinOp.H:66
int hidden_direction
Definition AMReX_MLLinOp.H:47
LPInfo & setAgglomerationGridSize(int x) noexcept
Definition AMReX_MLLinOp.H:52
LPInfo & setMaxCoarseningLevel(int n) noexcept
Definition AMReX_MLLinOp.H:57
LPInfo & setAgglomeration(bool x) noexcept
Definition AMReX_MLLinOp.H:49
int semicoarsening_direction
Definition AMReX_MLLinOp.H:46
Definition AMReX_MLLinOp.H:84
StateMode
Definition AMReX_MLLinOp.H:86
BCMode
Definition AMReX_MLLinOp.H:85
Location
Definition AMReX_MLLinOp.H:87
FabArray memory allocation information.
Definition AMReX_FabArray.H:66
Definition AMReX_MLLinOp.H:616
void operator=(const CommContainer &)=delete
~CommContainer()
Definition AMReX_MLLinOp.H:623
MPI_Comm comm
Definition AMReX_MLLinOp.H:617
CommContainer(const CommContainer &)=delete
CommContainer(CommContainer &&)=delete
CommContainer(MPI_Comm m) noexcept
Definition AMReX_MLLinOp.H:618