Skip to content

Commit

Permalink
Activate load-balancing
Browse files Browse the repository at this point in the history
  • Loading branch information
RemiLehe committed Jan 12, 2025
1 parent e6d93f9 commit eb43087
Show file tree
Hide file tree
Showing 3 changed files with 14 additions and 4 deletions.
10 changes: 10 additions & 0 deletions Source/Parallelization/WarpXRegrid.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -174,6 +174,14 @@ WarpX::RemakeLevel (int lev, Real /*time*/, const BoxArray& ba, const Distributi
using ablastr::fields::Direction;
using warpx::fields::FieldType;

const auto RemakeMultiFab = [&](auto& mf){
if (mf == nullptr) { return; }
const IntVect& ng = mf->nGrowVect();
auto pmf = std::remove_reference_t<decltype(mf)>{};
AllocInitMultiFab(pmf, mf->boxArray(), dm, mf->nComp(), ng, lev, mf->tags()[0]);
mf = std::move(pmf);
};

bool const eb_enabled = EB::enabled();
if (ba == boxArray(lev))
{
Expand All @@ -187,6 +195,8 @@ WarpX::RemakeLevel (int lev, Real /*time*/, const BoxArray& ba, const Distributi
{
if (eb_enabled) {
if (WarpX::electromagnetic_solver_id != ElectromagneticSolverAlgo::PSATD) {
RemakeMultiFab( m_eb_update_E[lev][idim] );
RemakeMultiFab( m_eb_update_B[lev][idim] );
if (WarpX::electromagnetic_solver_id == ElectromagneticSolverAlgo::ECT) {
m_borrowing[lev][idim] = std::make_unique<amrex::LayoutData<FaceInfoBox>>(amrex::convert(ba, Bfield_fp[lev][idim]->ixType().toIntVect()), dm);
}
Expand Down
4 changes: 2 additions & 2 deletions Source/WarpX.H
Original file line number Diff line number Diff line change
Expand Up @@ -1015,7 +1015,7 @@ public:
/** Set a flag to indicate on which grid points the field `field`
* should be updated.
*
* More specifically, this function fill the iMultiFabs in `eb_update`
* More specifically, this function fills the iMultiFabs in `eb_update`
* (which have the same indexType as the MultiFabs in `field`) with 1
* or 0, depending on whether the grid point is in the valid region
* or in the embedded boundary.
Expand All @@ -1034,7 +1034,7 @@ public:

// TODO documentation
void MarkECTUpdateECells (
std::array< std::unique_ptr<amrex::iMultiFab>,3> & eb_update_B,
std::array< std::unique_ptr<amrex::iMultiFab>,3> & eb_update_E,
ablastr::fields::VectorField const& edge_lengths );

// TODO documentation
Expand Down
4 changes: 2 additions & 2 deletions Source/WarpX.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -2304,14 +2304,14 @@ WarpX::AllocLevelMFs (int lev, const BoxArray& ba, const DistributionMapping& dm
if (lev == maxLevel()) {
if (WarpX::electromagnetic_solver_id != ElectromagneticSolverAlgo::PSATD) {

AllocInitMultiFab(m_eb_update_E[lev][0], amrex::convert(ba, Ex_nodal_flag), dm, ncomps,
AllocInitMultiFab(m_eb_update_E[lev][0], amrex::convert(ba, Ex_nodal_flag), dm, ncomps,
guard_cells.ng_FieldSolver, lev, "m_eb_update_E[x]");
AllocInitMultiFab(m_eb_update_E[lev][1], amrex::convert(ba, Ey_nodal_flag), dm, ncomps,
guard_cells.ng_FieldSolver, lev, "m_eb_update_E[y]");
AllocInitMultiFab(m_eb_update_E[lev][2], amrex::convert(ba, Ez_nodal_flag), dm, ncomps,
guard_cells.ng_FieldSolver, lev, "m_eb_update_E[z]");

AllocInitMultiFab(m_eb_update_B[lev][0], amrex::convert(ba, Bx_nodal_flag), dm, ncomps,
AllocInitMultiFab(m_eb_update_B[lev][0], amrex::convert(ba, Bx_nodal_flag), dm, ncomps,
guard_cells.ng_FieldSolver, lev, "m_eb_update_B[x]");
AllocInitMultiFab(m_eb_update_B[lev][1], amrex::convert(ba, By_nodal_flag), dm, ncomps,
guard_cells.ng_FieldSolver, lev, "m_eb_update_B[y]");
Expand Down

0 comments on commit eb43087

Please sign in to comment.