PIC Core

dimension.param

The spatial dimensionality of the simulation.

Defines

SIMDIM

Possible values: DIM3 for 3D3V and DIM2 for 2D3V.

namespace picongpu

Variables

constexpr uint32_t simDim = SIMDIM

grid.param

Definition of cell sizes and time step.

Our cells are defining a regular, cartesian grid. Our explicit FDTD field solvers define an upper bound for the time step value in relation to the cell size for convergence. Make sure to resolve important wavelengths of your simulation, e.g. shortest plasma wavelength and central laser wavelength both spatially and temporarily.

Units in reduced dimensions

In 2D3V simulations, the CELL_DEPTH_SI (Z) cell length is still used for normalization of densities, etc..

A 2D3V simulation in a cartesian PIC simulation such as ours only changes the degrees of freedom in motion for (macro) particles and all (field) information in z travels instantaneous, making the 2D3V simulation behave like the interaction of infinite “wire particles” in fields with perfect symmetry in Z.

namespace picongpu

Variables

constexpr uint32_t picongpu::ABSORBER_CELLS[3][2]= { {32, 32}, {32, 32}, {32, 32} }

Defines the size of the absorbing zone (in cells)

unit: none

constexpr float_X picongpu::ABSORBER_STRENGTH[3][2]= { {1.0e-3, 1.0e-3}, {1.0e-3, 1.0e-3}, {1.0e-3, 1.0e-3} }

Define the strength of the absorber for any direction.

unit: none

constexpr float_64 movePoint = 0.9

When to move the co-moving window.

An initial pseudo particle, flying with the speed of light, is fired at the begin of the simulation. When it reaches movePoint % of the absolute(*) simulation area, the co-moving window starts to move with the speed of light.

(*) Note: beware, that there is one “hidden” row of gpus at the y-front, when you use the co-moving window 0.75 means only 75% of simulation area is used for real simulation

Warning: this variable is deprecated, but currently still required for building purposes. Please keep the variable here. In case a moving window is enabled in your .cfg file, please set the move point using the ‘windowMovePoint’ parameter in that file, its default value is movePoint.

namespace SI

Variables

constexpr float_64 DELTA_T_SI = 0.8e-16

Duration of one timestep unit: seconds.

constexpr float_64 CELL_WIDTH_SI = 0.1772e-6

equals X unit: meter

constexpr float_64 CELL_HEIGHT_SI = 0.4430e-7

equals Y - the laser & moving window propagation direction unit: meter

constexpr float_64 CELL_DEPTH_SI = CELL_WIDTH_SI

equals Z unit: meter

components.param

Select a user-defined simulation class here, e.g.

with strongly modified initialization and/or PIC loop beyond the parametrization given in other .param files.

namespace simulation_starter

Simulation Starter Selection: This value does usually not need to be changed.

Change only if you want to implement your own SimulationHelper (e.g. MySimulation) class.

  • defaultPIConGPU : default PIConGPU configuration

fieldSolver.param

Configure the field solver.

Select the numerical Maxwell solver (e.g. Yee’s method).

Also allows to configure ad hoc mitigations for high frequency noise in some setups via current smoothing.

namespace picongpu
namespace fields

Typedefs

using CurrentInterpolation = currentInterpolation::None

Current Interpolation.

CurrentInterpolation is used to set a method performing the interpolate/assign operation from the generated currents of particle species to the electro-magnetic fields.

Allowed values are:

  • None:

    • default for staggered grids/Yee-scheme

    • updates E

  • Binomial: 2nd order Binomial filter

    • smooths the current before assignment in staggered grid

    • updates E & breaks local charge conservation slightly

  • NoneDS:

    • experimental assignment for all-centered/directional splitting

    • updates E & B at the same time

using Solver = maxwellSolver::Yee<CurrentInterpolation>

FieldSolver.

Field Solver Selection:

  • Yee< CurrentInterpolation > : standard Yee solver

  • YeePML< CurrentInterpolation >: standard Yee solver with PML absorber

  • Lehe< CurrentInterpolation >: Num. Cherenkov free field solver in a chosen direction

  • DirSplitting< CurrentInterpolation >: Sentoku’s Directional Splitting Method

  • None< CurrentInterpolation >: disable the vacuum update of E and B

laser.param

Configure laser profiles.

All laser propagate in y direction.

Available profiles:

  • None : no laser init

  • GaussianBeam : Gaussian beam (focusing)

  • PulseFrontTilt : Gaussian beam with a tilted pulse envelope in ‘x’ direction

  • PlaneWave : a plane wave (Gaussian in time)

  • Wavepacket : wavepacket (Gaussian in time and space, not focusing)

  • Polynom : a polynomial laser envelope

  • ExpRampWithPrepulse : wavepacket with exponential upramps and prepulse

In the end, this file needs to define a Selected class in namespace picongpu::fields::laserProfiles. A typical profile consists of a laser profile class and its parameters. For example:

using Selected = GaussianBeam< GaussianBeamParam >;

namespace picongpu
namespace fields
namespace laserProfiles

Typedefs

using Selected = None<>

currently selected laser profile

struct ExpRampWithPrepulseParam

Based on a wavepacket with Gaussian spatial envelope.

and the following temporal shape: A Gaussian peak (optionally lengthened by a plateau) is preceded by two pieces of exponential preramps, defined by 3 (time, intensity)- -points. The first two points get connected by an exponential, the 2nd and 3rd point are connected by another exponential, which is then extrapolated to the peak. The Gaussian is added everywhere, but typically contributes significantly only near the peak. It is advisable to set the third point far enough from the plateau (approx 3*FWHM), then the contribution from the Gaussian is negligible there, and the intensity can be set as measured from the laser profile. Optionally a Gaussian prepulse can be added, given by the parameters of the relative intensity and time point. The time of the prepulse and the three preramp points are given in SI, the intensities are given as multiples of the peak intensity.

Public Types

enum PolarisationType

Available polarisation types.

Values:

LINEAR_X = 1u
LINEAR_Z = 2u
CIRCULAR = 4u

Public Static Attributes

constexpr float_X INT_RATIO_PREPULSE = 0.
constexpr float_X INT_RATIO_POINT_1 = 1.e-8
constexpr float_X INT_RATIO_POINT_2 = 1.e-4
constexpr float_X INT_RATIO_POINT_3 = 1.e-4
constexpr float_64 TIME_PREPULSE_SI = -950.0e-15
constexpr float_64 TIME_PEAKPULSE_SI = 0.0e-15
constexpr float_64 TIME_POINT_1_SI = -1000.0e-15
constexpr float_64 TIME_POINT_2_SI = -300.0e-15
constexpr float_64 TIME_POINT_3_SI = -100.0e-15
constexpr float_64 WAVE_LENGTH_SI = 0.8e-6

