1#ifndef AMREX_FILLPATCHER_H_
2#define AMREX_FILLPATCHER_H_
3#include <AMReX_Config.H>
71template <
class MF = MultiFab>
124 template <
typename BC,
130 int scomp,
int dcomp,
int ncomp,
131 BC& cbc,
int cbccomp, BC& fbc,
int fbccomp,
133 PreInterpHook
const& pre_interp = {},
134 PostInterpHook
const& post_interp = {});
155 template <
typename BC,
156 typename PreInterpHook=NullInterpHook<MF>,
157 typename PostInterpHook=NullInterpHook<MF> >
161 int scomp,
int dcomp,
int ncomp,
162 BC& cbc,
int cbccomp,
164 PreInterpHook
const& pre_interp = {},
165 PostInterpHook
const& post_interp = {});
176 template <std::
size_t order>
196 template <
typename BC>
197 void fillRK (
int stage,
int iteration,
int ncycle, MF& mf,
Real time,
215 std::unique_ptr<MF> m_cf_crse_data_tmp;
216 std::unique_ptr<MF> m_cf_fine_data;
217 Real m_dt_coarse = std::numeric_limits<Real>::lowest();
238 m_eb_index_space(eb_index_space),
239 m_sfine(fba, fdm, 1, nghost,
MFInfo().SetAlloc(false))
242 "FillPatcher<MF>: MF must be FabArray type");
245 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
252template <
typename BC,
typename PreInterpHook,
typename PostInterpHook>
257 int scomp,
int dcomp,
int ncomp,
258 BC& cbc,
int cbccomp,
259 BC& fbc,
int fbccomp,
261 PreInterpHook
const& pre_interp,
262 PostInterpHook
const& post_interp)
269 fillCoarseFineBoundary(mf, nghost, time, cmf, ct, scomp, dcomp, ncomp,
270 cbc, cbccomp, bcs, bcscomp, pre_interp, post_interp);
273 m_fgeom, fbc, fbccomp);
282 m_fgeom, m_cgeom, m_eb_index_space);
286template <
typename BC,
typename PreInterpHook,
typename PostInterpHook>
291 int scomp,
int dcomp,
int ncomp,
292 BC& cbc,
int cbccomp,
294 PreInterpHook
const& pre_interp,
295 PostInterpHook
const& post_interp)
300 m_fba == mf.boxArray() &&
301 m_fdm == mf.DistributionMap() &&
302 m_cba == cmf[0]->boxArray() &&
303 m_cdm == cmf[0]->DistributionMap() &&
305 m_ncomp == cmf[0]->nComp());
307 auto const& fpc = getFPinfo();
309 if ( ! fpc.ba_crse_patch.empty())
311 if (m_cf_fine_data ==
nullptr) {
312 m_cf_fine_data = std::make_unique<MF>
313 (detail::make_mf_fine_patch<MF>(fpc, m_ncomp));
316 int ncmfs = cmf.
size();
317 for (
int icmf = 0; icmf < ncmfs; ++icmf) {
319 auto it = std::find_if(m_cf_crse_data.begin(), m_cf_crse_data.end(),
320 [=] (
auto const&
x) {
321 return amrex::almostEqual(x.first,t,5);
324 if (it == std::end(m_cf_crse_data)) {
325 MF mf_crse_patch = detail::make_mf_crse_patch<MF>(fpc, m_ncomp);
326 mf_crse_patch.ParallelCopy(*cmf[icmf], m_cgeom.periodicity());
328 std::pair<Real,std::unique_ptr<MF>> tmp;
330 tmp.second = std::make_unique<MF>(std::move(mf_crse_patch));
331 m_cf_crse_data.push_back(std::move(tmp));
335 if (m_cf_crse_data_tmp ==
nullptr) {
336 m_cf_crse_data_tmp = std::make_unique<MF>
337 (detail::make_mf_crse_patch<MF>(fpc, m_ncomp));
340 int const ng_space_interp = 8;
341 Box domain = m_cgeom.growPeriodicDomain(ng_space_interp);
345 if (m_cf_crse_data.size() == 1) {
347 }
else if (m_cf_crse_data.size() == 2) {
348 Real const teps = std::abs(m_cf_crse_data[1].first -
349 m_cf_crse_data[0].first) * 1.e-3_rt;
350 if (time > m_cf_crse_data[0].first - teps &&
351 time < m_cf_crse_data[0].first + teps) {
353 }
else if (time > m_cf_crse_data[1].first - teps &&
354 time < m_cf_crse_data[1].first + teps) {
361 if (idata == 0 || idata == 1) {
362 auto const& dst = m_cf_crse_data_tmp->arrays();
363 auto const& src = m_cf_crse_data[idata].second->const_arrays();
368 dst[bi](i,j,k,n) = src[bi](i,j,k,n+scomp);
371 }
else if (idata == 2) {
372 Real t0 = m_cf_crse_data[0].first;
373 Real t1 = m_cf_crse_data[1].first;
374 Real alpha = (t1-time)/(t1-t0);
375 Real beta = (time-t0)/(t1-t0);
376 auto const& a = m_cf_crse_data_tmp->arrays();
377 auto const& a0 = m_cf_crse_data[0].second->const_arrays();
378 auto const& a1 = m_cf_crse_data[1].second->const_arrays();
384 = alpha*a0[bi](i,j,k,scomp+n)
385 + beta*a1[bi](i,j,k,scomp+n);
391 amrex::Abort(
"FillPatcher: High order interpolation in time not supported. Or FillPatcher was not properly deleted.");
395 cbc(*m_cf_crse_data_tmp, 0, ncomp, m_cf_crse_data_tmp->nGrowVect(), time, cbccomp);
397 detail::call_interp_hook(pre_interp, *m_cf_crse_data_tmp, 0, ncomp);
400 ncomp,
IntVect(0), m_cgeom, m_fgeom,
402 mf.ixType()),nghost),
403 m_ratio, m_interp, bcs, bcscomp);
405 detail::call_interp_hook(post_interp, *m_cf_fine_data, scomp, ncomp);
407 mf.ParallelCopy(*m_cf_fine_data, scomp, dcomp, ncomp,
IntVect{0}, nghost);
411template <
typename MF>
412template <std::
size_t order>
416 BL_PROFILE(
"FillPatcher::storeRKCoarseData()");
418 m_cf_crse_data.resize(order+1);
420 auto const& fpc = getFPinfo();
422 for (
auto& tmf : m_cf_crse_data) {
423 tmf.first = std::numeric_limits<Real>::lowest();
424 tmf.second = std::make_unique<MF>(detail::make_mf_crse_patch<MF>(fpc, m_ncomp));
426 m_cf_crse_data[0].second->ParallelCopy(S_old, m_cgeom.periodicity());
427 for (std::size_t i = 0; i < order; ++i) {
428 m_cf_crse_data[i+1].second->ParallelCopy(RK_k[i], m_cgeom.periodicity());
432template <
typename MF>
433template <
typename BC>
435 MF& mf,
Real time, BC& cbc, BC& fbc,
439 int rk_order = m_cf_crse_data.size()-1;
440 if (rk_order != 3 && rk_order != 4) {
441 amrex::Abort(
"FillPatcher: unsupported RK order "+std::to_string(rk_order));
446 auto const& fpc = getFPinfo();
447 if (m_cf_crse_data_tmp ==
nullptr) {
448 m_cf_crse_data_tmp = std::make_unique<MF>
449 (detail::make_mf_crse_patch<MF>(fpc, m_ncomp));
452 auto const& u = m_cf_crse_data_tmp->arrays();
453 auto const& u0 = m_cf_crse_data[0].second->const_arrays();
454 auto const& k1 = m_cf_crse_data[1].second->const_arrays();
455 auto const& k2 = m_cf_crse_data[2].second->const_arrays();
456 auto const& k3 = m_cf_crse_data[3].second->const_arrays();
458 Real dtc = m_dt_coarse;
462 int const ng_space_interp = 8;
463 Box cdomain = m_cgeom.growPeriodicDomain(ng_space_interp);
464 cdomain.
convert(m_cf_crse_data_tmp->ixType());
468 Real b1 = xsi -
Real(5./6.)*xsi*xsi;
484 Real kk1 = k1[bi](i,j,k,n);
485 Real kk2 = k2[bi](i,j,k,n);
486 Real kk3 = k3[bi](i,j,k,n);
487 Real uu = b1*kk1 + b2*kk2 + b3*kk3;
488 u[bi](i,j,k,n) = u0[bi](i,j,k,n) + dtc*uu;
491 }
else if (stage == 2) {
496 Real kk1 = k1[bi](i,j,k,n);
497 Real kk2 = k2[bi](i,j,k,n);
498 Real kk3 = k3[bi](i,j,k,n);
499 Real uu = b1*kk1 + b2*kk2 + b3*kk3;
500 Real ut = c1*kk1 + c2*kk2 + c3*kk3;
501 u[bi](i,j,k,n) = u0[bi](i,j,k,n) + dtc*(uu + r*ut);
504 }
else if (stage == 3) {
509 Real kk1 = k1[bi](i,j,k,n);
510 Real kk2 = k2[bi](i,j,k,n);
511 Real kk3 = k3[bi](i,j,k,n);
512 Real uu = b1*kk1 + b2*kk2 + b3*kk3;
513 Real ut = c1*kk1 + c2*kk2 + c3*kk3;
514 Real utt = d1*kk1 + d2*kk2 + d3*kk3;
515 u[bi](i,j,k,n) = u0[bi](i,j,k,n) + dtc*
516 (uu + Real(0.5)*r*ut + Real(0.25)*r*r*utt);
520 }
else if (rk_order == 4) {
521 auto const& k4 = m_cf_crse_data[4].second->const_arrays();
523 Real xsi3 = xsi2*xsi;
549 Real kk1 = k1[bi](i,j,k,n);
550 Real kk2 = k2[bi](i,j,k,n);
551 Real kk3 = k3[bi](i,j,k,n);
552 Real kk4 = k4[bi](i,j,k,n);
553 Real uu = b1*kk1 + b2*kk2 + b3*kk3 + b4*kk4;
554 u[bi](i,j,k,n) = u0[bi](i,j,k,n) + dtc*uu;
557 }
else if (stage == 2) {
562 Real kk1 = k1[bi](i,j,k,n);
563 Real kk2 = k2[bi](i,j,k,n);
564 Real kk3 = k3[bi](i,j,k,n);
565 Real kk4 = k4[bi](i,j,k,n);
566 Real uu = b1*kk1 + b2*kk2 + b3*kk3 + b4*kk4;
567 Real ut = c1*kk1 + c2*kk2 + c3*kk3 + c4*kk4;
568 u[bi](i,j,k,n) = u0[bi](i,j,k,n) + dtc*(uu + Real(0.5)*r*ut);
571 }
else if (stage == 3 || stage == 4) {
574 Real at = (stage == 3) ?
Real(0.5)*r : r;
576 Real attt = (stage == 3) ?
Real(0.0625)*r3 :
Real(0.125)*r3;
582 Real kk1 = k1[bi](i,j,k,n);
583 Real kk2 = k2[bi](i,j,k,n);
584 Real kk3 = k3[bi](i,j,k,n);
585 Real kk4 = k4[bi](i,j,k,n);
586 Real uu = b1*kk1 + b2*kk2 + b3*kk3 + b4*kk4;
587 Real ut = c1*kk1 + c2*kk2 + c3*kk3 + c4*kk4;
588 Real utt = d1*kk1 + d2*kk2 + d3*kk3 + d4*kk4;
589 Real uttt = e1*kk1 + e2*kk2 + e3*kk3 + e4*kk4;
590 u[bi](i,j,k,n) = u0[bi](i,j,k,n) + dtc *
591 (uu + at*ut + att*utt + attt*(uttt+akk*(kk3-kk2)));
598 cbc(*m_cf_crse_data_tmp, 0, m_ncomp, m_cf_crse_data_tmp->nGrowVect(), time, 0);
600 if (m_cf_fine_data ==
nullptr) {
601 m_cf_fine_data = std::make_unique<MF>(detail::make_mf_fine_patch<MF>(fpc, m_ncomp));
605 m_ncomp,
IntVect(0), m_cgeom, m_fgeom,
607 mf.ixType()),m_nghost),
608 m_ratio, m_interp, bcs, 0);
611 mf.ParallelCopy(*m_cf_fine_data, 0, 0, m_ncomp,
IntVect(0), m_nghost);
613 mf.FillBoundary(m_fgeom.periodicity());
614 fbc(mf, 0, m_ncomp, m_nghost, time, 0);
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:567
IndexType ixType() const noexcept
Return index type of this BoxArray.
Definition AMReX_BoxArray.H:857
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:154
__host__ __device__ BoxND & convert(IndexTypeND< dim > typ) noexcept
Convert the BoxND from the current type into the argument type. This may change the BoxND coordinates...
Definition AMReX_Box.H:974
__host__ __device__ bool contains(const IntVectND< dim > &p) const noexcept
Return true if argument is contained within BoxND.
Definition AMReX_Box.H:212
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
Definition AMReX_EB2.H:28
static const FPinfo & TheFPinfo(const FabArrayBase &srcfa, const FabArrayBase &dstfa, const IntVect &dstng, const BoxConverter &coarsener, const Geometry &fgeom, const Geometry &cgeom, const EB2::IndexSpace *)
Definition AMReX_FabArrayBase.cpp:2068
FillPatcher is for filling a fine level MultiFab/FabArray.
Definition AMReX_FillPatcher.H:73
void fillCoarseFineBoundary(MF &mf, IntVect const &nghost, Real time, Vector< MF * > const &cmf, Vector< Real > const &ct, int scomp, int dcomp, int ncomp, BC &cbc, int cbccomp, Vector< BCRec > const &bcs, int bcscomp, PreInterpHook const &pre_interp={}, PostInterpHook const &post_interp={})
Function to fill data at coarse/fine boundary only.
Definition AMReX_FillPatcher.H:288
void fill(MF &mf, IntVect const &nghost, Real time, Vector< MF * > const &cmf, Vector< Real > const &ct, Vector< MF * > const &fmf, Vector< Real > const &ft, int scomp, int dcomp, int ncomp, BC &cbc, int cbccomp, BC &fbc, int fbccomp, Vector< BCRec > const &bcs, int bcscomp, PreInterpHook const &pre_interp={}, PostInterpHook const &post_interp={})
Function to fill data.
Definition AMReX_FillPatcher.H:254
void fillRK(int stage, int iteration, int ncycle, MF &mf, Real time, BC &cbc, BC &fbc, Vector< BCRec > const &bcs)
Fill ghost cells of fine AMR level for RK3 and RK4.
Definition AMReX_FillPatcher.H:434
void storeRKCoarseData(Real time, Real dt, MF const &S_old, Array< MF, order > const &RK_k)
Store coarse AMR level data for RK3 and RK4.
Definition AMReX_FillPatcher.H:413
FillPatcher(BoxArray const &fba, DistributionMapping const &fdm, Geometry const &fgeom, BoxArray const &cba, DistributionMapping const &cdm, Geometry const &cgeom, IntVect const &nghost, int ncomp, InterpBase *interp, EB2::IndexSpace const *eb_index_space=EB2::TopIndexSpaceIfPresent())
Constructor of FillPatcher.
Definition AMReX_FillPatcher.H:223
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:211
__host__ __device__ bool cellCentered() const noexcept
True if the IndexTypeND is CELL based in all directions.
Definition AMReX_IndexType.H:104
__host__ __device__ bool nodeCentered() const noexcept
True if the IndexTypeND is NODE based in all directions.
Definition AMReX_IndexType.H:110
__host__ __device__ constexpr bool allLE(const IntVectND< dim > &rhs) const noexcept
Returns true if this is less than or equal to argument for all components. NOTE: This is NOT a strict...
Definition AMReX_IntVect.H:400
Definition AMReX_InterpBase.H:26
Definition AMReX_InterpBase.H:15
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
Long size() const noexcept
Definition AMReX_Vector.H:53
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
__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
std::array< T, N > Array
Definition AMReX_Array.H:26
const IndexSpace * TopIndexSpaceIfPresent() noexcept
Definition AMReX_EB2.cpp:93
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:263
Definition AMReX_Amr.cpp:49
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1558
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
DistributionMapping const & DistributionMap(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2869
std::enable_if_t< IsFabArray< MF >::value > FillPatchSingleLevel(MF &mf, IntVect const &nghost, Real time, const Vector< MF * > &smf, const Vector< Real > &stime, int scomp, int dcomp, int ncomp, const Geometry &geom, BC &physbcf, int bcfcomp)
FillPatch with data from the current level.
Definition AMReX_FillPatchUtil_I.H:75
void FillPatchInterp(MultiFab &mf_fine_patch, int fcomp, MultiFab const &mf_crse_patch, int ccomp, int ncomp, IntVect const &ng, const Geometry &cgeom, const Geometry &fgeom, Box const &dest_domain, const IntVect &ratio, MFInterpolater *mapper, const Vector< BCRec > &bcs, int bcscomp)
Definition AMReX_FillPatchUtil.cpp:140
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230
BoxArray const & boxArray(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2864
Definition AMReX_FabArrayBase.H:305
Definition AMReX_TypeTraits.H:29
FabArray memory allocation information.
Definition AMReX_FabArray.H:66
Definition AMReX_FillPatchUtil.H:34