Particles

Particles are defined in modular steps. First, species need to be generally defined in speciesDefinition.param. Second, species are initialized with particles in speciesInitialization.param.

The following operations can be applied in the picongpu::particles::InitPipeline of the latter:

Initialization

CreateDensity

template<typename T_DensityFunctor, typename T_PositionFunctor, typename T_SpeciesType = boost::mpl::_1>
struct CreateDensity

Sample macroparticles according to the given spatial density profile.

Create macroparticles inside a species.

This function only concerns the number of macroparticles, positions, and weighting. So it basically performs sampling in the coordinate space, while not initializing other attributes. When needed, those should be set (for then-existing macroparticles) by subsequently calling Manipulate.

User input to this functor is two-fold. T_DensityFunctor represents spatial density of real particles, normalized according to our requirements. It describes the physical setup being simulated and only deals with real, not macro-, particles. T_PositionFunctor is more of a PIC simulation parameter. It defines how real particles in each cell will be represented with macroparticles. This concerns the count, weighting, and in-cell positions of the created macroparticles.

The sampling process operates independently for each cell, as follows:

  • Evaluate the amount of real particles in the cell, Nr, using T_DensityFunctor.

  • If Nr > 0, decide how to represent it with macroparticles using T_PositionFunctor:

    • (For simplicity we describe how all currently used functors but RandomPositionAndWeightingImpl and RandomBinomialImpl operate, see below for customization)

    • Try to have exactly T_PositionFunctor::numParticlesPerCell macroparticles with same weighting w = Nr / T_PositionFunctor::numParticlesPerCell.

    • If such w < MIN_WEIGHTING, instead use fewer macroparticles and higher weighting.

    • In any case the combined weighting of all new macroparticles will match Nr.

  • Create the selected number of macroparticles with selected weighting.

  • Set in-cell positions according to T_PositionFunctor.

In principle, one could override the logic inside the (If Nr > 0) block by implementing a custom functor. Then one could have an arbitrary number of macroparticles and weight distribution between them. The only requirement is that together it matches Nr. For an example of non-uniform weight distribution

See also

startPosition::RandomPositionAndWeightingImpl. Note that in this scheme almost all non-vacuum cells will start with the same number of macroparticles. Having a higher density in a cell would mean larger weighting, but not more macroparticles.

Note

FillAllGaps is automatically called after creation.

Template Parameters:
  • T_DensityFunctor – unary lambda functor with profile description, see density.param, example: picongpu::particles::densityProfiles::Homogenous

  • T_PositionFunctor – unary lambda functor with position description and number of macroparticles per cell, see particle.param, examples: picongpu::particles::startPosition::Quiet, picongpu::particles::startPosition::Random

  • T_SpeciesType – type or name as PMACC_CSTRING of the used species, see speciesDefinition.param

Derive

template<typename T_SrcSpeciesType, typename T_DestSpeciesType = boost::mpl::_1, typename T_Filter = filter::All>
struct Derive : public picongpu::particles::ManipulateDerive<manipulators::generic::None, T_SrcSpeciesType, boost::mpl::_1, filter::All>

Generate particles in a species by deriving from another species’ particles.

Create particles in T_DestSpeciesType by deriving (copying) all particles and their matching attributes (except particleId) from T_SrcSpeciesType.

Note

FillAllGaps is called on on T_DestSpeciesType after the derivation is finished.

Template Parameters:
  • T_SrcSpeciesType – type or name as PMACC_CSTRING of the source species

  • T_DestSpeciesType – type or name as PMACC_CSTRING of the destination species

  • T_Filter – picongpu::particles::filter, particle filter type to select source particles to derive

Manipulate

template<typename T_Manipulator, typename T_Species = boost::mpl::_1, typename T_Filter = filter::All, typename T_Area = std::integral_constant<uint32_t, CORE + BORDER>>
struct Manipulate : public pmacc::particles::algorithm::CallForEach<pmacc::particles::meta::FindByNameOrType<VectorAllSpecies, boost::mpl::_1>, detail::MakeUnaryFilteredFunctor<T_Manipulator, boost::mpl::_1, filter::All>, T_Area::value>

Run a user defined manipulation for each particle of a species in an area.

Allows to manipulate attributes of existing particles in a species with arbitrary unary functors (“manipulators”).

Provides two versions of operator() to either operate on T_Area or a custom area,

