Skip to content

Commit

Permalink
dont use tiling in the plasma pusher (#949)
Browse files Browse the repository at this point in the history
  • Loading branch information
AlexanderSinn authored May 22, 2023
1 parent a5ae726 commit bc98e67
Show file tree
Hide file tree
Showing 4 changed files with 16 additions and 23 deletions.
4 changes: 0 additions & 4 deletions src/Hipace.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -535,10 +535,6 @@ Hipace::SolveOneSlice (int islice, const int islice_local, int step)

m_multi_plasma.DoFieldIonization(lev, m_3D_geom[lev], m_fields);

if (m_multi_plasma.IonizationOn() && m_do_tiling) {
m_multi_plasma.TileSort(m_slice_geom[lev].Domain(), m_slice_geom[lev]);
}

// Push plasma particles
m_multi_plasma.AdvanceParticles(m_fields, m_multi_laser, m_3D_geom, false, lev);

Expand Down
2 changes: 1 addition & 1 deletion src/particles/plasma/MultiPlasma.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -105,7 +105,7 @@ MultiPlasma::AdvanceParticles (
{
for (int i=0; i<m_nplasmas; i++) {
AdvancePlasmaParticles(m_all_plasmas[i], fields, gm, temp_slice,
lev, m_all_bins[i], multi_laser);
lev, multi_laser);
}
}

Expand Down
3 changes: 1 addition & 2 deletions src/particles/pusher/PlasmaParticleAdvance.H
Original file line number Diff line number Diff line change
Expand Up @@ -20,12 +20,11 @@
* \param[in] gm Geometry of the simulation, to get the cell size etc.
* \param[in] temp_slice if true, the temporary data (x_temp, ...) will be used
* \param[in] lev MR level
* \param[in] bins objects containing indices of plasma particles in each tile
* \param[in] multi_laser Laser pulses, which affects the plasma via the ponderomotive force
*/
void
AdvancePlasmaParticles (PlasmaParticleContainer& plasma, const Fields & fields,
amrex::Vector<amrex::Geometry> const& gm, const bool temp_slice, int const lev,
PlasmaBins& bins, const MultiLaser& multi_laser);
const MultiLaser& multi_laser);

#endif // PLASMAPARTICLEADVANCE_H_
30 changes: 14 additions & 16 deletions src/particles/pusher/PlasmaParticleAdvance.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -28,13 +28,11 @@ template struct PlasmaMomentumDerivative<DualNumber>;
void
AdvancePlasmaParticles (PlasmaParticleContainer& plasma, const Fields & fields,
amrex::Vector<amrex::Geometry> const& gm, const bool temp_slice, int const lev,
PlasmaBins& bins, const MultiLaser& multi_laser)
const MultiLaser& multi_laser)
{
HIPACE_PROFILE("AdvancePlasmaParticles()");
using namespace amrex::literals;

const bool do_tiling = Hipace::m_do_tiling;

const PhysConst phys_const = get_phys_const();

// Loop over particle boxes
Expand Down Expand Up @@ -90,25 +88,25 @@ AdvancePlasmaParticles (PlasmaParticleContainer& plasma, const Fields & fields,
const amrex::Real clight = phys_const.c;
const amrex::Real clight_inv = 1._rt/phys_const.c;
const amrex::Real charge_mass_clight_ratio = plasma.m_charge/(plasma.m_mass * phys_const.c);

const int ntiles = do_tiling ? bins.numBins() : 1;

#ifdef AMREX_USE_OMP
#pragma omp parallel for if (amrex::Gpu::notInLaunchRegion())
#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
#endif
for (int itile=0; itile<ntiles; itile++){
BeamBins::index_type const * const indices =
do_tiling ? bins.permutationPtr() : nullptr;
BeamBins::index_type const * const offsets =
do_tiling ? bins.offsetsPtr() : nullptr;
int const num_particles =
do_tiling ? offsets[itile+1]-offsets[itile] : pti.numParticles();
{
int const num_particles = pti.numParticles();
#ifdef AMREX_USE_OMP
int const idx_begin = (num_particles * omp_get_thread_num()) / omp_get_num_threads();
int const idx_end = (num_particles * (omp_get_thread_num()+1)) / omp_get_num_threads();
#else
int constexpr idx_begin = 0;
int const idx_end = num_particles;
#endif

amrex::ParallelFor(
amrex::TypeList<amrex::CompileTimeOptions<0, 1, 2, 3>>{},
{Hipace::m_depos_order_xy},
num_particles,
idx_end - idx_begin,
[=] AMREX_GPU_DEVICE (long idx, auto depos_order) {
const int ip = do_tiling ? indices[offsets[itile]+idx] : idx;
const int ip = idx + idx_begin;
// only push plasma particles on their according MR level
if (setPositionEnforceBC.m_structs[ip].id() < 0 ||
setPositionEnforceBC.m_structs[ip].cpu() != lev) return;
Expand Down

0 comments on commit bc98e67

Please sign in to comment.