Binning

This binning plugin is a flexible binner for particles and field properties. Users can

  • Define their own axes

  • Define their own quantity which is binned

  • Choose which species are used for particle binning

  • Choose which fields are used for field binning

  • Choose how frequently they want the binning to be executed

  • Choose if the binning should be time averaging

  • Write custom output to file, for example other quantities related to the simulation which the user is interested in

  • Execute multiple binnings at the same time

  • Pass extra parameters as a tuple, if additional information is required by the kernels to do binning.

  • Define a host-side hook to execute code before binning at each notify step, for example to fill FieldTmp.

User Input

Users can set up their binning in the binningSetup.param file. After setup, PIConGPU needs to be recompiled.

Attention

Unlike other plugins, the binning plugin doesn’t provide any runtime configuration. To set up binning, users need to define it in the param file and then recompile.

A binner is created using the addParticleBinner() or addFieldBinner() functions, which describe the configuration options available to the user to set up particle and field binning respectively. Multiple binnings can be run at the same time by simply calling addParticleBinner() or addFieldBinner() multiple times with different parameters.

class BinningCreator

An object of this class is provided to the user to add their binning setups.

Public Functions

template<typename T_AccumulateOp = alpaka::AtomicAdd, typename T_Extras = std::tuple<>>
inline auto &addParticleBinner(std::string const &binnerOutputName, auto const &axisTupleObject, auto const &speciesTupleObject, auto const &depositionData, T_Extras const &extraData = {})

Creates a particle binner from user input and adds it to the vector of all binners.

The Particle binner will bin particle quantities to create histograms on a grid defined by axes. The results will be written to openPMD files.

Parameters:
  • binnerOutputName – filename and dataset name for openPMD output. Must be unique to avoid overwrites during data dumps and undefined behaviour during restarts.

  • axisTupleObject – tuple holding the axes. Each element in the tuple describes an axis for binning.

  • speciesTupleObject – tuple holding the species to bin. Each element specifies a species to be included in the binning.

  • depositionData – functor description of the deposited quantity. Defines which particle property to bin.

  • extraData – tuple holding extra data to be passed to the binner. Can be used for optional configurations.

Returns:

ParticleBinningData& reference to the created ParticleBinningData object. This can be used to further configure the binning setup if needed.

template<typename T_AccumulateOp = alpaka::AtomicAdd, typename T_Extras = std::tuple<>>
inline auto addFieldBinner(std::string const &binnerOutputName, auto const &axisTupleObject, auto const &fieldsTupleObject, auto const &depositionData, T_Extras const &extraData = {})

Creates a field binner from user input and adds it to the vector of all binners.

The Field binner will bin field quantities to create histograms on a grid defined by axes. The results will be written to openPMD files.

Parameters:
  • binnerOutputName – filename and dataset name for openPMD output. Must be unique to avoid overwrites during data dumps and undefined behaviour during restarts.

  • axisTupleObject – tuple holding the axes. Each element in the tuple describes an axis for binning.

  • fieldsTupleObject – tuple holding the fields to bin. Each element specifies a field to be included in the binning.

  • depositionData – functor description of the deposited quantity. Defines which field property to bin.

  • extraData – tuple holding extra data to be passed to the binner. Can be used for optional configurations.

Returns:

FieldBinningData& reference to the created FieldBinningData object. This can be used to further configure the binning setup if needed.

The most important parts of defining a binning are the axes (the axes of the histogram which define the bins) and the deposited quantity (the quantity to be binned). Both of these are described using the “Functor Description”.

Functor Description

The basic building block for the binning plugin is the Functor Description object, and it is used to describe both axes and the deposited quantity. It describes the properties which we find interesting and how we can calculate/get this property from the particle or field. A functor description is created using createFunctorDescription.

template<typename QuantityType, typename FunctorType>
inline auto picongpu::plugins::binning::createFunctorDescription(FunctorType functor, std::string name, std::array<double, numUnits> units = std::array<double, numUnits>({0., 0., 0., 0., 0., 0., 0.}))

Describes the functors, units and names for the axes and the deposited quantity.

Todo:

infer T_Quantity from T_Functor, needs particle type also, different particles may have different return types

Template Parameters:
  • QuantityType – The type returned by the functor

  • FunctorType – Automatically deduced type of the functor

Parameters:
  • functor – Functor which access the particle property

  • name – Name for the functor/axis written out with the openPMD output.

  • units – The dimensionality of the quantity returned by the functor in the 7D format. Defaults to unitless.

