aboutsummaryrefslogtreecommitdiff
path: root/Source/Utils/CoarsenIO.cpp
diff options
context:
space:
mode:
authorGravatar Edoardo Zoni <59625522+EZoni@users.noreply.github.com> 2020-05-11 08:01:49 -0700
committerGravatar GitHub <noreply@github.com> 2020-05-11 08:01:49 -0700
commite2eb322f5856994102237c8903781a5b3f77dd58 (patch)
tree2744f0ae6735ce2c427f344f72a90724c4f7f8af /Source/Utils/CoarsenIO.cpp
parentec3c2aab43130d96c599e38474e3467b47a4bfcd (diff)
downloadWarpX-e2eb322f5856994102237c8903781a5b3f77dd58.tar.gz
WarpX-e2eb322f5856994102237c8903781a5b3f77dd58.tar.zst
WarpX-e2eb322f5856994102237c8903781a5b3f77dd58.zip
Generalize coarsening for MR (#945)
* Move interpolation functions for MR to new folder * Preparatory clean-up of old namespace Average for future MR functions * Add interpolation for MR in new namespace Coarsen * Change file names Average.H/.cpp to Coarsen.H/.cpp * Remove Source/Parallelization/WarpXComm.H (not necessary anymore) * Coarsening for IO and MR in separate files * Clean up IO and MR Coarsen namespaces * Remove old interpolation functions (charge and current) * Void commit: trigger Travis CI build * Fix GPU build * Void commit: trigger Travis CI build * Add Python script to test interpolation points/weights in 1D * Move using-directives inside namespaces * Add Doxygen documentation and comments * Minor style changes * Improve new Python script
Diffstat (limited to 'Source/Utils/CoarsenIO.cpp')
-rw-r--r--Source/Utils/CoarsenIO.cpp106
1 files changed, 106 insertions, 0 deletions
diff --git a/Source/Utils/CoarsenIO.cpp b/Source/Utils/CoarsenIO.cpp
new file mode 100644
index 000000000..b75d054b4
--- /dev/null
+++ b/Source/Utils/CoarsenIO.cpp
@@ -0,0 +1,106 @@
+#include "CoarsenIO.H"
+
+using namespace amrex;
+
+void
+CoarsenIO::Loop ( MultiFab& mf_dst,
+ const MultiFab& mf_src,
+ const int dcomp,
+ const int scomp,
+ const int ncomp,
+ const int ngrow,
+ 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,
+ "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),
+ "source fine MultiFab does not have enough guard cells for this interpolation" );
+
+ // Auxiliary integer arrays (always 3D)
+ Gpu::ManagedVector<int> sf_gpuarr, sc_gpuarr, cr_gpuarr;
+ sf_gpuarr.resize( 3 ); // staggering of source fine MultiFab
+ sc_gpuarr.resize( 3 ); // staggering of destination coarse MultiFab
+ cr_gpuarr.resize( 3 ); // coarsening ratio
+
+ sf_gpuarr[0] = stag_src[0];
+ sf_gpuarr[1] = stag_src[1];
+#if (AMREX_SPACEDIM == 2)
+ sf_gpuarr[2] = 0;
+#elif (AMREX_SPACEDIM == 3)
+ sf_gpuarr[2] = stag_src[2];
+#endif
+
+ sc_gpuarr[0] = stag_dst[0];
+ sc_gpuarr[1] = stag_dst[1];
+#if (AMREX_SPACEDIM == 2)
+ sc_gpuarr[2] = 0;
+#elif (AMREX_SPACEDIM == 3)
+ sc_gpuarr[2] = stag_dst[2];
+#endif
+
+ cr_gpuarr[0] = crse_ratio[0];
+ cr_gpuarr[1] = crse_ratio[1];
+#if (AMREX_SPACEDIM == 2)
+ cr_gpuarr[2] = 1;
+#elif (AMREX_SPACEDIM == 3)
+ cr_gpuarr[2] = crse_ratio[2];
+#endif
+
+ int const* const AMREX_RESTRICT sf = sf_gpuarr.data();
+ int const* const AMREX_RESTRICT sc = sc_gpuarr.data();
+ int const* const AMREX_RESTRICT cr = cr_gpuarr.data();
+
+#ifdef _OPENMP
+#pragma omp parallel if (Gpu::notInLaunchRegion())
+#endif
+ // Loop over boxes (or tiles if not on GPU)
+ for (MFIter mfi( mf_dst, TilingIfNotGPU() ); mfi.isValid(); ++mfi)
+ {
+ // Tiles defined at the coarse level
+ const Box& bx = mfi.growntilebox( ngrow );
+ Array4<Real> const& arr_dst = mf_dst.array( mfi );
+ Array4<Real const> const& arr_src = mf_src.const_array( mfi );
+ ParallelFor( bx, ncomp,
+ [=] AMREX_GPU_DEVICE( int i, int j, int k, int n )
+ {
+ arr_dst(i,j,k,n+dcomp) = CoarsenIO::Interp(
+ arr_src, sf, sc, cr, i, j, k, n+scomp );
+ } );
+ }
+}
+
+void
+CoarsenIO::Coarsen ( MultiFab& mf_dst,
+ const MultiFab& mf_src,
+ const int dcomp,
+ const int scomp,
+ const int ncomp,
+ const int ngrow,
+ const IntVect crse_ratio )
+{
+ BL_PROFILE( "CoarsenIO::Coarsen" );
+
+ // Convert BoxArray of source MultiFab to staggering of destination MultiFab and coarsen it
+ BoxArray ba_tmp = amrex::convert( mf_src.boxArray(), mf_dst.ixType().toIntVect() );
+ AMREX_ALWAYS_ASSERT_WITH_MESSAGE( ba_tmp.coarsenable( crse_ratio ),
+ "source MultiFab converted to staggering of destination MultiFab is not coarsenable" );
+ 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 );
+ 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 );
+ // 3) copy from mf_tmp to mf_dst (with different BoxArray or DistributionMapping)
+ mf_dst.copy( mf_tmp, 0, dcomp, ncomp );
+ }
+}