unit: meter

constexpr float_64 UNITCONV_A0_to_Amplitude_SI = -2.0 * PI / WAVE_LENGTH_SI * picongpu::SI::ELECTRON_MASS_SI * picongpu::SI::SPEED_OF_LIGHT_SI * picongpu::SI::SPEED_OF_LIGHT_SI / picongpu::SI::ELECTRON_CHARGE_SI

UNITCONV.

constexpr float_64 _A0 = 20.

unit: W / m^2

unit: none

constexpr float_64 AMPLITUDE_SI = _A0 * UNITCONV_A0_to_Amplitude_SI

unit: Volt /meter

constexpr float_64 LASER_NOFOCUS_CONSTANT_SI = 0.0 * WAVE_LENGTH_SI / picongpu::SI::SPEED_OF_LIGHT_SI

unit: Volt /meter

The profile of the test Lasers 0 and 2 can be stretched by a constant area between the up and downramp unit: seconds

constexpr float_64 PULSE_LENGTH_SI = 3.0e-14 / 2.35482

Pulse length: sigma of std.

gauss for intensity (E^2) PULSE_LENGTH_SI = FWHM_of_Intensity / [ 2*sqrt{ 2* ln(2) } ] [ 2.354820045 ] Info: FWHM_of_Intensity = FWHM_Illumination = what a experimentalist calls “pulse duration” unit: seconds (1 sigma)

constexpr float_64 W0_X_SI = 2.5 * WAVE_LENGTH_SI

beam waist: distance from the axis where the pulse intensity (E^2) decreases to its 1/e^2-th part, WO_X_SI is this distance in x-direction W0_Z_SI is this distance in z-direction if both values are equal, the laser has a circular shape in x-z W0_SI = FWHM_of_Intensity / sqrt{ 2* ln(2) } [ 1.17741 ] unit: meter

constexpr float_64 W0_Z_SI = W0_X_SI
constexpr float_64 RAMP_INIT = 16.0

The laser pulse will be initialized half of PULSE_INIT times of the PULSE_LENGTH before plateau and half at the end of the plateau unit: none.

constexpr uint32_t initPlaneY = 0

cell from top where the laser is initialized

if initPlaneY == 0 than the absorber are disabled. if initPlaneY > absorbercells negative Y the negative absorber in y direction is enabled

valid ranges:

  • initPlaneY == 0

  • absorber cells negative Y < initPlaneY < cells in y direction of the top gpu

constexpr float_X LASER_PHASE = 0.0

laser phase shift (no shift: 0.0)

sin(omega*time + laser_phase): starts with phase=0 at center > E-field=0 at center

unit: rad, periodic in 2*pi

constexpr PolarisationType Polarisation = LINEAR_X

Polarization selection.

struct GaussianBeamParam

Public Types

enum PolarisationType

Available polarisation types.

Values:

LINEAR_X = 1u
LINEAR_Z = 2u
CIRCULAR = 4u
using LAGUERREMODES_t = gaussianBeam::LAGUERREMODES_t

Public Static Attributes

constexpr float_64 WAVE_LENGTH_SI = 0.8e-6

unit: meter

constexpr float_64 UNITCONV_A0_to_Amplitude_SI = -2.0 * PI / WAVE_LENGTH_SI * picongpu::SI::ELECTRON_MASS_SI * picongpu::SI::SPEED_OF_LIGHT_SI * picongpu::SI::SPEED_OF_LIGHT_SI / picongpu::SI::ELECTRON_CHARGE_SI

Convert the normalized laser strength parameter a0 to Volt per meter.

constexpr float_64 AMPLITUDE_SI = 1.738e13

unit: W / m^2

unit: none unit: Volt / meter unit: Volt / meter

constexpr float_64 PULSE_LENGTH_SI = 10.615e-15 / 4.0

Pulse length: sigma of std.

gauss for intensity (E^2) PULSE_LENGTH_SI = FWHM_of_Intensity / [ 2*sqrt{ 2* ln(2) } ] [ 2.354820045 ] Info: FWHM_of_Intensity = FWHM_Illumination = what a experimentalist calls “pulse duration”

unit: seconds (1 sigma)

constexpr float_64 W0_SI = 5.0e-6 / 1.17741

beam waist: distance from the axis where the pulse intensity (E^2) decreases to its 1/e^2-th part, at the focus position of the laser W0_SI = FWHM_of_Intensity / sqrt{ 2* ln(2) } [ 1.17741 ]

unit: meter

constexpr float_64 FOCUS_POS_SI = 4.62e-5

the distance to the laser focus in y-direction unit: meter

constexpr float_64 PULSE_INIT = 20.0

The laser pulse will be initialized PULSE_INIT times of the PULSE_LENGTH.

unit: none

constexpr uint32_t initPlaneY = 0

cell from top where the laser is initialized

if initPlaneY == 0 than the absorber are disabled. if initPlaneY > absorbercells negative Y the negative absorber in y direction is enabled

valid ranges:

  • initPlaneY == 0

  • absorber cells negative Y < initPlaneY < cells in y direction of the top gpu

constexpr float_X LASER_PHASE = 0.0

laser phase shift (no shift: 0.0)

sin(omega*time + laser_phase): starts with phase=0 at center > E-field=0 at center

unit: rad, periodic in 2*pi

constexpr uint32_t MODENUMBER = gaussianBeam::MODENUMBER
constexpr PolarisationType Polarisation = CIRCULAR

Polarization selection.

struct PlaneWaveParam

Public Types

enum PolarisationType

Available polarization types.

Values:

LINEAR_X = 1u
LINEAR_Z = 2u
CIRCULAR = 4u

Public Static Attributes

constexpr float_64 WAVE_LENGTH_SI = 0.8e-6

unit: meter

constexpr float_64 UNITCONV_A0_to_Amplitude_SI = -2.0 * PI / WAVE_LENGTH_SI * picongpu::SI::ELECTRON_MASS_SI * picongpu::SI::SPEED_OF_LIGHT_SI * picongpu::SI::SPEED_OF_LIGHT_SI / picongpu::SI::ELECTRON_CHARGE_SI

Convert the normalized laser strength parameter a0 to Volt per meter.

constexpr float_64 _A0 = 1.5