Returns:

FunctorDescription object

Functor

The functor is run on the GPU and constitutes the means by which the user can execute code on device while binning. This may be used to calculate the deposited quantity or to calculate the axis values.

The functor needs to follow the signature shown below. This provides the user access to the particle or field object and with information about the domain.

The return type is defined by the user.

The first two arguments are always the worker and the domainInfo. For particle binning this is followed by the particle object, and for field binning this is followed by the field object (if any). If extra data is passed, this is also passed as arguments. It is the responsibility of the user to ensure that their functors have an appropriate number and types of arguments to match the provided tuples. A few examples are shown below.

For particles:

auto myParticleFunctor = [] ALPAKA_FN_ACC(auto const& worker, auto const& domainInfo, auto const& particle) -> returnType
{
        // fn body
        return myParameter;
};

For particles with two extra data paramters:

auto myParticleWithExtraDataFunctor = [] ALPAKA_FN_ACC(auto const& worker, auto const& domainInfo, auto const& particle, auto const& extraData1, auto const& extraData2) -> returnType
{
        // fn body
        return myParameter;
};

For one field:

auto myFieldFunctor = [] ALPAKA_FN_ACC(auto const& worker, auto const& domainInfo, auto const& field) -> returnType
{
        // fn body
        return myParameter;
};

For two fields with one extra data parameter:

auto myFieldWithExtraDataFunctor = [] ALPAKA_FN_ACC(auto const& worker, auto const& domainInfo, auto const& field1, auto const& field2, auto const& extraData) -> returnType
{
        // fn body
        return myParameter;
};

Note

Fields and extra data may be tuples with multiple elements which are unpacked and passed to the user-defined functors. Therefore, it is the responsibility of the user to ensure that their functors have an appropriate number of arguments to match the provided tuples.

Domain Info

Enables the user to find the location of the particle or field in the simulation domain. Contains

For particle binning, the DomainInfo class contains:

class DomainInfoBase

Provides knowledge of the simulation domain to the user Names and concept are described at https://github.com/ComputationalRadiationPhysics/picongpu/wiki/PIConGPU-domain-definitions.

Subclassed by picongpu::plugins::binning::DomainInfo< BinningType::Field >, picongpu::plugins::binning::DomainInfo< BinningType::Particle >

Public Functions

inline HDINLINE DomainInfoBase(uint32_t simStep, pmacc::DataSpace<simDim> gOffset, pmacc::DataSpace<simDim> lOffset, pmacc::DataSpace<simDim> physicalSuperCellIdx)
Parameters:

physicalSuperCellIdx – supercell index relative to the border origin

Public Members

uint32_t currentStep

Current simulation timestep.

pmacc::DataSpace<simDim> globalOffset

Offset of the global domain on all GPUs.

pmacc::DataSpace<simDim> localOffset

Offset of the domain simulated on current GPU.

pmacc::DataSpace<simDim> blockCellOffset

Offset of domain simulated by current block wrt the border.

The global and local offsets can be understood by looking at the PIConGPU domain definitions.

For particle binning, the particle position is obtained at cell precision by default. To get sub-cell precision or SI units, use optional template parameters with getParticlePosition<DomainOrigin, PositionPrecision, PositionUnits>.

For field binning, the field DomainInfo additionally holds the localCellIndex in the supercell and has a method to getCellIndex<DomainOrigin, PositionUnits> to get the current cell index being binned relative to an origin (global, total, local). To get the exact position of the fields inside the cell, relative to the cell index, use the FieldPosition trait.

Dimensionality and units

Users can specify the units of their functor output using a 7 dimensional array. Each element of the array corresponds to an SI base unit, and the value stored in that index is the exponent of the unit. The dimensional base quantities are defined in SIBaseUnits_t following the international system of quantities (ISQ). If no units are given, the quantity is assumed to be dimensionless.

// Define units
std::array<double, numUnits> momentumDimension{};
momentumDimension[SIBaseUnits::length] = 1.0;
momentumDimension[SIBaseUnits::mass] = 1.0;
momentumDimension[SIBaseUnits::time] = -1.0;
enum picongpu::SIBaseUnits::SIBaseUnits_t

Values:

enumerator length
enumerator mass
enumerator time
enumerator electricCurrent
enumerator thermodynamicTemperature
enumerator amountOfSubstance
enumerator luminousIntensity

Axis

Axis is a combination of a functor description and an axis splitting. These are brought together by createAxis functions, depending on what kind of an axis you want. The name used in the functor description is used as the name of the axis for openPMD.

