1#ifndef AMREX_FillPatchUtil_I_H_
2#define AMREX_FillPatchUtil_I_H_
3#include <AMReX_Config.H>
10template <
typename F,
typename MF>
11auto call_interp_hook (
F const& f, MF& mf,
int icomp,
int ncomp)
12 ->
decltype(f(mf[0],
Box(),icomp,ncomp))
15#pragma omp parallel if (Gpu::notInLaunchRegion())
17 for (MFIter mfi(mf); mfi.isValid(); ++mfi) {
19 const Box& dbx = dfab.box();
20 f(dfab, dbx, icomp, ncomp);
24template <
typename F,
typename MF>
25auto call_interp_hook (
F const& f, MF& mf,
int icomp,
int ncomp)
26 ->
decltype(f(mf,icomp,ncomp))
34template <
typename Interp>
38 int ratio_max = ratio[0];
39#if (AMREX_SPACEDIM > 1)
40 ratio_max = std::max(ratio_max, ratio[1]);
42#if (AMREX_SPACEDIM == 3)
43 ratio_max = std::max(ratio_max, ratio[2]);
47 const IntVect& nbuf = blocking_factor / ratio_max;
55 fine_box.
refine(ratio_max);
57 const Box& fine_box_coarsened = mapper->CoarseBox(fine_box, ratio_max);
58 return crse_box.
contains(fine_box_coarsened);
61template <FabArrayType MF,
typename BC>
65 int scomp,
int dcomp,
int ncomp,
67 BC& physbcf,
int bcfcomp)
70 geom, physbcf, bcfcomp);
73template <FabArrayType MF,
typename BC>
77 int scomp,
int dcomp,
int ncomp,
79 BC& physbcf,
int bcfcomp)
84 (scomp+ncomp <= smf[0]->nComp()) &&
85 (dcomp+ncomp <= mf.nComp()) &&
87 nghost.
allLE(mf.nGrowVect()));
90 if constexpr (std::is_same_v<BC,PhysBCFunctUseCoarseGhost>) {
91 src_ghost = physbcf.fp1_src_ghost;
96 if (&mf == smf[0] && scomp == dcomp) {
97 mf.FillBoundary(dcomp, ncomp, nghost, geom.
periodicity());
99 mf.ParallelCopy(*smf[0], scomp, dcomp, ncomp, src_ghost, nghost, geom.
periodicity());
102 else if (smf.
size() == 2)
109 if (mf.boxArray() == smf[0]->boxArray() &&
110 mf.DistributionMap() == smf[0]->DistributionMap())
117 MFInfo(), smf[0]->Factory());
124 if ((dmf != smf[0] && dmf != smf[1]) || scomp != dcomp)
127 if constexpr (std::is_same_v<BC,PhysBCFunctUseCoarseGhost>) {
128 interp_ghost = physbcf.fp1_src_ghost;
130 interp_ghost.
min(nghost);
134#pragma omp parallel if (Gpu::notInLaunchRegion())
138 const Box& bx = mfi.growntilebox(interp_ghost);
139 const Real t0 = stime[0];
140 const Real t1 = stime[1];
141 auto const sfab0 = smf[0]->array(mfi);
142 auto const sfab1 = smf[1]->array(mfi);
143 auto dfab = dmf->array(mfi);
149 dfab(i,j,k,n+destcomp) = sfab0(i,j,k,n+scomp);
156 dfab(i,j,k,n+destcomp) = sfab1(i,j,k,n+scomp);
161 Real alpha = (t1-time)/(t1-t0);
162 Real beta = (time-t0)/(t1-t0);
165 dfab(i,j,k,n+destcomp) = alpha*sfab0(i,j,k,n+scomp)
166 + beta*sfab1(i,j,k,n+scomp);
173 dfab(i,j,k,n+destcomp) = sfab0(i,j,k,n+scomp);
183 mf.FillBoundary(dcomp, ncomp, nghost, geom.
periodicity());
187 mf.ParallelCopy(*dmf, 0, dcomp, ncomp, src_ghost, nghost, geom.
periodicity());
191 amrex::Abort(
"FillPatchSingleLevel: high-order interpolation in time not implemented yet");
194 physbcf(mf, dcomp, ncomp, nghost, time, bcfcomp);
214void FillPatchInterp (MultiFab& mf_fine_patch,
int fcomp, MultiFab
const& mf_crse_patch,
int ccomp,
215 int ncomp,
IntVect const& ng,
const Geometry& cgeom,
const Geometry& fgeom,
217 MFInterpolater* mapper,
const Vector<BCRec>& bcs,
int bcscomp);
219template <FabArrayType MF,
typename Interp>
220requires (!std::is_same_v<Interp, MFInterpolater>)
232#pragma omp parallel if (Gpu::notInLaunchRegion())
238 auto& sfab = mf_crse_patch[mfi];
239 const Box& sbx = sfab.box();
241 auto& dfab = mf_fine_patch[mfi];
245 mapper->interp(sfab, ccomp, dfab, fcomp, ncomp, dbx, ratio,
246 cgeom, fgeom, bcr, idummy, idummy,
RunOn::Gpu);
251template <FabArrayType MF>
260 ncomp, ng, cgeom, fgeom, dest_domain, ratio,
264 ncomp, ng, cgeom, fgeom, dest_domain, ratio,
271template <FabArrayType MF,
typename iMF,
typename Interp>
272requires (!std::is_same_v<Interp, MFInterpolater>)
275 MF
const& mf_crse_patch,
int crse_comp,
276 MF& mf_refined_patch,
int fine_comp,
277 int ncomp,
const IntVect& ratio,
278 const iMF& solve_mask,
const Geometry& crse_geom,
const Geometry& fine_geom,
279 int bcscomp,
RunOn gpu_or_cpu,
286 auto& sfab = mf_crse_patch[mfi];
287 const Box& sbx = sfab.box();
288 auto& dfab = mf_refined_patch[mfi];
289 Box const& dbx = dfab.box();
290 auto& ifab = solve_mask[mfi];
292 interp->interp_face(sfab,crse_comp,dfab,fine_comp,ncomp,
293 dbx, ratio, ifab, crse_geom, fine_geom,
294 bcr, bcscomp, gpu_or_cpu);
298template <FabArrayType MF,
typename iMF>
301 MF
const& mf_crse_patch,
int crse_comp,
302 MF& mf_refined_patch,
int fine_comp,
303 int ncomp,
const IntVect& ratio,
304 const iMF& solve_mask,
const Geometry& crse_geom,
const Geometry& fine_geom,
305 int bccomp,
RunOn gpu_or_cpu,
310 mf_crse_patch, crse_comp,mf_refined_patch, fine_comp,
311 ncomp, ratio, solve_mask, crse_geom, fine_geom, bccomp,
316 mf_crse_patch, crse_comp,mf_refined_patch, fine_comp,
317 ncomp, ratio, solve_mask, crse_geom, fine_geom, bccomp,
331 template <MultiFabLike MF>
332 requires (std::same_as<typename MF::FABType::value_type, FArrayBox>)
333 MF make_mf_crse_patch (FabArrayBase::FPinfo
const& fpc,
int ncomp)
335 MF mf_crse_patch(fpc.ba_crse_patch, fpc.dm_patch, ncomp, 0, MFInfo(),
336 *fpc.fact_crse_patch);
337 return mf_crse_patch;
340 template <MultiFabLike MF>
341 requires (std::same_as<typename MF::FABType::value_type, FArrayBox>)
342 MF make_mf_crse_patch (FabArrayBase::FPinfo
const& fpc,
int ncomp,
IndexType idx_type)
344 MF mf_crse_patch(
amrex::convert(fpc.ba_crse_patch, idx_type), fpc.dm_patch,
345 ncomp, 0, MFInfo(), *fpc.fact_crse_patch);
346 return mf_crse_patch;
349 template <MultiFabLike MF>
350 requires (std::same_as<typename MF::FABType::value_type, FArrayBox>)
351 MF make_mf_fine_patch (FabArrayBase::FPinfo
const& fpc,
int ncomp)
353 MF mf_fine_patch(fpc.ba_fine_patch, fpc.dm_patch, ncomp, 0, MFInfo(),
354 *fpc.fact_fine_patch);
355 return mf_fine_patch;
358 template <MultiFabLike MF>
359 requires (std::same_as<typename MF::FABType::value_type, FArrayBox>)
360 MF make_mf_fine_patch (FabArrayBase::FPinfo
const& fpc,
int ncomp,
IndexType idx_type)
362 MF mf_fine_patch(
amrex::convert(fpc.ba_fine_patch, idx_type), fpc.dm_patch,
363 ncomp, 0, MFInfo(), *fpc.fact_fine_patch);
364 return mf_fine_patch;
367 template <MultiFabLike MF>
368 requires (std::same_as<typename MF::FABType::value_type, FArrayBox>)
369 MF make_mf_refined_patch (FabArrayBase::FPinfo
const& fpc,
int ncomp,
IndexType idx_type,
IntVect ratio)
372 fpc.dm_patch, ncomp, 0, MFInfo(), *fpc.fact_fine_patch);
373 return mf_refined_patch;
376 template <MultiFabLike MF>
377 requires (std::same_as<typename MF::FABType::value_type, FArrayBox>)
378 MF make_mf_crse_mask (FabArrayBase::FPinfo
const& fpc,
int ncomp,
IndexType idx_type,
IntVect ratio)
381 fpc.dm_patch, ncomp, 0, MFInfo(), *fpc.fact_fine_patch);
385 template <MultiFabLike MF>
386 requires (std::same_as<typename MF::FABType::value_type, FArrayBox>)
387 void mf_set_domain_bndry (MF &mf, Geometry
const & geom)
389 mf.setDomainBndry(std::numeric_limits<Real>::quiet_NaN(), geom);
395 template <MultiFabLike MF>
396 requires (!std::same_as<typename MF::FABType::value_type, FArrayBox>)
397 MF make_mf_crse_patch (FabArrayBase::FPinfo
const& fpc,
int ncomp)
399 return MF(fpc.ba_crse_patch, fpc.dm_patch, ncomp, 0);
402 template <MultiFabLike MF>
403 requires (!std::same_as<typename MF::FABType::value_type, FArrayBox>)
404 MF make_mf_crse_patch (FabArrayBase::FPinfo
const& fpc,
int ncomp,
IndexType idx_type)
406 return MF(
amrex::convert(fpc.ba_crse_patch, idx_type), fpc.dm_patch, ncomp, 0);
409 template <MultiFabLike MF>
410 requires (!std::same_as<typename MF::FABType::value_type, FArrayBox>)
411 MF make_mf_fine_patch (FabArrayBase::FPinfo
const& fpc,
int ncomp)
413 return MF(fpc.ba_fine_patch, fpc.dm_patch, ncomp, 0);
416 template <MultiFabLike MF>
417 requires (!std::same_as<typename MF::FABType::value_type, FArrayBox>)
418 MF make_mf_fine_patch (FabArrayBase::FPinfo
const& fpc,
int ncomp,
IndexType idx_type)
420 return MF(
amrex::convert(fpc.ba_fine_patch, idx_type), fpc.dm_patch, ncomp, 0);
423 template <MultiFabLike MF>
424 requires (!std::same_as<typename MF::FABType::value_type, FArrayBox>)
425 MF make_mf_refined_patch (FabArrayBase::FPinfo
const& fpc,
int ncomp,
IndexType idx_type,
IntVect ratio)
430 template <MultiFabLike MF>
431 requires (!std::same_as<typename MF::FABType::value_type, FArrayBox>)
432 MF make_mf_crse_mask (FabArrayBase::FPinfo
const& fpc,
int ncomp,
IndexType idx_type,
IntVect ratio)
437 template <MultiFabLike MF>
438 requires (!std::same_as<typename MF::FABType::value_type, FArrayBox>)
439 void mf_set_domain_bndry (MF &, Geometry
const & )
444 template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
446 FillPatchTwoLevels_doit (MF& mf,
IntVect const& nghost,
Real time,
447 const Vector<MF*>& cmf,
const Vector<Real>& ct,
448 const Vector<MF*>& fmf,
const Vector<Real>& ft,
449 int scomp,
int dcomp,
int ncomp,
450 const Geometry& cgeom,
const Geometry& fgeom,
451 BC& cbc,
int cbccomp,
452 BC& fbc,
int fbccomp,
455 const Vector<BCRec>& bcs,
int bcscomp,
456 const PreInterpHook& pre_interp,
457 const PostInterpHook& post_interp,
458 EB2::IndexSpace
const* index_space,
459 bool return_error_code =
false)
463 int success_code = return_error_code ? 0 : -1;
464 int failure_code = 1;
466 if (nghost.max() > 0 || mf.getBDKey() != fmf[0]->getBDKey())
468 const InterpolaterBoxCoarsener& coarsener = mapper->BoxCoarsener(ratio);
472 + mf.ixType().nodeCentered(1),
473 + mf.ixType().nodeCentered(2) ) == 1 )
475 if ( !
dynamic_cast<Interpolater*
>(mapper) ){
476 amrex::Abort(
"This interpolater has not yet implemented a version for face-based data");
481 mf.DistributionMap(), ncomp, nghost, MFInfo().SetAlloc(
false) );
492 if ( ! fpc.ba_crse_patch.empty())
494 if (return_error_code) {
496 if (!cba.contains(fpc.ba_crse_patch,cgeom.periodicity())) {
501 MF mf_crse_patch = make_mf_crse_patch<MF> (fpc, ncomp, mf.boxArray().ixType());
506 MF mf_refined_patch = make_mf_refined_patch<MF> (fpc, ncomp, mf.boxArray().ixType(), ratio);
507 auto solve_mask = make_mf_crse_mask<iMultiFab>(fpc, ncomp, mf.boxArray().ixType(), ratio);
509 mf_set_domain_bndry(mf_crse_patch, cgeom);
510 if constexpr (std::is_same_v<BC,PhysBCFunctUseCoarseGhost>) {
511 cbc.fp1_src_ghost = cbc.cghost;
514 cgeom, cbc, cbccomp);
516 mf_set_domain_bndry(mf_refined_patch, fgeom);
517 if constexpr (std::is_same_v<BC,PhysBCFunctUseCoarseGhost>) {
518 fbc.fp1_src_ghost =
IntVect(0);
521 fgeom, fbc, fbccomp);
525 ncomp, nghost, MFInfo().SetAlloc(
false) );
526 MF mf_solution(
amrex::coarsen(mf_refined_patch.boxArray(), ratio), mf_refined_patch.DistributionMap(),
527 ncomp, 0, MFInfo().SetAlloc(
false) );
531 cgeom.periodicity());
533 solve_mask.setVal(1);
534 solve_mask.setVal(0, mask_cpc, 0, 1);
536 detail::call_interp_hook(pre_interp, mf_crse_patch, 0, ncomp);
538 InterpFace(mapper, mf_crse_patch, 0, mf_refined_patch, 0, ncomp,
539 ratio, solve_mask, cgeom, fgeom, bcscomp,
RunOn::Gpu, bcs);
541 detail::call_interp_hook(post_interp, mf_refined_patch, 0, ncomp);
543 bool aliasing =
false;
544 for (
auto const& fmf_a : fmf) {
545 aliasing = aliasing || (&mf == fmf_a);
548 mf.ParallelCopyToGhost(mf_refined_patch, 0, dcomp, ncomp,
551 mf.ParallelCopy(mf_refined_patch, 0, dcomp, ncomp,
565 if ( ! fpc.ba_crse_patch.empty())
567 if (return_error_code) {
568 BoxArray
const& cba = cmf[0]->boxArray();
569 if (!cba.contains(fpc.ba_crse_patch,cgeom.periodicity())) {
574 MF mf_crse_patch = make_mf_crse_patch<MF>(fpc, ncomp);
575 mf_set_domain_bndry (mf_crse_patch, cgeom);
577 if constexpr (std::is_same_v<BC,PhysBCFunctUseCoarseGhost>) {
578 cbc.fp1_src_ghost = cbc.cghost;
582 MF mf_fine_patch = make_mf_fine_patch<MF>(fpc, ncomp);
584 detail::call_interp_hook(pre_interp, mf_crse_patch, 0, ncomp);
587 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
588 if (fgeom.isPeriodic(i)) {
589 fdomain_g.grow(i, nghost[i]);
591 if constexpr (std::is_same_v
592 <BC, PhysBCFunctUseCoarseGhost>) {
593 fdomain_g.grow(i, fbc.nghost_outside_domain[i]);
598 ncomp,
IntVect(0), cgeom, fgeom,
599 fdomain_g, ratio, mapper, bcs, bcscomp);
601 detail::call_interp_hook(post_interp, mf_fine_patch, 0, ncomp);
603 mf.ParallelCopy(mf_fine_patch, 0, dcomp, ncomp,
IntVect{0}, nghost);
608 if constexpr(std::is_same_v<BC, PhysBCFunctUseCoarseGhost>) {
609 fbc.fp1_src_ghost =
IntVect(0);
612 fgeom, fbc, fbccomp);
617 template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
619 FillPatchTwoLevels_doit (Array<MF*, AMREX_SPACEDIM>
const& mf,
IntVect const& nghost,
Real time,
620 const Vector<Array<MF*, AMREX_SPACEDIM> >& cmf,
const Vector<Real>& ct,
621 const Vector<Array<MF*, AMREX_SPACEDIM> >& fmf,
const Vector<Real>& ft,
622 int scomp,
int dcomp,
int ncomp,
623 const Geometry& cgeom,
const Geometry& fgeom,
624 Array<BC, AMREX_SPACEDIM>& cbc,
const Array<int, AMREX_SPACEDIM>& cbccomp,
625 Array<BC, AMREX_SPACEDIM>& fbc,
const Array<int, AMREX_SPACEDIM>& fbccomp,
628 const Array<Vector<BCRec>, AMREX_SPACEDIM>& bcs,
const Array<int, AMREX_SPACEDIM>& bcscomp,
629 const PreInterpHook& pre_interp,
630 const PostInterpHook& post_interp,
631 EB2::IndexSpace
const* index_space)
633 BL_PROFILE(
"FillPatchTwoLevels (Array<MF*>)");
635 using FAB =
typename MF::FABType::value_type;
636 using iFAB =
typename iMultiFab::FABType::value_type;
639 && mf[1]->ixType().nodeCentered(1),
640 && mf[2]->ixType().nodeCentered(2)));
655 if (nghost.max() > 0 || mf[0]->getBDKey() != fmf[0][0]->getBDKey())
657 const InterpolaterBoxCoarsener& coarsener = mapper->BoxCoarsener(ratio);
663 fmf[0][0]->
DistributionMap(), ncomp, nghost, MFInfo().SetAlloc(
false) );
672 if ( !fpc.ba_crse_patch.empty() )
674 Array<MF, AMREX_SPACEDIM> mf_crse_patch;
675 Array<MF, AMREX_SPACEDIM> mf_refined_patch;
676 Array<iMultiFab, AMREX_SPACEDIM> solve_mask;
678 for (
int d=0; d<AMREX_SPACEDIM; ++d)
680 mf_crse_patch[d] = make_mf_crse_patch<MF> (fpc, ncomp, mf[d]->
boxArray().ixType());
681 mf_refined_patch[d] = make_mf_refined_patch<MF> (fpc, ncomp, mf[d]->
boxArray().ixType(), ratio);
682 solve_mask[d] = make_mf_crse_mask<iMultiFab>(fpc, ncomp, mf[d]->
boxArray().ixType(), ratio);
684 mf_set_domain_bndry(mf_crse_patch[d], cgeom);
685 Vector<MF*> cmf_time;
686 for (
const auto & mfab : cmf)
687 { cmf_time.push_back(mfab[d]); }
690 cgeom, cbc[d], cbccomp[d]);
692 mf_set_domain_bndry(mf_refined_patch[d], fgeom);
693 Vector<MF*> fmf_time;
694 for (
const auto & mfab : fmf)
695 { fmf_time.push_back(mfab[d]); }
698 fgeom, fbc[d], fbccomp[d]);
703 ncomp, nghost, MFInfo().SetAlloc(
false) );
705 ncomp, 0, MFInfo().SetAlloc(
false) );
709 cgeom.periodicity() );
711 solve_mask[d].setVal(1);
712 solve_mask[d].setVal(0, mask_cpc, 0, 1);
719#pragma omp parallel if (cc && Gpu::notInLaunchRegion() )
722 Vector<Array<BCRec, AMREX_SPACEDIM> > bcr(ncomp);
723 for (MFIter mfi(mf_refined_patch[0]); mfi.isValid(); ++mfi)
725 Array<FAB*, AMREX_SPACEDIM> sfab{
AMREX_D_DECL( &(mf_crse_patch[0][mfi]),
726 &(mf_crse_patch[1][mfi]),
727 &(mf_crse_patch[2][mfi]) )};
728 Array<FAB*, AMREX_SPACEDIM> dfab{
AMREX_D_DECL( &(mf_refined_patch[0][mfi]),
729 &(mf_refined_patch[1][mfi]),
730 &(mf_refined_patch[2][mfi]) )};
731 Array<iFAB*, AMREX_SPACEDIM> mfab{
AMREX_D_DECL( &(solve_mask[0][mfi]),
732 &(solve_mask[1][mfi]),
733 &(solve_mask[2][mfi]) )};
738 for (
int d=0; d<AMREX_SPACEDIM; ++d)
740 Vector<BCRec> bcr_d(ncomp);
743 bcscomp[d],0,ncomp,bcs[d],bcr_d);
745 for (
int n=0; n<ncomp; ++n)
746 { bcr[n][d] = bcr_d[n]; }
749 pre_interp(sfab, sbx_cc, 0, ncomp);
751 mapper->interp_arr(sfab, 0, dfab, 0, ncomp, dbx_cc, ratio, mfab,
752 cgeom, fgeom, bcr, idummy, idummy,
RunOn::Gpu);
754 post_interp(dfab, dbx_cc, 0, ncomp);
758 for (
int d=0; d<AMREX_SPACEDIM; ++d)
760 bool aliasing =
false;
761 for (
auto const& fmf_a : fmf) {
762 aliasing = aliasing || (mf[d] == fmf_a[d]);
765 mf[d]->ParallelCopyToGhost(mf_refined_patch[d], 0, dcomp, ncomp,
768 mf[d]->ParallelCopy(mf_refined_patch[d], 0, dcomp, ncomp,
775 for (
int d=0; d<AMREX_SPACEDIM; ++d)
777 Vector<MF*> fmf_time;
778 for (
auto const& ffab : fmf)
779 { fmf_time.push_back(ffab[d]); }
782 fgeom, fbc[d], fbccomp[d]);
789template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
794 int scomp,
int dcomp,
int ncomp,
796 BC& cbc,
int cbccomp,
797 BC& fbc,
int fbccomp,
801 const PreInterpHook& pre_interp,
802 const PostInterpHook& post_interp)
809 detail::FillPatchTwoLevels_doit(mf,nghost,time,cmf,ct,fmf,ft,
810 scomp,dcomp,ncomp,cgeom,fgeom,
811 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
812 pre_interp,post_interp,index_space);
815template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
820 int scomp,
int dcomp,
int ncomp,
822 BC& cbc,
int cbccomp,
823 BC& fbc,
int fbccomp,
827 const PreInterpHook& pre_interp,
828 const PostInterpHook& post_interp)
836 detail::FillPatchTwoLevels_doit(mf,mf.nGrowVect(),time,cmf,ct,fmf,ft,
837 scomp,dcomp,ncomp,cgeom,fgeom,
838 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
839 pre_interp,post_interp,index_space);
842template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
847 int scomp,
int dcomp,
int ncomp,
854 const PreInterpHook& pre_interp,
855 const PostInterpHook& post_interp)
863 detail::FillPatchTwoLevels_doit(mf,nghost,time,cmf,ct,fmf,ft,
864 scomp,dcomp,ncomp,cgeom,fgeom,
865 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
866 pre_interp,post_interp,index_space);
869template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
874 int scomp,
int dcomp,
int ncomp,
881 const PreInterpHook& pre_interp,
882 const PostInterpHook& post_interp)
894 detail::FillPatchTwoLevels_doit(mf,nghost,time,cmf,ct,fmf,ft,
895 scomp,dcomp,ncomp,cgeom,fgeom,
896 cbc,cbccomp_arr,fbc,fbccomp_arr,ratio,mapper,bcs,bcscomp_arr,
897 pre_interp,post_interp,index_space);
900template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
905 int scomp,
int dcomp,
int ncomp,
912 const PreInterpHook& pre_interp,
913 const PostInterpHook& post_interp)
925 detail::FillPatchTwoLevels_doit(mf,mf[0]->
nGrowVect(),time,cmf,ct,fmf,ft,
926 scomp,dcomp,ncomp,cgeom,fgeom,
927 cbc,cbccomp_arr,fbc,fbccomp_arr,ratio,mapper,bcs,bcscomp_arr,
928 pre_interp,post_interp,index_space);
932template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
938 int scomp,
int dcomp,
int ncomp,
940 BC& cbc,
int cbccomp,
941 BC& fbc,
int fbccomp,
945 const PreInterpHook& pre_interp,
946 const PostInterpHook& post_interp)
948 detail::FillPatchTwoLevels_doit(mf,nghost,time,cmf,ct,fmf,ft,
949 scomp,dcomp,ncomp,cgeom,fgeom,
950 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
951 pre_interp,post_interp,&index_space);
954template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
960 int scomp,
int dcomp,
int ncomp,
962 BC& cbc,
int cbccomp,
963 BC& fbc,
int fbccomp,
967 const PreInterpHook& pre_interp,
968 const PostInterpHook& post_interp)
970 detail::FillPatchTwoLevels_doit(mf,mf.nGrowVect(),time,cmf,ct,fmf,ft,
971 scomp,dcomp,ncomp,cgeom,fgeom,
972 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
973 pre_interp,post_interp,&index_space);
977template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
980 const MF& cmf,
int scomp,
int dcomp,
int ncomp,
982 BC& cbc,
int cbccomp,
983 BC& fbc,
int fbccomp,
987 const PreInterpHook& pre_interp,
988 const PostInterpHook& post_interp)
996 InterpFromCoarseLevel(mf,mf.nGrowVect(),time,index_space,cmf,scomp,dcomp,ncomp,cgeom,fgeom,
997 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
998 pre_interp,post_interp);
1001template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
1011 const PreInterpHook& pre_interp,
1012 const PostInterpHook& post_interp)
1015 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
1016 pre_interp,post_interp);
1019template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
1022 const MF& cmf,
int scomp,
int dcomp,
int ncomp,
1024 BC& cbc,
int cbccomp,
1025 BC& fbc,
int fbccomp,
1029 const PreInterpHook& pre_interp,
1030 const PostInterpHook& post_interp)
1038 InterpFromCoarseLevel(mf,nghost,time,index_space,cmf,scomp,dcomp,ncomp,cgeom,fgeom,
1039 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
1040 pre_interp,post_interp);
1043template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
1047 const MF& cmf,
int scomp,
int dcomp,
int ncomp,
1049 BC& cbc,
int cbccomp,
1050 BC& fbc,
int fbccomp,
1054 const PreInterpHook& pre_interp,
1055 const PostInterpHook& post_interp)
1059 using FAB =
typename MF::FABType::value_type;
1063 const BoxArray& ba = mf.boxArray();
1071 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
1073 fdomain_g.
grow(i, nghost[i]);
1075 if constexpr (std::is_same_v<BC, PhysBCFunctUseCoarseGhost>) {
1076 fdomain_g.
grow(i, fbc.nghost_outside_domain[i]);
1082 IntVect send_ghost(0), recv_ghost(0);
1083 if constexpr (std::is_same_v<BC, PhysBCFunctUseCoarseGhost>) {
1084 mf_crse_patch.define(
amrex::coarsen(ba,ratio), dm, ncomp, fbc.src_ghost);
1085 send_ghost = fbc.cghost;
1086 recv_ghost = fbc.src_ghost;
1090 for (
int i = 0, N = ba.
size(); i < N; ++i)
1094 ba_crse_patch.
set(i, coarsener.
doit(bx));
1104 mf_crse_patch.define(ba_crse_patch, dm, ncomp, 0,
MFInfo(), *factory);
1108 mf_crse_patch.define(ba_crse_patch, dm, ncomp, 0);
1110 detail::mf_set_domain_bndry (mf_crse_patch, cgeom);
1113 mf_crse_patch.ParallelCopy(cmf, scomp, 0, ncomp, send_ghost, recv_ghost,
1116 cbc(mf_crse_patch, 0, ncomp, mf_crse_patch.nGrowVect(), time, cbccomp);
1118 detail::call_interp_hook(pre_interp, mf_crse_patch, 0, ncomp);
1120 FillPatchInterp(mf, dcomp, mf_crse_patch, 0, ncomp, nghost, cgeom, fgeom, fdomain_g,
1121 ratio, mapper, bcs, bcscomp);
1124#pragma omp parallel if (Gpu::notInLaunchRegion())
1128 FAB& dfab = mf[mfi];
1129 Box dfab_bx = dfab.box();
1130 dfab_bx.
grow(nghost-mf.nGrowVect());
1131 const Box& dbx = dfab_bx & fdomain_g;
1133 post_interp(dfab, dbx, dcomp, ncomp);
1136 fbc(mf, dcomp, ncomp, nghost, time, fbccomp);
1139template <FabArrayType MF,
typename BC,
typename Interp,
typename PreInterpHook,
typename PostInterpHook>
1149 const PreInterpHook& pre_interp,
1150 const PostInterpHook& post_interp)
1154 using FAB =
typename MF::FABType::value_type;
1155 using iFAB =
typename iMultiFab::FABType::value_type;
1158 const BoxArray& ba = mf[0]->boxArray();
1162 && mf[1]->ixType().nodeCentered(1),
1163 && mf[2]->ixType().nodeCentered(2)));
1178 for (
int d=0; d<AMREX_SPACEDIM; ++d) {
1179 if (nghost[d] % ratio[d] != 0) {
1180 nghost_adj[d] += ratio[d] - (nghost[d] % ratio[d]);
1185 int dcomp_adj = dcomp;
1187 if (! nghost.
allGE(nghost_adj)) {
1188 for (
int d=0; d<AMREX_SPACEDIM; ++d) {
1189 mf_temp[d] = std::make_unique<MF>(mf[d]->
boxArray(),
1191 mf_local[d] = mf_temp[d].get();
1199 Box fdomain_g(fdomain);
1200 for (
int d = 0; d < AMREX_SPACEDIM; ++d) {
1202 fdomain_g.
grow(d,nghost_adj[d]);
1209 for (
int i = 0, N = ba.
size(); i < N; ++i)
1213 ba_crse_patch.
set(i, coarsener.
doit(bx));
1218 for (
int d = 0; d<AMREX_SPACEDIM; ++d)
1225 mf_crse_patch[d].define(ba_crse_idxed, dm, ncomp, 0,
MFInfo(), *crse_factory);
1227 mf_crse_patch[d].define(ba_crse_idxed, dm, ncomp, 0);
1229 detail::mf_set_domain_bndry(mf_crse_patch[d], cgeom);
1231 mf_crse_patch[d].ParallelCopy(*(cmf[d]), scomp, 0, ncomp, cgeom.
periodicity());
1232 cbc[d](mf_crse_patch[d], 0, ncomp, mf_crse_patch[d].nGrowVect(), time, cbccomp);
1235 int idummy1=0, idummy2=0;
1237#pragma omp parallel if (Gpu::notInLaunchRegion())
1248 &(mf_crse_patch[1][mfi]),
1249 &(mf_crse_patch[2][mfi]) )};
1251 &(*mf_local[1])[mfi],
1252 &(*mf_local[2])[mfi] )};
1256 const Box& dbx_cc = dfab_cc & fdomain_g;
1258 for (
int d=0; d<AMREX_SPACEDIM; ++d)
1264 bcscomp,0,ncomp,bcs[d],bcr_d);
1266 for (
int n=0; n<ncomp; ++n)
1267 { bcr[n][d] = bcr_d[n]; }
1270 pre_interp(sfab, sbx_cc, 0, ncomp);
1272 mapper->interp_arr(sfab, 0, dfab, 0, ncomp, dbx_cc, ratio, mfab,
1273 cgeom, fgeom, bcr, idummy1, idummy2,
RunOn::Gpu);
1275 post_interp(dfab, dbx_cc, 0, ncomp);
1279 for (
int d=0; d<AMREX_SPACEDIM; ++d)
1281 if (mf[d] != mf_local[d]) {
1282 amrex::Copy(*mf[d], *mf_local[d], 0, dcomp_adj, ncomp, nghost);
1285 fbc[d](*mf[d], dcomp, ncomp, nghost, time, fbccomp);
1289template <FabArrayType MF,
typename Interp>
1292 IntVect const& nghost_outside_domain,
1293 const MF& cmf,
int scomp,
int dcomp,
int ncomp,
1295 const IntVect& ratio, Interp* mapper,
1300 cgeom, fgeom, erfbc, 0, erfbc, 0, ratio, mapper,
1304template <FabArrayType MF>
1308 const Vector<Real>& stime,
int scomp,
int dcomp,
int ncomp,
1316template <FabArrayType MF,
typename Interp>
1322 int scomp,
int dcomp,
int ncomp,
1324 const IntVect& ratio, Interp* mapper,
1330 FillPatchTwoLevels(mf, nghost, time, cmf, ct, fmf, ft, scomp, dcomp, ncomp,
1331 cgeom, fgeom, erfbc, 0, erfbc, 0, ratio, mapper,
1335template <FabArrayType MF,
typename BC,
typename Interp>
1339 int scomp,
int dcomp,
int ncomp,
1354 auto get_clayout = [&] () -> std::tuple<BoxArray,BoxArray,DistributionMapping>
1359 BoxArray const& ba = mf.boxArray();
1360 auto const& typ = ba.
ixType();
1361 std::map<int,Vector<Box>> extra_boxes_map;
1364 for (
int i = 0, N =
int(ba.
size()); i < N; ++i) {
1365 Box const& cbox = mapper->CoarseBox(
amrex::grow(ba[i],nghost),ratio[level-1]);
1367 Box gdomain = geom[level-1].growNonPeriodicDomain(cbox.
length());
1370 auto& extra_boxes = extra_boxes_map[i];
1371 auto const& pshift = geom[level-1].periodicity().shiftIntVect();
1372 for (
auto const& piv : pshift) {
1375 extra_boxes.push_back(ibox);
1383 if (!extra_boxes_map.empty()) {
1385 auto& lbox = cbl2.
data();
1388 for (
auto const& [i, vb] : extra_boxes_map) {
1390 for (
int j = 1, nj =
int(vb.size()); j < nj; ++j) {
1391 lbox.push_back(vb[j]);
1392 procmap2.push_back(dm[i]);
1399 return std::make_tuple(
BoxArray(std::move(cbl)), cba2, dm2);
1410 level <
int(bc.
size()) &&
1411 level <
int(ratio.
size()+1));
1415 }
else if (level >=
int(smf.size()))
1417 auto const& [ba1, ba2, dm2] = get_clayout();
1422 mf.DistributionMap(), {0,0,0},
1424 cmf1.define(ba1, mf.DistributionMap(), ncomp, 0,
MFInfo(), *factory);
1429 cmf2.define(ba2, dm2, ncomp, 0,
MFInfo(), *factory2);
1434 cmf1.define(ba1, mf.DistributionMap(), ncomp, 0);
1436 cmf2.define(ba2, dm2, ncomp, 0);
1440 MF* p_mf_inside = (ba2.empty()) ? &cmf1 : &cmf2;
1442 geom, bc, bccomp, ratio, mapper, bcr, bcrcomp);
1443 if (&cmf1 != p_mf_inside) {
1444 cmf1.ParallelCopy(*p_mf_inside, geom[level-1].periodicity());
1446 Box domain_g = geom[level].growPeriodicDomain(nghost);
1447 domain_g.
convert(mf.ixType());
1448 FillPatchInterp(mf, dcomp, cmf1, 0, ncomp, nghost, geom[level-1], geom[level],
1449 domain_g, ratio[level-1], mapper, bcr, bcrcomp);
1452 int error_code = detail::FillPatchTwoLevels_doit(mf, nghost, time,
1453 smf[level-1], st[level-1],
1454 smf[level ], st[level ],
1455 scomp, dcomp, ncomp,
1456 geom[level-1], geom[level],
1457 bc[level-1], bccomp,
1459 ratio[level-1], mapper, bcr, bcrcomp,
1460 hook, hook, index_space,
true);
1461 if (error_code == 0) {
return; }
1463 auto const& [ba1, ba2, dm2] = get_clayout();
1469 mf.DistributionMap(), {0,0,0},
1471 cmf_tmp.define(ba1, mf.DistributionMap(), ncomp, 0,
MFInfo(), *factory);
1476 cmf_tmp.define(ba2, dm2, ncomp, 0,
MFInfo(), *factory);
1482 cmf_tmp.define(ba1, mf.DistributionMap(), ncomp, 0);
1484 cmf_tmp.define(ba2, dm2, ncomp, 0);
1489 geom, bc, bccomp, ratio, mapper, bcr, bcrcomp);
1495 for (
auto const* p : fmf) {
1499 for (
auto& a : fmf_raii) {
1504 detail::FillPatchTwoLevels_doit(mf, nghost, time,
1508 geom[level-1], geom[level],
1509 bc[level-1], bccomp,
1511 ratio[level-1], mapper, bcr, bcrcomp,
1512 hook, hook, index_space);
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:49
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:111
#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:564
IndexType ixType() const noexcept
Return index type of this BoxArray.
Definition AMReX_BoxArray.H:854
static bool SameRefs(const BoxArray &lhs, const BoxArray &rhs)
whether two BoxArrays share the same data
Definition AMReX_BoxArray.H:837
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition AMReX_BoxArray.H:611
void set(int i, const Box &ibox)
Set element i in this BoxArray to Box ibox.
Definition AMReX_BoxArray.cpp:878
A class for managing a List of Boxes that share a common IndexType. This class implements operations ...
Definition AMReX_BoxList.H:52
Vector< Box > & data() noexcept
Returns a reference to the Vector<Box>.
Definition AMReX_BoxList.H:215
void reserve(std::size_t n)
Definition AMReX_BoxList.H:90
void push_back(const Box &bn)
Append a Box to this BoxList.
Definition AMReX_BoxList.H:93
__host__ __device__ BoxND & grow(int i) noexcept
Definition AMReX_Box.H:649
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:155
__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:982
__host__ __device__ bool contains(const IntVectND< dim > &p) const noexcept
Return true if argument is contained within BoxND.
Definition AMReX_Box.H:216
__host__ __device__ IndexTypeND< dim > ixType() const noexcept
Return the indexing type.
Definition AMReX_Box.H:136
__host__ __device__ BoxND & refine(int ref_ratio) noexcept
Refine BoxND by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo*ratio...
Definition AMReX_Box.H:706
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
const Vector< int > & ProcessorMap() const noexcept
Returns a constant reference to the mapping of boxes in the underlying BoxArray to the CPU that holds...
Definition AMReX_DistributionMapping.cpp:49
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:2066
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:75
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:216
Periodicity periodicity() const noexcept
Definition AMReX_Geometry.H:361
bool isPeriodic(int dir) const noexcept
Is the domain periodic in the specified direction?
Definition AMReX_Geometry.H:337
__host__ __device__ constexpr CellIndex ixType(int dir) const noexcept
Returns the CellIndex in direction dir.
Definition AMReX_IndexType.H:117
__host__ __device__ constexpr int min() const noexcept
minimum (no absolute values) value
Definition AMReX_IntVect.H:324
__host__ __device__ constexpr 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:542
__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:771
__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:492
Definition AMReX_InterpBase.H:34
Definition AMReX_InterpBase.H:20
Box doit(const Box &fine) const override
Apply the coarse-box logic to the supplied fine box.
Definition AMReX_InterpBase.cpp:10
Virtual base class for interpolaters.
Definition AMReX_Interpolater.H:32
Definition AMReX_MFInterpolater.H:20
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
Definition AMReX_PhysBCFunct.H:127
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
std::unique_ptr< EBFArrayBoxFactory > makeEBFabFactory(const Geometry &a_geom, const BoxArray &a_ba, const DistributionMapping &a_dm, const Vector< int > &a_ngrow, EBSupport a_support)
Definition AMReX_EBFabFactory.cpp:217
__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
__host__ __device__ BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1289
__host__ __device__ BoxND< dim > refine(const BoxND< dim > &b, int ref_ratio) noexcept
Definition AMReX_Box.H:1468
std::array< T, N > Array
Definition AMReX_Array.H:26
const IndexSpace * TopIndexSpaceIfPresent() noexcept
Return the top IndexSpace if one has been built (nullptr otherwise).
Definition AMReX_EB2.cpp:93
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 Copy(FabArray< DFAB > &dst, FabArray< SFAB > const &src, int srccomp, int dstcomp, int numcomp, int nghost)
Definition AMReX_FabArray.H:180
DistributionMapping const & DistributionMap(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2867
IntVect nGrowVect(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2857
RunOn
Definition AMReX_GpuControl.H:65
bool ProperlyNested(const IntVect &ratio, const IntVect &blocking_factor, int ngrow, const IndexType &boxType, Interp *mapper)
Test if AMR grids are properly nested.
Definition AMReX_FillPatchUtil_I.H:35
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
__host__ __device__ BoxND< dim > shift(const BoxND< dim > &b, int dir, int nzones) noexcept
Return a BoxND with indices shifted by nzones in dir direction.
Definition AMReX_Box.H:1504
IndexTypeND< 3 > IndexType
IndexType is an alias for amrex::IndexTypeND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:36
void InterpFace(Interp *interp, MF const &mf_crse_patch, int crse_comp, MF &mf_refined_patch, int fine_comp, int ncomp, const IntVect &ratio, const iMF &solve_mask, const Geometry &crse_geom, const Geometry &fine_geom, int bcscomp, RunOn gpu_or_cpu, const Vector< BCRec > &bcs)
Definition AMReX_FillPatchUtil_I.H:274
void 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
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
__host__ __device__ bool almostEqual(T x, T y, int ulp=2)
Definition AMReX_Algorithm.H:123
bool TilingIfNotGPU() noexcept
Definition AMReX_MFIter.H:12
void InterpFromCoarseLevel(MF &mf, Real time, const MF &cmf, int scomp, int dcomp, int ncomp, const Geometry &cgeom, const Geometry &fgeom, BC &cbc, int cbccomp, BC &fbc, int fbccomp, const IntVect &ratio, Interp *mapper, const Vector< BCRec > &bcs, int bcscomp, const PreInterpHook &pre_interp={}, const PostInterpHook &post_interp={})
Fill with interpolation of coarse level data.
Definition AMReX_FillPatchUtil_I.H:979
void setBC(const Box &bx, const Box &domain, int src_comp, int dest_comp, int ncomp, const Vector< BCRec > &bc_dom, Vector< BCRec > &bcr) noexcept
Function for setting array of BCs.
Definition AMReX_BCRec.cpp:8
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:241
void FillPatchTwoLevels(MF &mf, IntVect const &nghost, Real time, const Vector< MF * > &cmf, const Vector< Real > &ct, const Vector< MF * > &fmf, const Vector< Real > &ft, int scomp, int dcomp, int ncomp, const Geometry &cgeom, const Geometry &fgeom, BC &cbc, int cbccomp, BC &fbc, int fbccomp, const IntVect &ratio, Interp *mapper, const Vector< BCRec > &bcs, int bcscomp, const PreInterpHook &pre_interp={}, const PostInterpHook &post_interp={})
FillPatch with data from the current level and the level below.
Definition AMReX_FillPatchUtil_I.H:791
BoxArray const & boxArray(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2862
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)
Helper that applies a MFInterpolater to fill a fine patch from a coarse patch.
Definition AMReX_FillPatchUtil.cpp:139
void FillPatchNLevels(MF &mf, int level, const IntVect &nghost, Real time, const Vector< Vector< MF * > > &smf, const Vector< Vector< Real > > &st, int scomp, int dcomp, int ncomp, const Vector< Geometry > &geom, Vector< BC > &bc, int bccomp, const Vector< IntVect > &ratio, Interp *mapper, const Vector< BCRec > &bcr, int bcrcomp)
FillPatch with data from AMR levels.
Definition AMReX_FillPatchUtil_I.H:1337
FabArray memory allocation information.
Definition AMReX_FabArray.H:68
Definition AMReX_FillPatchUtil.H:39