unit: W / m^2

unit: none

constexpr float_64 AMPLITUDE_SI = _A0 * UNITCONV_A0_to_Amplitude_SI

unit: Volt / meter

constexpr float_64 LASER_NOFOCUS_CONSTANT_SI = 13.34e-15

unit: Volt / meter

The profile of the test Lasers 0 and 2 can be stretched by a constant area between the up and downramp unit: seconds

constexpr float_64 PULSE_LENGTH_SI = 10.615e-15 / 4.0

Pulse length: sigma of std.

gauss for intensity (E^2) PULSE_LENGTH_SI = FWHM_of_Intensity / [ 2*sqrt{ 2* ln(2) } ] [ 2.354820045 ] Info: FWHM_of_Intensity = FWHM_Illumination = what a experimentalist calls “pulse duration” unit: seconds (1 sigma)

constexpr uint32_t initPlaneY = 0

cell from top where the laser is initialized

if initPlaneY == 0 than the absorber are disabled. if initPlaneY > absorbercells negative Y the negative absorber in y direction is enabled

valid ranges:

  • initPlaneY == 0

  • absorber cells negative Y < initPlaneY < cells in y direction of the top gpu

constexpr float_64 RAMP_INIT = 20.6146

The laser pulse will be initialized half of PULSE_INIT times of the PULSE_LENGTH before and after the plateau unit: none.

constexpr float_X LASER_PHASE = 0.0

laser phase shift (no shift: 0.0)

sin(omega*time + laser_phase): starts with phase=0 at center > E-field=0 at center

unit: rad, periodic in 2*pi

constexpr PolarisationType Polarisation = LINEAR_X

Polarization selection.

struct PolynomParam

Based on a wavepacket with Gaussian spatial envelope.

Wavepacket with a polynomial temporal intensity shape.

Public Types

enum PolarisationType

Available polarization types.

Values:

LINEAR_X = 1u
LINEAR_Z = 2u
CIRCULAR = 4u

Public Static Attributes

constexpr float_64 WAVE_LENGTH_SI = 0.8e-6

unit: meter

constexpr float_64 UNITCONV_A0_to_Amplitude_SI = -2.0 * PI / WAVE_LENGTH_SI * picongpu::SI::ELECTRON_MASS_SI * picongpu::SI::SPEED_OF_LIGHT_SI * picongpu::SI::SPEED_OF_LIGHT_SI / picongpu::SI::ELECTRON_CHARGE_SI

Convert the normalized laser strength parameter a0 to Volt per meter.

constexpr float_64 AMPLITUDE_SI = 1.738e13

unit: W / m^2

unit: none unit: Volt / meter unit: Volt / meter

constexpr float_64 LASER_NOFOCUS_CONSTANT_SI = 13.34e-15

The profile of the test Lasers 0 and 2 can be stretched by a constant area between the up and downramp unit: seconds.

constexpr float_64 PULSE_LENGTH_SI = 10.615e-15 / 4.0

Pulse length: sigma of std.

gauss for intensity (E^2) PULSE_LENGTH_SI = FWHM_of_Intensity / [ 2*sqrt{ 2* ln(2) } ] [ 2.354820045 ] Info: FWHM_of_Intensity = FWHM_Illumination = what a experimentalist calls “pulse duration” unit: seconds (1 sigma)

constexpr float_64 W0_X_SI = 4.246e-6

beam waist: distance from the axis where the pulse intensity (E^2) decreases to its 1/e^2-th part, at the focus position of the laser unit: meter

constexpr float_64 W0_Z_SI = W0_X_SI
constexpr uint32_t initPlaneY = 0

cell from top where the laser is initialized

if initPlaneY == 0 than the absorber are disabled. if initPlaneY > absorbercells negative Y the negative absorber in y direction is enabled

valid ranges:

  • initPlaneY == 0

  • absorber cells negative Y < initPlaneY < cells in y direction of the top gpu

constexpr float_64 PULSE_INIT = 20.0

The laser pulse will be initialized PULSE_INIT times of the PULSE_LENGTH.

unit: none

constexpr float_X LASER_PHASE = 0.0

laser phase shift (no shift: 0.0)

sin(omega*time + laser_phase): starts with phase=0 at center > E-field=0 at center

unit: rad, periodic in 2*pi

constexpr PolarisationType Polarisation = LINEAR_X

Polarization selection.

struct PulseFrontTiltParam

Public Types

enum PolarisationType

Available polarisation types.

Values:

LINEAR_X = 1u
LINEAR_Z = 2u
CIRCULAR = 4u

Public Static Attributes

constexpr float_64 WAVE_LENGTH_SI = 0.8e-6

unit: meter

constexpr float_64 UNITCONV_A0_to_Amplitude_SI = -2.0 * PI / WAVE_LENGTH_SI * picongpu::SI::ELECTRON_MASS_SI * picongpu::SI::SPEED_OF_LIGHT_SI * picongpu::SI::SPEED_OF_LIGHT_SI / picongpu::SI::ELECTRON_CHARGE_SI

Convert the normalized laser strength parameter a0 to Volt per meter.

constexpr float_64 AMPLITUDE_SI = 1.738e13

unit: W / m^2

unit: none unit: Volt / meter unit: Volt / meter

constexpr float_64 PULSE_LENGTH_SI = 10.615e-15 / 4.0

Pulse length: sigma of std.

gauss for intensity (E^2) PULSE_LENGTH_SI = FWHM_of_Intensity / [ 2*sqrt{ 2* ln(2) } ] [ 2.354820045 ] Info: FWHM_of_Intensity = FWHM_Illumination = what a experimentalist calls “pulse duration”

unit: seconds (1 sigma)

constexpr float_64 W0_SI = 5.0e-6 / 1.17741

beam waist: distance from the axis where the pulse intensity (E^2) decreases to its 1/e^2-th part, at the focus position of the laser W0_SI = FWHM_of_Intensity / sqrt{ 2* ln(2) } [ 1.17741 ]

unit: meter

constexpr float_64 FOCUS_POS_SI = 4.62e-5

the distance to the laser focus in y-direction unit: meter

constexpr float_64 TILT_X_SI = 0.0

the tilt angle between laser propagation in y-direction and laser axis in x-direction (0 degree == no tilt) unit: degree

constexpr float_64 PULSE_INIT = 20.0

The laser pulse will be initialized PULSE_INIT times of the PULSE_LENGTH.

unit: none

constexpr uint32_t initPlaneY = 0

cell from top where the laser is initialized