See also

pmacc::particles::algorithm::CallForEach.

See also

picongpu::particles::manipulators

Warning

Does NOT call FillAllGaps after manipulation! If the manipulation deactivates particles or creates “gaps” in any other way, FillAllGaps needs to be called for the T_Species manually in the next step!

Template Parameters:
  • T_Manipulator – unary lambda functor accepting one particle species,

  • T_Species – type or name as PMACC_CSTRING of the used species

  • T_Filter – picongpu::particles::filter, particle filter type to select particles in T_Species to manipulate

  • T_Area – area to process particles in operator()(currentStep), wrapped into std::integral_constant for boost::mpl::apply to work; does not affect operator()(currentStep, areaMapperFactory)

ManipulateDerive

template<typename T_Manipulator, typename T_SrcSpeciesType, typename T_DestSpeciesType = boost::mpl::_1, typename T_SrcFilter = filter::All>
struct ManipulateDerive

Generate particles in a species by deriving and manipulating from another species’ particles.

Create particles in T_DestSpeciesType by deriving (copying) all particles and their matching attributes (except particleId) from T_SrcSpeciesType. During the derivation, the particle attributes in can be manipulated with T_ManipulateFunctor.

See also

picongpu::particles::manipulators

Note

FillAllGaps is called on on T_DestSpeciesType after the derivation is finished. If the derivation also manipulates the T_SrcSpeciesType, e.g. in order to deactivate some particles for a move, FillAllGaps needs to be called for the T_SrcSpeciesType manually in the next step!

Template Parameters:
  • T_Manipulator – a pseudo-binary functor accepting two particle species: destination and source,

  • T_SrcSpeciesType – type or name as PMACC_CSTRING of the source species

  • T_DestSpeciesType – type or name as PMACC_CSTRING of the destination species

  • T_SrcFilter – picongpu::particles::filter, particle filter type to select particles in T_SrcSpeciesType to derive into T_DestSpeciesType

FillAllGaps

template<typename T_SpeciesType = boost::mpl::_1>
struct FillAllGaps

Generate a valid, contiguous list of particle frames.

Some operations, such as deactivating or adding particles to a particle species can generate “gaps” in our internal particle storage, a list of frames.

This operation copies all particles from the end of the frame list to “gaps” in the beginning of the frame list. After execution, the requirement that all particle frames must be filled contiguously with valid particles and that all frames but the last are full is fulfilled.

Template Parameters:

T_SpeciesType – type or name as PMACC_CSTRING of the particle species to fill gaps in memory

Manipulation Functors

Some of the particle operations above can take the following functors as arguments to manipulate attributes of particle species. A particle filter (see following section) is used to only manipulated selected particles of a species with a functor.

Free

template<typename T_Functor>
struct Free : protected picongpu::particles::functor::User<T_Functor>

call simple free user defined manipulators

example for particle.param: set in cell position to zero

struct FunctorInCellPositionZero
{
    template< typename T_Particle >
    HDINLINE void operator()( T_Particle & particle )
    {
        particle[ position_ ] = floatD_X::create( 0.0 );
    }
    static constexpr char const * name = "inCellPositionZero";
};

using InCellPositionZero = generic::Free<
   FunctorInCellPositionZero
>;

Template Parameters:

T_Functor – user defined manipulators optional: can implement one host side constructor T_Functor() or T_Functor(uint32_t currentTimeStep)

FreeRng

template<typename T_Functor, typename T_Distribution>
struct FreeRng : protected picongpu::particles::functor::User<T_Functor>, private picongpu::particles::functor::misc::Rng<T_Distribution>

call simple free user defined functor and provide a random number generator

example for particle.param: add

#include <pmacc/random/distributions/Uniform.hpp>

struct FunctorRandomX
{
    template< typename T_Rng, typename T_Particle >
    HDINLINE void operator()( T_Rng& rng, T_Particle& particle )
    {
        particle[ position_ ].x() = rng();
    }
    static constexpr char const * name = "randomXPos";
};

using RandomXPos = generic::FreeRng<
   FunctorRandomX,
   pmacc::random::distributions::Uniform< float_X >
>;

and to InitPipeline in speciesInitialization.param:

Manipulate< manipulators::RandomXPos, SPECIES_NAME >

Template Parameters:
  • T_Functor – user defined unary functor

  • T_Distribution – pmacc::random::distributions, random number distribution

