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)
94template <
typename T>
class MLMGT;
95template <
typename T>
class MLCGSolverT;
96template <
typename T>
class MLPoissonT;
97template <
typename T>
class MLABecLaplacianT;
98template <
typename T>
class GMRESMLMGT;
101template <
typename MF>
106 template <
typename T>
friend class MLMGT;
134 bool eb_limit_coarsening =
true);
136 [[nodiscard]]
virtual std::string
name ()
const {
return std::string(
"Unspecified"); }
211 template <
typename AMF>
212 requires (!std::same_as<MF,AMF>)
216 template <
typename AMF>
217 requires (!std::same_as<MF,AMF>)
243 const MF* =
nullptr) = 0;
245 template <MultiFabLike AMF>
246 requires (!std::same_as<MF,AMF>)
248 const AMF* robinbc_a =
nullptr,
249 const AMF* robinbc_b =
nullptr,
250 const AMF* robinbc_f =
nullptr);
270 [[nodiscard]]
virtual int getNComp ()
const {
return 1; }
272 [[nodiscard]]
virtual int getNGrow (
int = 0,
int = 0)
const {
return 0; }
311 amrex::Abort(
"MLLinOpT::interpAssign: Must be implemented for FMG cycle");
326 amrex::Abort(
"MLLinOpT::interpolationAmr: Must be implemented for composite solves across multiple AMR levels");
339 const MF& fine_sol,
const MF& fine_rhs)
342 amrex::Abort(
"MLLinOpT::averageDownSolutionRHS: Must be implemented for composite solves across multiple AMR levels");
356 virtual void apply (
int amrlev,
int mglev, MF& out, MF& in,
BCMode bc_mode,
369 virtual void smooth (
int amrlev,
int mglev, MF& sol,
const MF& rhs,
370 bool skip_fillboundary,
int niter)
const = 0;
373 virtual void normalize (
int amrlev,
int mglev, MF& mf)
const {
387 const MF* crse_bcdata=
nullptr) = 0;
403 BCMode bc_mode,
const MF* crse_bcdata=
nullptr) = 0;
417 MF& res,
const MF& crse_sol,
const MF& crse_rhs,
418 MF& fine_res, MF& fine_sol,
const MF& fine_rhs)
const
422 amrex::Abort(
"MLLinOpT::reflux: Must be implemented for composite solves across multiple AMR levels");
437 amrex::Abort(
"AMReX_MLLinOp::compFlux::How did we get here?");
452 amrex::Abort(
"AMReX_MLLinOp::compGrad::How did we get here?");
470 [[nodiscard]]
virtual bool scaleRHS (
int , MF* )
const {
476 MF
const& )
const {
return {}; }
491 amrex::Warning(
"This function might need to be implemented for GMRES to work with this LinOp.");
495 [[nodiscard]]
virtual bool isSingular (
int amrlev)
const = 0;
500 virtual RT xdoty (
int amrlev,
int mglev,
const MF&
x,
const MF&
y,
bool local)
const = 0;
506 virtual std::unique_ptr<MLLinOpT<MF>>
makeNLinOp (
int )
const
508 amrex::Abort(
"MLLinOp::makeNLinOp: N-Solve not supported");
515 amrex::Abort(
"MLLinOp::getFluxes: How did we get here?");
519 amrex::Abort(
"MLLinOp::getFluxes: How did we get here?");
525 amrex::Abort(
"MLLinOp::getEBFluxes: How did we get here?");
529#if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
531 amrex::Abort(
"MLLinOp::makeHypre: How did we get here?");
534 [[nodiscard]]
virtual std::unique_ptr<HypreNodeLap> makeHypreNodeLap(
536 const std::string& )
const
538 amrex::Abort(
"MLLinOp::makeHypreNodeLap: How did we get here?");
543#if defined(AMREX_USE_PETSC) && (AMREX_SPACEDIM > 1)
544 [[nodiscard]]
virtual std::unique_ptr<PETScABecLap> makePETSc ()
const {
545 amrex::Abort(
"MLLinOp::makePETSc: How did we get here?");
556 [[nodiscard]]
virtual RT normInf (
int amrlev, MF
const& mf,
bool local)
const = 0;
563 amrex::Abort(
"MLLinOpT::avgDownResAmr: Must be implemented for composite solves across multiple AMR levels");
572 [[nodiscard]]
bool isMFIterSafe (
int amrlev,
int mglev1,
int mglev2)
const;
580 [[nodiscard]]
const Geometry&
Geom (
int amr_lev,
int mglev=0) const noexcept {
return m_geom[amr_lev][mglev]; }
628 struct CommContainer {
630 CommContainer (
MPI_Comm m) noexcept : comm(m) {}
631 CommContainer (
const CommContainer&) =
delete;
632 CommContainer (CommContainer&&) =
delete;
633 void operator= (
const CommContainer&) =
delete;
634 void operator= (CommContainer&&) =
delete;
637 if (comm != MPI_COMM_NULL) { MPI_Comm_free(&comm); }
713 [[nodiscard]]
virtual MF
make (
int amrlev,
int mglev,
IntVect const& ng)
const;
715 [[nodiscard]]
virtual MF
makeAlias (MF
const& mf)
const;
721 [[nodiscard]]
virtual std::unique_ptr<FabFactory<FAB> >
makeFactory (
int ,
int )
const {
722 return std::make_unique<DefaultFabFactory<FAB>>();
731 template <
typename T>
735 return Array4<T>(a.dataPtr(), {a.begin[1],a.begin[2],0}, {a.end[1],a.end[2],1}, a.nComp());
737 return Array4<T>(a.dataPtr(), {a.begin[0],a.begin[2],0}, {a.end[0],a.end[2],1}, a.nComp());
739 return Array4<T>(a.dataPtr(), {a.begin[0],a.begin[1],0}, {a.end[0],a.end[1],1}, a.nComp());
745 template <
typename T>
746 [[nodiscard]] T
get_d0 (T
const& dx, T
const& dy, T
const&)
const noexcept
755 template <
typename T>
756 [[nodiscard]] T
get_d1 (T
const&, T
const& dy, T
const& dz)
const noexcept
774 int ratio,
int strategy);
777 virtual void checkPoint (std::string
const& )
const {
781 Vector<std::unique_ptr<MF>> levelbc_raii;
782 Vector<std::unique_ptr<MF>> robin_a_raii;
783 Vector<std::unique_ptr<MF>> robin_b_raii;
784 Vector<std::unique_ptr<MF>> robin_f_raii;
787template <
typename MF>
794 [[maybe_unused]]
bool eb_limit_coarsening)
803 if (info.con_grid_size <= 0) { info.con_grid_size =
AMREX_D_PICK(32, 16, 8); }
813 if (!a_factory.empty() && eb_limit_coarsening) {
816 info.max_coarsening_level = std::min(info.max_coarsening_level,
817 f->maxCoarseningLevel());
821 defineGrids(a_geom, a_grids, a_dmap, a_factory);
825template <
typename MF>
835 if ( ! a_factory.empty() ) {
837 if (ebf && !(ebf->isAllRegular())) {
838 mg_domain_min_width = 4;
843 m_num_amr_levels = 0;
844 for (
int amrlev = 0; amrlev < std::ssize(a_geom); amrlev++) {
845 if (!a_grids[amrlev].empty()) {
850 m_amr_ref_ratio.resize(m_num_amr_levels);
851 m_num_mg_levels.resize(m_num_amr_levels);
853 m_geom.resize(m_num_amr_levels);
854 m_grids.resize(m_num_amr_levels);
855 m_dmap.resize(m_num_amr_levels);
856 m_factory.resize(m_num_amr_levels);
860 const RealBox& rb = a_geom[0].ProbDomain();
861 const int coord = a_geom[0].Coord();
862 const Array<int,AMREX_SPACEDIM>& is_per = a_geom[0].isPeriodic();
864 IntVect mg_coarsen_ratio_v(mg_coarsen_ratio);
865 IntVect mg_box_min_width_v(mg_box_min_width);
866 IntVect mg_domain_min_width_v(mg_domain_min_width);
867 if (hasHiddenDimension()) {
869 "Hidden direction only supported for 3d");
870 mg_coarsen_ratio_v[info.hidden_direction] = 1;
871 mg_box_min_width_v[info.hidden_direction] = 0;
872 mg_domain_min_width_v[info.hidden_direction] = 0;
876 for (
int amrlev = m_num_amr_levels-1; amrlev > 0; --amrlev)
878 m_num_mg_levels[amrlev] = 1;
879 m_geom[amrlev].push_back(a_geom[amrlev]);
880 m_grids[amrlev].push_back(a_grids[amrlev]);
881 m_dmap[amrlev].push_back(a_dmap[amrlev]);
882 if (amrlev < std::ssize(a_factory)) {
883 m_factory[amrlev].emplace_back(a_factory[amrlev]->clone());
885 m_factory[amrlev].push_back(std::make_unique<DefaultFabFactory<FAB>>());
888 IntVect rr = mg_coarsen_ratio_v;
889 const Box& dom = a_geom[amrlev].Domain();
890 for (
int i = 0; i < 2; ++i)
892 if (!dom.coarsenable(rr)) {
amrex::Abort(
"MLLinOp: Uncoarsenable domain"); }
895 if (cdom == a_geom[amrlev-1].Domain()) {
break; }
897 ++(m_num_mg_levels[amrlev]);
899 m_geom[amrlev].emplace_back(cdom, rb, coord, is_per);
901 m_grids[amrlev].push_back(a_grids[amrlev]);
903 m_grids[amrlev].back().coarsen(rr);
905 m_dmap[amrlev].push_back(a_dmap[amrlev]);
907 rr *= mg_coarsen_ratio_v;
910#if (AMREX_SPACEDIM > 1)
911 if (hasHiddenDimension()) {
912 m_amr_ref_ratio[amrlev-1] = rr[(info.hidden_direction+1) % AMREX_SPACEDIM];
916 m_amr_ref_ratio[amrlev-1] = rr[0];
921 m_num_mg_levels[0] = 1;
922 m_geom[0].push_back(a_geom[0]);
923 m_grids[0].push_back(a_grids[0]);
924 m_dmap[0].push_back(a_dmap[0]);
925 if (!a_factory.empty()) {
926 m_factory[0].emplace_back(a_factory[0]->clone());
928 m_factory[0].push_back(std::make_unique<DefaultFabFactory<FAB>>());
931 m_domain_covered.resize(m_num_amr_levels,
false);
932 auto npts0 = m_grids[0][0].numPts();
933 m_domain_covered[0] = (npts0 == compactify(m_geom[0][0].Domain()).numPts());
934 for (
int amrlev = 1; amrlev < m_num_amr_levels; ++amrlev)
936 if (!m_domain_covered[amrlev-1]) {
break; }
937 m_domain_covered[amrlev] = (m_grids[amrlev][0].numPts() ==
938 compactify(m_geom[amrlev][0].Domain()).numPts());
942 bool aggable =
false;
944 if (m_grids[0][0].size() > 1 && info.do_agglomeration)
946 if (m_domain_covered[0])
948 aggbox = m_geom[0][0].Domain();
949 if (hasHiddenDimension()) {
950 aggbox.
makeSlab(hiddenDirection(), m_grids[0][0][0].smallEnd(hiddenDirection()));
956 aggbox = m_grids[0][0].minimalBox();
957 aggable = (aggbox.numPts() == npts0);
963 int agg_lev = 0, con_lev = 0;
966 && info.semicoarsening_direction >= -1
967 && info.semicoarsening_direction < AMREX_SPACEDIM );
969 if (info.do_agglomeration && aggable)
971 Box dbx = m_geom[0][0].Domain();
973 Real const nbxs =
static_cast<Real>(m_grids[0][0].size());
976 *info.agg_grid_size));
977 Vector<Box> domainboxes{dbx};
978 Vector<Box> boundboxes{bbx};
979 Vector<int> agg_flag{
false};
980 Vector<IntVect> accum_coarsen_ratio{
IntVect(1)};
983 for (
int lev = 0; lev < info.max_coarsening_level; ++lev)
985 IntVect rr_level = mg_coarsen_ratio_v;
986 bool const do_semicoarsening_level = info.do_semicoarsening
987 && numsclevs < info.max_semicoarsening_level;
988 if (do_semicoarsening_level
989 && info.semicoarsening_direction != -1)
991 rr_level[info.semicoarsening_direction] = 1;
994 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
996 rr_dir[idim] = rr_level[idim];
997 is_coarsenable[idim] = dbx.coarsenable(rr_dir, mg_domain_min_width_v)
998 && bbx.coarsenable(rr_dir, mg_box_min_width_v);
999 if (!is_coarsenable[idim] && do_semicoarsening_level
1000 && info.semicoarsening_direction == -1)
1002 is_coarsenable[idim] =
true;
1009 if (do_semicoarsening_level && info.semicoarsening_direction == -1) {
1011 int n_ones =
AMREX_D_TERM(
static_cast<int>(rr_level[0] == 1),
1012 +
static_cast<int>(rr_level[1] == 1),
1013 +
static_cast<int>(rr_level[2] == 1));
1014 if (n_ones > 1) {
break; }
1016 if (rr_level != mg_coarsen_ratio_v) {
1020 accum_coarsen_ratio.push_back(accum_coarsen_ratio.back()*rr_level);
1021 domainboxes.push_back(dbx.coarsen(rr_level));
1022 boundboxes.push_back(bbx.coarsen(rr_level));
1023 bool to_agg = (bbx.d_numPts() / nbxs) < 0.999*threshold_npts;
1024 agg_flag.push_back(to_agg);
1027 for (
int lev = 1, nlevs =
static_cast<int>(domainboxes.size()); lev < nlevs; ++lev) {
1028 if (!agged && !agg_flag[lev] &&
1029 a_grids[0].coarsenable(accum_coarsen_ratio[lev], mg_box_min_width_v))
1031 m_grids[0].push_back(
amrex::coarsen(a_grids[0], accum_coarsen_ratio[lev]));
1032 m_dmap[0].push_back(a_dmap[0]);
1034 IntVect cr = domainboxes[lev-1].length() / domainboxes[lev].length();
1035 if (!m_grids[0].back().coarsenable(cr)) {
1038 m_grids[0].emplace_back(boundboxes[lev]);
1039 IntVect max_grid_size(info.agg_grid_size);
1040 if (info.do_semicoarsening && info.max_semicoarsening_level >= lev
1041 && info.semicoarsening_direction != -1)
1044 AMREX_D_TERM(
int mgs_0 = (max_grid_size[0]+blen[0]-1) / blen[0];,
1045 int mgs_1 = (max_grid_size[1]+blen[1]-1) / blen[1];,
1046 int mgs_2 = (max_grid_size[2]+blen[2]-1) / blen[2]);
1047 max_grid_size[info.semicoarsening_direction]
1050 m_grids[0].back().maxSize(max_grid_size);
1051 m_dmap[0].push_back(DistributionMapping());
1057 m_geom[0].emplace_back(domainboxes[lev],rb,coord,is_per);
1062 Long consolidation_threshold = 0;
1063 Real avg_npts = 0.0;
1064 if (info.do_consolidation) {
1067 *info.con_grid_size,
1068 *info.con_grid_size);
1071 Box const& dom0 = a_geom[0].Domain();
1074 for (
int lev = 0; lev < info.max_coarsening_level; ++lev)
1076 IntVect rr_level = mg_coarsen_ratio_v;
1077 bool do_semicoarsening_level = info.do_semicoarsening
1078 && numsclevs < info.max_semicoarsening_level;
1079 if (do_semicoarsening_level
1080 && info.semicoarsening_direction != -1)
1082 rr_level[info.semicoarsening_direction] = 1;
1085 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
1087 rr_dir[idim] = rr_vec[idim] * rr_level[idim];
1088 is_coarsenable[idim] = dom0.coarsenable(rr_dir, mg_domain_min_width_v)
1089 && a_grids[0].coarsenable(rr_dir, mg_box_min_width_v);
1090 if (!is_coarsenable[idim] && do_semicoarsening_level
1091 && info.semicoarsening_direction == -1)
1093 is_coarsenable[idim] =
true;
1100 if (do_semicoarsening_level && info.semicoarsening_direction == -1) {
1102 int n_ones =
AMREX_D_TERM(
static_cast<int>(rr_level[0] == 1),
1103 +
static_cast<int>(rr_level[1] == 1),
1104 +
static_cast<int>(rr_level[2] == 1));
1105 if (n_ones > 1) {
break; }
1107 if (rr_level != mg_coarsen_ratio_v) {
1112 m_geom[0].emplace_back(
amrex::coarsen(dom0, rr_vec), rb, coord, is_per);
1115 if (info.do_consolidation)
1117 if (avg_npts/
static_cast<Real>(
AMREX_D_TERM(rr_vec[0], *rr_vec[1], *rr_vec[2]))
1118 <
Real(0.999)*
static_cast<Real>(consolidation_threshold))
1121 con_lev = m_dmap[0].size();
1122 m_dmap[0].push_back(DistributionMapping());
1126 m_dmap[0].push_back(m_dmap[0].back());
1131 m_dmap[0].push_back(a_dmap[0]);
1136 m_num_mg_levels[0] = m_grids[0].size();
1138 for (
int mglev = 0; mglev < m_num_mg_levels[0] - 1; mglev++){
1139 const Box& fine_domain = m_geom[0][mglev].Domain();
1140 const Box& crse_domain = m_geom[0][mglev+1].Domain();
1141 mg_coarsen_ratio_vec.push_back(fine_domain.length()/crse_domain.length());
1144 for (
int amrlev = 0; amrlev < m_num_amr_levels; ++amrlev) {
1145 if (AMRRefRatio(amrlev) == 4 && mg_coarsen_ratio_vec.empty()) {
1146 mg_coarsen_ratio_vec.push_back(
IntVect(2));
1152 makeAgglomeratedDMap(m_grids[0], m_dmap[0]);
1156 makeConsolidatedDMap(m_grids[0], m_dmap[0], info.con_ratio, info.con_strategy);
1161 m_bottom_comm = makeSubCommunicator(m_dmap[0].back());
1165 m_bottom_comm = m_default_comm;
1168 m_do_agglomeration = agged;
1169 m_do_consolidation = coned;
1173 Print() <<
"MLLinOp::defineGrids(): agglomerated AMR level 0 starting at MG level "
1174 << agg_lev <<
" of " << m_num_mg_levels[0] <<
"\n";
1176 Print() <<
"MLLinOp::defineGrids(): consolidated AMR level 0 starting at MG level "
1177 << con_lev <<
" of " << m_num_mg_levels[0]
1178 <<
" (ratio = " << info.con_ratio <<
")" <<
"\n";
1180 Print() <<
"MLLinOp::defineGrids(): no agglomeration or consolidation of AMR level 0\n";
1184 for (
int amrlev = 0; amrlev < m_num_amr_levels; ++amrlev)
1186 for (
int mglev = 1; mglev < m_num_mg_levels[amrlev]; ++mglev)
1188 m_factory[amrlev].emplace_back(makeFactory(amrlev,mglev));
1192 for (
int amrlev = 1; amrlev < m_num_amr_levels; ++amrlev)
1195 "MLLinOp: grids not coarsenable between AMR levels");
1199template <
typename MF>
1201MLLinOpT<MF>::defineBC ()
1203 m_needs_coarse_data_for_bc = !m_domain_covered[0];
1205 levelbc_raii.resize(m_num_amr_levels);
1206 robin_a_raii.resize(m_num_amr_levels);
1207 robin_b_raii.resize(m_num_amr_levels);
1208 robin_f_raii.resize(m_num_amr_levels);
1211template <
typename MF>
1216 const int ncomp = getNComp();
1221template <
typename MF>
1226 const int ncomp = getNComp();
1228 "MLLinOp::setDomainBC: wrong size");
1231 m_lobc_orig = m_lobc;
1232 m_hibc_orig = m_hibc;
1233 for (
int icomp = 0; icomp < ncomp; ++icomp) {
1234 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
1235 if (m_geom[0][0].isPeriodic(idim)) {
1237 m_hibc[icomp][idim] == BCType::Periodic);
1240 m_hibc[icomp][idim] != BCType::Periodic);
1257 if (hasHiddenDimension()) {
1258 const int hd = hiddenDirection();
1259 for (
int n = 0; n < ncomp; ++n) {
1265 if (hasInhomogNeumannBC() && !supportInhomogNeumannBC()) {
1266 amrex::Abort(
"Inhomogeneous Neumann BC not supported");
1268 if (hasRobinBC() && !supportRobinBC()) {
1273template <
typename MF>
1277 int ncomp = m_lobc_orig.size();
1278 for (
int n = 0; n < ncomp; ++n) {
1279 for (
int idim = 0; idim <AMREX_SPACEDIM; ++idim) {
1280 if (m_lobc_orig[n][idim] == bct || m_hibc_orig[n][idim] == bct) {
1288template <
typename MF>
1292 return hasBC(BCType::inhomogNeumann);
1295template <
typename MF>
1299 return hasBC(BCType::Robin);
1302template <
typename MF>
1306#if (AMREX_SPACEDIM == 3)
1307 if (info.hasHiddenDimension()) {
1309 const auto& hi = b.bigEnd();
1310 if (info.hidden_direction == 0) {
1312 }
else if (info.hidden_direction == 1) {
1324template <
typename MF>
1332 for (
int i = 1, N=
static_cast<int>(ba.
size()); i < N; ++i)
1342 for (
int iproc = 0; iproc < nprocs; ++iproc) {
1344 for (
int ibox : sfc[iproc]) {
1348 dm[i].define(std::move(pmap));
1353template <
typename MF>
1355MLLinOpT<MF>::makeConsolidatedDMap (
const Vector<BoxArray>& ba,
1356 Vector<DistributionMapping>& dm,
1357 int ratio,
int strategy)
1359 BL_PROFILE(
"MLLinOp::makeConsolidatedDMap()");
1363 for (
int i = 1, N=
static_cast<int>(ba.size()); i < N; ++i)
1370 const auto& pmap_fine = dm[i-1].ProcessorMap();
1371 Vector<int> pmap(pmap_fine.size());
1373 if (strategy == 1) {
1374 for (
auto&
x: pmap) {
1377 }
else if (strategy == 2) {
1378 int nprocs_con =
static_cast<int>(std::ceil(
static_cast<Real>(nprocs)
1379 /
static_cast<Real>(factor)));
1380 for (
auto&
x: pmap) {
1381 auto d = std::div(
x,nprocs_con);
1384 }
else if (strategy == 3) {
1385 if (factor == ratio) {
1387 for (
int iproc = 0; iproc < nprocs; ++iproc) {
1388 for (
int ibox : sfc[iproc]) {
1393 for (
auto&
x: pmap) {
1399 dm[i].define(std::move(pmap));
1401 Vector<int> pmap_g(pmap.size());
1403 dm[i].define(std::move(pmap_g));
1409template <
typename MF>
1411MLLinOpT<MF>::makeSubCommunicator (
const DistributionMapping& dm)
1413 BL_PROFILE(
"MLLinOp::makeSubCommunicator()");
1417 Vector<int> newgrp_ranks = dm.ProcessorMap();
1418 std::sort(newgrp_ranks.begin(), newgrp_ranks.end());
1419 auto last = std::unique(newgrp_ranks.begin(), newgrp_ranks.end());
1420 newgrp_ranks.erase(last, newgrp_ranks.end());
1424 MPI_Comm_group(m_default_comm, &defgrp);
1426 MPI_Group_incl(defgrp,
static_cast<int>(newgrp_ranks.size()), newgrp_ranks.data(), &newgrp);
1428 Vector<int> local_newgrp_ranks(newgrp_ranks.size());
1430 newgrp_ranks.data(),
static_cast<int>(newgrp_ranks.size()));
1431 MPI_Group_incl(defgrp,
static_cast<int>(local_newgrp_ranks.size()), local_newgrp_ranks.data(), &newgrp);
1434 MPI_Comm_create(m_default_comm, newgrp, &newcomm);
1436 m_raii_comm = std::make_unique<CommContainer>(newcomm);
1438 MPI_Group_free(&defgrp);
1439 MPI_Group_free(&newgrp);
1444 return m_default_comm;
1448template <
typename MF>
1453 m_domain_bloc_lo = lo_bcloc;
1454 m_domain_bloc_hi = hi_bcloc;
1457template <
typename MF>
1462 setCoarseFineBC(
crse,
IntVect(crse_ratio), bc_type);
1465template <
typename MF>
1470 m_coarse_data_for_bc =
crse;
1471 m_coarse_data_crse_ratio = crse_ratio;
1472 m_coarse_fine_bc_type = bc_type;
1475template <
typename MF>
1476template <
typename AMF>
1477requires (!std::same_as<MF,AMF>)
1482 setCoarseFineBC(
crse,
IntVect(crse_ratio), bc_type);
1485template <
typename MF>
1486template <
typename AMF>
1487requires (!std::same_as<MF,AMF>)
1492 m_coarse_data_for_bc_raii = MF(
crse->boxArray(),
crse->DistributionMap(),
1494 m_coarse_data_for_bc_raii.LocalCopy(*
crse, 0, 0,
crse->nComp(),
1496 m_coarse_data_for_bc = &m_coarse_data_for_bc_raii;
1497 m_coarse_data_crse_ratio = crse_ratio;
1498 m_coarse_fine_bc_type = bc_type;
1501template <
typename MF>
1506 mf.
resize(m_num_amr_levels);
1507 for (
int alev = 0; alev < m_num_amr_levels; ++alev) {
1508 mf[alev].resize(m_num_mg_levels[alev]);
1509 for (
int mlev = 0; mlev < m_num_mg_levels[alev]; ++mlev) {
1510 mf[alev][mlev] = make(alev, mlev, ng);
1515template <
typename MF>
1519 if constexpr (IsMultiFabLike_v<MF>) {
1521 m_dmap[amrlev][mglev], getNComp(), ng,
MFInfo(),
1522 *m_factory[amrlev][mglev]);
1530template <
typename MF>
1534 if constexpr (IsMultiFabLike_v<MF>) {
1538 amrex::Abort(
"MLLinOpT::makeAlias: how did we get here?");
1543template <
typename MF>
1547 if constexpr (IsMultiFabLike_v<MF>) {
1548 BoxArray cba = m_grids[amrlev][mglev];
1549 IntVect ratio = (amrlev > 0) ?
IntVect(2) : mg_coarsen_ratio_vec[mglev];
1552 return MF(cba, m_dmap[amrlev][mglev], getNComp(), ng);
1555 amrex::Abort(
"MLLinOpT::makeCoarseMG: how did we get here?");
1560template <
typename MF>
1564 if constexpr (IsMultiFabLike_v<MF>) {
1565 BoxArray cba = m_grids[famrlev][0];
1566 IntVect ratio(AMRRefRatioVect(famrlev-1));
1569 return MF(cba, m_dmap[famrlev][0], getNComp(), ng);
1572 amrex::Abort(
"MLLinOpT::makeCoarseAmr: how did we get here?");
1577template <
typename MF>
1581 if (new_size <= 0 || new_size >= m_num_mg_levels[0]) {
return; }
1583 m_num_mg_levels[0] = new_size;
1585 m_geom[0].resize(new_size);
1586 m_grids[0].resize(new_size);
1587 m_dmap[0].resize(new_size);
1588 m_factory[0].resize(new_size);
1590 if (m_bottom_comm != m_default_comm) {
1591 m_bottom_comm = makeSubCommunicator(m_dmap[0].back());
1595template <
typename MF>
1601 const int ncomp = this->getNComp();
1603 if (!fres.isAllRegular()) {
1604 if constexpr (std::is_same<MF,MultiFab>()) {
1606 mg_coarsen_ratio_vec[clev-1]);
1608 amrex::Abort(
"EB_average_down only works with MultiFab");
1616 amrex::Abort(
"For non-FabArray, MLLinOpT<MF>::avgDownResMG should be overridden.");
1620template <
typename MF>
1624 return m_dmap[amrlev][mglev1] == m_dmap[amrlev][mglev2]
1628template <
typename MF>
1629template <MultiFabLike AMF>
1630requires (!std::same_as<MF,AMF>)
1633 const AMF* robinbc_a,
const AMF* robinbc_b,
1634 const AMF* robinbc_f)
1636 const int ncomp = this->getNComp();
1638 levelbc_raii[amrlev] = std::make_unique<MF>(levelbcdata->boxArray(),
1639 levelbcdata->DistributionMap(),
1640 ncomp, levelbcdata->nGrowVect());
1641 levelbc_raii[amrlev]->LocalCopy(*levelbcdata, 0, 0, ncomp,
1642 levelbcdata->nGrowVect());
1644 levelbc_raii[amrlev].reset();
1648 robin_a_raii[amrlev] = std::make_unique<MF>(robinbc_a->boxArray(),
1649 robinbc_a->DistributionMap(),
1650 ncomp, robinbc_a->nGrowVect());
1651 robin_a_raii[amrlev]->LocalCopy(*robinbc_a, 0, 0, ncomp,
1652 robinbc_a->nGrowVect());
1654 robin_a_raii[amrlev].reset();
1658 robin_b_raii[amrlev] = std::make_unique<MF>(robinbc_b->boxArray(),
1659 robinbc_b->DistributionMap(),
1660 ncomp, robinbc_b->nGrowVect());
1661 robin_b_raii[amrlev]->LocalCopy(*robinbc_b, 0, 0, ncomp,
1662 robinbc_b->nGrowVect());
1664 robin_b_raii[amrlev].reset();
1668 robin_f_raii[amrlev] = std::make_unique<MF>(robinbc_f->boxArray(),
1669 robinbc_f->DistributionMap(),
1670 ncomp, robinbc_f->nGrowVect());
1671 robin_f_raii[amrlev]->LocalCopy(*robinbc_f, 0, 0, ncomp,
1672 robinbc_f->nGrowVect());
1674 robin_f_raii[amrlev].reset();
1677 this->setLevelBC(amrlev, levelbc_raii[amrlev].
get(), robin_a_raii[amrlev].
get(),
1678 robin_b_raii[amrlev].
get(), robin_f_raii[amrlev].
get());
1681template <
typename MF>
1686 return xdoty(0,0,*
x[0],*
y[0],
false);
1689template <
typename MF>
1694 auto r = xdoty(0,0,*
x[0],*
x[0],
false);
1695 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
Infrastructure for storing per-face boundary data in FabSets.
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:172
#define AMREX_D_PICK(a, b, c)
Definition AMReX_SPACE.H:173
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:564
static bool SameRefs(const BoxArray &lhs, const BoxArray &rhs)
whether two BoxArrays share the same data
Definition AMReX_BoxArray.H:837
BoxArray & coarsen(int refinement_ratio)
Coarsen each Box in the BoxArray to the specified ratio.
Definition AMReX_BoxArray.cpp:672
BoxArray & convert(IndexType typ)
Apply Box::convert(IndexType) to each Box in the BoxArray.
Definition AMReX_BoxArray.cpp:813
__host__ __device__ BoxND & makeSlab(int direction, int slab_index) noexcept
Definition AMReX_Box.H:834
__host__ __device__ const IntVectND< dim > & smallEnd() const &noexcept
Return the inclusive lower bound of the box.
Definition AMReX_Box.H:112
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
static DistributionMapping makeSFC(const MultiFab &weight, bool sort=true)
Definition AMReX_DistributionMapping.cpp:1766
Definition AMReX_EBFabFactory.H:32
Definition AMReX_FabFactory.H:50
Solve using GMRES with multigrid as preconditioner.
Definition AMReX_GMRES_MLMG.H:28
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:75
Interface
HYPRE interface modes supported.
Definition AMReX_Hypre.H:37
__host__ static __device__ constexpr std::size_t size() noexcept
Definition AMReX_IntVect.H:824
__host__ __device__ constexpr 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:867
Definition AMReX_MLABecLaplacian.H:15
Definition AMReX_MLCGSolver.H:12
Definition AMReX_MLLinOp.H:103
const MF * m_coarse_data_for_bc
Definition AMReX_MLLinOp.H:651
Vector< Vector< std::unique_ptr< FabFactory< FAB > > > > m_factory
Definition AMReX_MLLinOp.H:622
virtual void avgDownResMG(int clev, MF &cres, MF const &fres) const
Definition AMReX_MLLinOp.H:1597
int NAMRLevels() const noexcept
Return the number of AMR levels.
Definition AMReX_MLLinOp.H:575
bool m_do_consolidation
Definition AMReX_MLLinOp.H:613
bool isCellCentered() const noexcept
Definition AMReX_MLLinOp.H:705
IntVect m_ixtype
Definition AMReX_MLLinOp.H:610
void setVerbose(int v) noexcept
Set verbosity.
Definition AMReX_MLLinOp.H:253
bool isMFIterSafe(int amrlev, int mglev1, int mglev2) const
Definition AMReX_MLLinOp.H:1622
RealVect m_coarse_bc_loc
Definition AMReX_MLLinOp.H:650
virtual bool needsUpdate() const
Does it need update if it's reused?
Definition AMReX_MLLinOp.H:275
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:1517
FabFactory< FAB > const * Factory(int amr_lev, int mglev=0) const noexcept
Definition AMReX_MLLinOp.H:669
void setDomainBC(const Vector< Array< BCType, 3 > > &lobc, const Vector< Array< BCType, 3 > > &hibc)
Boundary of the whole domain.
Definition AMReX_MLLinOp.H:1223
Array< Real, 3 > m_domain_bloc_hi
Definition AMReX_MLLinOp.H:645
T get_d0(T const &dx, T const &dy, T const &) const noexcept
Definition AMReX_MLLinOp.H:746
MPI_Comm BottomCommunicator() const noexcept
Definition AMReX_MLLinOp.H:696
void setEnforceSingularSolvable(bool o) noexcept
Definition AMReX_MLLinOp.H:262
MPI_Comm Communicator() const noexcept
Definition AMReX_MLLinOp.H:697
int mg_domain_min_width
Definition AMReX_MLLinOp.H:594
void setMaxOrder(int o) noexcept
Set order of interpolation at coarse/fine boundary.
Definition AMReX_MLLinOp.H:256
virtual void compGrad(int amrlev, const Array< MF *, 3 > &grad, MF &sol, Location loc) const
Compute gradients of the solution.
Definition AMReX_MLLinOp.H:448
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:308
virtual std::string name() const
Definition AMReX_MLLinOp.H:136
GpuArray< BCType, 3 > LoBC(int icomp=0) const noexcept
Definition AMReX_MLLinOp.H:673
bool doAgglomeration() const noexcept
Definition AMReX_MLLinOp.H:701
MF m_coarse_data_for_bc_raii
Definition AMReX_MLLinOp.H:652
virtual void setDirichletNodesToZero(int, int, MF &) const
Definition AMReX_MLLinOp.H:488
MLLinOpT< MF > & operator=(const MLLinOpT< MF > &)=delete
std::unique_ptr< CommContainer > m_raii_comm
Definition AMReX_MLLinOp.H:642
bool m_do_semicoarsening
Definition AMReX_MLLinOp.H:615
bool hasRobinBC() const noexcept
Definition AMReX_MLLinOp.H:1297
Vector< Array< BCType, 3 > > m_hibc
Definition AMReX_MLLinOp.H:584
virtual void resizeMultiGrid(int new_size)
Definition AMReX_MLLinOp.H:1579
Vector< Vector< BoxArray > > m_grids
Definition AMReX_MLLinOp.H:620
virtual MF makeCoarseAmr(int famrlev, IntVect const &ng) const
Definition AMReX_MLLinOp.H:1562
bool m_do_agglomeration
Definition AMReX_MLLinOp.H:612
virtual MF makeAlias(MF const &mf) const
Definition AMReX_MLLinOp.H:1532
Array4< T > compactify(Array4< T > const &a) const noexcept
Definition AMReX_MLLinOp.H:732
static constexpr int mg_coarsen_ratio
Definition AMReX_MLLinOp.H:592
virtual void copyNSolveSolution(MF &, MF const &) const
Definition AMReX_MLLinOp.H:552
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:258
virtual int getNComp() const
Return number of components.
Definition AMReX_MLLinOp.H:270
virtual void getFluxes(const Vector< MF * > &, const Vector< MF * > &) const
Definition AMReX_MLLinOp.H:517
void setCoarseFineBCLocation(const RealVect &cloc) noexcept
Definition AMReX_MLLinOp.H:699
Vector< int > m_amr_ref_ratio
Definition AMReX_MLLinOp.H:605
MPI_Comm m_default_comm
Definition AMReX_MLLinOp.H:625
bool isBottomActive() const noexcept
Definition AMReX_MLLinOp.H:694
virtual Vector< RT > getSolvabilityOffset(int, int, MF const &) const
get offset for fixing solvability
Definition AMReX_MLLinOp.H:475
virtual void getFluxes(const Vector< Array< MF *, 3 > > &, const Vector< MF * > &, Location) const
Definition AMReX_MLLinOp.H:512
virtual BottomSolver getDefaultBottomSolver() const
Definition AMReX_MLLinOp.H:267
typename FabDataType< MF >::fab_type FAB
Definition AMReX_MLLinOp.H:113
virtual RT normInf(int amrlev, MF const &mf, bool local) const =0
Vector< int > m_num_mg_levels
Definition AMReX_MLLinOp.H:607
bool hasBC(BCType bct) const noexcept
Definition AMReX_MLLinOp.H:1275
Vector< Vector< DistributionMapping > > m_dmap
Definition AMReX_MLLinOp.H:621
int verbose
Definition AMReX_MLLinOp.H:598
IntVect m_coarse_data_crse_ratio
Definition AMReX_MLLinOp.H:649
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)
const Vector< int > & AMRRefRatio() const noexcept
Return AMR refinement ratios.
Definition AMReX_MLLinOp.H:657
MLLinOpT(MLLinOpT< MF > &&)=delete
Vector< Array< BCType, 3 > > m_hibc_orig
Definition AMReX_MLLinOp.H:588
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:1459
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)
bool needsCoarseDataForBC() const noexcept
Needs coarse data for bc?
Definition AMReX_MLLinOp.H:182
typename FabDataType< MF >::value_type RT
Definition AMReX_MLLinOp.H:114
virtual void update()
Update for reuse.
Definition AMReX_MLLinOp.H:277
Vector< Array< BCType, 3 > > m_lobc_orig
Definition AMReX_MLLinOp.H:587
bool m_precond_mode
Definition AMReX_MLLinOp.H:654
virtual std::unique_ptr< FabFactory< FAB > > makeFactory(int, int) const
Definition AMReX_MLLinOp.H:721
virtual void applyMetricTerm(int, int, MF &) const
apply metric terms if there are any
Definition AMReX_MLLinOp.H:456
virtual void unapplyMetricTerm(int, int, MF &) const
unapply metric terms if there are any
Definition AMReX_MLLinOp.H:458
bool hasHiddenDimension() const noexcept
Definition AMReX_MLLinOp.H:727
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:416
virtual IntVect getNGrowVectRestriction() const
Definition AMReX_MLLinOp.H:707
virtual void compFlux(int amrlev, const Array< MF *, 3 > &fluxes, MF &sol, Location loc) const
Compute fluxes.
Definition AMReX_MLLinOp.H:433
static constexpr int mg_box_min_width
Definition AMReX_MLLinOp.H:593
virtual RT dotProductPrecond(Vector< MF const * > const &x, Vector< MF const * > const &y) const
Definition AMReX_MLLinOp.H:1683
virtual MF makeCoarseMG(int amrlev, int mglev, IntVect const &ng) const
Definition AMReX_MLLinOp.H:1545
virtual void fixSolvabilityByOffset(int, int, MF &, Vector< RT > const &) const
fix solvability by subtracting offset from RHS
Definition AMReX_MLLinOp.H:479
const Geometry & Geom(int amr_lev, int mglev=0) const noexcept
Definition AMReX_MLLinOp.H:580
virtual void endPrecondBC()
Definition AMReX_MLLinOp.H:570
int hiddenDirection() const noexcept
Definition AMReX_MLLinOp.H:728
void setDomainBCLoc(const Array< Real, 3 > &lo_bcloc, const Array< Real, 3 > &hi_bcloc) noexcept
Set location of domain boundaries.
Definition AMReX_MLLinOp.H:1450
Vector< Array< BCType, 3 > > m_lobc
Definition AMReX_MLLinOp.H:583
Vector< int > m_domain_covered
Definition AMReX_MLLinOp.H:623
const MLLinOpT< MF > * m_parent
Definition AMReX_MLLinOp.H:608
bool doSemicoarsening() const noexcept
Definition AMReX_MLLinOp.H:703
virtual bool supportNSolve() const
Definition AMReX_MLLinOp.H:550
virtual bool supportRobinBC() const noexcept
Definition AMReX_MLLinOp.H:688
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:373
virtual void applyInhomogNeumannTerm(int, MF &) const
Extra terms introduced when we treat inhomogeneous Nuemann BC as homogeneous.
Definition AMReX_MLLinOp.H:464
virtual void avgDownResAmr(int clev, MF &cres, MF const &fres) const
Definition AMReX_MLLinOp.H:560
Vector< Vector< Geometry > > m_geom
first Vector is for amr level and second is mg level
Definition AMReX_MLLinOp.H:619
virtual RT norm2Precond(Vector< MF const * > const &x) const
Definition AMReX_MLLinOp.H:1691
MLLinOpT(const MLLinOpT< MF > &)=delete
virtual bool scaleRHS(int, MF *) const
scale RHS to fix solvability
Definition AMReX_MLLinOp.H:470
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:1467
Box compactify(Box const &b) const noexcept
Definition AMReX_MLLinOp.H:1304
bool m_needs_coarse_data_for_bc
Definition AMReX_MLLinOp.H:647
GpuArray< BCType, 3 > HiBC(int icomp=0) const noexcept
Definition AMReX_MLLinOp.H:678
int maxorder
Definition AMReX_MLLinOp.H:600
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:789
Vector< IntVect > mg_coarsen_ratio_vec
Definition AMReX_MLLinOp.H:616
MF MFType
Definition AMReX_MLLinOp.H:112
virtual void preparePrecond()
Definition AMReX_MLLinOp.H:484
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:338
LPInfo info
Definition AMReX_MLLinOp.H:596
int NMGLevels(int amrlev) const noexcept
Return the number of MG levels at given AMR level.
Definition AMReX_MLLinOp.H:578
T get_d1(T const &, T const &dy, T const &dz) const noexcept
Definition AMReX_MLLinOp.H:756
bool enforceSingularSolvable
Definition AMReX_MLLinOp.H:602
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:322
bool doConsolidation() const noexcept
Definition AMReX_MLLinOp.H:702
virtual std::unique_ptr< MLLinOpT< MF > > makeNLinOp(int) const
Definition AMReX_MLLinOp.H:506
virtual bool supportInhomogNeumannBC() const noexcept
Definition AMReX_MLLinOp.H:689
virtual void prepareForFluxes(int, const MF *=nullptr)
Definition AMReX_MLLinOp.H:389
LinOpBCType m_coarse_fine_bc_type
Definition AMReX_MLLinOp.H:648
Array< Real, 3 > m_domain_bloc_lo
Definition AMReX_MLLinOp.H:644
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:554
int AMRRefRatio(int amr_lev) const noexcept
Return AMR refinement ratio at given AMR level.
Definition AMReX_MLLinOp.H:660
MPI_Comm m_bottom_comm
Definition AMReX_MLLinOp.H:626
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:569
virtual void applyOverset(int, MF &) const
for overset solver only
Definition AMReX_MLLinOp.H:467
virtual void getEBFluxes(const Vector< MF * > &, const Vector< MF * > &) const
Definition AMReX_MLLinOp.H:523
virtual void unimposeNeumannBC(int, MF &) const
This is needed for our nodal projection solver.
Definition AMReX_MLLinOp.H:461
bool getEnforceSingularSolvable() const noexcept
Definition AMReX_MLLinOp.H:265
virtual void prepareForSolve()=0
IntVect AMRRefRatioVect(int amr_lev) const noexcept
Return AMR refinement ratio as IntVect (1 in hidden direction)
Definition AMReX_MLLinOp.H:663
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:1503
bool hasInhomogNeumannBC() const noexcept
Definition AMReX_MLLinOp.H:1290
void setDomainBC(const Array< BCType, 3 > &lobc, const Array< BCType, 3 > &hibc) noexcept
Boundary of the whole domain.
Definition AMReX_MLLinOp.H:1213
virtual int getNGrow(int=0, int=0) const
Definition AMReX_MLLinOp.H:272
int m_num_amr_levels
Definition AMReX_MLLinOp.H:604
Definition AMReX_MLMGBndry.H:12
Definition AMReX_MLMG.H:17
Definition AMReX_MLPoisson.H:19
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:29
Long size() const noexcept
Definition AMReX_Vector.H:54
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
amrex_long Long
Definition AMReX_INT.H:30
__host__ __device__ BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) coarsening ratio.
Definition AMReX_Box.H:1418
std::array< T, N > Array
Definition AMReX_Array.H:26
bool notInLaunchRegion() noexcept
Definition AMReX_GpuControl.H:89
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:223
int MPI_Comm
Definition AMReX_ccse-mpi.H:51
int MPI_Group
Definition AMReX_ccse-mpi.H:52
static constexpr int MPI_COMM_NULL
Definition AMReX_ccse-mpi.H:59
Definition AMReX_Amr.cpp:50
@ make_alias
Definition AMReX_MakeType.H:7
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1567
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:359
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
LinOpBCType
Definition AMReX_LO_BCTYPES.H:27
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)
Volume-weighted average-down from fine to coarse using EB volume fractions.
Definition AMReX_EBMultiFabUtil.cpp:336
__host__ __device__ BoxND< dim > enclosedCells(const BoxND< dim > &b, int dir) noexcept
Return a BoxND with CELL based coordinates in direction dir that is enclosed by b....
Definition AMReX_Box.H:1595
BottomSolver
Definition AMReX_MLLinOp.H:32
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
std::unique_ptr< Hypre > makeHypre(const BoxArray &grids, const DistributionMapping &dmap, const Geometry &geom, MPI_Comm comm_, Hypre::Interface interface, const iMultiFab *overset_mask)
Factory that instantiates the requested HYPRE interface.
Definition AMReX_Hypre.cpp:12
void Warning(const std::string &msg)
Print out warning message to cerr.
Definition AMReX.cpp:247
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:241
__host__ __device__ constexpr int get(IntVectND< dim > const &iv) noexcept
Get I'th element of IntVectND<dim>
Definition AMReX_IntVect.H:1334
A multidimensional array accessor.
Definition AMReX_Array4.H:285
Definition AMReX_FabDataType.H:10
Fixed-size array that can be used on GPU.
Definition AMReX_Array.H:43
Definition AMReX_TypeTraits.H:27
Definition AMReX_MLLinOp.H:37
LPInfo & setConsolidationRatio(int x) noexcept
Definition AMReX_MLLinOp.H:57
int con_strategy
Definition AMReX_MLLinOp.H:44
bool do_semicoarsening
Definition AMReX_MLLinOp.H:40
bool has_metric_term
Definition AMReX_MLLinOp.H:45
LPInfo & setConsolidationGridSize(int x) noexcept
Definition AMReX_MLLinOp.H:56
int max_semicoarsening_level
Definition AMReX_MLLinOp.H:47
bool hasHiddenDimension() const noexcept
Definition AMReX_MLLinOp.H:66
int con_ratio
Definition AMReX_MLLinOp.H:43
bool do_consolidation
Definition AMReX_MLLinOp.H:39
int con_grid_size
Definition AMReX_MLLinOp.H:42
LPInfo & setSemicoarsening(bool x) noexcept
Definition AMReX_MLLinOp.H:54
LPInfo & setHiddenDirection(int n) noexcept
Definition AMReX_MLLinOp.H:63
LPInfo & setConsolidation(bool x) noexcept
Definition AMReX_MLLinOp.H:53
LPInfo & setSemicoarseningDirection(int n) noexcept
Definition AMReX_MLLinOp.H:62
LPInfo & setMaxSemicoarseningLevel(int n) noexcept
Definition AMReX_MLLinOp.H:61
bool do_agglomeration
Definition AMReX_MLLinOp.H:38
LPInfo & setMetricTerm(bool x) noexcept
Definition AMReX_MLLinOp.H:59
bool deterministic
Enable deterministic mode for GPU operations.
Definition AMReX_MLLinOp.H:50
static constexpr int getDefaultConsolidationGridSize()
Definition AMReX_MLLinOp.H:78
LPInfo & setConsolidationStrategy(int x) noexcept
Definition AMReX_MLLinOp.H:58
int max_coarsening_level
Definition AMReX_MLLinOp.H:46
int agg_grid_size
Definition AMReX_MLLinOp.H:41
static constexpr int getDefaultAgglomerationGridSize()
Definition AMReX_MLLinOp.H:70
int hidden_direction
Definition AMReX_MLLinOp.H:49
LPInfo & setAgglomerationGridSize(int x) noexcept
Definition AMReX_MLLinOp.H:55
LPInfo & setMaxCoarseningLevel(int n) noexcept
Definition AMReX_MLLinOp.H:60
LPInfo & setDeterministic(bool x) noexcept
Definition AMReX_MLLinOp.H:64
LPInfo & setAgglomeration(bool x) noexcept
Definition AMReX_MLLinOp.H:52
int semicoarsening_direction
Definition AMReX_MLLinOp.H:48
Definition AMReX_MLLinOp.H:88
StateMode
Definition AMReX_MLLinOp.H:90
BCMode
Definition AMReX_MLLinOp.H:89
Location
Definition AMReX_MLLinOp.H:91
FabArray memory allocation information.
Definition AMReX_FabArray.H:68