Skip to content

Commit

Permalink
Fix MFInterpolater
Browse files Browse the repository at this point in the history
  • Loading branch information
WeiqunZhang committed Jul 16, 2023
1 parent e182e35 commit 04ad7f3
Show file tree
Hide file tree
Showing 2 changed files with 20 additions and 16 deletions.
8 changes: 1 addition & 7 deletions Src/AmrCore/AMReX_Interpolater.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -850,13 +850,7 @@ CellConservativeQuartic::interp (const FArrayBox& crse,
RunOn runon)
{
BL_PROFILE("CellConservativeQuartic::interp()");
AMREX_ASSERT(ratio[0] == 2);
#if (AMREX_SPACEDIM >= 2)
AMREX_ASSERT(ratio[0] == ratio[1]);
#endif
#if (AMREX_SPACEDIM == 3)
AMREX_ASSERT(ratio[1] == ratio[2]);
#endif
AMREX_ASSERT(ratio == 2);
amrex::ignore_unused(ratio);

//
Expand Down
28 changes: 19 additions & 9 deletions Src/AmrCore/AMReX_MFInterpolater.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -99,6 +99,11 @@ MFCellConsLinInterp::interp (MultiFab const& crsemf, int ccomp, MultiFab& finemf

Box const& cdomain = cgeom.Domain();

IntVect minus1;
for (int dim = 0; dim < AMREX_SPACEDIM; dim++) {
minus1[dim] = (ratio[dim] > 1) ? -1 : 0;
}

#ifdef AMREX_USE_GPU
if (Gpu::inLaunchRegion()) {
MultiFab crse_tmp(crsemf.boxArray(), crsemf.DistributionMap(), AMREX_SPACEDIM*nc, 0);
Expand All @@ -116,14 +121,14 @@ MFCellConsLinInterp::interp (MultiFab const& crsemf, int ccomp, MultiFab& finemf
Real drf = fgeom.CellSize(0);
Real rlo = fgeom.Offset(0);
if (do_linear_limiting) {
ParallelFor(crsemf, IntVect(-1),
ParallelFor(crsemf, minus1,
[=] AMREX_GPU_DEVICE (int box_no, int i, int, int) noexcept
{
mf_cell_cons_lin_interp_llslope(i,0,0, tmp[box_no], crse[box_no], ccomp, nc,
cdomain, ratio, pbc);
});
} else {
ParallelFor(crsemf, IntVect(-1), nc,
ParallelFor(crsemf, minus1, nc,
[=] AMREX_GPU_DEVICE (int box_no, int i, int, int, int n) noexcept
{
mf_cell_cons_lin_interp_mcslope_sph(i, n, tmp[box_no], crse[box_no], ccomp, nc,
Expand All @@ -145,14 +150,14 @@ MFCellConsLinInterp::interp (MultiFab const& crsemf, int ccomp, MultiFab& finemf
Real drf = fgeom.CellSize(0);
Real rlo = fgeom.Offset(0);
if (do_linear_limiting) {
ParallelFor(crsemf, IntVect(-1),
ParallelFor(crsemf, minus1,
[=] AMREX_GPU_DEVICE (int box_no, int i, int j, int) noexcept
{
mf_cell_cons_lin_interp_llslope(i,j,0, tmp[box_no], crse[box_no], ccomp, nc,
cdomain, ratio, pbc);
});
} else {
ParallelFor(crsemf, IntVect(-1), nc,
ParallelFor(crsemf, minus1, nc,
[=] AMREX_GPU_DEVICE (int box_no, int i, int j, int, int n) noexcept
{
mf_cell_cons_lin_interp_mcslope_rz(i,j,n, tmp[box_no], crse[box_no], ccomp, nc,
Expand All @@ -172,14 +177,14 @@ MFCellConsLinInterp::interp (MultiFab const& crsemf, int ccomp, MultiFab& finemf
#endif
{
if (do_linear_limiting) {
ParallelFor(crsemf, IntVect(-1),
ParallelFor(crsemf, minus1,
[=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept
{
mf_cell_cons_lin_interp_llslope(i,j,k, tmp[box_no], crse[box_no], ccomp, nc,
cdomain, ratio, pbc);
});
} else {
ParallelFor(crsemf, IntVect(-1), nc,
ParallelFor(crsemf, minus1, nc,
[=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
{
mf_cell_cons_lin_interp_mcslope(i,j,k,n, tmp[box_no], crse[box_no], ccomp, nc,
Expand Down Expand Up @@ -212,7 +217,7 @@ MFCellConsLinInterp::interp (MultiFab const& crsemf, int ccomp, MultiFab& finemf
auto const& fine = finemf.array(mfi);
auto const& crse = crsemf.const_array(mfi);

Box const& cbox = amrex::grow(crsemf[mfi].box(), -1);
Box const& cbox = amrex::grow(crsemf[mfi].box(), minus1);
tmpfab.resize(cbox, AMREX_SPACEDIM*nc);
auto const& tmp = tmpfab.array();
auto const& ctmp = tmpfab.const_array();
Expand Down Expand Up @@ -335,6 +340,11 @@ MFCellConsLinMinmaxLimitInterp::interp (MultiFab const& crsemf, int ccomp, Multi

Box const& cdomain = cgeom.Domain();

IntVect minus1;
for (int dim = 0; dim < AMREX_SPACEDIM; dim++) {
minus1[dim] = (ratio[dim] > 1) ? -1 : 0;
}

#ifdef AMREX_USE_GPU
if (Gpu::inLaunchRegion()) {
MultiFab crse_tmp(crsemf.boxArray(), crsemf.DistributionMap(), AMREX_SPACEDIM*nc, 0);
Expand All @@ -347,7 +357,7 @@ MFCellConsLinMinmaxLimitInterp::interp (MultiFab const& crsemf, int ccomp, Multi
BCRec const* pbc = d_bc.data();
Gpu::copyAsync(Gpu::hostToDevice, bcs.begin()+bcomp, bcs.begin()+bcomp+nc, d_bc.begin());

ParallelFor(crsemf, IntVect(-1),
ParallelFor(crsemf, minus1,
[=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept
{
mf_cell_cons_lin_interp_limit_minmax_llslope(i,j,k, tmp[box_no], crse[box_no], ccomp, nc,
Expand Down Expand Up @@ -378,7 +388,7 @@ MFCellConsLinMinmaxLimitInterp::interp (MultiFab const& crsemf, int ccomp, Multi
auto const& fine = finemf.array(mfi);
auto const& crse = crsemf.const_array(mfi);

Box const& cbox = amrex::grow(crsemf[mfi].box(), -1);
Box const& cbox = amrex::grow(crsemf[mfi].box(), minus1);
tmpfab.resize(cbox, AMREX_SPACEDIM*nc);
auto const& tmp = tmpfab.array();
auto const& ctmp = tmpfab.const_array();
Expand Down

0 comments on commit 04ad7f3

Please sign in to comment.