FreeTotalCellOffset

template<typename T_Functor>
struct FreeTotalCellOffset : protected picongpu::particles::functor::User<T_Functor>, private picongpu::particles::functor::misc::TotalCellOffset

call simple free user defined manipulators and provide the cell information

The functor passes the cell offset of the particle relative to the total domain origin into the functor.

example for particle.param: set a user-defined species attribute y0 (type: uint32_t) to the current total y-cell index

struct FunctorSaveYcell
{
    template< typename T_Particle >
    HDINLINE void operator()(
       DataSpace< simDim > const & particleOffsetToTotalOrigin,
       T_Particle & particle
    )
    {
        particle[ y0_ ] = particleOffsetToTotalOrigin.y();
    }
    static constexpr char const * name = "saveYcell";
};

using SaveYcell = unary::FreeTotalCellOffset<
   FunctorSaveYcell
>;

Template Parameters:

T_Functor – user defined unary functor

CopyAttribute

template<typename T_DestAttribute, typename T_SrcAttribute>
using picongpu::particles::manipulators::unary::CopyAttribute = generic::Free<acc::CopyAttribute<T_DestAttribute, T_SrcAttribute>>

copy a particle source attribute to a destination attribute

This is an unary functor and operates on one particle.

Template Parameters:
  • T_DestAttribute – type of the destination attribute e.g. momentumPrev1

  • T_SrcAttribute – type of the source attribute e.g. momentum

Drift

template<typename T_ParamClass = param::DriftCfg, typename T_ValueFunctor = pmacc::math::operation::Add>
using picongpu::particles::manipulators::unary::Drift = generic::Free<acc::Drift<T_ParamClass, T_ValueFunctor>>

change particle’s momentum based on speed

allow to manipulate a speed to a particle

Template Parameters:
  • T_ParamClass – param::DriftCfg, configuration parameter

  • T_ValueFunctor – pmacc::math::operation::*, binary functor type to manipulate the momentum attribute

RandomPosition

using picongpu::particles::manipulators::unary::RandomPosition = generic::FreeRng<acc::RandomPosition, pmacc::random::distributions::Uniform<float_X>>

Change the in cell position.

This functor changes the in-cell position of a particle. The new in-cell position is uniformly distributed position between [0.0;1.0).

example: add

particles::Manipulate<RandomPosition,SPECIES_NAME>
to InitPipeline in speciesInitialization.param

Temperature

template<typename T_ParamClass = param::TemperatureCfg, typename T_ValueFunctor = pmacc::math::operation::Add>
using picongpu::particles::manipulators::unary::Temperature = generic::FreeRng<acc::Temperature<T_ParamClass, T_ValueFunctor>, pmacc::random::distributions::Normal<float_X>>

Modify particle momentum based on temperature.

Sample a random momentum value distributed according to the given temperature and add it to the existing particle momentum.

Note: initial electron temperature should generally be chosen so that the resulting Debye length is resolved by the grid.

Template Parameters:
  • T_ValueFunctor – pmacc::math::operation::*, binary functor type to add a new momentum to an old one Version with fixed temperature given via parameter struct

  • T_ParamClass – configuration parameter, follows requirements of param::TemperatureCfg

template<typename T_TemperatureFunctor = param::TemperatureFunctor, typename T_ValueFunctor = pmacc::math::operation::Add>
using picongpu::particles::manipulators::unary::FreeTemperature = FreeTotalCellOffsetRng<acc::FreeTemperature<T_TemperatureFunctor, T_ValueFunctor>, pmacc::random::distributions::Normal<float_X>>

Version with user-provided temperature functor.

Template Parameters:

T_TemperatureFunctor – temperature functor, follows requirements of param::TemperatureFunctor

Assign

using picongpu::particles::manipulators::binary::Assign = generic::Free<acc::Assign>

assign attributes of one particle to another

Can be used as binary and higher order operator but only the first two particles are used for the assign operation.

Assign all matching attributes of a source particle to the destination particle. Attributes that only exist in the destination species are initialized with the default value. Attributes that only exists in the source particle will be ignored.

DensityWeighting

using picongpu::particles::manipulators::binary::DensityWeighting = generic::Free<acc::DensityWeighting>

Re-scale the weighting of a cloned species by densityRatio.

