aboutsummaryrefslogtreecommitdiff
path: root/Source/Utils/CoarsenIO.cpp
diff options
context:
space:
mode:
authorGravatar David Grote <grote1@llnl.gov> 2021-03-17 15:41:53 -0700
committerGravatar GitHub <noreply@github.com> 2021-03-17 15:41:53 -0700
commit93cdd24708e8e975409daebce758ddbcaaac9325 (patch)
tree8a10ebd185306f613c4509126e211fd5a6e1ff93 /Source/Utils/CoarsenIO.cpp
parent161a0f9af3251aa2691ec90753de70632c5226bc (diff)
downloadWarpX-93cdd24708e8e975409daebce758ddbcaaac9325.tar.gz
WarpX-93cdd24708e8e975409daebce758ddbcaaac9325.tar.zst
WarpX-93cdd24708e8e975409daebce758ddbcaaac9325.zip
Replaced almost all nGrow with nGrowVect (#1801)
Diffstat (limited to 'Source/Utils/CoarsenIO.cpp')
-rw-r--r--Source/Utils/CoarsenIO.cpp31
1 files changed, 25 insertions, 6 deletions
diff --git a/Source/Utils/CoarsenIO.cpp b/Source/Utils/CoarsenIO.cpp
index c98fb902e..b51419bbd 100644
--- a/Source/Utils/CoarsenIO.cpp
+++ b/Source/Utils/CoarsenIO.cpp
@@ -8,17 +8,17 @@ CoarsenIO::Loop ( MultiFab& mf_dst,
const int dcomp,
const int scomp,
const int ncomp,
- const int ngrow,
+ const IntVect ngrowvect,
const IntVect crse_ratio )
{
// Staggering of source fine MultiFab and destination coarse MultiFab
const IntVect stag_src = mf_src.boxArray().ixType().toIntVect();
const IntVect stag_dst = mf_dst.boxArray().ixType().toIntVect();
- if ( crse_ratio > IntVect(1) ) AMREX_ALWAYS_ASSERT_WITH_MESSAGE( ngrow == 0,
+ if ( crse_ratio > IntVect(1) ) AMREX_ALWAYS_ASSERT_WITH_MESSAGE( ngrowvect == IntVect(0),
"option of filling guard cells of destination MultiFab with coarsening not supported for this interpolation" );
- AMREX_ALWAYS_ASSERT_WITH_MESSAGE( mf_src.nGrowVect() >= stag_dst-stag_src+IntVect(ngrow),
+ AMREX_ALWAYS_ASSERT_WITH_MESSAGE( mf_src.nGrowVect() >= stag_dst-stag_src+ngrowvect,
"source fine MultiFab does not have enough guard cells for this interpolation" );
// Auxiliary integer arrays (always 3D)
@@ -57,7 +57,7 @@ CoarsenIO::Loop ( MultiFab& mf_dst,
for (MFIter mfi( mf_dst, TilingIfNotGPU() ); mfi.isValid(); ++mfi)
{
// Tiles defined at the coarse level
- const Box& bx = mfi.growntilebox( ngrow );
+ const Box& bx = mfi.growntilebox( ngrowvect );
Array4<Real> const& arr_dst = mf_dst.array( mfi );
Array4<Real const> const& arr_src = mf_src.const_array( mfi );
ParallelFor( bx, ncomp,
@@ -78,6 +78,25 @@ CoarsenIO::Coarsen ( MultiFab& mf_dst,
const int ngrow,
const IntVect crse_ratio )
{
+ amrex::IntVect ngrowvect(ngrow);
+ Coarsen(mf_dst,
+ mf_src,
+ dcomp,
+ scomp,
+ ncomp,
+ ngrowvect,
+ crse_ratio);
+}
+
+void
+CoarsenIO::Coarsen ( MultiFab& mf_dst,
+ const MultiFab& mf_src,
+ const int dcomp,
+ const int scomp,
+ const int ncomp,
+ const IntVect ngrowvect,
+ const IntVect crse_ratio )
+{
BL_PROFILE("CoarsenIO::Coarsen()");
// Convert BoxArray of source MultiFab to staggering of destination MultiFab and coarsen it
@@ -87,14 +106,14 @@ CoarsenIO::Coarsen ( MultiFab& mf_dst,
ba_tmp.coarsen( crse_ratio );
if ( ba_tmp == mf_dst.boxArray() and mf_src.DistributionMap() == mf_dst.DistributionMap() )
- CoarsenIO::Loop( mf_dst, mf_src, dcomp, scomp, ncomp, ngrow, crse_ratio );
+ CoarsenIO::Loop( mf_dst, mf_src, dcomp, scomp, ncomp, ngrowvect, crse_ratio );
else
{
// Cannot coarsen into MultiFab with different BoxArray or DistributionMapping:
// 1) create temporary MultiFab on coarsened version of source BoxArray with same DistributionMapping
MultiFab mf_tmp( ba_tmp, mf_src.DistributionMap(), ncomp, 0, MFInfo(), FArrayBoxFactory() );
// 2) interpolate from mf_src to mf_tmp (start writing into component 0)
- CoarsenIO::Loop( mf_tmp, mf_src, 0, scomp, ncomp, ngrow, crse_ratio );
+ CoarsenIO::Loop( mf_tmp, mf_src, 0, scomp, ncomp, ngrowvect, crse_ratio );
// 3) copy from mf_tmp to mf_dst (with different BoxArray or DistributionMapping)
mf_dst.copy( mf_tmp, 0, dcomp, ncomp );
}