aboutsummaryrefslogtreecommitdiff
path: root/Source/Particles/WarpXParticleContainer.cpp
diff options
context:
space:
mode:
authorGravatar Andrew Myers <atmyers@lbl.gov> 2021-10-18 14:06:08 -0700
committerGravatar GitHub <noreply@github.com> 2021-10-18 14:06:08 -0700
commit2bc1322aa76f7026e5b42639bb3d12125da2407c (patch)
tree7b8dbd4285fea5d47bf2b2b18a24af767342f34a /Source/Particles/WarpXParticleContainer.cpp
parent3e9903c07561c24bb8309a7bfa8feeb76ba2eebd (diff)
downloadWarpX-2bc1322aa76f7026e5b42639bb3d12125da2407c.tar.gz
WarpX-2bc1322aa76f7026e5b42639bb3d12125da2407c.tar.zst
WarpX-2bc1322aa76f7026e5b42639bb3d12125da2407c.zip
Option to do single precision mesh communication. (#2294)
* option to use single precision guard cell exhanges * add missing files * fix namespace issue * change precision of comms to float * ParmParse the single_precision_comms flag * set back to real * test * make sure dst is filled * nGrow -> nGrowVect * restore float * don't override valid cells * single precision mesh * whitespace * wrap SumBoundary * Wrap additional uses of FillBoundary. * catch missing copies of ParallelCopy * missing OverrideSyncs * add wrapper for iMultifab * fix typo * moar typos * typo * strip single_precision_mesh option * fix original copy * update fusible syntax Co-authored-by: Weiqun Zhang <WeiqunZhang@lbl.gov> * Fix: Single Precision Builds Should not copy around data for `do_single_precision_comms` * Docs: warpx.do_single_precision_comms * initialize this tmp multifab to 0.0 * fix tiny profile label * remove orig copies, they are only correct for FillBoundary * loosen tolerance for space charge tests for single precision * missing some setVal * another missing setVal * missing setVal * add wrapper for new version of sumboundary * add explicit cast to silence compiler warning * add a test for single precision comms * revert change to test precision * add benchmark for single precision comms test * restore tolerance I removed by mistake * tolerance * copyright headers * drop tolerance for single precision tests in default analysis script * missing python module * bump tol again * fix bad merge * Apply suggestions from code review Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> Co-authored-by: Remi Lehe <remi.lehe@normalesup.org> Co-authored-by: Axel Huebl <axel.huebl@plasma.ninja> Co-authored-by: Weiqun Zhang <WeiqunZhang@lbl.gov>
Diffstat (limited to 'Source/Particles/WarpXParticleContainer.cpp')
-rw-r--r--Source/Particles/WarpXParticleContainer.cpp21
1 files changed, 14 insertions, 7 deletions
diff --git a/Source/Particles/WarpXParticleContainer.cpp b/Source/Particles/WarpXParticleContainer.cpp
index 5a76a24df..8610238db 100644
--- a/Source/Particles/WarpXParticleContainer.cpp
+++ b/Source/Particles/WarpXParticleContainer.cpp
@@ -13,6 +13,7 @@
#include "Deposition/CurrentDeposition.H"
#include "Pusher/GetAndSetPosition.H"
#include "Pusher/UpdatePosition.H"
+#include "Parallelization/WarpXCommUtil.H"
#include "ParticleBoundaries_K.H"
#include "Utils/CoarsenMR.H"
#include "Utils/WarpXAlgorithmSelection.H"
@@ -783,23 +784,29 @@ WarpXParticleContainer::DepositCharge (amrex::Vector<std::unique_ptr<amrex::Mult
#endif
// Exchange guard cells
- if (local == false) rho[lev]->SumBoundary(m_gdb->Geom(lev).periodicity());
+ if (local == false) {
+ WarpXCommUtil::SumBoundary(*rho[lev], m_gdb->Geom(lev).periodicity());
+ }
}
// Now that the charge has been deposited at each level,
// we average down from fine to crse
if (interpolate_across_levels)
{
- for (int lev = finest_level - 1; lev >= 0; --lev)
- {
+ for (int lev = finest_level - 1; lev >= 0; --lev) {
const DistributionMapping& fine_dm = rho[lev+1]->DistributionMap();
BoxArray coarsened_fine_BA = rho[lev+1]->boxArray();
coarsened_fine_BA.coarsen(m_gdb->refRatio(lev));
MultiFab coarsened_fine_data(coarsened_fine_BA, fine_dm, rho[lev+1]->nComp(), 0);
coarsened_fine_data.setVal(0.0);
- const int refinement_ratio = 2;
- CoarsenMR::Coarsen(coarsened_fine_data, *rho[lev+1], IntVect(refinement_ratio));
- rho[lev]->ParallelAdd(coarsened_fine_data, m_gdb->Geom(lev).periodicity());
+
+ int const refinement_ratio = 2;
+
+ CoarsenMR::Coarsen( coarsened_fine_data, *rho[lev+1], IntVect(refinement_ratio) );
+ WarpXCommUtil::ParallelAdd(*rho[lev], coarsened_fine_data, 0, 0, rho[lev]->nComp(),
+ amrex::IntVect::TheZeroVector(),
+ amrex::IntVect::TheZeroVector(),
+ m_gdb->Geom(lev).periodicity());
}
}
}
@@ -860,7 +867,7 @@ WarpXParticleContainer::GetChargeDensity (int lev, bool local)
WarpX::GetInstance().ApplyInverseVolumeScalingToChargeDensity(rho.get(), lev);
#endif
- if (local == false) rho->SumBoundary(gm.periodicity());
+ if (local == false) { WarpXCommUtil::SumBoundary(*rho, gm.periodicity()); }
return rho;
}