Particles
Particle containers
Particle structures and functions are defined in Source/Particles/
. WarpX uses the Particle
class from AMReX for single particles. An ensemble of particles (e.g., a plasma species, or laser particles) is stored as a WarpXParticleContainer
(see description below) in a per-box (and even per-tile on CPU) basis.
-
class WarpXParticleContainer : public NamedComponentParticleContainer<amrex::DefaultAllocator>
WarpXParticleContainer is the base polymorphic class from which all concrete particle container classes (that store a collection of particles) derive. Derived classes can be used for plasma particles, photon particles, or non-physical particles (e.g., for the laser antenna). It derives from amrex::ParticleContainer<0,0,PIdx::nattribs>, where the template arguments stand for the number of int and amrex::Real SoA and AoS data in amrex::Particle.
AoS amrex::Real: x, y, z (default), 0 additional (first template parameter)
AoS int: id, cpu (default), 0 additional (second template parameter)
SoA amrex::Real: PIdx::nattribs (third template parameter), see PIdx for the list.
WarpXParticleContainer contains the main functions for initialization, interaction with the grid (field gather and current deposition) and particle push.
Note: many functions are pure virtual (meaning they MUST be defined in derived classes, e.g., Evolve) or actual functions (e.g. CurrentDeposition).
Subclassed by LaserParticleContainer, PhysicalParticleContainer
Physical species are stored in PhysicalParticleContainer
, that derives from WarpXParticleContainer
. In particular, the main function to advance all particles in a physical species is PhysicalParticleContainer::Evolve
(see below).
-
virtual void PhysicalParticleContainer::Evolve(int lev, const amrex::MultiFab &Ex, const amrex::MultiFab &Ey, const amrex::MultiFab &Ez, const amrex::MultiFab &Bx, const amrex::MultiFab &By, const amrex::MultiFab &Bz, amrex::MultiFab &jx, amrex::MultiFab &jy, amrex::MultiFab &jz, amrex::MultiFab *cjx, amrex::MultiFab *cjy, amrex::MultiFab *cjz, amrex::MultiFab *rho, amrex::MultiFab *crho, const amrex::MultiFab *cEx, const amrex::MultiFab *cEy, const amrex::MultiFab *cEz, const amrex::MultiFab *cBx, const amrex::MultiFab *cBy, const amrex::MultiFab *cBz, amrex::Real t, amrex::Real dt, DtType a_dt_type = DtType::Full, bool skip_deposition = false) override
Finally, all particle species (physical plasma species PhysicalParticleContainer
, photon species PhotonParticleContainer
or non-physical species LaserParticleContainer
) are stored in MultiParticleContainer
. The class WarpX
holds one instance of MultiParticleContainer
as a member variable, called WarpX::mypc
(where mypc stands for “my particle containers”):
-
class MultiParticleContainer
The class MultiParticleContainer holds multiple instances of the polymorphic class WarpXParticleContainer, stored in its member variable “allcontainers”. The class WarpX typically has a single (pointer to an) instance of MultiParticleContainer.
MultiParticleContainer typically has two types of functions:
Functions that loop over all instances of WarpXParticleContainer in allcontainers and calls the corresponding function (for instance, MultiParticleContainer::Evolve loops over all particles containers and calls the corresponding WarpXParticleContainer::Evolve function).
Functions that specifically handle multiple species (for instance ReadParameters or mapSpeciesProduct).
Loop over particles
A typical loop over particles reads:
// pc is a std::unique_ptr<WarpXParticleContainer>
// Loop over MR levels
for (int lev = 0; lev <= finest_level; ++lev) {
// Loop over particles, box by box
for (WarpXParIter pti(*this, lev); pti.isValid(); ++pti) {
// Do something on particles
// [MY INNER LOOP]
}
}
The innermost step [MY INNER LOOP]
typically calls amrex::ParallelFor
to perform operations on all particles in a portable way. For this reasons, the particle data needs to be converted in plain-old-data structures. The innermost loop in the code snippet above could look like:
// Get Array-Of-Struct particle data, also called data
// (x, y, z, id, cpu)
const auto& particles = pti.GetArrayOfStructs();
// Get Struct-Of-Array particle data, also called attribs
// (ux, uy, uz, w, Exp, Ey, Ez, Bx, By, Bz)
auto& attribs = pti.GetAttribs();
auto& Exp = attribs[PIdx::Ex];
// [...]
// Number of particles in this box
const long np = pti.numParticles();
Link fields and particles?
In WarpX, the loop over boxes through a MultiFab
iterator MFIter
and the loop over boxes through a ParticleContainer
iterator WarpXParIter
are consistent.
On a loop over boxes in a MultiFab
(MFIter
), it can be useful to access particle data on a GPU-friendly way. This can be done by:
// Index of grid (= box)
const int grid_id = mfi.index();
// Index of tile within the grid
const int tile_id = mfi.LocalTileIndex();
// Get GPU-friendly arrays of particle data
auto& ptile = GetParticles(lev)[std::make_pair(grid_id,tile_id)];
ParticleType* pp = particle_tile.GetArrayOfStructs()().data();
// Only need attribs (i.e., SoA data)
auto& soa = ptile.GetStructOfArrays();
// As an example, let's get the ux momentum
const ParticleReal * const AMREX_RESTRICT ux = soa.GetRealData(PIdx::ux).data();
On a loop over particles it can be useful to access the fields on the box we are looping over (typically when we use both field and particle data on the same box, for field gather or current deposition for instance). This is done for instance by adding this snippet in [MY INNER LOOP]
:
// E is a reference to, say, WarpX::Efield_aux
// Get the Ex field on the grid
const FArrayBox& exfab = (*E[lev][0])[pti];
// Let's be generous and also get the underlying box (i.e., index info)
const Box& box = pti.validbox();
Main functions
-
virtual void PhysicalParticleContainer::PushPX(WarpXParIter &pti, amrex::FArrayBox const *exfab, amrex::FArrayBox const *eyfab, amrex::FArrayBox const *ezfab, amrex::FArrayBox const *bxfab, amrex::FArrayBox const *byfab, amrex::FArrayBox const *bzfab, const amrex::IntVect ngEB, const int, const long offset, const long np_to_push, int lev, int gather_lev, amrex::Real dt, ScaleFields scaleFields, DtType a_dt_type = DtType::Full)
-
void WarpXParticleContainer::DepositCurrent(amrex::Vector<std::array<std::unique_ptr<amrex::MultiFab>, 3>> &J, const amrex::Real dt, const amrex::Real relative_time)
Deposit current density.
- Parameters:
J – [inout] vector of current densities (one three-dimensional array of pointers to MultiFabs per mesh refinement level)
dt – [in] Time step for particle level
relative_time – [in] Time at which to deposit J, relative to the time of the current positions of the particles. When different than 0, the particle position will be temporarily modified to match the time of the deposition.
Note
The current deposition is used both by PhysicalParticleContainer
and LaserParticleContainer
, so it is in the parent class WarpXParticleContainer
.
Buffers
To reduce numerical artifacts at the boundary of a mesh-refinement patch, WarpX has an option to use buffers: When particles evolve on the fine level, they gather from the coarse level (e.g., Efield_cax
, a copy of the aux
data from the level below) if they are located on the fine level but fewer than WarpX::n_field_gather_buffer
cells away from the coarse-patch boundary. Similarly, when particles evolve on the fine level, they deposit on the coarse level (e.g., Efield_cp
) if they are located on the fine level but fewer than WarpX::n_current_deposition_buffer
cells away from the coarse-patch boundary.
WarpX::gather_buffer_masks
and WarpX::current_buffer_masks
contain masks indicating if a cell is in the interior of the fine-resolution patch or in the buffers. Then, particles depending on this mask in
-
void PhysicalParticleContainer::PartitionParticlesInBuffers(long &nfine_current, long &nfine_gather, long const np, WarpXParIter &pti, int const lev, amrex::iMultiFab const *current_masks, amrex::iMultiFab const *gather_masks)
Note
Buffers are complex!
Particle attributes
WarpX adds the following particle attributes by default to WarpX particles. These attributes are either stored in an Array-of-Struct (AoS) or Struct-of-Array (SoA) location of the AMReX particle containers. The data structures for those are either pre-described at compile-time (CT) or runtime (RT).
Attribute name |
|
Description |
Where |
When |
Notes |
---|---|---|---|---|---|
|
|
Particle position. |
AoS |
CT |
|
|
|
CPU index where the particle was created. |
AoS |
CT |
|
|
|
CPU-local particle index where the particle was created. |
AoS |
CT |
|
|
|
Ion ionization level |
SoA |
RT |
Added when ionization physics is used. |
|
|
QED: optical depth of the Quantum- Synchrotron process |
SoA |
RT |
Added when PICSAR QED physics is used. |
|
|
QED: optical depth of the Breit- Wheeler process |
SoA |
RT |
Added when PICSAR QED physics is used. |
WarpX allows extra runtime attributes to be added to particle containers (through AddRealComp("attrname")
or AddIntComp("attrname")
).
The attribute name can then be used to access the values of that attribute.
For example, using a particle iterator, pti
, to loop over the particles the command pti.GetAttribs(particle_comps["attrname"]).dataPtr();
will return the values of the "attrname"
attribute.
User-defined integer or real attributes are initialized when particles are generated in AddPlasma()
.
The attribute is initialized with a required user-defined parser function.
Please see the input options addIntegerAttributes
and addRealAttributes
for a user-facing documentation.
Commonly used runtime attributes are described in the table below and are all part of SoA particle storage:
Attribute name |
|
Description |
Default value |
---|---|---|---|
|
|
The coordinates of the particles at the previous timestep. |
user-defined |
|
|
The coordinates of the particles when they were created. |
user-defined |
A Python example that adds runtime options can be found in Examples/Tests/particle_data_python
Note
Only use _
to separate components of vectors!