if initPlaneY == 0 than the absorber are disabled. if initPlaneY > absorbercells negative Y the negative absorber in y direction is enabled

valid ranges:

  • initPlaneY == 0

  • absorber cells negative Y < initPlaneY < cells in y direction of the top gpu

constexpr float_X LASER_PHASE = 0.0

laser phase shift (no shift: 0.0)

sin(omega*time + laser_phase): starts with phase=0 at center > E-field=0 at center

unit: rad, periodic in 2*pi

constexpr PolarisationType Polarisation = CIRCULAR

Polarization selection.

struct WavepacketParam

Public Types

enum PolarisationType

Available polarisation types.

Values:

LINEAR_X = 1u
LINEAR_Z = 2u
CIRCULAR = 4u

Public Static Attributes

constexpr float_64 WAVE_LENGTH_SI = 0.8e-6

unit: meter

constexpr float_64 UNITCONV_A0_to_Amplitude_SI = -2.0 * PI / WAVE_LENGTH_SI * picongpu::SI::ELECTRON_MASS_SI * picongpu::SI::SPEED_OF_LIGHT_SI * picongpu::SI::SPEED_OF_LIGHT_SI / picongpu::SI::ELECTRON_CHARGE_SI

Convert the normalized laser strength parameter a0 to Volt per meter.

constexpr float_64 AMPLITUDE_SI = 1.738e13

unit: W / m^2

unit: none unit: Volt / meter unit: Volt / meter

constexpr float_64 LASER_NOFOCUS_CONSTANT_SI = 7.0 * WAVE_LENGTH_SI / picongpu::SI::SPEED_OF_LIGHT_SI

The profile of the test Lasers 0 and 2 can be stretched by a constant area between the up and downramp unit: seconds.

constexpr float_64 PULSE_LENGTH_SI = 10.615e-15 / 4.0

Pulse length: sigma of std.

gauss for intensity (E^2) PULSE_LENGTH_SI = FWHM_of_Intensity / [ 2*sqrt{ 2* ln(2) } ] [ 2.354820045 ] Info: FWHM_of_Intensity = FWHM_Illumination = what a experimentalist calls “pulse duration”

unit: seconds (1 sigma)

constexpr float_64 W0_X_SI = 4.246e-6

beam waist: distance from the axis where the pulse intensity (E^2) decreases to its 1/e^2-th part, at the focus position of the laser W0_SI = FWHM_of_Intensity / sqrt{ 2* ln(2) } [ 1.17741 ]

unit: meter

constexpr float_64 W0_Z_SI = W0_X_SI
constexpr float_64 PULSE_INIT = 20.0

The laser pulse will be initialized PULSE_INIT times of the PULSE_LENGTH.

unit: none

constexpr uint32_t initPlaneY = 0

cell from top where the laser is initialized

if initPlaneY == 0 than the absorber are disabled. if initPlaneY > absorbercells negative Y the negative absorber in y direction is enabled

valid ranges:

  • initPlaneY == 0

  • absorber cells negative Y < initPlaneY < cells in y direction of the top gpu

constexpr float_X LASER_PHASE = 0.0

laser phase shift (no shift: 0.0)

sin(omega*time + laser_phase): starts with phase=0 at center > E-field=0 at center

unit: rad, periodic in 2*pi

constexpr PolarisationType Polarisation = LINEAR_X

Polarization selection.

namespace gaussianBeam

Functions

picongpu::fields::laserProfiles::gaussianBeam::PMACC_CONST_VECTOR(float_X, MODENUMBER+ 1, LAGUERREMODES, 1. 0)

Variables

constexpr uint32_t MODENUMBER = 0

Use only the 0th Laguerremode for a standard Gaussian.

List of available laser profiles.

pml.param

Configure the perfectly matched layer (PML).

To enable PML use YeePML field solver.

namespace picongpu
namespace fields
namespace maxwellSolver
namespace yeePML

Variables

constexpr uint32_t THICKNESS = 8
constexpr uint32_t picongpu::fields::maxwellSolver::yeePML::NUM_CELLS[3][2]= { { THICKNESS, THICKNESS }, { THICKNESS, THICKNESS }, { THICKNESS, THICKNESS } }

Thickness of the absorbing layer, in number of cells.

PML is located inside the global simulation area, near the outer borders. Setting size to 0 results in disabling absorption at the corresponding boundary. Normally thickness is between 6 and 16 cells, with larger values providing less reflections. 8 cells should be good enough for most simulations. There are no requirements on thickness being a multiple of the supercell size. It is only required that PML is small enough to be fully contained in a single layer of local domains near the global simulation area boundary (Note that the domains of this layer might be changing, e.g. due to moving window.) Unit: number of cells.

constexpr float_64 SIGMA_KAPPA_GRADING_ORDER = 4.0

Order of polynomial grading for artificial electric conductivity and stretching coefficient.

The conductivity (sigma) is polynomially scaling from 0 at the internal border of PML to the maximum value (defined below) at the external border. The stretching coefficient (kappa) scales from 1 to the corresponding maximum value (defined below) with the same polynomial. The grading is given in [Taflove, Hagness], eq. (7.60a, b), with the order denoted ‘m’. Must be >= 0. Normally between 3 and 4, not required to be integer. Unitless.

constexpr float_64 SIGMA_OPT_SI[3] = {0.8 * (SIGMA_KAPPA_GRADING_ORDER + 1.0) / (SI::Z0_SI * SI::CELL_WIDTH_SI), , }
constexpr float_64 SIGMA_OPT_MULTIPLIER = 1.0
constexpr float_64 SIGMA_MAX_SI[3] = {SIGMA_OPT_SI[0] * SIGMA_OPT_MULTIPLIER, , }

Max value of artificial electric conductivity in PML.

Components correspond to directions: element 0 corresponds to absorption along x direction, 1 = y, 2 = z. Grading is described in comments for SIGMA_KAPPA_GRADING_ORDER. Too small values lead to significant reflections from the external border, too large - to reflections due to discretization errors. Artificial magnetic permeability will be chosen to perfectly match this. Must be >= 0. Normally between 0.7 * SIGMA_OPT_SI and 1.1 * SIGMA_OPT_SI. Unit: siemens / m.

constexpr float_64 KAPPA_MAX[3] = {1.0, , }

Max value of coordinate stretching coefficient in PML.

