aboutsummaryrefslogtreecommitdiff
path: root/Source
diff options
context:
space:
mode:
Diffstat (limited to 'Source')
-rw-r--r--Source/Diagnostics/ReducedDiags/LoadBalanceCosts.cpp14
-rw-r--r--Source/Evolve/WarpXEvolve.cpp4
-rw-r--r--Source/FieldSolver/WarpX_QED_Field_Pushers.cpp2
-rw-r--r--Source/Initialization/WarpXInitData.cpp6
-rw-r--r--Source/Laser/LaserParticleContainer.cpp2
-rw-r--r--Source/Parallelization/WarpXRegrid.cpp25
-rw-r--r--Source/Particles/PhysicalParticleContainer.cpp6
-rw-r--r--Source/Particles/WarpXParticleContainer.cpp2
-rw-r--r--Source/WarpX.H8
-rw-r--r--Source/WarpX.cpp4
10 files changed, 31 insertions, 42 deletions
diff --git a/Source/Diagnostics/ReducedDiags/LoadBalanceCosts.cpp b/Source/Diagnostics/ReducedDiags/LoadBalanceCosts.cpp
index 3d893977a..c65a38abf 100644
--- a/Source/Diagnostics/ReducedDiags/LoadBalanceCosts.cpp
+++ b/Source/Diagnostics/ReducedDiags/LoadBalanceCosts.cpp
@@ -25,7 +25,7 @@ void LoadBalanceCosts::ComputeDiags (int step)
// get WarpX class object
auto& warpx = WarpX::GetInstance();
- const amrex::Vector<amrex::Real>* cost = warpx.getCosts(0);
+ const amrex::LayoutData<amrex::Real>* cost = warpx.getCosts(0);
// judge if the diags should be done
// costs is initialized only if we're doing load balance
@@ -51,20 +51,12 @@ void LoadBalanceCosts::ComputeDiags (int step)
m_data.assign(dataSize, 0.0);
// read in WarpX costs to local copy; compute if using `Heuristic` update
- amrex::Vector<std::unique_ptr<amrex::Vector<amrex::Real> > > costs;
+ amrex::Vector<std::unique_ptr<amrex::LayoutData<amrex::Real> > > costs;
costs.resize(nLevels);
for (int lev = 0; lev < nLevels; ++lev)
{
- costs[lev].reset(new amrex::Vector<Real>);
- const int nBoxesLev = warpx.getCosts(lev)->size();
- costs[lev]->resize(nBoxesLev);
- for (int i = 0; i < nBoxesLev; ++i)
- {
- // If `Heuristic` update, this fills with zeros;
- // if `Timers` update, this fills with timer-based costs
- (*costs[lev])[i] = (*warpx.getCosts(lev))[i];
- }
+ costs[lev].reset(new amrex::LayoutData<Real>(*warpx.getCosts(lev)));
}
if (warpx.load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Heuristic)
diff --git a/Source/Evolve/WarpXEvolve.cpp b/Source/Evolve/WarpXEvolve.cpp
index bc1e31806..8ca059184 100644
--- a/Source/Evolve/WarpXEvolve.cpp
+++ b/Source/Evolve/WarpXEvolve.cpp
@@ -61,7 +61,7 @@ WarpX::Evolve (int numsteps)
if (warpx_py_beforestep) warpx_py_beforestep();
#endif
- amrex::Vector<amrex::Real>* cost = WarpX::getCosts(0);
+ amrex::LayoutData<amrex::Real>* cost = WarpX::getCosts(0);
if (cost) {
#ifdef WARPX_USE_PSATD
amrex::Abort("LoadBalance for PSATD: TODO");
@@ -82,7 +82,7 @@ WarpX::Evolve (int numsteps)
// (Giving more importance to most recent costs; only needed
// for timers update, heuristic load balance considers the
// instantaneous costs)
- for (int i=0; i<cost->size(); ++i)
+ for (int i : cost->IndexArray())
{
(*cost)[i] *= (1. - 2./load_balance_intervals.localPeriod(step+1));
}
diff --git a/Source/FieldSolver/WarpX_QED_Field_Pushers.cpp b/Source/FieldSolver/WarpX_QED_Field_Pushers.cpp
index f9db7b39a..ddb0359e5 100644
--- a/Source/FieldSolver/WarpX_QED_Field_Pushers.cpp
+++ b/Source/FieldSolver/WarpX_QED_Field_Pushers.cpp
@@ -90,7 +90,7 @@ WarpX::Hybrid_QED_Push (int lev, PatchType patch_type, Real a_dt)
Jz = current_cp[lev][2].get();
}
- amrex::Vector<amrex::Real>* cost = WarpX::getCosts(lev);
+ amrex::LayoutData<amrex::Real>* cost = WarpX::getCosts(lev);
// Loop through the grids, and over the tiles within each grid
#ifdef _OPENMP
diff --git a/Source/Initialization/WarpXInitData.cpp b/Source/Initialization/WarpXInitData.cpp
index 3afd16667..0943d549e 100644
--- a/Source/Initialization/WarpXInitData.cpp
+++ b/Source/Initialization/WarpXInitData.cpp
@@ -420,9 +420,9 @@ WarpX::InitLevelData (int lev, Real /*time*/)
}
if (costs[lev]) {
- std::fill((*costs[lev]).begin(),
- (*costs[lev]).end(),
- 0.0);
+ for (int i : costs[lev]->IndexArray()) {
+ (*costs[lev])[i] = 0.0;
+ }
}
}
diff --git a/Source/Laser/LaserParticleContainer.cpp b/Source/Laser/LaserParticleContainer.cpp
index a8847c79a..6100a4408 100644
--- a/Source/Laser/LaserParticleContainer.cpp
+++ b/Source/Laser/LaserParticleContainer.cpp
@@ -407,7 +407,7 @@ LaserParticleContainer::Evolve (int lev,
BL_ASSERT(OnSameGrids(lev,jx));
- amrex::Vector<amrex::Real>* cost = WarpX::getCosts(lev);
+ amrex::LayoutData<amrex::Real>* cost = WarpX::getCosts(lev);
#ifdef _OPENMP
#pragma omp parallel
diff --git a/Source/Parallelization/WarpXRegrid.cpp b/Source/Parallelization/WarpXRegrid.cpp
index f04c59c78..6c6710121 100644
--- a/Source/Parallelization/WarpXRegrid.cpp
+++ b/Source/Parallelization/WarpXRegrid.cpp
@@ -35,12 +35,6 @@ WarpX::LoadBalance ()
const int nLevels = finestLevel();
for (int lev = 0; lev <= nLevels; ++lev)
{
- LayoutData<Real> costsLayoutData(boxArray(lev), DistributionMap(lev));
- for (auto i : costsLayoutData.IndexArray())
- {
- costsLayoutData[i] = (*costs[lev])[i];
- }
-
// Compute the new distribution mapping
DistributionMapping newdm;
const amrex::Real nboxes = costs[lev]->size();
@@ -52,11 +46,11 @@ WarpX::LoadBalance ()
amrex::Real proposedEfficiency = 0.0;
newdm = (load_balance_with_sfc)
- ? DistributionMapping::makeSFC(costsLayoutData,
+ ? DistributionMapping::makeSFC(*costs[lev],
currentEfficiency, proposedEfficiency,
false,
ParallelDescriptor::IOProcessorNumber())
- : DistributionMapping::makeKnapSack(costsLayoutData,
+ : DistributionMapping::makeKnapSack(*costs[lev],
currentEfficiency, proposedEfficiency,
nmax,
false,
@@ -283,9 +277,11 @@ WarpX::RemakeLevel (int lev, Real /*time*/, const BoxArray& ba, const Distributi
if (costs[lev] != nullptr)
{
- costs[lev].reset(new amrex::Vector<Real>);
- const int nboxes = Efield_fp[lev][0].get()->size();
- costs[lev]->resize(nboxes, 0.0);
+ costs[lev].reset(new amrex::LayoutData<Real>(ba, dm));
+ for (int i : costs[lev]->IndexArray())
+ {
+ (*costs[lev])[i] = 0.0;
+ }
}
SetDistributionMap(lev, dm);
@@ -299,7 +295,7 @@ WarpX::RemakeLevel (int lev, Real /*time*/, const BoxArray& ba, const Distributi
}
void
-WarpX::ComputeCostsHeuristic (amrex::Vector<std::unique_ptr<amrex::Vector<amrex::Real> > >& a_costs)
+WarpX::ComputeCostsHeuristic (amrex::Vector<std::unique_ptr<amrex::LayoutData<amrex::Real> > >& a_costs)
{
for (int lev = 0; lev <= finest_level; ++lev)
{
@@ -333,6 +329,9 @@ WarpX::ResetCosts ()
{
for (int lev = 0; lev <= finest_level; ++lev)
{
- costs[lev]->assign((*costs[lev]).size(), 0.0);
+ for (int i : costs[lev]->IndexArray())
+ {
+ (*costs[lev])[i] = 0.0;
+ }
}
}
diff --git a/Source/Particles/PhysicalParticleContainer.cpp b/Source/Particles/PhysicalParticleContainer.cpp
index ce2871b3f..ebfea28c7 100644
--- a/Source/Particles/PhysicalParticleContainer.cpp
+++ b/Source/Particles/PhysicalParticleContainer.cpp
@@ -478,7 +478,7 @@ PhysicalParticleContainer::AddPlasma (int lev, RealBox part_realbox)
defineAllParticleTiles();
- amrex::Vector<amrex::Real>* cost = WarpX::getCosts(lev);
+ amrex::LayoutData<amrex::Real>* cost = WarpX::getCosts(lev);
const int nlevs = numLevels();
static bool refine_injection = false;
@@ -950,7 +950,7 @@ PhysicalParticleContainer::FieldGather (int lev,
{
BL_ASSERT(OnSameGrids(lev,Ex));
- amrex::Vector<amrex::Real>* cost = WarpX::getCosts(lev);
+ amrex::LayoutData<amrex::Real>* cost = WarpX::getCosts(lev);
#ifdef _OPENMP
#pragma omp parallel
@@ -1021,7 +1021,7 @@ PhysicalParticleContainer::Evolve (int lev,
BL_ASSERT(OnSameGrids(lev,jx));
- amrex::Vector<amrex::Real>* cost = WarpX::getCosts(lev);
+ amrex::LayoutData<amrex::Real>* cost = WarpX::getCosts(lev);
const iMultiFab* current_masks = WarpX::CurrentBufferMasks(lev);
const iMultiFab* gather_masks = WarpX::GatherBufferMasks(lev);
diff --git a/Source/Particles/WarpXParticleContainer.cpp b/Source/Particles/WarpXParticleContainer.cpp
index 2344d7985..49c8264e9 100644
--- a/Source/Particles/WarpXParticleContainer.cpp
+++ b/Source/Particles/WarpXParticleContainer.cpp
@@ -768,7 +768,7 @@ WarpXParticleContainer::PushX (int lev, amrex::Real dt)
if (do_not_push) return;
- amrex::Vector<amrex::Real>* cost = WarpX::getCosts(lev);
+ amrex::LayoutData<amrex::Real>* cost = WarpX::getCosts(lev);
#ifdef _OPENMP
#pragma omp parallel
diff --git a/Source/WarpX.H b/Source/WarpX.H
index abca491f4..a19d9bf81 100644
--- a/Source/WarpX.H
+++ b/Source/WarpX.H
@@ -242,7 +242,7 @@ public:
/** get low-high-low-high-... vector for each direction indicating if mother grid PMLs are enabled */
std::vector<bool> getPMLdirections() const;
- static amrex::Vector<amrex::Real>* getCosts (int lev) {
+ static amrex::LayoutData<amrex::Real>* getCosts (int lev) {
if (m_instance) {
return m_instance->costs[lev].get();
} else
@@ -467,7 +467,7 @@ public:
* @param[in] costs vector of (`unique_ptr` to) vectors; expected to be initialized
* to correct number of boxes and boxes per level
*/
- void ComputeCostsHeuristic (amrex::Vector<std::unique_ptr<amrex::Vector<amrex::Real> > >& costs);
+ void ComputeCostsHeuristic (amrex::Vector<std::unique_ptr<amrex::LayoutData<amrex::Real> > >& costs);
protected:
@@ -671,9 +671,9 @@ private:
/** Load balancing intervals that reads the "load_balance_int" string int the input file
* for getting steps at which load balancing is performed */
IntervalsParser load_balance_intervals;
- /** Collection of vectors to keep track of weights used in load balancing
+ /** Collection of LayoutData to keep track of weights used in load balancing
* routines. Contains timer-based or heuristic-based costs depending on input option */
- amrex::Vector<std::unique_ptr<amrex::Vector<amrex::Real> > > costs;
+ amrex::Vector<std::unique_ptr<amrex::LayoutData<amrex::Real> > > costs;
/** Load balance with 'space filling curve' strategy. */
int load_balance_with_sfc = 0;
/** Controls the maximum number of boxes that can be assigned to a rank during
diff --git a/Source/WarpX.cpp b/Source/WarpX.cpp
index 72adbe889..17a47caaa 100644
--- a/Source/WarpX.cpp
+++ b/Source/WarpX.cpp
@@ -1062,9 +1062,7 @@ WarpX::AllocLevelMFs (int lev, const BoxArray& ba, const DistributionMapping& dm
if (load_balance_intervals.isActivated())
{
- costs[lev].reset(new amrex::Vector<Real>);
- const int nboxes = Efield_fp[lev][0].get()->size();
- costs[lev]->resize(nboxes);
+ costs[lev].reset(new amrex::LayoutData<Real>(ba, dm));
}
}