Attention

The return type of the functor as specified in the functor description is required to be the same as the type of the range (min, max).

Currently implemented axis types
  • Linear Axis

  • Log Axis

template<typename T_Attribute, typename T_AttrFunctor>
class LinearAxis

Linear axis with contiguous fixed sized bins.

Axis splitting is defined with min, max and n_bins. Bin size = (max-min)/n_bins. Bins are closed open [) intervals [min, min + size), [min + size, min + 2*size) ,…, [max-size, max). Allocates 2 extra bins, for under and overflow. These are bin index 0 and (n_bins+2)-1 T_Attribute is a Numeric/Arithmetic type

template<typename T_Attribute, typename T_AttrFunctor>
class LogAxis

Log axis with logarithmically sized bins.

Axis splitting is defined with min, max and n_bins. A bin with edges a and b is a closed-open interval [a,b) If overflow bins are enabled, allocates 2 extra bins, for under and overflow. These are bin index 0 and (n_bins+2)-1

Binning can be done over an arbitrary number of axes, by creating a tuple of all the axes. Limited by memory depending on number of bins in each axis.

Axis Splitting

Defines the axis range and how it is split into bins. Bins are defined as closed open intervals from the lower edge to the upper edge. In the future, this plugin will support other ways to split the domain, eg. using the binWidth or by auto-selecting the parameters.

template<typename T_Data>
class AxisSplitting

Holds the axis range in SI units and information on how this range is split into bins.

Public Members

Range<T_Data> m_range

Range object in SI units.

uint32_t nBins

Number of bins in range.

bool enableOverflowBins

Enable or Disable overflow bins.

Number of overflow bis is the responsibility of the axis implementaiton Defaults to true

Range
template<typename T_Data>
class Range

Holds the range in SI Units in axis space over which the binning will be done.

Public Members

T_Data min

Minimum of binning range in SI Units.

T_Data max

Maximum of binning range in SI Units.

Note

Axes are passed to addParticleBinner or addFieldBinner grouped in a tuple. This is just a collection of axis objects and is of arbitrary size. Users can make a tuple for axes by using the createTuple() function and passing in the axis objects as arguments.

Species

PIConGPU species which should be used in particle binning. Species can be instances of a species type or a particle species name as a PMACC_CSTRING. For example,

auto electronsObj = PMACC_CSTRING("e"){};

Optionally, users can specify a filter to be used with the species. This is a predicate functor, i.e. it is a functor with a signature as described above and which returns a boolean. If the filter returns true it means the particle is included in the binning. They can then create a FilteredSpecies object which contains the species and the filter.

/**
 * Optionally pass a particle filter to be used for this
 * Simply a particle predicate functor
 */
auto allParticles
    = [] ALPAKA_FN_ACC(auto const& worker, auto const& domainInfo, auto const& particle) -> bool
{ return true; };

auto filteredElectrons = FilteredSpecies{electronsObj, allParticles};

Note

Species are passed to addParticleBinner in the form of a tuple. This is just a collection of Species and FilteredSpecies objects (the tuple can be a mixture of both) and is of arbitrary size. Users can make a species tuple by using the createSpeciesTuple() function and passing in the objects as arguments.

Fields

PIConGPU fields which should be used in field binning. Fields can be instances of a field type. For example,

// Bring the fields together in a tuple
DataConnector& dc = Environment<>::get().DataConnector();
auto fieldsTuple = createTuple(
    dc.get<FieldB>(FieldB::getName())->getDeviceDataBox(),
    dc.get<FieldTmp>(FieldTmp::getUniqueId(0))->getDeviceDataBox());

Fields are passed to addFieldBinner in the form of a tuple. This is just a collection of field objects and is of arbitrary size. Users can make a fields tuple by using the createTuple() function and passing in the objects as arguments.

Note

It is possible to have an empty tuple for fields when doing field binning, in which case the functor will be called with no fields. This may be useful if you are passing in extra data and want field traversal over it.

Deposited Quantity

Quantity to be deposited is simply a functor description. The signature of quantity functors is

auto myQuantityFunctor = [] ALPAKA_FN_ACC(auto const& worker, auto const& domainInfo, ...) -> returnType

