1#ifndef AMREX_YAFLUXREGISTER_H_
2#define AMREX_YAFLUXREGISTER_H_
3#include <AMReX_Config.H>
8#include <AMReX_YAFluxRegister_K.H>
36 using T =
typename MF::value_type;
37 using FAB =
typename MF::fab_type;
58 const IntVect& ref_ratio,
int fine_lev,
int nvar);
64 const IntVect& ref_ratio,
int fine_lev,
int nvar);
79 const std::array<FAB const*, AMREX_SPACEDIM>& flux,
95 const std::array<FAB const*, AMREX_SPACEDIM>& flux,
96 const Real* dx,
Real dt,
int srccomp,
int destcomp,
97 int numcomp,
RunOn runon)
noexcept;
109 const std::array<FAB const*, AMREX_SPACEDIM>& flux,
125 const std::array<FAB const*, AMREX_SPACEDIM>& a_flux,
126 const Real* dx,
Real dt,
int srccomp,
int destcomp,
127 int numcomp,
RunOn runon)
noexcept;
135 void Reflux (MF& state,
int dc = 0);
138 void Reflux (MF& state,
int srccomp,
int destcomp,
int numcomp);
147 return !(
m_cfp_fab[mfi.LocalIndex()].empty());
192template <
typename MF>
196 const IntVect& ref_ratio,
int fine_lev,
int nvar)
198 define(fba, cba, fdm, cdm, fgeom, cgeom, ref_ratio, fine_lev, nvar);
201template <
typename MF>
206 const IntVect& ref_ratio,
int fine_lev,
int nvar)
211 m_fine_level = fine_lev;
214 m_crse_data.define(cba, cdm, nvar, 0);
216 m_crse_flag.define(cba, cdm, 1, 1);
218 const auto& cperiod = m_crse_geom.periodicity();
219 const std::vector<IntVect>& pshifts = cperiod.shiftIntVect();
224 Box cdomain = m_crse_geom.Domain();
225 for (
int idim=0; idim < AMREX_SPACEDIM; ++idim) {
226 if (m_crse_geom.isPeriodic(idim)) {
227 cdomain.
grow(idim,1);
231 m_crse_fab_flag.
resize(m_crse_flag.local_size(), crse_cell);
233 m_crse_flag.setVal(crse_cell);
237 m_crse_flag.setVal(crse_fine_boundary_cell, cpc1, 0, 1);
239 m_crse_flag.setVal(fine_cell, cpc0, 0, 1);
240 auto recv_layout_mask = m_crse_flag.RecvLayoutMask(cpc0);
242#pragma omp parallel if (Gpu::notInLaunchRegion())
245 if (recv_layout_mask[mfi]) {
246 m_crse_fab_flag[mfi.LocalIndex()] = fine_cell;
255 const auto n_cfba =
static_cast<int>(cfba.
size());
260 const int nthreads = omp_get_max_threads();
267 const int tid = omp_get_thread_num();
272 for (
int i = 0; i < n_cfba; ++i)
278 const auto ntmp =
static_cast<int>(bl_tmp.
size());
282 for (
int j = 0; j < ntmp; ++j) {
286 if (proc == myproc) {
292 for (
auto const& bl : bl_priv) {
296 for (
auto const& pmp : procmap_priv) {
297 cfp_procmap.insert(std::end(cfp_procmap), std::begin(pmp), std::end(pmp));
300 for (
auto& lid : localindex_priv) {
302 for (
int j = 0; j < nl; ++j) {
303 m_cfp_localindex.push_back(nlocal);
312 for (
int i = 0; i < n_cfba; ++i)
318 const auto ntmp =
static_cast<int>(bl_tmp.
size());
322 for (
int j = 0; j < ntmp; ++j) {
323 cfp_procmap.push_back(proc);
326 if (proc == myproc) {
327 for (
int j = 0; j < ntmp; ++j) {
328 m_cfp_localindex.push_back(nlocal);
340 m_cfpatch.define(cfp_ba, cfp_dm, nvar, 0);
342 m_cfp_fab.resize(nlocal);
345 const int li = mfi.LocalIndex();
346 const int flgi = m_cfp_localindex[li];
347 FAB& fab = m_cfpatch[mfi];
348 m_cfp_fab[flgi].push_back(&fab);
351 bool is_periodic = m_fine_geom.isAnyPeriodic();
353 m_cfp_mask.define(cfp_ba, cfp_dm, 1, 0);
354 m_cfp_mask.setVal(
T(1.0));
361 const Box& domainbox = m_crse_geom.Domain();
364#pragma omp parallel if (!run_on_gpu)
367 std::vector< std::pair<int,Box> > isects;
371 const Box& bx = mfi.fabbox();
374 FAB& fab = m_cfp_mask[mfi];
376 auto const& arr = m_cfp_mask.array(mfi);
378 for (
const auto& iv : pshifts)
383 for (
const auto& is : isects)
385 const Box& ibx = is.second - iv;
388 tags.push_back({arr,ibx});
392 fab.template setVal<RunOn::Host>(
T(0.0), ibx);
405 tag.
dfab(i,j,k,n) =
T(0);
411template <
typename MF>
415 m_crse_data.setVal(
T(0.0));
416 m_cfpatch.setVal(
T(0.0));
419template <
typename MF>
422 const std::array<FAB const*, AMREX_SPACEDIM>& flux,
425 BL_ASSERT(m_crse_data.nComp() == flux[0]->nComp());
428 int numcomp = m_crse_data.nComp();
429 CrseAdd(mfi, flux, dx, dt, srccomp, destcomp, numcomp, runon);
432template <
typename MF>
435 const std::array<FAB const*, AMREX_SPACEDIM>& flux,
436 const Real* dx,
Real dt,
int srccomp,
int destcomp,
437 int numcomp,
RunOn runon)
noexcept
439 BL_ASSERT(m_crse_data.nComp() >= destcomp+numcomp &&
440 flux[0]->nComp() >= srccomp+numcomp);
447 if (m_crse_fab_flag[mfi.LocalIndex()] == crse_cell) {
451 const Box& bx = mfi.tilebox();
453 auto dtdy =
static_cast<T>(dt/dx[1]);,
454 auto dtdz =
static_cast<T>(dt/dx[2]););
456 FAB const* fy = flux[1];,
457 FAB const* fz = flux[2];);
463 auto dest_arr = m_crse_data.array(mfi,destcomp);
464 auto const flag = m_crse_flag.const_array(mfi);
472 yafluxreg_crseadd(tbx, dest_arr, flag,
AMREX_D_DECL(fxarr,fyarr,fzarr),
477template <
typename MF>
480 const std::array<FAB const*, AMREX_SPACEDIM>& flux,
483 BL_ASSERT(m_crse_data.nComp() == flux[0]->nComp());
486 int numcomp = m_crse_data.nComp();
487 FineAdd(mfi, flux, dx, dt, srccomp, destcomp, numcomp, runon);
490template <
typename MF>
493 const std::array<FAB const*, AMREX_SPACEDIM>& a_flux,
494 const Real* dx,
Real dt,
int srccomp,
int destcomp,
495 int numcomp,
RunOn runon)
noexcept
497 BL_ASSERT(m_cfpatch.nComp() >= destcomp+numcomp &&
498 a_flux[0]->nComp() >= srccomp+numcomp);
504 const int li = mfi.LocalIndex();
506 if (cfp_fabs.empty()) {
return; }
508 const Box& tbx = mfi.tilebox();
512 const T ratio =
static_cast<T>(
AMREX_D_TERM(m_ratio[0],*m_ratio[1],*m_ratio[2]));
513 std::array<T,AMREX_SPACEDIM> dtdx{{
AMREX_D_DECL(
static_cast<T>(dt/(dx[0]*ratio)),
514 static_cast<T>(dt/(dx[1]*ratio)),
515 static_cast<T>(dt/(dx[2]*ratio)))}};
516 const Dim3 rr = m_ratio.dim3();
519 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
524 int fluxcomp = srccomp;
525 std::array<FAB const*,AMREX_SPACEDIM> flux{{
AMREX_D_DECL(a_flux[0],a_flux[1],a_flux[2])}};
528 std::array<FAB,AMREX_SPACEDIM> ftmp;
531 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
533 ftmp[idim].resize(b,numcomp);
534 ftmp[idim].template setVal<RunOn::Host>(
T(0.0));
535 ftmp[idim].template copy<RunOn::Host>(*a_flux[idim], srccomp, 0, numcomp);
536 flux[idim] = &ftmp[idim];
543 for (
int idim=0; idim < AMREX_SPACEDIM; ++idim)
547 FAB const* f = flux[idim];
548 for (
FAB* cfp : cfp_fabs)
551 const Box& lobx_is = lobx & cfp->box();
555 auto d = cfp->array(destcomp);
556 auto dtdxs = dtdx[idim];
557 int dirside = idim*2+side;
561 yafluxreg_fineadd(tmpbox, d, farr, dtdxs, numcomp, dirside, rr);
566 const Box& hibx_is = hibx & cfp->box();
570 auto d = cfp->array(destcomp);
571 auto dtdxs = dtdx[idim];
572 int dirside = idim*2+side;
576 yafluxreg_fineadd(tmpbox, d, farr, dtdxs, numcomp, dirside, rr);
584template <
typename MF>
590 int numcomp = m_ncomp;
591 Reflux(state, srccomp, destcomp, numcomp);
594template <
typename MF>
602 if (!m_cfp_mask.empty())
605#pragma omp parallel if (Gpu::notInLaunchRegion())
609 const Box& bx = m_cfpatch[mfi].box();
610 auto const maskfab = m_cfp_mask.array(mfi);
611 auto cfptfab = m_cfpatch.array(mfi,srccomp);
614 cfptfab(i,j,k,n) *= maskfab(i,j,k);
619 m_crse_data.ParallelCopy(m_cfpatch, srccomp, srccomp, numcomp,
623 BL_ASSERT(state.nComp() >= destcomp + numcomp);
625 auto const& dst = state.arrays();
626 auto const& src = m_crse_data.const_arrays();
627 auto const& vol = m_cvol->const_arrays();
631 dst[bno](i,j,k,destcomp+n) += src[bno](i,j,k,srccomp+n) / vol[bno](i,j,k);
634 amrex::Add(state, m_crse_data, srccomp, destcomp, numcomp, 0);
638template <
typename MF>
645template <
typename MF>
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:111
#define AMREX_LAUNCH_HOST_DEVICE_LAMBDA_FLAG(where_to_run, box, tbox, block)
Definition AMReX_GpuLaunch.nolint.H:136
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
std::vector< std::pair< int, Box > > intersections(const Box &bx) const
Return intersections of Box and BoxArray.
Definition AMReX_BoxArray.cpp:1186
BoxList complementIn(const Box &b) const
Return box - boxarray.
Definition AMReX_BoxArray.cpp:1314
BoxArray & coarsen(int refinement_ratio)
Coarsen each Box in the BoxArray to the specified ratio.
Definition AMReX_BoxArray.cpp:672
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition AMReX_BoxArray.H:615
void uniqify()
Make ourselves unique.
Definition AMReX_BoxArray.cpp:1617
A class for managing a List of Boxes that share a common IndexType. This class implements operations ...
Definition AMReX_BoxList.H:52
Long size() const noexcept
The number of Boxes in this BoxList.
Definition AMReX_BoxList.H:113
void join(const BoxList &blist)
Join the BoxList to ourselves.
Definition AMReX_BoxList.cpp:71
__host__ __device__ BoxND & grow(int i) noexcept
Definition AMReX_Box.H:641
__host__ __device__ BoxND< new_dim > resize() const noexcept
Return a new BoxND of size new_dim by either shrinking or expanding this BoxND.
Definition AMReX_Box.H:893
__host__ __device__ bool cellCentered() const noexcept
Return true if BoxND is cell-centered in all indexing directions.
Definition AMReX_Box.H:327
__host__ __device__ bool contains(const IntVectND< dim > &p) const noexcept
Return true if argument is contained within BoxND.
Definition AMReX_Box.H:212
__host__ __device__ bool ok() const noexcept
Checks if it is a proper BoxND (including a valid type).
Definition AMReX_Box.H:208
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
@ ADD
Definition AMReX_FabArrayBase.H:394
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
__host__ static __device__ 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:679
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:85
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:169
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
Definition AMReX_YAFluxRegister.H:33
MF & getFineData()
Mutable access to fine-side accumulation data.
Definition AMReX_YAFluxRegister.H:640
Vector< Vector< FAB * > > m_cfp_fab
The size of this is (# of local fine grids (# of crse/fine patches for that grid))
Definition AMReX_YAFluxRegister.H:178
Vector< int > m_crse_fab_flag
Definition AMReX_YAFluxRegister.H:174
MF m_cfpatch
This is built on crse/fine patches.
Definition AMReX_YAFluxRegister.H:176
void Reflux(MF &state, int dc=0)
Apply the accumulated flux divergence to the coarse MultiFab state.
Definition AMReX_YAFluxRegister.H:586
Geometry m_fine_geom
Definition AMReX_YAFluxRegister.H:181
IntVect m_ratio
Definition AMReX_YAFluxRegister.H:184
void FineAdd(const MFIter &mfi, const std::array< FAB const *, 3 > &flux, const Real *dx, Real dt, RunOn runon) noexcept
Add fine-level fluxes for the tile identified by mfi.
Definition AMReX_YAFluxRegister.H:479
void define(const BoxArray &fba, const BoxArray &cba, const DistributionMapping &fdm, const DistributionMapping &cdm, const Geometry &fgeom, const Geometry &cgeom, const IntVect &ref_ratio, int fine_lev, int nvar)
Define the register using the same arguments as the constructor.
Definition AMReX_YAFluxRegister.H:203
bool CrseHasWork(const MFIter &mfi) const noexcept
Return true if coarse tile mfi abuts a coarse/fine boundary.
Definition AMReX_YAFluxRegister.H:141
YAFluxRegisterT()=default
Construct an empty register; call define() before use.
MF const * m_cvol
Definition AMReX_YAFluxRegister.H:188
iMultiFab m_crse_flag
Definition AMReX_YAFluxRegister.H:173
void setCrseVolume(MF const *cvol)
Definition AMReX_YAFluxRegister.H:164
MF m_crse_data
Definition AMReX_YAFluxRegister.H:172
Geometry m_crse_geom
Definition AMReX_YAFluxRegister.H:182
int m_fine_level
Definition AMReX_YAFluxRegister.H:185
void setDeterministic(bool flag)
Enable deterministic mode for GPU operations via flag (uses deterministic reductions).
Definition AMReX_YAFluxRegister.H:167
typename MF::fab_type FAB
Definition AMReX_YAFluxRegister.H:37
CellType
Definition AMReX_YAFluxRegister.H:156
@ fine_cell
Definition AMReX_YAFluxRegister.H:158
@ crse_cell
Definition AMReX_YAFluxRegister.H:158
@ crse_fine_boundary_cell
Definition AMReX_YAFluxRegister.H:158
bool FineHasWork(const MFIter &mfi) const noexcept
Return true if fine tile mfi has flux registers to update.
Definition AMReX_YAFluxRegister.H:146
MF & getCrseData()
Mutable access to coarse-side accumulation data.
Definition AMReX_YAFluxRegister.H:647
MF m_cfp_mask
Definition AMReX_YAFluxRegister.H:177
Vector< int > m_cfp_localindex
Definition AMReX_YAFluxRegister.H:179
typename MF::value_type T
Definition AMReX_YAFluxRegister.H:36
bool m_deterministic
Definition AMReX_YAFluxRegister.H:189
int m_ncomp
Definition AMReX_YAFluxRegister.H:186
bool getDeterministic() const
Definition AMReX_YAFluxRegister.H:168
void reset()
Reset stored fluxes and flags before a coarse advance.
Definition AMReX_YAFluxRegister.H:413
void CrseAdd(const MFIter &mfi, const std::array< FAB const *, 3 > &flux, const Real *dx, Real dt, RunOn runon) noexcept
Add coarse-level fluxes for the tile indicated by mfi.
Definition AMReX_YAFluxRegister.H:421
A Collection of IArrayBoxes.
Definition AMReX_iMultiFab.H:34
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
__host__ __device__ BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) coarsening ratio.
Definition AMReX_Box.H:1409
__host__ __device__ BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1280
__host__ __device__ BoxND< dim > refine(const BoxND< dim > &b, int ref_ratio) noexcept
Definition AMReX_Box.H:1459
int MyProc() noexcept
Definition AMReX_ParallelDescriptor.H:128
bool inLaunchRegion() noexcept
Definition AMReX_GpuControl.H:92
Definition AMReX_Amr.cpp:49
__host__ __device__ BoxND< dim > adjCellHi(const BoxND< dim > &b, int dir, int len=1) noexcept
Similar to adjCellLo but builds an adjacent BoxND on the high end.
Definition AMReX_Box.H:1735
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
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:193
__host__ __device__ BoxND< dim > adjCellLo(const BoxND< dim > &b, int dir, int len=1) noexcept
Return the cell centered BoxND of length len adjacent to b on the low end along the coordinate direct...
Definition AMReX_Box.H:1714
__host__ __device__ BoxND< dim > surroundingNodes(const BoxND< dim > &b, int dir) noexcept
Return a BoxND with NODE based coordinates in direction dir that encloses BoxND b....
Definition AMReX_Box.H:1522
RunOn
Definition AMReX_GpuControl.H:69
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
void Add(FabArray< FAB > &dst, FabArray< FAB > const &src, int srccomp, int dstcomp, int numcomp, int nghost)
Definition AMReX_FabArray.H:243
Definition AMReX_TagParallelFor.H:58
Array4< T > dfab
Definition AMReX_TagParallelFor.H:59
A multidimensional array accessor.
Definition AMReX_Array4.H:283
Definition AMReX_Dim3.H:12
parallel copy or add
Definition AMReX_FabArrayBase.H:538
FabArray memory allocation information.
Definition AMReX_FabArray.H:66