Components correspond to directions: element 0 corresponds to absorption along x direction, 1 = y, 2 = z. Grading is described in comments for SIGMA_KAPPA_GRADING_ORDER. Must be >= 1. For relatively homogeneous domains 1.0 is a reasonable value. Highly elongated domains can have better absorption with values between 7.0 and 20.0, for example, see section 7.11.2 in [Taflove, Hagness]. Unitless.

constexpr float_64 ALPHA_GRADING_ORDER = 1.0

Order of polynomial grading for complex frequency shift.

The complex frequency shift (alpha) is polynomially downscaling from the maximum value (defined below) at the internal border of PML to 0 at the external border. The grading is given in [Taflove, Hagness], eq. (7.79), with the order denoted ‘m_a’. Must be >= 0. Normally values are around 1.0. Unitless.

constexpr float_64 ALPHA_MAX_SI[3] = {0.2, , }

Complex frequency shift in PML.

Components correspond to directions: element 0 corresponds to absorption along x direction, 1 = y, 2 = z. Setting it to 0 will make PML behave as uniaxial PML. Setting it to a positive value helps to attenuate evanescent modes, but can degrade absorption of propagating modes, as described in section 7.7 and 7.11.3 in [Taflove, Hagness]. Must be >= 0. Normally values are 0 or between 0.15 and 0.3. Unit: siemens / m.

pusher.param

Configure particle pushers.

Those pushers can then be selected by a particle species in species.param and speciesDefinition.param

namespace picongpu
struct particlePusherAccelerationParam

Subclassed by picongpu::particlePusherAcceleration::UnitlessParam

Public Static Attributes

constexpr float_64 AMPLITUDEx_SI = 0.0

Define strength of constant and homogeneous accelerating electric field in SI per dimension.

unit: Volt / meter

constexpr float_64 AMPLITUDEy_SI = -1.e11

The moving window propagation direction unit: Volt / meter (1e11 V/m = 1 GV/cm)

constexpr float_64 AMPLITUDEz_SI = 0.0

unit: Volt / meter

constexpr float_64 ACCELERATION_TIME_SI = 10000.0 * picongpu::SI::DELTA_T_SI

Acceleration duration unit: second.

namespace particlePusherAxel

Enums

enum TrajectoryInterpolationType

Values:

LINEAR = 1u
NONLINEAR = 2u

Variables

constexpr TrajectoryInterpolationType TrajectoryInterpolation = LINEAR
namespace particlePusherProbe

Typedefs

using ActualPusher = void

Also push the probe particles?

In many cases, probe particles are static throughout the simulation. This option allows to set an “actual” pusher that shall be used to also change the probe particle positions.

Examples:

  • particles::pusher::Boris

  • particles::pusher::[all others from above]

  • void (no push)

density.param

Configure existing or define new normalized density profiles here.

During particle species creation in speciesInitialization.param, those profiles can be translated to spatial particle distributions.

namespace picongpu
namespace densityProfiles

Typedefs

using Gaussian = GaussianImpl<GaussianParam>
using Homogenous = HomogenousImpl
using LinearExponential = LinearExponentialImpl<LinearExponentialParam>
using GaussianCloud = GaussianCloudImpl<GaussianCloudParam>
using SphereFlanks = SphereFlanksImpl<SphereFlanksParam>
using FromHDF5 = FromHDF5Impl<FromHDF5Param>
using FreeFormula = FreeFormulaImpl<FreeFormulaFunctor>

Functions

picongpu::densityProfiles::PMACC_STRUCT(GaussianParam, ( PMACC_C_VALUE (float_X, gasFactor, -1.0))( PMACC_C_VALUE (float_X, gasPower, 4.0))( PMACC_C_VALUE (uint32_t, vacuumCellsY, 50))( PMACC_C_VALUE (float_64, gasCenterLeft_SI, 4.62e-5))(PMACC_C_VALUE(float_64, gasCenterRight_SI, 4.62e-5))(PMACC_C_VALUE(float_64, gasSigmaLeft_SI, 4.62e-5))(PMACC_C_VALUE(float_64, gasSigmaRight_SI, 4.62e-5)))

Profile Formula: const float_X exponent = abs((y - gasCenter_SI) / gasSigma_SI); const float_X density = exp(gasFactor * pow(exponent, gasPower));

takes gasCenterLeft_SI for y < gasCenterLeft_SI, gasCenterRight_SI for y > gasCenterRight_SI, and exponent = 0.0 for gasCenterLeft_SI < y < gasCenterRight_SI

picongpu::densityProfiles::PMACC_STRUCT(LinearExponentialParam, ( PMACC_C_VALUE (uint32_t, vacuumCellsY, 50))( PMACC_C_VALUE (float_64, gasYMax_SI, 1.0e-3))(PMACC_C_VALUE(float_64, gasA_SI, 1.0e-3))(PMACC_C_VALUE(float_64, gasD_SI, 1.0e-3))(PMACC_C_VALUE(float_64, gasB, 0.0)))

parameter for LinearExponential profile

* Density Profile: /\
*                 /  -,_
*   linear       /      -,_    exponential
*   slope       /  |       -,_ slope
*                  MAX
*

picongpu::densityProfiles::PMACC_STRUCT(GaussianCloudParam, ( PMACC_C_VALUE (float_X, gasFactor, -0.5))( PMACC_C_VALUE (float_X, gasPower, 2.0))( PMACC_C_VALUE (uint32_t, vacuumCellsY, 50))( PMACC_C_VECTOR_DIM (float_64, simDim, center_SI, 1.134e-5, 1.134e-5, 1.134e-5))(PMACC_C_VECTOR_DIM(float_64, simDim, sigma_SI, 7.0e-6, 7.0e-6, 7.0e-6)))
picongpu::densityProfiles::PMACC_STRUCT(SphereFlanksParam, ( PMACC_C_VALUE (uint32_t, vacuumCellsY, 50))( PMACC_C_VALUE (float_64, r_SI, 1.0e-3))(PMACC_C_VALUE(float_64, ri_SI, 0.0))(PMACC_C_VECTOR_DIM(float_64, simDim, center_SI, 8.0e-3, 8.0e-3, 8.0e-3))(PMACC_C_VALUE(float_64, exponent_SI, 1.0e3)))

The profile consists out of the composition of 3 1D profiles with the scheme: exponential increasing flank, constant sphere, exponential decreasing flank.

*           ___
*  1D:  _,./   \.,_   rho(r)
*
*  2D:  ..,x,..   density: . low
*       .,xxx,.            , middle
*       ..,x,..            x high (constant)
*