This option makes it evident that the binning is more than just about creating histograms. While histograms are a common use case, the binning plugin allows for the accumulation of various quantities within bins and not just noting the frequencies of occurrences. This means that users can define custom quantities to be accumulated in each bin, such as charge, energy, momentum, or any other property of interest. The flexibility of the functor description enables users to specify exactly what and how they want to accumulate data in the bins. For example, you might want to accumulate the total charge of particles within each bin, or the average kinetic energy of particles in a specific region. The deposited quantity functor provides the mechanism to calculate and return these values, which are then accumulated in the corresponding bins during the simulation. By default the deposited quantity is added for each bin, but this is configurable by the user by setting an accumulate operation.

Extra Data

Users can pass extra data to the functors if additional information is required by the kernels to do binning. Extra data is passed to binningCreator’s addParticleBinner or addFieldBinner functions in the form of a tuple. This tuple is just a collection of objects and is of arbitrary size. Users can make a tuple of the extra data by using the createTuple() function and passing in the objects as arguments.

Note

Ensure that the functors have an appropriate number of arguments to match the provided tuples. The extra data may be tuples with multiple elements which are unpacked and passed to the user-defined functors.

Host-side hook

Users can define a host-side hook to execute code before binning starts for a notify step. This hook is set using the setHostSideHook method of the BinningCreator class and takes a callable (lambda/functor/std::function) which returns void.

[]() -> void{
        // host-side code to be executed before binning
}

For example, this can be used to pre-process data on the host side, such as filling temporary fields before they are used in field binning.

Notify period

Set the periodicity of the output. Follows the period syntax defined here. Notifies every time step by default.

Dump Period

Defines the number of notify steps to accumulate over. Note that this is not accumulating over actual PIC iterations, but over the notify periods. If time averaging is enabled, this is also the period to do time averaging over. For example a value of 10 means that after every 10 notifies, an accumulated file will be written out. If PIConGPU exits before executing 10 notifies, then there will be no output. The plugin dumps on every notify if this is set to either 0 or 1. This is the default behaviour.

Time Averaging

When dumping the accumulated output, whether or not to divide by the dump period, i.e. do a time averaging. Enabled by default.

Attention

The user needs to set a dump period to enable time averaging.

Binning Particles Leaving the Simulation Volume

enum picongpu::plugins::binning::ParticleRegion

Bin particles in enabled region.

All regions must be represented by a unique bit

Values:

enumerator Bounded

Bounded - Particles inside the global simulation volume, 01 in binary, corresponds to the first bit.

enumerator Leaving

Leaving - Particles that have left the global simulation volume in this timestep, 10 in binary, corresponds to the second bit.

By default, only particles within the simulation volume are binned by the particle binner. However, users can modify this behavior to include or exclusively bin particles that are leaving the global simulation volume. This can be configured using the enableRegion and disableRegion options with the regions defined by the ParticleRegion enum.

Users must carefully configure the notify period when using the binning plugin for leaving particles. The plugin bins particles leaving the global simulation volume at every timestep (except 0) after particles are pushed, regardless of the notify period. If the plugin is not notified at every timestep, this can cause discrepancies between the binning process and time-averaged data or histogram dumps, which follow the notify period. Additionally, the binning plugin is first notified at timestep 0, allowing users to bin initial conditions. However, leaving particles are first binned at timestep 1, after the initial particle push. Therefore, users should consider setting the notify period’s start at timestep 1, depending on their specific needs.

writeOpenPMDFunctor

Users can also write out custom output to file, for example other quantities related to the simulation which the user is interested in. This is a lambda with the following signature, and is set using the setWriteOpenPMDFunctor method.

[=](::openPMD::Series& series, ::openPMD::Iteration& iteration, ::openPMD::Mesh& mesh) -> void

Note

Make sure to capture by copy only, as the objects defined in the param file are not kept alive

OpenPMD JSON Configuration

Users can set a json configuration string for the OpenPMD output format using the setOpenPMDJsonCfg setter method.

OpenPMD Output

The binning outputs are stored in HDF5 files in simOutput/binningOpenPMD/ directory.

The files are named as <binnerOutputName>_<timestep>.h5.

The OpenPMD mesh is call “Binning”.

The outputs in written in SI units.

The output histogram has 2 bins more in each dimension than the user-defined nBins in that dimension, to deal with under and overflow.

The number of bin edges written out for an axis is one more than the user-defined nBins. These represent the bins in [min,max]. Since there are actually nBins + 2 bins, two edges are not written out. These are the first and last edge, corresponding to the overflow bins, and they have the value of -inf and + inf.

Attribute

Description

unitSI

Scaling factor for the deposited quantity to convert to SI

<axisName>_bin_edges

The edges of the bins of an axis in SI units

<axisName>_units

The units of an axis