When deriving species from each other, the new species “inherits” the macro-particle weighting of the first one. This functor can be used to manipulate the weighting of the new species’ macro particles to satisfy the input densityRatio of it.

note: needs the densityRatio flag on both species, used by the GetDensityRatio trait.

ProtonTimesWeighting

using picongpu::particles::manipulators::binary::ProtonTimesWeighting = generic::Free<acc::ProtonTimesWeighting>

Re-scale the weighting of a cloned species by numberOfProtons.

When deriving species from each other, the new species “inherits” the macro-particle weighting of the first one. This functor can be used to manipulate the weighting of the new species’ macro particles to be a multiplied by the number of protons of the initial species.

As an example, this is useful when initializing a quasi-neutral, pre-ionized plasma of ions and electrons. Electrons can be created from ions via deriving and increasing their weight to avoid simulating multiple macro electrons per macro ion (with Z>1).

note: needs the atomicNumbers flag on the initial species, used by the GetAtomicNumbers trait.

Manipulation Filters

Most of the particle functors shall operate on all valid particles, where filter::All is the default assumption. One can limit the domain or subset of particles with filters such as the ones below (or define new ones).

All

struct All

RelativeGlobalDomainPosition

template<typename T_Params>
struct RelativeGlobalDomainPosition

filter particle dependent on the global position

Check if a particle is within a relative area in one direction of the global domain.

Template Parameters:

T_Params – picongpu::particles::filter::param::RelativeGlobalDomainPosition, parameter to configure the functor

Free

template<typename T_Functor>
struct Free : protected picongpu::particles::functor::User<T_Functor>

call simple free user defined filter

example for particleFilters.param: each particle with in-cell position greater than 0.5

struct FunctorEachParticleAboveMiddleOfTheCell
{
    template< typename T_Particle >
    HDINLINE bool operator()( T_Particle const & particle )
    {
        bool result = false;
        if( particle[ position_ ].y() >= float_X( 0.5 ) )
            result = true;
        return result;
    }
    static constexpr char const * name = "eachParticleAboveMiddleOfTheCell";

    static constexpr bool isDeterministic = true;
};

using EachParticleAboveMiddleOfTheCell = generic::Free<
   FunctorEachParticleAboveMiddleOfTheCell
>;

Template Parameters:

T_Functor – user defined filter optional: can implement one host side constructor T_Functor() or T_Functor(uint32_t currentTimeStep)

FreeRng

template<typename T_Functor, typename T_Distribution>
struct FreeRng : protected picongpu::particles::functor::User<T_Functor>, private picongpu::particles::functor::misc::Rng<T_Distribution>

call simple free user defined functor and provide a random number generator

example for particleFilters.param: get every second particle (random sample of 50%)

struct FunctorEachSecondParticle
{
    template< typename T_Rng, typename T_Particle >
    HDINLINE bool operator()(
        T_Rng & rng,
        T_Particle const & particle
    )
    {
        bool result = false;
        if( rng() >= float_X( 0.5 ) )
            result = true;
        return result;
    }
    static constexpr char const * name = "eachSecondParticle";

    static constexpr bool isDeterministic = false;
};

using EachSecondParticle = generic::FreeRng<
   FunctorEachSecondParticle,
   pmacc::random::distributions::Uniform< float_X >
>;

Template Parameters:
  • T_Functor – user defined unary functor

  • T_Distribution – pmacc::random::distributions, random number distribution

FreeTotalCellOffset

template<typename T_Functor>
struct FreeTotalCellOffset : protected picongpu::particles::functor::User<T_Functor>, private picongpu::particles::functor::misc::TotalCellOffset

call simple free user defined functor and provide the cell information

The functor passes the cell offset of the particle relative to the total domain origin into the functor.

example for particleFilters.param: each particle with a cell offset of 5 in X direction

struct FunctorEachParticleInXCell5
{
    template< typename T_Particle >
    HDINLINE bool operator()(
        DataSpace< simDim > const & particleOffsetToTotalOrigin,
        T_Particle const & particle
    )
    {
        bool result = false;
        if( particleOffsetToTotalOrigin.x() == 5 )
            result = true;
        return result;
    }
    static constexpr char const * name = "eachParticleInXCell5";

    static constexpr bool isDeterministic = true;
};

using EachParticleInXCell5 = generic::FreeTotalCellOffset<
   FunctorEachParticleInXCell5
>;

Template Parameters:

T_Functor – user defined unary functor