Skip to content

Commit

Permalink
Add a new InterFromCoarseLevel for ERF (#4150)
Browse files Browse the repository at this point in the history
  • Loading branch information
WeiqunZhang committed Sep 17, 2024
1 parent dc4b76e commit eecb9ea
Show file tree
Hide file tree
Showing 2 changed files with 73 additions and 0 deletions.
34 changes: 34 additions & 0 deletions Src/AmrCore/AMReX_FillPatchUtil.H
Original file line number Diff line number Diff line change
Expand Up @@ -694,6 +694,40 @@ namespace amrex
const PreInterpHook& pre_interp = {},
const PostInterpHook& post_interp = {});

/**
* \brief Fill with interpolation of coarse level data
*
* It's the CALLER's responsibility to make sure all ghost cells of the
* coarse MF needed for interpolation are filled already before calling
* this function. It's assumed that the fine level MultiFab mf's
* BoxArray is coarsenable by the refinement ratio. There is no support
* for EB.
*
* \tparam MF the MultiFab/FabArray type
* \tparam Interp spatial interpolater
*
* \param mf destination MF on the fine level
* \param nghost number of ghost cells of mf inside domain needed to be filled
* \param nghost_outside_domain number of ghost cells of mf outside domain needed to be filled
* \param cmf source MF on the coarse level
* \param scomp starting component of the source MF
* \param dcomp starting component of the destination MF
* \param ncomp number of components
* \param cgeom Geometry for the coarse level
* \param fgeom Geometry for the fine level
* \param ratio refinement ratio
* \param mapper spatial interpolater
* \param bcs boundar types for each component
*/
template <typename MF, typename Interp>
std::enable_if_t<IsFabArray<MF>::value>
InterpFromCoarseLevel (MF& mf, IntVect const& nghost,
IntVect const& nghost_outside_domain,
const MF& cmf, int scomp, int dcomp, int ncomp,
const Geometry& cgeom, const Geometry& fgeom,
const IntVect& ratio, Interp* mapper,
const Vector<BCRec>& bcs);

#ifndef BL_NO_FORT
enum InterpEM_t { InterpE, InterpB};

Expand Down
39 changes: 39 additions & 0 deletions Src/AmrCore/AMReX_FillPatchUtil_I.H
Original file line number Diff line number Diff line change
Expand Up @@ -1210,6 +1210,45 @@ InterpFromCoarseLevel (Array<MF*, AMREX_SPACEDIM> const& mf, IntVect const& ngho
}
}

template <typename MF, typename Interp>
std::enable_if_t<IsFabArray<MF>::value>
InterpFromCoarseLevel (MF& mf, IntVect const& nghost,
IntVect const& nghost_outside_domain,
const MF& cmf, int scomp, int dcomp, int ncomp,
const Geometry& cgeom, const Geometry& fgeom,
const IntVect& ratio, Interp* mapper,
const Vector<BCRec>& bcs)
{
const BoxArray& ba = mf.boxArray();
const DistributionMapping& dm = mf.DistributionMap();

const IndexType& typ = ba.ixType();
BL_ASSERT(typ == cmf.boxArray().ixType());

const InterpolaterBoxCoarsener& coarsener = mapper->BoxCoarsener(ratio);
Box tmp(-nghost, IntVect(32), typ);
Box tmp2 = coarsener.doit(tmp);
IntVect src_ghost = -tmp2.smallEnd();

tmp = Box(-nghost_outside_domain, IntVect(32), typ);
tmp2 = coarsener.doit(tmp);
IntVect src_ghost_outside_domain = -tmp2.smallEnd();

IntVect cghost = cmf.nGrowVect();
cghost.min(src_ghost);

AMREX_ALWAYS_ASSERT(cghost.allGE(src_ghost_outside_domain));

MF mf_crse_patch(amrex::coarsen(ba,ratio), dm, ncomp, src_ghost);
mf_crse_patch.ParallelCopy(cmf, scomp, 0, ncomp, cghost, src_ghost,
cgeom.periodicity());

Box fdomain_g = amrex::convert(fgeom.Domain(),typ);
fdomain_g.grow(nghost_outside_domain);
FillPatchInterp(mf, dcomp, mf_crse_patch, 0, ncomp, nghost, cgeom, fgeom,
fdomain_g, ratio, mapper, bcs, 0);
}

template <typename MF, typename BC, typename Interp>
std::enable_if_t<IsFabArray<MF>::value>
FillPatchNLevels (MF& mf, int level, const IntVect& nghost, Real time,
Expand Down

0 comments on commit eecb9ea

Please sign in to comment.