Accumulation

The binning plugin provides flexible options for accumulating data within bins. Instead of simply adding values to accumulate in a bin, users can utilize any alpaka atomic operation. This includes operations such as subtraction, minimum, maximum, AND, OR, XOR, and more. This flexibility enables users to tailor the accumulation process to their specific needs. For instance, you might want to track the minimum or maximum value of a certain property within each bin. A list of available alpaka atomic operations can be found in the alpaka documentation.

Note

Only binary atomics on arithmetic types are supported. In particular CAS operations are not supported.

The accumulate operation is passed as a template parameter to the createBinner functions. For example, to use the maximum operation for accumulation, you would pass the corresponding Alpaka atomic operation as a template parameter:

binningCreator
    .addParticleBinner<alpaka::AtomicMax>(
        "particleBinning",
        createTuple(ax_timeStep, ax_py),
        speciesTuple,
        getCounts)
    .setDumpPeriod(2000);

Similarly, you can use other alpaka atomic operations such as alpaka::AtomicMin, alpaka::AtomicSub, alpaka::AtomicAnd, alpaka::AtomicOr, alpaka::AtomicXor, etc.

Example binning plugin usage: Laser Wakefield electron spectrometer

The LWFA example contains a sample binning plugin setup to calculate an in-situ electron spectrometer. The kinetic energy of the electrons \(E = (\gamma - 1) m_o c^2\) is plotted on axis 1 and the direction of the electrons \(\theta = \mathrm{atan}(p_x/p_y)\) is plotted on axis 2. The charge \(Q\) moving in the bin direction \(\theta\) at the bin energy \(E\) is calculated for each bin. Such spectrometers are a common tool in plasma based electron acceleration experiments [Kurz2018].

Note

Please note that if you specify the SI units of an axis, e.g. via energyDimension in the LWFA example, PIConGPU automatically converts to the internal unit system, but then also expects SI-compliant values for the axis range (in the case of energy Joules).

To read the electron spectrometer data in python, one could load and plot it like this:

import numpy as np
import matplotlib.pyplot as plt
from matplotlib.colors import LogNorm
import openpmd_api as io
import scipy.constants as const

# access openPMD series of eSpec
series = io.Series("./LWFA/simOutput/binningOpenPMD/eSpec_%T.h5", access=io.Access_Type.read_only)

last_iter = list(series.iterations)[-1]
it = series.iterations[last_iter]
espec_h = it.meshes['Binning'][io.Mesh_Record_Component.SCALAR]

# load data
espec = espec_h[:,:]
series.flush()

# convert to SI units and make positve (electrons have a negative charge)
espec *= espec_h.get_attribute('unitSI') * -1

# get axes (they are already in the correct SI unit)
E_bins = espec_h.get_attribute('Energy_bin_edges')
theta_bins = espec_h.get_attribute('pointingXY_bin_edges')

# convert C/J/rad -> C/MeV/mrad
convert_C_per_Joule_per_rad_to_pC_per_MeV_per_mrad = 1./1e-12 * const.elementary_charge/1e6 * 1/1e3

# plot
plt.pcolormesh(np.array(E_bins) / const.elementary_charge / 1e6,
               np.array(theta_bins) / 0.001,
               espec[1:-1, 1:-1] * convert_C_per_Joule_per_rad_to_pC_per_MeV_per_mrad,
               norm=LogNorm(), cmap=plt.cm.inferno)
cb = plt.colorbar()

plt.xlabel(r"$E \, \mathrm{[MeV]}$", fontsize=18)
plt.xticks(fontsize=14)

plt.ylabel(r"$\theta \, \mathrm{[mrad]}$", fontsize=18)
plt.yticks(fontsize=14)

cb.set_label(r"$\frac{\mathrm{d}^2 Q}{\mathrm{d} E \mathrm{d}\theta} \, \mathrm{[pC/MeV/mrad]}$", fontsize=20)
for i in cb.ax.get_yticklabels():
    i.set_fontsize(14)

plt.tight_layout()
plt.show()

References

[Kurz2018]

T. Kurz, J.P. Couperus, J.M. Krämer, et al. Calibration and cross-laboratory implementation of scintillating screens for electron bunch charge determination, Review of Scientific Instruments (2018), https://doi.org/10.1063/1.5041755

More examples of binning plugin usage

Further examples of the binning plugin can be found in the atomic physics example at share/picongpu/examples/AtomicPhysics/include/picongpu/param/binningSetup.param and the share/picongpu/tests/compile2 test.