1#ifndef AMREX_MultiFabUtil_H_
2#define AMREX_MultiFabUtil_H_
3#include <AMReX_Config.H>
19 const MultiFab& nd,
int scomp,
20 int ncomp,
int ngrow = 0);
29 const Vector<const MultiFab*>& edge,
33 const Vector<const MultiFab*>& fc,
36 template <
typename CMF,
typename FMF,
37 std::enable_if_t<IsFabArray_v<CMF> && IsFabArray_v<FMF>,
int> = 0>
39 const Array<const FMF*,AMREX_SPACEDIM>& fc,
43 const Vector<const MultiFab*>& fc,
44 const Geometry& geom);
47 const Array<const MultiFab*,AMREX_SPACEDIM>& fc,
48 const Geometry& geom);
54 bool use_harmonic_averaging =
false);
60 bool use_harmonic_averaging =
false);
63 template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int> = 0>
65 const Vector<MF*>&
crse,
69 template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int> = 0>
71 const Vector<MF*>&
crse,
75 template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int> = 0>
77 const Array<MF*,AMREX_SPACEDIM>&
crse,
81 template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int> = 0>
83 const Array<MF*,AMREX_SPACEDIM>&
crse,
92 template <
typename FAB>
94 const IntVect& ratio,
int ngcrse=0);
97 template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int> = 0>
99 const Array<MF*,AMREX_SPACEDIM>&
crse,
100 const IntVect& ratio,
const Geometry& crse_geom);
102 template <
typename FAB>
104 const IntVect& ratio,
const Geometry& crse_geom);
108 const Vector<MultiFab*>&
crse,
112 const Array<MultiFab*,AMREX_SPACEDIM>&
crse,
119 const IntVect& ratio,
int ngcrse=0);
122 template <
typename FAB>
124 FabArray<FAB>& S_crse,
127 bool mfiter_is_definitely_safe=
false);
135 void average_down (
const MultiFab& S_fine, MultiFab& S_crse,
136 const Geometry& fgeom,
const Geometry& cgeom,
137 int scomp,
int ncomp,
const IntVect& ratio);
138 void average_down (
const MultiFab& S_fine, MultiFab& S_crse,
139 const Geometry& fgeom,
const Geometry& cgeom,
140 int scomp,
int ncomp,
int rr);
145 template<
typename FAB>
146 void average_down (
const FabArray<FAB>& S_fine, FabArray<FAB>& S_crse,
147 int scomp,
int ncomp,
const IntVect& ratio);
148 template<
typename FAB>
149 void average_down (
const FabArray<FAB>& S_fine, FabArray<FAB>& S_crse,
150 int scomp,
int ncomp,
int rr);
155 int scomp,
int ncomp,
157 const Geometry& cgeom,
const Geometry& fgeom);
164 void writeFabs (
const MultiFab& mf,
const std::string& name);
165 void writeFabs (
const MultiFab& mf,
int comp,
int ncomp,
const std::string& name);
171 const Geometry& geom,
int start_comp,
int ncomp,
172 bool interpolate=
false,
173 RealBox
const& bnd_rbx = RealBox());
182 template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int> FOO = 0>
191 template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int> FOO = 0>
197 template <
typename FAB>
198 iMultiFab
makeFineMask (
const FabArray<FAB>& cmf,
const BoxArray& fba,
const IntVect& ratio,
199 int crse_value = 0,
int fine_value = 1);
200 iMultiFab
makeFineMask (
const BoxArray& cba,
const DistributionMapping& cdm,
201 const BoxArray& fba,
const IntVect& ratio,
202 int crse_value = 0,
int fine_value = 1);
203 template <
typename FAB>
204 iMultiFab
makeFineMask (
const FabArray<FAB>& cmf,
const BoxArray& fba,
const IntVect& ratio,
205 Periodicity
const& period,
int crse_value,
int fine_value);
206 iMultiFab
makeFineMask (
const BoxArray& cba,
const DistributionMapping& cdm,
207 const IntVect& cnghost,
const BoxArray& fba,
const IntVect& ratio,
208 Periodicity
const& period,
int crse_value,
int fine_value);
209 template <
typename FAB>
210 iMultiFab
makeFineMask (
const FabArray<FAB>& cmf,
const FabArray<FAB>& fmf,
212 Periodicity
const& period,
int crse_value,
int fine_value);
213 template <
typename FAB>
214 iMultiFab
makeFineMask (
const FabArray<FAB>& cmf,
const FabArray<FAB>& fmf,
216 Periodicity
const& period,
int crse_value,
int fine_value,
217 LayoutData<int>& has_cf);
219 MultiFab
makeFineMask (
const BoxArray& cba,
const DistributionMapping& cdm,
220 const BoxArray& fba,
const IntVect& ratio,
221 Real crse_value, Real fine_value);
224 void computeDivergence (MultiFab& divu,
const Array<MultiFab const*,AMREX_SPACEDIM>& umac,
225 const Geometry& geom);
228 void computeGradient (MultiFab& grad,
const Array<MultiFab const*,AMREX_SPACEDIM>& umac,
229 const Geometry& geom);
238 Periodicity
const& period);
241 template <
typename T,
typename U>
244 T mf_out(mf_in.boxArray(), mf_in.DistributionMap(), mf_in.nComp(), mf_in.nGrowVect());
247#pragma omp parallel if (Gpu::notInLaunchRegion())
251 const Long n = mfi.fabbox().numPts() * mf_in.nComp();
252 auto *
pdst = mf_out[mfi].dataPtr();
253 auto const* psrc = mf_in [mfi].dataPtr();
256 pdst[i] =
static_cast<typename U::value_type
>(psrc[i]);
322 template <
typename Op,
typename T,
typename FAB,
typename F,
323 std::enable_if_t<IsBaseFab<FAB>::value
324#ifndef AMREX_USE_CUDA
325 && IsCallableR<T,F,int,int,int,int>::value
329 ReduceToPlane (
int direction,
Box const& domain, FabArray<FAB>
const& mf,
F const& f);
375 template <
typename Op,
typename FA,
typename F,
376 std::enable_if_t<IsMultiFabLike_v<FA>
377#ifndef AMREX_USE_CUDA
378 && IsCallableR<
typename FA::value_type,
426 template <
typename Op,
typename FA,
typename F,
427 std::enable_if_t<IsMultiFabLike_v<FA>
428#ifndef AMREX_USE_CUDA
429 && IsCallableR<
typename FA::value_type,
451 Gpu::HostVector<Real>
sumToLine (MultiFab
const& mf,
int icomp,
int ncomp,
452 Box const& domain,
int direction,
bool local =
false);
462 Vector<Geometry>
const& geom,
463 Vector<IntVect>
const& ratio,
493 void FillRandom (MultiFab& mf,
int scomp,
int ncomp);
506 void FillRandomNormal (MultiFab& mf,
int scomp,
int ncomp, Real mean, Real stddev);
519 [[nodiscard]] Vector<MultiFab>
convexify (Vector<MultiFab const*>
const& mf,
520 Vector<IntVect>
const& refinement_ratio);
525template <
typename FAB>
528 int crse_value,
int fine_value)
534template <
typename FAB>
537 Periodicity const& period,
int crse_value,
int fine_value)
540 fba, ratio, period, crse_value, fine_value);
543template <
typename FAB>
547 Periodicity const& period,
int crse_value,
int fine_value)
550 mask.setVal(crse_value);
553 1, 0,
MFInfo().SetAlloc(
false));
555 mask.setVal(fine_value, cpc, 0, 1);
560template <
typename FAB>
564 Periodicity const& period,
int crse_value,
int fine_value,
568 mask.setVal(crse_value);
571 1, 0,
MFInfo().SetAlloc(
false));
573 mask.setVal(fine_value, cpc, 0, 1);
575 has_cf =
mask.RecvLayoutMask(cpc);
582template <
typename FAB>
584 const IntVect& ratio,
int ngcrse,
bool mfiter_is_definitely_safe)
590 int ncomp =
crse.nComp();
591 using value_type =
typename FAB::value_type;
596#pragma omp parallel if (Gpu::notInLaunchRegion())
600 const Box& bx = mfi.growntilebox(ngcrse);
615 crse.ParallelCopy(ctmp,0,0,ncomp,ngcrse,ngcrse);
624template<
typename FAB>
630template<
typename FAB>
632 int scomp,
int ncomp,
const IntVect& ratio)
639 using value_type =
typename FAB::value_type;
652 auto const& crsema = S_crse.
arrays();
654 if (is_cell_centered) {
658 amrex_avgdown(i,j,k,n,crsema[box_no],finema[box_no],scomp,scomp,ratio);
674#pragma omp parallel if (Gpu::notInLaunchRegion())
679 const Box& bx = mfi.tilebox();
683 if (is_cell_centered) {
703 auto const& crsema = crse_S_fine.
arrays();
705 if (is_cell_centered) {
709 amrex_avgdown(i,j,k,n,crsema[box_no],finema[box_no],0,scomp,ratio);
725#pragma omp parallel if (Gpu::notInLaunchRegion())
730 const Box& bx = mfi.tilebox();
738 if (is_cell_centered) {
772 int numcomp,
IntVect nghost,
bool local)
781 auto const& xma =
x.const_arrays();
787 auto const& xfab = xma[box_no];
788 auto const& yfab = yma[box_no];
789 for (
int n = 0; n < numcomp; ++n) {
790 t += f(xfab(i,j,k,xcomp+n) , yfab(i,j,k,ycomp+n));
798#pragma omp parallel if (!system::regtest_reduction) reduction(+:sm)
802 Box const& bx = mfi.growntilebox(nghost);
807 sm += f(xfab(i,j,k,xcomp+n) , yfab(i,j,k,ycomp+n));
827template <
typename MMF,
typename Pred,
typename F>
834 int numcomp,
IntVect nghost,
bool local)
846 auto const& xma =
x.const_arrays();
848 auto const& mma =
mask.const_arrays();
853 if (pf(mma[box_no](i,j,k))) {
854 auto const& xfab = xma[box_no];
855 auto const& yfab = yma[box_no];
856 for (
int n = 0; n < numcomp; ++n) {
857 t += f(xfab(i,j,k,xcomp+n) , yfab(i,j,k,ycomp+n));
866#pragma omp parallel if (!system::regtest_reduction) reduction(+:sm)
870 Box const& bx = mfi.growntilebox(nghost);
873 auto const& mfab =
mask.const_array(mfi);
876 if (pf(mfab(i,j,k))) {
877 sm += f(xfab(i,j,k,xcomp+n) , yfab(i,j,k,ycomp+n));
890template <
typename CMF,
typename FMF,
891 std::enable_if_t<IsFabArray_v<CMF> && IsFabArray_v<FMF>,
int> FOO>
901 auto const& ccma = cc.arrays();
903 auto const& fyma = fc[1]->const_arrays();,
904 auto const& fzma = fc[2]->const_arrays(););
908#if (AMREX_SPACEDIM == 1)
916#
if (AMREX_SPACEDIM == 1)
928#pragma omp parallel if (Gpu::notInLaunchRegion())
932 const Box bx = mfi.growntilebox(ngrow);
933 auto const& ccarr = cc.array(mfi);
934 AMREX_D_TERM(
auto const& fxarr = fc[0]->const_array(mfi);,
935 auto const& fyarr = fc[1]->const_array(mfi);,
936 auto const& fzarr = fc[2]->const_array(mfi););
938#if (AMREX_SPACEDIM == 1)
955template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int>>
958 const IntVect& ratio,
int ngcrse)
968template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int>>
975template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int>>
978 int ratio,
int ngcrse)
983template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int>>
986 const IntVect& ratio,
int ngcrse)
988 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
993template <
typename FAB>
995 const IntVect& ratio,
int ngcrse)
1001 const auto type =
fine.ixType();
1003 for (dir = 0; dir < AMREX_SPACEDIM; ++dir) {
1004 if (type.nodeCentered(dir)) {
break; }
1006 auto tmptype = type;
1008 if (dir >= AMREX_SPACEDIM || !tmptype.cellCentered()) {
1009 amrex::Abort(
"average_down_faces: not face index type");
1011 const int ncomp =
crse.nComp();
1016 auto const& crsema =
crse.arrays();
1017 auto const& finema =
fine.const_arrays();
1030#pragma omp parallel if (Gpu::notInLaunchRegion())
1034 const Box& bx = mfi.growntilebox(ngcrse);
1035 auto const& crsearr =
crse.array(mfi);
1036 auto const& finearr =
fine.const_array(mfi);
1049 crse.ParallelCopy(ctmp,0,0,ncomp,ngcrse,ngcrse);
1053template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int>>
1058 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
1063template <
typename FAB>
1073template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int> FOO>
1076 using T =
typename MF::value_type;
1077 const int ncomp = mf.nComp();
1079 auto* dp = dv.
data();
1081 auto loc = cell.
dim3();
1084 Box const& box = mfi.validbox();
1087 auto const& fab = mf.const_array(mfi);
1090 for (
int n = 0; n < ncomp; ++n) {
1091 dp[n] = fab(loc.x,loc.y,loc.z,n);
1104template <typename MF, std::enable_if_t<IsFabArray<MF>::value,
int> FOO>
1107 bool do_bnd = (!bnd_bx.
isEmpty());
1109 BoxArray const& ba = mf.boxArray();
1111 const auto nboxes =
static_cast<int>(ba.
size());
1117 for (
int i = 0; i < nboxes; ++i) {
1118 Box const&
b = ba[i];
1120 lo[dir] =
b.smallEnd(dir);
1121 if (
b.contains(lo)) {
1123 hi[dir] =
b.bigEnd(dir);
1124 Box b1d(lo,hi,
b.ixType());
1126 procmap.push_back(dm[i]);
1127 index_map.push_back(i);
1131 for (
int i = 0; i < nboxes; ++i) {
1132 Box const&
b = ba[i];
1133 Box const& b1d = bnd_bx &
b;
1136 procmap.push_back(dm[i]);
1137 index_map.push_back(i);
1147 MF rmf(rba, rdm, mf.nComp(),
IntVect(0),
1148 MFInfo().SetArena(mf.arena()));
1150#pragma omp parallel if (Gpu::notInLaunchRegion())
1153 Box const&
b = mfi.validbox();
1154 auto const& dfab = rmf.array(mfi);
1155 auto const& sfab = mf.const_array(index_map[mfi.index()]);
1159 dfab(i,j,k,n) = sfab(i,j,k,n);
1167template <
typename Op,
typename T,
typename F>
1171#if defined(AMREX_USE_GPU)
1176 constexpr int nthreads = 128;
1177 auto nblocks =
static_cast<int>(b2d.
numPts());
1178#ifdef AMREX_USE_SYCL
1180 amrex::launch<nthreads>(nblocks, shared_mem_bytes,
Gpu::gpuStream(),
1183 int bid = h.blockIdx();
1184 int tid = h.threadIdx();
1189 int bid = blockIdx.x;
1190 int tid = threadIdx.x;
1195 if (direction == 0) {
1196 int k = bid / len.y;
1197 int j = bid - k*len.y;
1200 for (
int i = blo.x + tid; i < blo.x+len.x; i += nthreads) {
1201 Op().local_update(tmp, f(box_no,i,j,k));
1204 }
else if (direction == 1) {
1205 int k = bid / len.x;
1206 int i = bid - k*len.x;
1209 for (
int j = blo.y + tid; j < blo.y+len.y; j += nthreads) {
1210 Op().local_update(tmp, f(box_no,i,j,k));
1214 int j = bid / len.x;
1215 int i = bid - j*len.x;
1218 for (
int k = blo.z + tid; k < blo.z+len.z; k += nthreads) {
1219 Op().local_update(tmp, f(box_no,i,j,k));
1223#ifdef AMREX_USE_SYCL
1224 Op().template parallel_update<T>(*p, tmp, h);
1226 Op().template parallel_update<T,nthreads>(*p, tmp);
1231 if (direction == 0) {
1234 Op().local_update(ar(0,j,k), f(box_no,i,j,k));
1236 }
else if (direction == 1) {
1239 Op().local_update(ar(i,0,k), f(box_no,i,j,k));
1244 Op().local_update(ar(i,j,0), f(box_no,i,j,k));
1251template <
typename Op,
typename T,
typename FAB,
typename F,
1252 std::enable_if_t<IsBaseFab<FAB>::value
1253#ifndef AMREX_USE_CUDA
1262 Box domain2d = domain;
1269 r.template setVal<RunOn::Device>(initval);
1270 auto const& ar =
r.array();
1272 for (
MFIter mfi(mf,
MFItInfo().UseDefaultStream().DisableDeviceSync());
1275 Box bx = mfi.validbox() & domain;
1277 int box_no = mfi.LocalIndex();
1278 detail::reduce_to_plane<Op, T>(ar, direction, bx, box_no, f);
1287template <
typename Op,
typename FA,
typename F>
1290 using T =
typename FA::value_type;
1295 if (npts != ndomain.
numPts()) {
1296 amrex::Abort(
"ReduceToPlaneMF: mf's BoxArray must have a rectangular domain.");
1302 BoxList bl = mf.boxArray().boxList();
1303 for (
auto&
b : bl) {
1304 b.setRange(direction, 0);
1307 FA tmpfa(ba, mf.DistributionMap(), 1, 0);
1308 tmpfa.setVal(initval);
1312 Box bx = mfi.validbox() & ndomain;
1314 int box_no = mfi.LocalIndex();
1315 detail::reduce_to_plane<Op, T>(tmpfa.array(mfi), direction, bx, box_no, f);
1323template <
typename Op,
typename FA,
typename F,
1324 std::enable_if_t<IsMultiFabLike_v<FA>
1325#ifndef AMREX_USE_CUDA
1326 && IsCallableR<
typename FA::value_type,
1332 auto [fa3, fa2] = ReduceToPlaneMF2<Op>(direction, domain, mf, f);
1333 fa3.ParallelCopy(fa2);
1334 return std::move(fa3);
1337template <
typename Op,
typename FA,
typename F,
1338 std::enable_if_t<IsMultiFabLike_v<FA>
1339#ifndef AMREX_USE_CUDA
1340 && IsCallableR<
typename FA::value_type,
1347 using T =
typename FA::value_type;
1352 auto tmpmf = detail::reduce_to_plane<Op>(direction, domain, mf, f);
1358 ba2d.
maxSize(mf.boxArray()[0].length());
1361 FA mf2d(ba2d, dm2d, 1, 0);
1362 mf2d.setVal(initval);
1364 static_assert(std::is_same_v<Op, ReduceOpSum>,
"Currently only ReduceOpSum is supported.");
1365 mf2d.ParallelAdd(tmpmf);
1367 return std::make_pair(std::move(tmpmf), std::move(mf2d));
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_HOST_DEVICE_PARALLEL_FOR_1D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:109
#define AMREX_HOST_DEVICE_PARALLEL_FOR_3D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:110
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:111
#define AMREX_LAUNCH_HOST_DEVICE_LAMBDA(...)
Definition AMReX_GpuLaunch.nolint.H:16
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1089
Real * pdst
Definition AMReX_HypreMLABecLap.cpp:1090
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< int const > mask
Definition AMReX_InterpFaceRegister.cpp:93
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
#define AMREX_LOOP_3D(bx, i, j, k, block)
Definition AMReX_Loop.nolint.H:4
#define AMREX_LOOP_4D(bx, ncomp, i, j, k, n, block)
Definition AMReX_Loop.nolint.H:16
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:129
A FortranArrayBox(FAB)-like object.
Definition AMReX_BaseFab.H:183
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:550
IndexType ixType() const noexcept
Return index type of this BoxArray.
Definition AMReX_BoxArray.H:840
BoxArray & coarsen(int refinement_ratio)
Coarsen each Box in the BoxArray to the specified ratio.
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition AMReX_BoxArray.H:597
BoxArray & maxSize(int block_size)
Forces each Box in BoxArray to have sides <= block_size.
A class for managing a List of Boxes that share a common IndexType. This class implements operations ...
Definition AMReX_BoxList.H:52
bool isEmpty() const noexcept
Is this BoxList empty?
Definition AMReX_BoxList.H:135
void push_back(const Box &bn)
Append a Box to this BoxList.
Definition AMReX_BoxList.H:93
AMREX_GPU_HOST_DEVICE IndexTypeND< dim > ixType() const noexcept
Returns the indexing type.
Definition AMReX_Box.H:127
AMREX_GPU_HOST_DEVICE bool isEmpty() const noexcept
Checks if it is an empty BoxND.
Definition AMReX_Box.H:196
AMREX_GPU_HOST_DEVICE bool ok() const noexcept
Checks if it is a proper BoxND (including a valid type).
Definition AMReX_Box.H:200
AMREX_GPU_HOST_DEVICE Long numPts() const noexcept
Returns the number of points contained in the BoxND.
Definition AMReX_Box.H:346
AMREX_GPU_HOST_DEVICE bool contains(const IntVectND< dim > &p) const noexcept
Returns true if argument is contained within BoxND.
Definition AMReX_Box.H:204
AMREX_GPU_HOST_DEVICE BoxND & setRange(int dir, int sm_index, int n_cells=1) noexcept
Set the entire range in a given direction, starting at sm_index with length n_cells....
Definition AMReX_Box.H:1046
Definition AMReX_FabFactory.H:76
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:41
IntVect nGrowVect() const noexcept
Definition AMReX_FabArrayBase.H:79
bool isFusingCandidate() const noexcept
Is this a good candidate for kernel fusing?
bool is_cell_centered() const noexcept
This tests on whether the FabArray is cell-centered.
bool is_nodal() const noexcept
This tests on whether the FabArray is fully nodal.
IndexType ixType() const noexcept
Return index type.
Definition AMReX_FabArrayBase.H:85
const DistributionMapping & DistributionMap() const noexcept
Return constant reference to associated DistributionMapping.
Definition AMReX_FabArrayBase.H:130
int nComp() const noexcept
Return number of variables (aka components) associated with each point.
Definition AMReX_FabArrayBase.H:82
const BoxArray & boxArray() const noexcept
Return a constant reference to the BoxArray that defines the valid region associated with this FabArr...
Definition AMReX_FabArrayBase.H:94
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:344
Array4< typename FabArray< FAB >::value_type const > const_array(const MFIter &mfi) const noexcept
Definition AMReX_FabArray.H:584
void ParallelCopy(const FabArray< FAB > &src, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
Definition AMReX_FabArray.H:840
Array4< typename FabArray< FAB >::value_type const > array(const MFIter &mfi) const noexcept
Definition AMReX_FabArray.H:560
MultiArray4< typename FabArray< FAB >::value_type > arrays() noexcept
Definition AMReX_FabArray.H:632
MultiArray4< typename FabArray< FAB >::value_type const > const_arrays() const noexcept
Definition AMReX_FabArray.H:646
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:73
Periodicity periodicity() const noexcept
Definition AMReX_Geometry.H:355
Definition AMReX_Tuple.H:93
static AMREX_EXPORT constexpr int warp_size
Definition AMReX_GpuDevice.H:173
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool allGE(const IntVectND< dim > &rhs) const noexcept
Returns true if this is greater than or equal to argument for all components. NOTE: This is NOT a str...
Definition AMReX_IntVect.H:441
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 dim3() const noexcept
Definition AMReX_IntVect.H:163
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE constexpr IntVectND< dim > TheUnitVector() noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition AMReX_IntVect.H:680
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE constexpr IntVectND< dim > TheZeroVector() noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition AMReX_IntVect.H:670
a one-thingy-per-box distributed object
Definition AMReX_LayoutData.H:13
Definition AMReX_MFIter.H:57
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:141
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:38
Definition AMReX_PODVector.H:262
iterator begin() noexcept
Definition AMReX_PODVector.H:617
iterator end() noexcept
Definition AMReX_PODVector.H:621
T * data() noexcept
Definition AMReX_PODVector.H:609
This provides length of period for periodic domains. 0 means it is not periodic in that direction....
Definition AMReX_Periodicity.H:17
static const Periodicity & NonPeriodic() noexcept
Definition AMReX_Periodicity.cpp:52
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:27
Definition AMReX_iMultiFab.H:32
void copy(HostToDevice, InIter begin, InIter end, OutIter result) noexcept
A host-to-device copy routine. Note this is just a wrapper around memcpy, so it assumes contiguous st...
Definition AMReX_GpuContainers.H:121
static constexpr DeviceToHost deviceToHost
Definition AMReX_GpuContainers.H:99
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:237
bool inLaunchRegion() noexcept
Definition AMReX_GpuControl.H:86
bool inNoSyncRegion() noexcept
Definition AMReX_GpuControl.H:146
gpuStream_t gpuStream() noexcept
Definition AMReX_GpuDevice.H:218
void Sum(T &v, MPI_Comm comm)
Definition AMReX_ParallelReduce.H:204
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition AMReX_ParallelContext.H:70
void reduce_to_plane(Array4< T > const &ar, int direction, Box const &bx, int box_no, F const &f)
Definition AMReX_MultiFabUtil.H:1168
Definition AMReX_Amr.cpp:49
void FillRandomNormal(MultiFab &mf, int scomp, int ncomp, Real mean, Real stddev)
Fill MultiFab with random numbers from normal distribution.
Definition AMReX_MultiFabUtil.cpp:1212
int nComp(FabArrayBase const &fa)
void FillRandom(MultiFab &mf, int scomp, int ncomp)
Fill MultiFab with random numbers from uniform distribution.
Definition AMReX_MultiFabUtil.cpp:1199
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:191
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 void amrex_avg_fc_to_cc(int i, int, int, Array4< CT > const &cc, Array4< FT const > const &fx, int cccomp, GeometryData const &gd) noexcept
Definition AMReX_MultiFabUtil_1D_C.H:33
std::unique_ptr< MultiFab > get_slice_data(int dir, Real coord, const MultiFab &cc, const Geometry &geom, int start_comp, int ncomp, bool interpolate, RealBox const &bnd_rbx)
Definition AMReX_MultiFabUtil.cpp:552
Vector< typename MF::value_type > get_cell_data(MF const &mf, IntVect const &cell)
Get data in a cell of MultiFab/FabArray.
Definition AMReX_MultiFabUtil.H:1074
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
ReduceData< Ts... >::Type ParReduce(TypeList< Ops... > operation_list, TypeList< Ts... > type_list, FabArray< FAB > const &fa, IntVect const &nghost, F &&f)
Parallel reduce for MultiFab/FabArray.
Definition AMReX_ParReduce.H:47
BaseFab< T > ReduceToPlane(int direction, Box const &domain, FabArray< FAB > const &mf, F const &f)
Reduce FabArray/MultiFab data to a plane Fab.
Definition AMReX_MultiFabUtil.H:1258
std::pair< FA, FA > ReduceToPlaneMF2(int direction, Box const &domain, FA const &mf, F const &f)
Reduce FabArray/MultiFab data to plane FabArray.
Definition AMReX_MultiFabUtil.H:1345
Vector< MultiFab > convexify(Vector< MultiFab const * > const &mf, Vector< IntVect > const &refinement_ratio)
Convexify AMR data.
Definition AMReX_MultiFabUtil.cpp:1225
FA ReduceToPlaneMF(int direction, Box const &domain, FA const &mf, F const &f)
Reduce FabArray/MultiFab data to plane FabArray.
Definition AMReX_MultiFabUtil.H:1330
void average_down_faces(const Vector< const MF * > &fine, const Vector< MF * > &crse, const IntVect &ratio, int ngcrse=0)
Average fine face-based FabArray onto crse face-based FabArray.
Definition AMReX_MultiFabUtil.H:956
MultiFab periodicShift(MultiFab const &mf, IntVect const &offset, Periodicity const &period)
Periodic shift MultiFab.
Definition AMReX_MultiFabUtil.cpp:792
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void amrex_avgdown_faces(Box const &bx, Array4< T > const &crse, Array4< T const > const &fine, int ccomp, int fcomp, int ncomp, IntVect const &ratio, int) noexcept
Definition AMReX_MultiFabUtil_1D_C.H:130
Real volumeWeightedSum(Vector< MultiFab const * > const &mf, int icomp, Vector< Geometry > const &geom, Vector< IntVect > const &ratio, bool local)
Volume weighted sum for a vector of MultiFabs.
Definition AMReX_MultiFabUtil.cpp:958
bool isMFIterSafe(const FabArrayBase &x, const FabArrayBase &y)
Definition AMReX_MFIter.H:219
MultiFab ToMultiFab(const iMultiFab &imf)
Convert iMultiFab to MultiFab.
Definition AMReX_MultiFabUtil.cpp:542
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 lbound(Array4< T > const &a) noexcept
Definition AMReX_Array4.H:312
void writeFabs(const MultiFab &mf, const std::string &name)
Write each fab individually.
Definition AMReX_MultiFabUtil.cpp:529
Real NormHelper(const MultiFab &x, int xcomp, const MultiFab &y, int ycomp, F const &f, int numcomp, IntVect nghost, bool local)
Returns part of a norm based on two MultiFabs.
Definition AMReX_MultiFabUtil.H:769
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 average_down_edges(const Vector< const MultiFab * > &fine, const Vector< MultiFab * > &crse, const IntVect &ratio, int ngcrse)
Average fine edge-based MultiFab onto crse edge-based MultiFab.
Definition AMReX_MultiFabUtil.cpp:447
IntVectND< AMREX_SPACEDIM > IntVect
Definition AMReX_BaseFwd.H:30
Gpu::HostVector< Real > sumToLine(MultiFab const &mf, int icomp, int ncomp, Box const &domain, int direction, bool local)
Sum MultiFab data to line.
Definition AMReX_MultiFabUtil.cpp:816
void average_face_to_cellcenter(MultiFab &cc, int dcomp, const Vector< const MultiFab * > &fc, int ngrow)
Average face-based MultiFab onto cell-centered MultiFab.
Definition AMReX_MultiFabUtil.cpp:141
AMREX_GPU_HOST_DEVICE void cast(BaseFab< Tto > &tofab, BaseFab< Tfrom > const &fromfab, Box const &bx, SrcComp scomp, DestComp dcomp, NumComps ncomp) noexcept
Definition AMReX_BaseFabUtility.H:13
void print_state(const MultiFab &mf, const IntVect &cell, const int n, const IntVect &ng)
Output state data for a single zone.
Definition AMReX_MultiFabUtil.cpp:524
void computeDivergence(MultiFab &divu, const Array< MultiFab const *, AMREX_SPACEDIM > &umac, const Geometry &geom)
Computes divergence of face-data stored in the umac MultiFab.
Definition AMReX_MultiFabUtil.cpp:691
bool TilingIfNotGPU() noexcept
Definition AMReX_MFIter.H:12
void sum_fine_to_coarse(const MultiFab &S_fine, MultiFab &S_crse, int scomp, int ncomp, const IntVect &ratio, const Geometry &cgeom, const Geometry &)
Definition AMReX_MultiFabUtil.cpp:393
void average_cellcenter_to_face(const Vector< MultiFab * > &fc, const MultiFab &cc, const Geometry &geom, int ncomp, bool use_harmonic_averaging)
Average cell-centered MultiFab onto face-based MultiFab with geometric weighting.
Definition AMReX_MultiFabUtil.cpp:218
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 length(Array4< T > const &a) noexcept
Definition AMReX_Array4.H:326
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void amrex_avgdown(Box const &bx, Array4< T > const &crse, Array4< T const > const &fine, int ccomp, int fcomp, int ncomp, IntVect const &ratio) noexcept
Definition AMReX_MultiFabUtil_1D_C.H:195
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230
const int[]
Definition AMReX_BLProfiler.cpp:1664
void average_node_to_cellcenter(MultiFab &cc, int dcomp, const MultiFab &nd, int scomp, int ncomp, int ngrow)
Average nodal-based MultiFab onto cell-centered MultiFab.
Definition AMReX_MultiFabUtil.cpp:62
void FourthOrderInterpFromFineToCoarse(MultiFab &cmf, int scomp, int ncomp, MultiFab const &fmf, IntVect const &ratio)
Fourth-order interpolation from fine to coarse level.
Definition AMReX_MultiFabUtil.cpp:1112
FabArray< BaseFab< Long > > ToLongMultiFab(const iMultiFab &imf)
Convert iMultiFab to Long.
Definition AMReX_MultiFabUtil.cpp:547
MF get_line_data(MF const &mf, int dir, IntVect const &cell, Box const &bnd_bx=Box())
Get data in a line of MultiFab/FabArray.
Definition AMReX_MultiFabUtil.H:1105
void average_edge_to_cellcenter(MultiFab &cc, int dcomp, const Vector< const MultiFab * > &edge, int ngrow)
Average edge-based MultiFab onto cell-centered MultiFab.
Definition AMReX_MultiFabUtil.cpp:97
void computeGradient(MultiFab &grad, const Array< MultiFab const *, AMREX_SPACEDIM > &umac, const Geometry &geom)
Computes gradient of face-data stored in the umac MultiFab.
Definition AMReX_MultiFabUtil.cpp:743
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void amrex_avgdown_nodes(Box const &bx, Array4< T > const &crse, Array4< T const > const &fine, int ccomp, int fcomp, int ncomp, IntVect const &ratio) noexcept
Definition AMReX_MultiFabUtil_1D_C.H:253
iMultiFab makeFineMask(const BoxArray &cba, const DistributionMapping &cdm, const BoxArray &fba, const IntVect &ratio, int crse_value, int fine_value)
Definition AMReX_MultiFabUtil.cpp:607
std::array< T, N > Array
Definition AMReX_Array.H:24
void average_down_nodal(const FabArray< FAB > &S_fine, FabArray< FAB > &S_crse, const IntVect &ratio, int ngcrse=0, bool mfiter_is_definitely_safe=false)
Average fine node-based MultiFab onto crse node-centered MultiFab.
Definition AMReX_MultiFabUtil.H:583
Definition AMReX_FabArrayCommI.H:896
Definition AMReX_Array4.H:61
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE T * ptr(int i, int j, int k) const noexcept
Definition AMReX_Array4.H:149
parallel copy or add
Definition AMReX_FabArrayBase.H:536
Definition AMReX_Geometry.H:25
int coord
Definition AMReX_Geometry.H:59
Definition AMReX_GpuTypes.H:86
Test if a given type T is callable with arguments of type Args...
Definition AMReX_TypeTraits.H:207
FabArray memory allocation information.
Definition AMReX_FabArray.H:66
Definition AMReX_MFIter.H:20
Struct for holding types.
Definition AMReX_TypeList.H:12