picongpu::densityProfiles::PMACC_STRUCT(FromHDF5Param, ( PMACC_C_STRING (filename,"gas"))(PMACC_C_STRING(datasetName,"fields/e_chargeDensity"))(PMACC_C_VALUE(uint32_t, iteration, 0))( PMACC_C_VALUE (float_X, defaultDensity, 0.0)))
struct FreeFormulaFunctor

Public Functions

HDINLINE float_X picongpu::densityProfiles::FreeFormulaFunctor::operator()(const floatD_64 & position_SI, const float3_64 & cellSize_SI)

This formula uses SI quantities only.

The profile will be multiplied by BASE_DENSITY_SI.

Return

float_X density [normalized to 1.0]

Parameters
  • position_SI: total offset including all slides [meter]

  • cellSize_SI: cell sizes [meter]

namespace SI

Variables

constexpr float_64 BASE_DENSITY_SI = 1.e25

Base density in particles per m^3 in the density profiles.

This is often taken as reference maximum density in normalized profiles. Individual particle species can define a densityRatio flag relative to this value.

unit: ELEMENTS/m^3

speciesAttributes.param

This file defines available attributes that can be stored with each particle of a particle species.

Each attribute defined here needs to implement furthermore the traits

  • Unit

  • UnitDimension

  • WeightingPower

  • MacroWeighted in speciesAttributes.unitless for further information about these traits see therein.

namespace picongpu

Functions

alias(position)

relative (to cell origin) in-cell position of a particle

With this definition we do not define any type like float3_X, float3_64, … This is only a name without a specialization.

value_identifier(uint64_t, particleId, IdProvider<simDim>::getNewId())

unique identifier for a particle

picongpu::value_identifier(floatD_X, position_pic, floatD_X::create (0.))

specialization for the relative in-cell position

picongpu::value_identifier(float3_X, momentum, float3_X::create (0.))

momentum at timestep t

picongpu::value_identifier(float3_X, momentumPrev1, float3_X::create (0._X))

momentum at (previous) timestep t-1

picongpu::value_identifier(float_X, weighting, 0. _X)

weighting of the macro particle

picongpu::value_identifier(int16_t, voronoiCellId, - 1)

Voronoi cell of the macro particle.

picongpu::value_identifier(float3_X, probeE, float3_X::create (0.))

interpolated electric field with respect to particle shape

picongpu::value_identifier(float3_X, probeB, float3_X::create (0.))

interpolated electric field with respect to particle shape

picongpu::value_identifier(bool, radiationMask, false)

masking a particle for radiation

The mask is used by the user defined filter RadiationParticleFilter in radiation.param to (de)select particles for the radiation calculation.

picongpu::value_identifier(bool, transitionRadiationMask, false)

masking a particle for transition radiation

The mask is used by the user defined filter TransitionRadiationParticleFilter in transitionRadiation.param to (de)select particles for the transition radiation calculation.

picongpu::value_identifier(float_X, boundElectrons, 0. _X)

number of electrons bound to the atom / ion

value type is float_X to avoid casts during the runtime

  • float_X instead of integer types are reasonable because effective charge numbers are possible

  • required for ion species if ionization is enabled

  • setting it requires atomicNumbers to also be set

picongpu::value_identifier(flylite::Superconfig, superconfig, flylite::Superconfig::create (0.))

atomic superconfiguration

atomic configuration of an ion for collisional-radiative modeling, see also flylite.param

value_identifier(DataSpace<simDim>, totalCellIdx, DataSpace<simDim>())

Total cell index of a particle.

The total cell index is a N-dimensional DataSpace given by a GPU’s globalDomain.offset + localDomain.offset added to the N-dimensional cell index the particle belongs to on that GPU.

alias(shape)

alias for particle shape, see also species.param

alias(particlePusher)

alias for particle pusher, see alsospecies.param

alias(ionizers)

alias for particle ionizers, see also ionizer.param

alias(ionizationEnergies)

alias for ionization energy container, see also ionizationEnergies.param

alias(synchrotronPhotons)

alias for synchrotronPhotons, see also speciesDefinition.param

alias for ion species used for bremsstrahlung

alias(bremsstrahlungPhotons)

alias for photon species used for bremsstrahlung

alias(interpolation)

alias for particle to field interpolation, see also species.param

alias(current)

alias for particle current solver, see also species.param

alias(atomicNumbers)

alias for particle flag: atomic numbers, see also ionizer.param

  • only reasonable for atoms / ions / nuclei

  • is required when boundElectrons is set

alias(effectiveNuclearCharge)

alias for particle flag: effective nuclear charge,

  • see also ionizer.param

  • only reasonable for atoms / ions / nuclei

alias(populationKinetics)

alias for particle population kinetics model (e.g.

FLYlite)

see also flylite.param

alias(massRatio)

alias for particle mass ratio

mass ratio between base particle, see also speciesConstants.param SI::BASE_MASS_SI and a user defined species

default value: 1.0 if unset

alias(chargeRatio)

alias for particle charge ratio

charge ratio between base particle, see also speciesConstants.param SI::BASE_CHARGE_SI and a user defined species

default value: 1.0 if unset

alias(densityRatio)

alias for particle density ratio

density ratio between default density, see also density.param SI::BASE_DENSITY_SI and a user defined species

default value: 1.0 if unset

alias(exchangeMemCfg)

alias to reserved bytes for each communication direction

This is an optional flag and overwrites the default species configuration in memory.param.

A memory config must be of the following form:

struct ExampleExchangeMemCfg
{
    static constexpr uint32_t BYTES_EXCHANGE_X = 5 * 1024 * 1024;
    static constexpr uint32_t BYTES_EXCHANGE_Y = 5 * 1024 * 1024;
    static constexpr uint32_t BYTES_EXCHANGE_Z = 5 * 1024 * 1024;
    static constexpr uint32_t BYTES_CORNER = 16 * 1024;
    static constexpr uint32_t BYTES_EDGES = 16 * 1024;
};

alias(boundaryCondition)

alias to specify the boundary condition for particles

The default behavior if this alias is not given to a species is that the particles which leave the global simulation box where deleted. This also notifies all plugins that can handle leaving particles.

Note: alias boundaryCondition will be ignored if the runtime parameter --periodic is set.

The following species attributes are defined by PMacc and always stored with a particle:

namespace pmacc

Functions

pmacc::value_identifier(lcellId_t, localCellIdx, 0)

cell of a particle inside a supercell

Value is a linear cell index inside the supercell

pmacc::value_identifier(uint8_t, multiMask, 0)

state of a particle

