1 #ifndef AMREX_FILLPATCHER_H_
2 #define AMREX_FILLPATCHER_H_
3 #include <AMReX_Config.H>
71 template <
class MF = MultiFab>
124 template <
typename BC,
127 void fill (MF& mf,
IntVect const& nghost, Real time,
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,
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>
194 template <
typename BC>
195 void fillRK (
int stage,
int iteration,
int ncycle, MF& mf, Real time,
236 m_eb_index_space(eb_index_space),
237 m_sfine(fba, fdm, 1, nghost,
MFInfo().SetAlloc(false))
240 "FillPatcher<MF>: MF must be FabArray type");
243 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
250 template <
typename BC,
typename PreInterpHook,
typename PostInterpHook>
255 int scomp,
int dcomp,
int ncomp,
256 BC& cbc,
int cbccomp,
257 BC& fbc,
int fbccomp,
259 PreInterpHook
const& pre_interp,
260 PostInterpHook
const& post_interp)
267 fillCoarseFineBoundary(mf, nghost, time, cmf, ct, scomp, dcomp, ncomp,
268 cbc, cbccomp, bcs, bcscomp, pre_interp, post_interp);
271 m_fgeom, fbc, fbccomp);
280 m_fgeom, m_cgeom, m_eb_index_space);
284 template <
typename BC,
typename PreInterpHook,
typename PostInterpHook>
289 int scomp,
int dcomp,
int ncomp,
290 BC& cbc,
int cbccomp,
292 PreInterpHook
const& pre_interp,
293 PostInterpHook
const& post_interp)
298 m_fba == mf.boxArray() &&
299 m_fdm == mf.DistributionMap() &&
300 m_cba == cmf[0]->boxArray() &&
301 m_cdm == cmf[0]->DistributionMap() &&
303 m_ncomp == cmf[0]->nComp());
305 auto const& fpc = getFPinfo();
307 if ( ! fpc.ba_crse_patch.empty())
309 if (m_cf_fine_data ==
nullptr) {
310 m_cf_fine_data = std::make_unique<MF>
311 (detail::make_mf_fine_patch<MF>(fpc, m_ncomp));
314 int ncmfs = cmf.
size();
315 for (
int icmf = 0; icmf < ncmfs; ++icmf) {
317 auto it = std::find_if(m_cf_crse_data.begin(), m_cf_crse_data.end(),
318 [=] (
auto const&
x) {
319 return amrex::almostEqual(x.first,t,5);
322 if (it ==
std::end(m_cf_crse_data)) {
323 MF mf_crse_patch = detail::make_mf_crse_patch<MF>(fpc, m_ncomp);
324 mf_crse_patch.ParallelCopy(*cmf[icmf], m_cgeom.periodicity());
326 std::pair<Real,std::unique_ptr<MF>> tmp;
328 tmp.second = std::make_unique<MF>(std::move(mf_crse_patch));
329 m_cf_crse_data.push_back(std::move(tmp));
333 if (m_cf_crse_data_tmp ==
nullptr) {
334 m_cf_crse_data_tmp = std::make_unique<MF>
335 (detail::make_mf_crse_patch<MF>(fpc, m_ncomp));
338 int const ng_space_interp = 8;
339 Box domain = m_cgeom.growPeriodicDomain(ng_space_interp);
343 if (m_cf_crse_data.size() == 1) {
345 }
else if (m_cf_crse_data.size() == 2) {
346 Real
const teps =
std::abs(m_cf_crse_data[1].first -
347 m_cf_crse_data[0].first) * 1.e-3_rt;
348 if (time > m_cf_crse_data[0].first - teps &&
349 time < m_cf_crse_data[0].first + teps) {
351 }
else if (time > m_cf_crse_data[1].first - teps &&
352 time < m_cf_crse_data[1].first + teps) {
359 if (idata == 0 || idata == 1) {
360 auto const& dst = m_cf_crse_data_tmp->arrays();
361 auto const& src = m_cf_crse_data[idata].second->const_arrays();
366 dst[bi](i,j,k,n) = src[bi](i,j,k,n+scomp);
369 }
else if (idata == 2) {
370 Real t0 = m_cf_crse_data[0].first;
371 Real t1 = m_cf_crse_data[1].first;
372 Real alpha = (t1-time)/(t1-t0);
373 Real beta = (time-t0)/(t1-t0);
374 auto const& a = m_cf_crse_data_tmp->arrays();
375 auto const& a0 = m_cf_crse_data[0].second->const_arrays();
376 auto const& a1 = m_cf_crse_data[1].second->const_arrays();
382 = alpha*a0[bi](i,j,k,scomp+n)
383 + beta*a1[bi](i,j,k,scomp+n);
389 amrex::Abort(
"FillPatcher: High order interpolation in time not supported. Or FillPatcher was not properly deleted.");
393 cbc(*m_cf_crse_data_tmp, 0, ncomp, m_cf_crse_data_tmp->nGrowVect(), time, cbccomp);
398 ncomp,
IntVect(0), m_cgeom, m_fgeom,
400 mf.ixType()),nghost),
401 m_ratio, m_interp, bcs, bcscomp);
405 mf.ParallelCopy(*m_cf_fine_data, scomp, dcomp, ncomp,
IntVect{0}, nghost);
409 template <
typename MF>
410 template <std::
size_t order>
414 BL_PROFILE(
"FillPatcher::storeRKCoarseData()");
416 m_cf_crse_data.resize(order+1);
418 auto const& fpc = getFPinfo();
420 for (
auto& tmf : m_cf_crse_data) {
421 tmf.first = std::numeric_limits<Real>::lowest();
422 tmf.second = std::make_unique<MF>(detail::make_mf_crse_patch<MF>(fpc, m_ncomp));
424 m_cf_crse_data[0].second->ParallelCopy(S_old, m_cgeom.periodicity());
425 for (std::size_t i = 0; i < order; ++i) {
426 m_cf_crse_data[i+1].second->ParallelCopy(RK_k[i], m_cgeom.periodicity());
430 template <
typename MF>
431 template <
typename BC>
433 MF& mf, Real time, BC& cbc, BC& fbc,
437 int rk_order = m_cf_crse_data.size()-1;
438 if (rk_order != 3 && rk_order != 4) {
439 amrex::Abort(
"FillPatcher: unsupported RK order "+std::to_string(rk_order));
444 auto const& fpc = getFPinfo();
445 if (m_cf_crse_data_tmp ==
nullptr) {
446 m_cf_crse_data_tmp = std::make_unique<MF>
447 (detail::make_mf_crse_patch<MF>(fpc, m_ncomp));
450 auto const& u = m_cf_crse_data_tmp->arrays();
451 auto const& u0 = m_cf_crse_data[0].second->const_arrays();
452 auto const& k1 = m_cf_crse_data[1].second->const_arrays();
453 auto const& k2 = m_cf_crse_data[2].second->const_arrays();
454 auto const& k3 = m_cf_crse_data[3].second->const_arrays();
456 Real dtc = m_dt_coarse;
457 Real
r = Real(1) / Real(ncycle);
458 Real xsi = Real(iteration-1) / Real(ncycle);
460 int const ng_space_interp = 8;
461 Box cdomain = m_cgeom.growPeriodicDomain(ng_space_interp);
462 cdomain.
convert(m_cf_crse_data_tmp->ixType());
466 Real b1 = xsi - Real(5./6.)*xsi*xsi;
467 Real b2 = Real(1./6.)*xsi*xsi;
468 Real b3 = Real(2./3)*xsi*xsi;
470 Real c1 = Real(1.) - Real(5./3.)*xsi;
471 Real c2 = Real(1./3.)*xsi;
472 Real c3 = Real(4./3.)*xsi;
474 constexpr Real d1 = Real(-5./3.);
475 constexpr Real d2 = Real(1./3.);
476 constexpr Real d3 = Real(4./3.);
482 Real kk1 = k1[bi](i,j,k,n);
483 Real kk2 = k2[bi](i,j,k,n);
484 Real kk3 = k3[bi](i,j,k,n);
485 Real uu = b1*kk1 + b2*kk2 + b3*kk3;
486 u[bi](i,j,k,n) = u0[bi](i,j,k,n) + dtc*uu;
489 }
else if (stage == 2) {
494 Real kk1 = k1[bi](i,j,k,n);
495 Real kk2 = k2[bi](i,j,k,n);
496 Real kk3 = k3[bi](i,j,k,n);
497 Real uu = b1*kk1 + b2*kk2 + b3*kk3;
498 Real ut = c1*kk1 + c2*kk2 + c3*kk3;
499 u[bi](i,j,k,n) = u0[bi](i,j,k,n) + dtc*(uu + r*ut);
502 }
else if (stage == 3) {
507 Real kk1 = k1[bi](i,j,k,n);
508 Real kk2 = k2[bi](i,j,k,n);
509 Real kk3 = k3[bi](i,j,k,n);
510 Real uu = b1*kk1 + b2*kk2 + b3*kk3;
511 Real ut = c1*kk1 + c2*kk2 + c3*kk3;
512 Real utt = d1*kk1 + d2*kk2 + d3*kk3;
513 u[bi](i,j,k,n) = u0[bi](i,j,k,n) + dtc*
514 (uu + Real(0.5)*r*ut + Real(0.25)*r*r*utt);
518 }
else if (rk_order == 4) {
519 auto const& k4 = m_cf_crse_data[4].second->const_arrays();
521 Real xsi3 = xsi2*xsi;
523 Real b1 = xsi - Real(1.5)*xsi2 + Real(2./3.)*xsi3;
524 Real b2 = xsi2 - Real(2./3.)*xsi3;
526 Real b4 = Real(-0.5)*xsi2 + Real(2./3.)*xsi3;
528 Real c1 = Real(1.) - Real(3.)*xsi + Real(2.)*xsi2;
529 Real c2 = Real(2.)*xsi - Real(2.)*xsi2;
531 Real c4 = -xsi + Real(2.)*xsi2;
533 Real d1 = Real(-3.) + Real(4.)*xsi;
534 Real d2 = Real( 2.) - Real(4.)*xsi;
536 Real d4 = Real(-1.) + Real(4.)*xsi;
538 constexpr Real e1 = Real( 4.);
539 constexpr Real e2 = Real(-4.);
540 constexpr Real e3 = Real(-4.);
541 constexpr Real e4 = Real( 4.);
547 Real kk1 = k1[bi](i,j,k,n);
548 Real kk2 = k2[bi](i,j,k,n);
549 Real kk3 = k3[bi](i,j,k,n);
550 Real kk4 = k4[bi](i,j,k,n);
551 Real uu = b1*kk1 + b2*kk2 + b3*kk3 + b4*kk4;
552 u[bi](i,j,k,n) = u0[bi](i,j,k,n) + dtc*uu;
555 }
else if (stage == 2) {
560 Real kk1 = k1[bi](i,j,k,n);
561 Real kk2 = k2[bi](i,j,k,n);
562 Real kk3 = k3[bi](i,j,k,n);
563 Real kk4 = k4[bi](i,j,k,n);
564 Real uu = b1*kk1 + b2*kk2 + b3*kk3 + b4*kk4;
565 Real ut = c1*kk1 + c2*kk2 + c3*kk3 + c4*kk4;
566 u[bi](i,j,k,n) = u0[bi](i,j,k,n) + dtc*(uu + Real(0.5)*r*ut);
569 }
else if (stage == 3 || stage == 4) {
572 Real
at = (stage == 3) ? Real(0.5)*
r :
r;
573 Real att = (stage == 3) ? Real(0.25)*r2 : Real(0.5)*r2;
574 Real attt = (stage == 3) ? Real(0.0625)*r3 : Real(0.125)*r3;
575 Real akk = (stage == 3) ? Real(-4.) : Real(4.);
580 Real kk1 = k1[bi](i,j,k,n);
581 Real kk2 = k2[bi](i,j,k,n);
582 Real kk3 = k3[bi](i,j,k,n);
583 Real kk4 = k4[bi](i,j,k,n);
584 Real uu = b1*kk1 + b2*kk2 + b3*kk3 + b4*kk4;
585 Real ut = c1*kk1 + c2*kk2 + c3*kk3 + c4*kk4;
586 Real utt = d1*kk1 + d2*kk2 + d3*kk3 + d4*kk4;
587 Real uttt = e1*kk1 + e2*kk2 + e3*kk3 + e4*kk4;
588 u[bi](i,j,k,n) = u0[bi](i,j,k,n) + dtc *
589 (uu + at*ut + att*utt + attt*(uttt+akk*(kk3-kk2)));
596 cbc(*m_cf_crse_data_tmp, 0, m_ncomp, m_cf_crse_data_tmp->nGrowVect(), time, 0);
598 if (m_cf_fine_data ==
nullptr) {
599 m_cf_fine_data = std::make_unique<MF>(detail::make_mf_fine_patch<MF>(fpc, m_ncomp));
603 m_ncomp,
IntVect(0), m_cgeom, m_fgeom,
605 mf.ixType()),m_nghost),
606 m_ratio, m_interp, bcs, 0);
609 mf.ParallelCopy(*m_cf_fine_data, 0, 0, m_ncomp,
IntVect(0), m_nghost);
611 mf.FillBoundary(m_fgeom.periodicity());
612 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:550
IndexType ixType() const noexcept
Return index type of this BoxArray.
Definition: AMReX_BoxArray.H:837
AMREX_GPU_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:912
AMREX_GPU_HOST_DEVICE IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition: AMReX_Box.H:146
AMREX_GPU_HOST_DEVICE bool contains(const IntVectND< dim > &p) const noexcept
Returns true if argument is contained within BoxND.
Definition: AMReX_Box.H:204
Calculates the distribution of FABs to MPI processes.
Definition: AMReX_DistributionMapping.H:41
Definition: AMReX_EB2.H:26
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 *)
FillPatcher is for filling a fine level MultiFab/FabArray.
Definition: AMReX_FillPatcher.H:73
Geometry m_cgeom
Definition: AMReX_FillPatcher.H:205
Vector< std::pair< Real, std::unique_ptr< MF > > > m_cf_crse_data
Definition: AMReX_FillPatcher.H:212
Geometry m_fgeom
Definition: AMReX_FillPatcher.H:204
DistributionMapping m_cdm
Definition: AMReX_FillPatcher.H:203
InterpBase * m_interp
Definition: AMReX_FillPatcher.H:208
BoxArray m_fba
Definition: AMReX_FillPatcher.H:200
std::unique_ptr< MF > m_cf_fine_data
Definition: AMReX_FillPatcher.H:214
std::unique_ptr< MF > m_cf_crse_data_tmp
Definition: AMReX_FillPatcher.H:213
BoxArray m_cba
Definition: AMReX_FillPatcher.H:201
IntVect m_nghost
Definition: AMReX_FillPatcher.H:206
int m_ncomp
Definition: AMReX_FillPatcher.H:207
EB2::IndexSpace const * m_eb_index_space
Definition: AMReX_FillPatcher.H:209
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:286
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:252
IntVect m_ratio
Definition: AMReX_FillPatcher.H:211
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:432
MF m_sfine
Definition: AMReX_FillPatcher.H:210
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:411
DistributionMapping m_fdm
Definition: AMReX_FillPatcher.H:202
FabArrayBase::FPinfo const & getFPinfo()
Definition: AMReX_FillPatcher.H:276
Real m_dt_coarse
Definition: AMReX_FillPatcher.H:215
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=nullptr)
Constructor of FillPatcher.
Definition: AMReX_FillPatcher.H:221
Rectangular problem domain geometry.
Definition: AMReX_Geometry.H:73
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition: AMReX_Geometry.H:210
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool cellCentered() const noexcept
True if the IndexTypeND is CELL based in all directions.
Definition: AMReX_IndexType.H:101
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool nodeCentered() const noexcept
True if the IndexTypeND is NODE based in all directions.
Definition: AMReX_IndexType.H:107
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE 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:393
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:27
Long size() const noexcept
Definition: AMReX_Vector.H:50
const IndexSpace * TopIndexSpaceIfPresent() noexcept
Definition: AMReX_EB2.cpp:76
AMREX_GPU_HOST_DEVICE Long at(T const &, Long offset) noexcept
Definition: AMReX_GpuRange.H:30
void streamSynchronize() noexcept
Definition: AMReX_GpuDevice.H:237
auto call_interp_hook(F const &f, MF &mf, int icomp, int ncomp) -> decltype(f(mf[0], Box(), icomp, ncomp))
Definition: AMReX_FillPatchUtil_I.H:10
Definition: AMReX_Amr.cpp:49
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:200
DistributionMapping const & DistributionMap(FabArrayBase const &fa)
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:73
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition: AMReX_Box.H:1211
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Returns a BoxND with different type.
Definition: AMReX_Box.H:1435
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE T abs(const GpuComplex< T > &a_z) noexcept
Return the absolute value of a complex number.
Definition: AMReX_GpuComplex.H:356
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 end(BoxND< dim > const &box) noexcept
Definition: AMReX_Box.H:1890
BoxArray const & boxArray(FabArrayBase const &fa)
IntVectND< AMREX_SPACEDIM > IntVect
Definition: AMReX_BaseFwd.H:30
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
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > refine(const BoxND< dim > &b, int ref_ratio) noexcept
Definition: AMReX_Box.H:1342
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition: AMReX.cpp:225
std::array< T, N > Array
Definition: AMReX_Array.H:24
Definition: AMReX_FabArrayBase.H:304
Definition: AMReX_TypeTraits.H:29
FabArray memory allocation information.
Definition: AMReX_FabArray.H:66
Definition: AMReX_FillPatchUtil.H:34