Skip to content

Commit

Permalink
Add ionization tracking framework to collision routine
Browse files Browse the repository at this point in the history
  • Loading branch information
budjensen committed Jan 9, 2025
1 parent 5e7eb29 commit 08a5c10
Showing 1 changed file with 41 additions and 0 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -520,6 +520,8 @@ void BackgroundMCCCollision::doBackgroundIonization
( int lev, amrex::LayoutData<amrex::Real>* cost,
WarpXParticleContainer& species1, WarpXParticleContainer& species2, amrex::Real t)
{
using warpx::fields::FieldType;

WARPX_PROFILE("BackgroundMCCCollision::doBackgroundIonization()");

const SmartCopyFactory copy_factory_elec(species1, species1);
Expand All @@ -535,6 +537,15 @@ void BackgroundMCCCollision::doBackgroundIonization

const amrex::ParticleReal sqrt_kb_m = std::sqrt(PhysConst::kb / m_background_mass);

if (WarpX::do_MCC_ionization_tracking)
{
const amrex::MultiFab &coll_ionization = *warpx.m_fields.get(FieldType::collision_ionization, lev);

// get parameters for getParticleCell
const auto plo = species1.Geom(lev).ProbLoArray();
const auto dxi = species1.Geom(lev).InvCellSizeArray();
}

#ifdef AMREX_USE_OMP
#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
#endif
Expand Down Expand Up @@ -565,6 +576,36 @@ void BackgroundMCCCollision::doBackgroundIonization
setNewParticleIDs(elec_tile, np_elec, num_added);
setNewParticleIDs(ion_tile, np_ion, num_added);

if (WarpX::do_MCC_ionization_tracking)
{
// get parameters for getParticleCell (method used here is similar to TemperatureFunctor.cpp)
// should these lines be removed and we just do ptd = elec_tile.getParticleTileData()?
auto& tile = pti.GetParticleTile();
auto ptd = tile.getParticleTileData();

// get multi-fab for ionization creation tracking
const auto& coll_ionization_arr = coll_ionization[pti].array();

// get Struct-Of-Array particle data
// as long as w is set after the copy, will it include newly created particles (be of length np_elec + num_added)?
auto& attribs = pti.GetAttribs();
amrex::ParticleReal* const AMREX_RESTRICT w = attribs[PIdx::w].dataPtr();

// add each created particle to tracking buffer
amrex::ParallelFor(num_added,
[=] AMREX_GPU_HOST_DEVICE (long ip)
{
// store the energy (add E_lost * w[ip] to the buffer of energy sent to the background)
// get the particle's cell index (method used here is similar to TemperatureFunctor.cpp)
const auto p = WarpXParticleContainer::ParticleType(ptd, ip + np_elec);
const auto [ii, jj, kk] = amrex::getParticleCell(p , plo, dxi).dim3();

// store energy lost in the buffer
coll_ionization_arr(ii, jj, kk, 0) += w[ip + np_elec];
}
);
}

if (cost && WarpX::load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::Timers)
{
amrex::Gpu::synchronize();
Expand Down

0 comments on commit 08a5c10

Please sign in to comment.