Particle might be valid or invalid in a particle frame. Valid particles can further be marked as candidates to leave a supercell. Possible multiMask values are:

  • 0 (zero): no particle (invalid)

  • 1: particle (valid)

  • 2 to 27: (valid) particle that is about to leave its supercell but is still stored in the current particle frame. Directions to leave the supercell are defined as follows. An ExchangeType = value - 1 (e.g. 27 - 1 = 26) means particle leaves supercell in the direction of FRONT(value=18) && TOP(value=6) && LEFT(value=2) which defines a diagonal movement over a supercell corner (18+6+2=26).

speciesConstants.param

Constants and thresholds for particle species.

Defines the reference mass and reference charge to express species with (default: electrons with negative charge).

namespace picongpu

Variables

constexpr float_X picongpu::GAMMA_THRESH = 1.005_X

Threshold between relativistic and non-relativistic regime.

Threshold used for calculations that want to separate between high-precision formulas for relativistic and non-relativistic use-cases, e.g. energy-binning algorithms.

constexpr float_X picongpu::GAMMA_INV_SQUARE_RAD_THRESH = 0.18_X

Threshold in radiation plugin between relativistic and non-relativistic regime.

This limit is used to decide between a pure 1-sqrt(1-x) calculation and a 5th order Taylor approximation of 1-sqrt(1-x) to avoid halving of significant digits due to the sqrt() evaluation at x = 1/gamma^2 near 0.0. With 0.18 the relative error between Taylor approximation and real value will be below 0.001% = 1e-5 * for x=1/gamma^2 < 0.18

namespace SI

Variables

constexpr float_64 BASE_MASS_SI = ELECTRON_MASS_SI

base particle mass

reference for massRatio in speciesDefinition.param

unit: kg

constexpr float_64 BASE_CHARGE_SI = ELECTRON_CHARGE_SI

base particle charge

reference for chargeRatio in speciesDefinition.param

unit: C

species.param

Forward declarations for speciesDefinition.param in case one wants to use the same particle shape, interpolation, current solver and particle pusher for all particle species.

namespace picongpu

Typedefs

using UsedParticleShape = particles::shapes::TSC

Particle Shape definitions.

  • particles::shapes::CIC : 1st order

  • particles::shapes::TSC : 2nd order

  • particles::shapes::PCS : 3rd order

  • particles::shapes::P4S : 4th order

example: using UsedParticleShape = particles::shapes::CIC;

using UsedField2Particle = FieldToParticleInterpolation<UsedParticleShape, AssignedTrilinearInterpolation>

define which interpolation method is used to interpolate fields to particles

using UsedParticleCurrentSolver = currentSolver::Esirkepov<UsedParticleShape>

select current solver method

  • currentSolver::Esirkepov< SHAPE > : particle shapes - CIC, TSC, PCS, P4S (1st to 4th order)

  • currentSolver::VillaBune<> : particle shapes - CIC (1st order) only

  • currentSolver::EmZ< SHAPE > : particle shapes - CIC, TSC, PCS, P4S (1st to 4th order)

For development purposes:

  • currentSolver::currentSolver::EsirkepovNative< SHAPE > : generic version of currentSolverEsirkepov without optimization (~4x slower and needs more shared memory)

using UsedParticlePusher = particles::pusher::Boris

particle pusher configuration

Defining a pusher is optional for particles

  • particles::pusher::Vay : better suited relativistic boris pusher

  • particles::pusher::Boris : standard boris pusher

  • particles::pusher::ReducedLandauLifshitz : 4th order RungeKutta pusher with classical radiation reaction

For diagnostics & modeling: ———————————————

  • particles::pusher::Acceleration : Accelerate particles by applying a constant electric field

  • particles::pusher::Free : free propagation, ignore fields (= free stream model)

  • particles::pusher::Photon : propagate with c in direction of normalized mom.

  • particles::pusher::Probe : Probe particles that interpolate E & B For development purposes: ———————————————–

  • particles::pusher::Axel : a pusher developed at HZDR during 2011 (testing)

speciesDefinition.param

Define particle species.

This file collects all previous declarations of base (reference) quantities and configured solvers for species and defines particle species. This includes “attributes” (lvalues to store with each species) and “flags” (rvalues & aliases for solvers to perform with the species for each timestep and ratios to base quantities). With those information, a Particles class is defined for each species and then collected in the list VectorAllSpecies.

namespace picongpu

Typedefs

using DefaultParticleAttributes = MakeSeq_t<position<position_pic>, momentum, weighting>

describe attributes of a particle

using ParticleFlagsPhotons = MakeSeq_t<particlePusher<particles::pusher::Photon>, shape<UsedParticleShape>, interpolation<UsedField2Particle>, massRatio<MassRatioPhotons>, chargeRatio<ChargeRatioPhotons>>
using PIC_Photons = Particles<PMACC_CSTRING("ph"), ParticleFlagsPhotons, DefaultParticleAttributes>
using ParticleFlagsElectrons = MakeSeq_t<particlePusher<UsedParticlePusher>, shape<UsedParticleShape>, interpolation<UsedField2Particle>, current<UsedParticleCurrentSolver>, massRatio<MassRatioElectrons>, chargeRatio<ChargeRatioElectrons>>
using PIC_Electrons = Particles<PMACC_CSTRING("e"), ParticleFlagsElectrons, DefaultParticleAttributes>
using ParticleFlagsIons = MakeSeq_t<particlePusher<UsedParticlePusher>, shape<UsedParticleShape>, interpolation<UsedField2Particle>, current<UsedParticleCurrentSolver>, massRatio<MassRatioIons>, chargeRatio<ChargeRatioIons>, densityRatio<DensityRatioIons>, atomicNumbers<ionization::atomicNumbers::Hydrogen_t>>
using PIC_Ions = Particles<PMACC_CSTRING("i"), ParticleFlagsIons, DefaultParticleAttributes>
using VectorAllSpecies = MakeSeq_t<PIC_Electrons, PIC_Ions>

All known particle species of the simulation.

List all defined particle species from above in this list to make them available to the PIC algorithm.

Functions

picongpu::value_identifier(float_X, MassRatioPhotons, 0. 0)
picongpu::value_identifier(float_X, ChargeRatioPhotons, 0. 0)
picongpu::value_identifier(float_X, MassRatioElectrons, 1. 0)
picongpu::value_identifier(float_X, ChargeRatioElectrons, 1. 0)
picongpu::value_identifier(float_X, MassRatioIons, 1836. 152672)
picongpu::value_identifier(float_X, ChargeRatioIons, -1. 0)
picongpu::value_identifier(float_X, DensityRatioIons, 1. 0)

