aboutsummaryrefslogtreecommitdiff
path: root/Source/Parallelization/WarpXComm.cpp
diff options
context:
space:
mode:
authorGravatar Axel Huebl <axel.huebl@plasma.ninja> 2019-09-27 11:14:00 -0700
committerGravatar Axel Huebl <axel.huebl@plasma.ninja> 2019-09-27 15:17:34 -0700
commitc9577f8d200d99c40d15d5ff0d2fdacbddc2026f (patch)
tree5606dba411201d581ae386bb9694b933c4b07bcd /Source/Parallelization/WarpXComm.cpp
parent9db3a53e959a6279c313a0d8d3e80288f6775b5a (diff)
downloadWarpX-c9577f8d200d99c40d15d5ff0d2fdacbddc2026f.tar.gz
WarpX-c9577f8d200d99c40d15d5ff0d2fdacbddc2026f.tar.zst
WarpX-c9577f8d200d99c40d15d5ff0d2fdacbddc2026f.zip
Rename, Profile & Remove TODO
Diffstat (limited to 'Source/Parallelization/WarpXComm.cpp')
-rw-r--r--Source/Parallelization/WarpXComm.cpp19
1 files changed, 10 insertions, 9 deletions
diff --git a/Source/Parallelization/WarpXComm.cpp b/Source/Parallelization/WarpXComm.cpp
index 40f8203a9..4f870e79c 100644
--- a/Source/Parallelization/WarpXComm.cpp
+++ b/Source/Parallelization/WarpXComm.cpp
@@ -2,7 +2,7 @@
#include <WarpX.H>
#include <WarpX_f.H>
#include <WarpXSumGuardCells.H>
-#include <Parallelization/CurrentSynchronize.H>
+#include <Parallelization/InterpolateCurrentFineToCoarse.H>
#include <algorithm>
#include <cstdlib>
@@ -337,7 +337,7 @@ WarpX::SyncCurrent ()
std::array< MultiFab*,3> crse { current_cp[lev][0].get(),
current_cp[lev][1].get(),
current_cp[lev][2].get() };
- SyncCurrent(fine, crse, refinement_ratio[0]);
+ interpolateCurrentFineToCoarse(fine, crse, refinement_ratio[0]);
}
// For each level
@@ -350,10 +350,11 @@ WarpX::SyncCurrent ()
}
void
-WarpX::SyncCurrent (const std::array<const amrex::MultiFab*,3>& fine,
- const std::array< amrex::MultiFab*,3>& coarse,
- int const refinement_ratio)
+WarpX::interpolateCurrentFineToCoarse ( std::array< amrex::MultiFab const *, 3 > const & fine,
+ std::array< amrex::MultiFab *, 3 > const & coarse,
+ int const refinement_ratio)
{
+ BL_PROFILE("InterpolateCurrentFineToCoarse()");
BL_ASSERT(refinement_ratio == 2);
const IntVect& ng = (fine[0]->nGrowVect() + 1) / refinement_ratio; // add equivalent no. of guards to coarse patch
@@ -372,11 +373,11 @@ WarpX::SyncCurrent (const std::array<const amrex::MultiFab*,3>& fine,
auto const & arrCoarse = coarse[idim]->array(mfi);
if( idim == 0 )
- amrex::ParallelFor( bx, WarpxSyncCurrent<0>(arrFine, arrCoarse, refinement_ratio) );
+ amrex::ParallelFor( bx, InterpolateCurrentFineToCoarse<0>(arrFine, arrCoarse, refinement_ratio) );
else if( idim == 1 )
- amrex::ParallelFor( bx, WarpxSyncCurrent<1>(arrFine, arrCoarse, refinement_ratio) );
+ amrex::ParallelFor( bx, InterpolateCurrentFineToCoarse<1>(arrFine, arrCoarse, refinement_ratio) );
else if( idim == 2 )
- amrex::ParallelFor( bx, WarpxSyncCurrent<2>(arrFine, arrCoarse, refinement_ratio) );
+ amrex::ParallelFor( bx, InterpolateCurrentFineToCoarse<2>(arrFine, arrCoarse, refinement_ratio) );
}
}
}
@@ -452,7 +453,7 @@ WarpX::RestrictCurrentFromFineToCoarsePatch (int lev)
std::array< MultiFab*,3> crse { current_cp[lev][0].get(),
current_cp[lev][1].get(),
current_cp[lev][2].get() };
- SyncCurrent(fine, crse, refinement_ratio[0]);
+ interpolateCurrentFineToCoarse(fine, crse, refinement_ratio[0]);
}
void