particle.param

Configurations for particle manipulators.

Set up and declare functors that can be used in speciesInitalization.param for particle species initialization and manipulation, such as temperature distributions, drifts, pre-ionization and in-cell position.

namespace picongpu
namespace particles

Variables

constexpr float_X MIN_WEIGHTING = 10.0

a particle with a weighting below MIN_WEIGHTING will not be created / will be deleted

unit: none

constexpr uint32_t TYPICAL_PARTICLES_PER_CELL = 2u

Number of maximum particles per cell during density profile evaluation.

Determines the weighting of a macro particle and with it, the number of particles “sampling” dynamics in phase space.

namespace manipulators

Typedefs

using AssignXDrift = unary::Drift<DriftParam, nvidia::functors::Assign>

definition of manipulator that assigns a drift in X

using AddTemperature = unary::Temperature<TemperatureParam>
using DoubleWeighting = generic::Free<DoubleWeightingFunctor>

definition of a free particle manipulator: double weighting

using RandomEnabledRadiation = generic::FreeRng<RandomEnabledRadiationFunctor, pmacc::random::distributions::Uniform<float_X>>
using RandomPosition = unary::RandomPosition

changes the in-cell position of each particle of a species

Functions

picongpu::particles::manipulators::CONST_VECTOR(float_X, 3, DriftParam_direction, 1. 0, 0. 0, 0. 0)

Parameter for DriftParam.

struct DoubleWeightingFunctor

Unary particle manipulator: double each weighting.

Public Functions

template<typename T_Particle>DINLINE void picongpu::particles::manipulators::DoubleWeightingFunctor::operator()(T_Particle & particle)
struct DriftParam

Parameter for a particle drift assignment.

Public Members

const DriftParam_direction_t direction

Public Static Attributes

constexpr float_64 gamma = 1.0
struct RandomEnabledRadiationFunctor

Public Functions

template<typename T_Rng, typename T_Particle>DINLINE void picongpu::particles::manipulators::RandomEnabledRadiationFunctor::operator()(T_Rng & rng, T_Particle & particle)
struct TemperatureParam

Parameter for a temperature assignment.

Public Static Attributes

constexpr float_64 temperature = 0.0
namespace startPosition

Typedefs

using Random = RandomImpl<RandomParameter>

definition of random particle start

using Quiet = QuietImpl<QuietParam>

definition of quiet particle start

using OnePosition = OnePositionImpl<OnePositionParameter>

definition of one specific position for particle start

Functions

picongpu::particles::startPosition::CONST_VECTOR(float_X, 3, InCellOffset, 0. 0, 0. 0, 0. 0)

sit directly in lower corner of the cell

struct OnePositionParameter

Public Members

const InCellOffset_t inCellOffset

Public Static Attributes

constexpr uint32_t numParticlesPerCell = TYPICAL_PARTICLES_PER_CELL

Count of particles per cell at initial state.

unit: none

struct QuietParam

Public Types

using numParticlesPerDimension = mCT::shrinkTo<mCT::Int<1, TYPICAL_PARTICLES_PER_CELL, 1>, simDim>::type

Count of particles per cell per direction at initial state.

unit: none

struct RandomParameter

Public Static Attributes

constexpr uint32_t numParticlesPerCell = TYPICAL_PARTICLES_PER_CELL

Count of particles per cell at initial state.

unit: none

More details on the order of initialization of particles inside a particle species can be found here.

List of all pre-defined particle manipulators.

unit.param

In this file we define typical scales for normalization of physical quantities aka “units”.

Usually, a user would not change this file but might use the defined constants in other input files.

namespace picongpu

Variables

constexpr float_64 UNIT_TIME = SI::DELTA_T_SI

Unit of time.

constexpr float_64 UNIT_LENGTH = UNIT_TIME * UNIT_SPEED

Unit of length.

constexpr float_64 UNIT_MASS = SI::BASE_MASS_SI * double(particles::TYPICAL_NUM_PARTICLES_PER_MACROPARTICLE)

Unit of mass.

constexpr float_64 UNIT_CHARGE = -1.0 * SI::BASE_CHARGE_SI * double(particles::TYPICAL_NUM_PARTICLES_PER_MACROPARTICLE)

Unit of charge.

constexpr float_64 UNIT_ENERGY = (UNIT_MASS * UNIT_LENGTH * UNIT_LENGTH / (UNIT_TIME * UNIT_TIME))

Unit of energy.

constexpr float_64 UNIT_EFIELD = 1.0 / (UNIT_TIME * UNIT_TIME / UNIT_MASS / UNIT_LENGTH * UNIT_CHARGE)

Unit of EField: V/m.

constexpr float_64 UNIT_BFIELD = (UNIT_MASS / (UNIT_TIME * UNIT_CHARGE))
namespace particles

Variables

constexpr float_X TYPICAL_NUM_PARTICLES_PER_MACROPARTICLE = float_64(SI::BASE_DENSITY_SI * SI::CELL_WIDTH_SI * SI::CELL_HEIGHT_SI * SI::CELL_DEPTH_SI) / float_64(particles::TYPICAL_PARTICLES_PER_CELL)

Number of particles per makro particle (= macro particle weighting) unit: none.

particleFilters.param

A common task in both modeling and in situ processing (output) is the selection of particles of a particle species by attributes.

Users can define such selections as particle filters in this file.

Particle filters are simple mappings assigning each particle of a species either true or false (ignore / filter out).

All active filters need to be listed in AllParticleFilters. They are then combined with VectorAllSpecies at compile-time, e.g. for plugins.

namespace picongpu
namespace particles
namespace filter

Typedefs

using AllParticleFilters = MakeSeq_t<All>

Plugins: collection of all available particle filters.

Create a list of all filters here that you want to use in plugins.

Note: filter All is defined in picongpu/particles/filter/filter.def

List of all pre-defined particle filters.

speciesInitialization.param

Initialize particles inside particle species.

This is the final step in setting up particles (defined in speciesDefinition.param) via density profiles (defined in density.param). One can then further derive particles from one species to another and manipulate attributes with “manipulators” and “filters” (defined in particle.param and particleFilters.param).

namespace picongpu
namespace particles

Typedefs

using InitPipeline = bmpl::vector<>

InitPipeline defines in which order species are initialized.

the functors are called in order (from first to last functor)

List of all initialization methods for particle species.