PIC Core

dimension.param

The spatial dimensionality of the simulation.

Defines

SIMDIM

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

namespace picongpu

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

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 require 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, Debye length 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 instantaneously, making the 2D3V simulation behave like the interaction of infinite “wire particles” in fields with perfect symmetry in Z.

namespace picongpu

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

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 picongpu

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

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. Simulation) class.

  • defaultPIConGPU : default PIConGPU configuration

iterationStart.param

Specify a sequence of functors to be called at start of each time iteration.

namespace picongpu

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

Typedefs

using IterationStartPipeline = pmacc::mp_list<>

IterationStartPipeline defines the functors called at each iteration start.

The functors will be called in the given order.

The functors must be default-constructible and take the current time iteration as the only parameter. These are the same requirements as for functors in particles::InitPipeline.

fieldSolver.param

Configure the field solver.

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

Attention

Currently, the laser initialization in PIConGPU is implemented to work with the standard Yee solver. Using a solver of higher order will result in a slightly increased laser amplitude and energy than expected.

namespace picongpu

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

namespace fields

Typedefs

using Solver = maxwellSolver::Yee

FieldSolver.

Field Solver Selection (note <> for some solvers), all in namespace maxwellSolver:

  • Yee: Standard Yee solver approximating derivatives with respect to time and space by second order finite differences.

  • CKC: Cole-Karkkainen-Cowan Solver, Dispersion free solver in the direction of the smallest grid size.

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

  • ArbitraryOrderFDTD<4>: Solver using 4 neighbors to each direction to approximate spatial derivatives by finite differences. The number of neighbors can be changed from 4 to any positive, integer number. The order of the solver will be twice the number of neighbors in each direction. Yee’s method is a special case of this using one neighbor to each direction.

  • Substepping<Solver, 4>: use the given Solver (Yee, etc.) and substep each time step by factor 4

  • None: disable the vacuum update of E and B, including no J contribution to E

fieldAbsorber.param

Configure the field absorber parameters.

Field absorber type is set by command-line option &#8212;fieldAbsorber.

namespace picongpu

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

namespace fields
namespace absorber

Variables

constexpr uint32_t THICKNESS = 12
constexpr uint32_t NUM_CELLS[3][2] = {{THICKNESS, THICKNESS}, {THICKNESS, THICKNESS}, {THICKNESS, THICKNESS}}

Thickness of the absorbing layer, in number of cells.

This setting applies to applies for all absorber kinds. The absorber layer is located inside the global simulation area, near the outer borders. Setting size to 0 results in disabling absorption at the corresponding boundary. Note that for non-absorbing boundaries the actual thickness will be 0 anyways. There are no requirements on thickness being a multiple of the supercell size.

For PML the recommended thickness is between 6 and 16 cells. For the exponential damping it is 32.

Unit: number of cells.

namespace exponential

Settings for the Exponential absorber.

Variables

constexpr float_X 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 all directions.

Elements corredponding to non-absorber borders will have no effect.

Unit: none

namespace pml

Settings for the Pml absorber.

These parameters can generally be left with default values. For more details on the meaning of the parameters, refer to the following references. J.A. Roden, S.D. Gedney. Convolution PML (CPML): An efficient FDTD implementation of the CFS - PML for arbitrary media. Microwave and optical technology letters. 27 (5), 334-339 (2000). A. Taflove, S.C. Hagness. Computational Electrodynamics. The Finite-Difference Time-Domain Method. Third Edition. Artech house, Boston (2005), referred to as [Taflove, Hagness].

Variables

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), 0.8 * (SIGMA_KAPPA_GRADING_ORDER + 1.0) / (SI::Z0_SI * SI::CELL_HEIGHT_SI), 0.8 * (SIGMA_KAPPA_GRADING_ORDER + 1.0) / (SI::Z0_SI * SI::CELL_DEPTH_SI)}
constexpr float_64 SIGMA_OPT_MULTIPLIER = 1.0
constexpr float_64 SIGMA_MAX_SI[3] = {SIGMA_OPT_SI[0] * SIGMA_OPT_MULTIPLIER, SIGMA_OPT_SI[1] * SIGMA_OPT_MULTIPLIER, SIGMA_OPT_SI[2] * 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, 1.0, 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, 0.2, 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.

incidentField.param

Configure incident field profile and offset of the Huygens surface for each boundary.

Available profiles:

  • profiles::DispersivePulse<> : Gaussian pulse allowing to set first-, second-, and third-order dispersion in focus. That is, SD, AD, GDD, and TOD, respectively.

  • profiles::ExpRampWithPrepulse<> : exponential ramp with prepulse wavepacket with given parameters

  • profiles::Free<> : custom profile with user-provided functors to calculate incident E and B

  • profiles::GaussianPulse<> : Pulse with Gaussian profile in all three dimensions with given parameters

  • profiles::None : no incident field

  • profiles::PlaneWave<> : plane wave profile with given parameters

  • profiles::Polynom<> : wavepacket with a polynomial temporal intensity shape profile with given parameters

  • profiles::PulseFrontTilt<> : GaussianPulse with tilted pulse front with given parameters

  • profiles::Wavepacket<> : wavepacket with Gaussian spatial and temporal envelope profile with given parameters

All profiles but Free<> and None are parametrized with a profile-specific structure. Their interfaces are defined in the corresponding .def files inside directory picongpu/fields/incidentField/profiles/. Note that all these parameter structures inherit common base structures from BaseParam.def. Thus, a user-provided structure must also define all members according to the base struct.

In the end, this file needs to define XMin, XMax, YMin, YMax, ZMin, ZMax (the latter two can be skipped in 2d) type aliases in namespace picongpu::fields::incidentField. Each of them could be a single profile or a typelist of profiles created with MakeSeq_t. In case a typelist is used, the resulting field is a sum of effects of all profiles in the list. This file also has to define constexpr array POSITION that controls positioning of the generating surface relative to total domain. For example:

using XMin = profiles::Free<UserFunctorIncidentE, UserFunctorIncidentB>;
using XMax = profiles::None;
using YMin = MakeSeq_t<profiles::PlaneWave<UserPlaneWaveParams>, profiles::Wavepacket<UserWavepacketParams>>;
using YMax = profiles::None;
using ZMin = profiles::Polynom<UserPolynomParams>;
using ZMax = profiles::GaussianPulse<UserGaussianPulseParams>;

constexpr int32_t POSITION[3][2] = { {16, -16}, {16, -16}, {16, -16} };
namespace picongpu

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

namespace fields
namespace incidentField

Unnamed Group

using XMin = profiles::None

Incident field profile types along each boundary, these 6 types (or aliases) are required.

using XMax = profiles::None
using YMin = profiles::None
using YMax = profiles::None
using ZMin = profiles::None
using ZMax = profiles::None

Typedefs

using MyProfile = profiles::Free<UserFunctorIncidentE, UserFunctorIncidentB>

Make a profile with the user-provided functors defined above.

The functors for incident field E and B must be consistent to each other. They should work for all boundaries the profile is applied to.

Variables

constexpr int32_t POSITION[3][2] = {{16, -16}, {16, -16}, {16, -16}}

Position in cells of the Huygens surface relative to start of the total domain.

The position is set as an offset, in cells, counted from the start of the total domain. For the max boundaries, negative position values are allowed. These negative values are treated as position at (global_domain_size[d] + POSITION[d][1]). It is also possible to specify the position explicitly as a positive number. Then it is on a user to make sure the position is correctly calculated wrt the grid size.

Except moving window simulations, the position must be inside the global domain. The distance between the Huygens surface and each global domain boundary must be at least absorber_thickness + (FDTD_spatial_order / 2 - 1). However beware of setting position = direction * (absorber_thickness + const), as then changing absorber parameters will affect laser positioning. When all used profiles are None, the check for POSITION validity is skipped.

For moving window simulations, POSITION for the YMax side can be located outside the initially simulated volume. In this case, parts of the generation surface outside of the currently simulated volume is are treated as if they had zero incident field and it is user’s responsibility to apply a source matching such a case.

class UserFunctorIncidentE

User-defined functor to set values of incident E field.

Public Functions

PMACC_ALIGN(m_unitField, const float3_64)
inline HDINLINE UserFunctorIncidentE(float_X const currentStep, const float3_64 unitField)

Create a functor.

Parameters:
  • currentStep – current time step index, note that it is fractional

  • unitField – conversion factor from SI to internal units, field_internal = field_SI / unitField

inline HDINLINE float3_X operator() (const floatD_X &totalCellIdx) const

Calculate incident field E_inc(r, t) for a source.

Parameters:

totalCellIdx – cell index in the total domain (including all moving window slides), note that it is fractional

Returns:

incident field value in internal units

class UserFunctorIncidentB

User-defined functor to set values of incident B field.

Public Functions

PMACC_ALIGN(m_unitField, const float3_64)
inline HDINLINE UserFunctorIncidentB(float_X const currentStep, const float3_64 unitField)

Create a functor.

Parameters:
  • currentStep – current time step index, note that it is fractional

  • unitField – conversion factor from SI to internal units, field_internal = field_SI / unitField

inline HDINLINE float3_X operator() (const floatD_X &totalCellIdx) const

Calculate incident field B_inc(r, t) for a source.

Parameters:

totalCellIdx – cell index in the total domain (including all moving window slides), note that it is fractional

Returns:

incident field value in internal units

pusher.param

Configure particle pushers.

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

namespace picongpu

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

struct particlePusherAccelerationParam

Subclassed by picongpu::particlePusherAcceleration::UnitlessParam

Public Static Attributes

static constexpr float_64 AMPLITUDEx_SI = 0.0

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

unit: Volt / meter

static constexpr float_64 AMPLITUDEy_SI = -1.e11

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

static constexpr float_64 AMPLITUDEz_SI = 0.0

unit: Volt / meter

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

Acceleration duration unit: second.

namespace particlePusherAxel

Enums

enum TrajectoryInterpolationType

Values:

enumerator LINEAR
enumerator NONLINEAR

Variables

constexpr TrajectoryInterpolationType TrajectoryInterpolation = LINEAR
namespace particlePusherHigueraCary
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)

namespace particlePusherVay
namespace particles
namespace pusher

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.

This profile is normalized to units of picongpu::SI::BASE_DENSITY_SI, also defined in this file. Note that it only operates with physical density, and does not concern macroparticles. The number and sampling of macroparticles per cell are defined independently of a density profile. Please refer to documentation of picongpu::particles::CreateDensity<> for further details on this interaction.

Available profiles:

  • HomogenousImpl : homogeneous density in whole simulation volume

  • GaussianImpl<> : Gaussian profile in ‘y’, optionally with preceeding vacuum

  • GaussianCloudImpl<> : Gaussian profile in all axes, optionally with preceeding vacuum in ‘y’

  • LinearExponentialImpl<> : linear ramping of density in ‘y’ into exponential slope after

  • SphereFlanksImpl<> : composition of 1D profiles, each in form of exponential increasing flank, constant sphere, exponential decreasing flank

  • EveryNthCellImpl<> : checkerboard profile matching the grid, particles are only present in cells with the given stride from one another in all directions

  • FreeFormulaImpl<> : apply user-defined functor for calculating density, refer to picongpu::densityProfiles::IProfile for interface requirements

  • FromOpenPMDImpl<> : load density values from a given file, requires openPMD API dependency

In the end, this file typically defines an alias for each density profile to be used. These aliases do not have to follow any naming convention, but serve as template parameters for invocations of picongpu::particles::CreateDensity<> in speciesInitialization.param.

namespace picongpu

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

namespace densityProfiles

Typedefs

using Gaussian = GaussianImpl<GaussianParam>
using Homogenous = HomogenousImpl
using LinearExponential = LinearExponentialImpl<LinearExponentialParam>
using GaussianCloud = GaussianCloudImpl<GaussianCloudParam>
using SphereFlanks = SphereFlanksImpl<SphereFlanksParam>
using EveryFourthCell = EveryNthCellImpl<mCT::UInt32<4, 4, 4>>
using FreeFormula = FreeFormulaImpl<FreeFormulaFunctor>
struct GaussianParam

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

Public Static Attributes

static constexpr float_X gasFactor = -1.0
static constexpr float_X gasPower = 4.0
static constexpr uint32_t vacuumCellsY = 50

height of vacuum area on top border

this vacuum is important because of the laser initialization, which is done in the first cells of the simulation and assumes a charge-free volume unit: cells

static constexpr float_64 gasCenterLeft_SI = 4.62e-5

The central position of the distribution unit: meter.

static constexpr float_64 gasCenterRight_SI = 4.62e-5
static constexpr float_64 gasSigmaLeft_SI = 4.62e-5

the distance from gasCenter_SI until the gas density decreases to its 1/e-th part unit: meter

static constexpr float_64 gasSigmaRight_SI = 4.62e-5
struct LinearExponentialParam

parameter for LinearExponential profile

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

Public Static Attributes

static constexpr uint32_t vacuumCellsY = 50

height of vacuum area on top border

this vacuum is important because of the laser initialization, which is done in the first cells of the simulation and assumes a charge-free volume unit: cells

static constexpr float_64 gasYMax_SI = 1.0e-3

Y-Position where the linear slope ends and the exponential slope begins unit: meter.

static constexpr float_64 gasA_SI = 1.0e-3

Parameters for the linear slope: For Y <= gasYMax_SI: \rho / BASE_DENSITY = A * Y + B = element [0.0; 1.0] unit for A: 1/m unit for B: none.

static constexpr float_64 gasD_SI = 1.0e-3

Parameters for the exponential slope For Y > gasYMax_SI: let Y’ = Y - gasYMax_SI \rho = exp[ - Y’ * D ] = element [0.0; 1.0] unit: 1/m.

static constexpr float_64 gasB = 0.0
struct GaussianCloudParam

Public Static Attributes

static constexpr float_X gasFactor = -0.5

Profile Formula: exponent = |globalCellPos - center| / sigma density = e^[ gasFactor * exponent^gasPower ].

static constexpr float_X gasPower = 2.0
static constexpr uint32_t vacuumCellsY = 50

height of vacuum area on top border

this vacuum is important because of the laser initialization, which is done in the first cells of the simulation and assumes a charge-free volume unit: cells

static constexpr floatD_64 center_SI = float3_64(1.134e-5, 1.134e-5, 1.134e-5).shrink<simDim>()

The central position of the gas distribution unit: meter.

static constexpr floatD_64 sigma_SI = float3_64(7.0e-6, 7.0e-6, 7.0e-6).shrink<simDim>()

the distance from gasCenter_SI until the gas density decreases to its 1/e-th part unit: meter

struct SphereFlanksParam

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)
*

Public Static Attributes

static constexpr uint32_t vacuumCellsY = 50

height of vacuum area on top border

this vacuum is important because of the laser initialization, which is done in the first cells of the simulation and assumes a charge-free volume unit: cells

static constexpr float_64 r_SI = 1.0e-3

Radius of the constant sphere unit: meter.

static constexpr float_64 ri_SI = 0.0

Inner radius if you want to build a shell/ring unit: meter.

static constexpr floatD_64 center_SI = float3_64(8.0e-3, 8.0e-3, 8.0e-3).shrink<simDim>()

Middle of the constant sphere unit: meter.

static constexpr float_64 exponent_SI = 1.0e3

Parameters for the exponential slope For distance > r_SI: let distance’ = distance - r \rho = exp[ - distance’ * exponent ] unit: 1/m.

struct FromOpenPMDParam

Density values taken from an openPMD file.

The density values must be a scalar dataset of type float_X, type mismatch would cause errors. This implementation would ignore all openPMD metadata but axisLabels. Each value in the dataset defines density in the cell with the corresponding total coordinate minus the given offset. When the functor is instantiated, it will load the part matching the current domain position. Density in points not present in the file would be set to the given default density. Dimensionality of the file indexing must match the simulation dimensionality. Density values are in BASE_DENSITY_SI units.

Public Static Attributes

static constexpr char const *filename = "density.h5"

Path to the openPMD input file.

This value can alternatively be controlled at runtime by setting it to “” here and providing command-line option <species>_runtimeDensityFile. Refer to description of this command-line option for details. Note that runtime option only exists for this parameter, and not others in this struct.

File-based iteration format is also supported, with the usual openPMD API naming scheme.

It is recommended to use a full path to make it independent of how PIConGPU is launched. Relative paths require consistency to the current directory when PIConGPU is started. With tbg and the standard .tpl files, relative to the resulting simOutput directory.

static constexpr char const *datasetName = "e_density"

Name of the openPMD dataset inside the file.

By default, this dataset indexing is assumed to be in (x, y, z) coordinates. This can be changed by setting openPMD attribute “axisLabels” of the correponding dataset. For example, PIConGPU output uses z-y-x via this attribute, and that automatically works. Note that only C dataOrder is supported.

Warning

it is only the dataset itself, a simple text name and not something like “/[directories]/density/[iteration]/fields/e_density”.

static constexpr uint32_t iteration = 0

Iteration inside the file (only file, not related to the current simulation time iteration)

static constexpr DataSpace<simDim> offset = DataSpace<DIM3>(0, 0, 0).shrink<simDim>()

Offset in the file in cells: each value is density at (total cell index - offset)

This offset is in (x, y, z) coordinates. Positive offset means the file values “start” from index == offset in the total coordinates. Negative offset is also supported.

static constexpr float_X defaultDensity = 0.0_X

Default value to be used for cells with no corresponding file value.

struct FreeFormulaFunctor

Public Functions

inline HDINLINE float_X 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.

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

  • cellSize_SI – cell sizes [meter]

Returns:

float_X density [normalized to 1.0]

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

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

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

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

specialization for the relative in-cell position

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

momentum at timestep t

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

momentum at (previous) timestep t-1

value_identifier (float_X, weighting, 0._X)

weighting of the macro particle

value_identifier (float_X, weightingDampingFactor, 1._X)

optional damping factor for weighting of the macro particle

This factor is only used for current deposition with absorbing particle boundaries and PML. In case this attribute is defined, affected particles would have their weightings damped only for that procedure. Otherwise, the damping would affect weighting directly and so potentially overall properties like density. It should be generally not a considerable issue, except lots of particles move from PML to internal area.

Thus, for efficiency reasons it is recommended to not enable this attribute by default.

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

interpolated electric field with respect to particle shape

Attribute can be added to any species.

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

interpolated magnetic field with respect to particle shape

Attribute can be added to any species.

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.

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.

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

    Todo:

    connect default to proton number

alias(atomicLevels)

alias for one of two representations of the atomic state, see also atomicPhysics.param

alias(atomicConfigNumber)

alias for atomic state of ion, see also atomicPhysics.param

alias(atomicPhysicsSolver)
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(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(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;
    using REF_LOCAL_DOM_SIZE = mCT::Int<0, 0, 0>;
    const std::array<float_X, 3> DIR_SCALING_FACTOR = {{0.0, 0.0, 0.0}};
};
alias(boundaryCondition)

alias to specify the internal pmacc boundary treatment for particles

It controls the internal behavior and intented for special cases only. To set physical boundary conditions for a species, instead use <species>_boundary command-line option.

The default behavior if this alias is not given to a species is to do nothing. The existing boundary implementations already take care of the particles leaving the global simulation volume.

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

namespace pmacc

this specifies how an object of the ConfigNumber class can be written to an external file for storage.

Functions

value_identifier (lcellId_t, localCellIdx, 0)

cell of a particle inside a supercell

Value is a linear cell index inside the supercell

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

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

Variables

constexpr float_X 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 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

Particle shape, field to particle interpolation, current solver, and particle pusher can be declared here for usage in speciesDefinition.param.

See also

MODELS / Hierarchy of Charge Assignment Schemes in the online documentation for information on particle shapes.

namespace picongpu

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

Typedefs

using UsedParticleShape = particles::shapes::TSC

select macroparticle shape

WARNING the shape names are redefined and diverge from PIConGPU versions before 0.6.0.

  • particles::shapes::CIC : Assignment function is a piecewise linear spline

  • particles::shapes::TSC : Assignment function is a piecewise quadratic spline

  • particles::shapes::PQS : Assignment function is a piecewise cubic spline

  • particles::shapes::PCS : Assignment function is a piecewise quartic spline

using UsedField2Particle = FieldToParticleInterpolation<UsedParticleShape, AssignedTrilinearInterpolation>

select interpolation method to be used for interpolation of grid-based field values to particle positions

using UsedParticleCurrentSolver = currentSolver::Esirkepov<UsedParticleShape>

select current solver method

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

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

  • currentSolver::EZ< SHAPE, STRATEGY > : same as EmZ (just an alias to match the currently used naming)

STRATEGY (optional):

using UsedParticlePusher = particles::pusher::Boris

particle pusher configuration

Defining a pusher is optional for particles

  • particles::pusher::HigueraCary : Higuera & Cary’s relativistic pusher preserving both volume and ExB velocity

  • particles::pusher::Vay : Vay’s relativistic pusher preserving ExB velocity

  • particles::pusher::Boris : Boris’ relativistic pusher preserving volume

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

  • particles::pusher::Composite : composite of two given pushers, switches between using one (or none) of those

For diagnostics & modeling: ———————————————&#8212;

  • 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: ———————————————–&#8212;

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

Current solver details.

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

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

Typedefs

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

describe attributes of a particle

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

value_identifier (float_X, MassRatioElectrons, 1.0)
value_identifier (float_X, ChargeRatioElectrons, 1.0)
value_identifier (float_X, MassRatioIons, 1836.152672)
value_identifier (float_X, ChargeRatioIons, -1.0)
value_identifier (float_X, DensityRatioIons, 1.0)

particle.param

Configurations for particle manipulators.

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

namespace picongpu

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

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

(Approximate) Number of maximum macro-particles per cell.

Used internally for unit normalization. And used in startPosition functors further below to set real maximum number of macro-particles per cell.

namespace manipulators

Typedefs

using AssignXDrift = unary::Drift<DriftParam, pmacc::math::operation::Assign>

Definition of manipulator that assigns a drift in X using parameters from struct DriftParam.

using AddTemperature = unary::Temperature<TemperatureParam>

Definition of manipulator assigning a temperature using parameters from struct TemperatureParam.

using AddFreeTemperature = unary::FreeTemperature<TemperatureFunctor>

Definition of manipulator that assigns a temperature according to the functor defined by struct TemperaturFunctor.

using DoubleWeighting = generic::Free<DoubleWeightingFunctor>

Definition of the free particle manipulator which doubles each weighting.

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

Definition of manipulator that selects macro-particles for Radiation plugin.

using RandomPosition = unary::RandomPosition

Definition of manipulator that changes the in-cell position of each particle of a species.

using SetNeutral = unary::ChargeState<0u>

definition of manipulator that sets the boundElectron attribute(charge state) of each particle of an ion of a species to neutral

struct DriftParam

Define Lorentz factor of initial particle drift.

Public Static Attributes

static constexpr float_64 gamma = 1.0
static constexpr auto driftDirection = float3_X(1.0, 0.0, 0.0)

Define initial particle drift direction vector.

struct TemperatureParam

Define initial particle temperature.

Public Static Attributes

static constexpr float_64 temperature = 0.0

Initial temperature unit: keV.

struct TemperatureFunctor

Define initial particle temperature as a function of position.

This is a functor which needs to follow the requirements of param::TemperatureFunctor.

Public Functions

inline TemperatureFunctor()

Constructor, can take currentStep or no parameters (can also be auto-generated by a compiler)

Parameters:

currentStep – current time iteration

inline HDINLINE float_X operator() (const DataSpace< simDim > &totalCellOffset)

Return the temperature in keV for the given position.

Return type may be float_X or float_64.

Parameters:

totalCellOffset – total offset including all slides [in cells]

struct DoubleWeightingFunctor

Unary particle manipulator: double each weighting.

Public Functions

template<typename T_Particle> inline DINLINE void operator() (T_Particle &particle)
struct RandomEnabledRadiationFunctor

Define mask which randomly marks macro-particles used by the radiation plugin to calculate far field radiation.

Public Functions

template<typename T_Rng, typename T_Particle> inline DINLINE void operator() (T_Rng &rng, T_Particle &particle)
namespace startPosition

Typedefs

using Random = RandomImpl<RandomParameter>

Definition of start position functor that randomly distributes macro-particles within a cell.

using RandomPositionAndWeighting = RandomPositionAndWeightingImpl<RandomParameter>

Definition of start position functor that randomly distributes macro-particles within a cell and randomly distributes weightings across macro-particles within a cell.

using Quiet = QuietImpl<QuietParam>

Definition of Quiet start position functor that positions macro-particles regularly on the grid.

No random number generator used.

using OnePosition = OnePositionImpl<OnePositionParameter>

Definition of OnePosition start position functor that places macro-particles at the initial in-cell position defined above.

Functions

CONST_VECTOR (float_X, 3, InCellOffset, 0.0, 0.0, 0.0)

Define initial in-cell particle position used as parameter in OnePosition functor.

Here, macro-particles sit directly in lower corner of the cell.

This defines the type InCellOffset_t where every instance of this type has the preset values defined here.

struct RandomParameter

Define target number for marco-particles per cell to be used in Random start position functor.

Public Static Attributes

static constexpr uint32_t numParticlesPerCell = TYPICAL_PARTICLES_PER_CELL

Maximum number of macro-particles per cell during density profile evaluation.

Determines the weighting of a macro particle as well as the number of macro-particles which sample the evolution of the particle distribution function in phase space.

unit: none

struct QuietParam

Define target number for marco-particles per cell along a direction.

To be used in Quiet start position functor.

Here, one macro-particle per cell along x, one macro-particle per cell along z, and TYPICAL_PARTICLES_PER_CELL macro-particles per cell along y.

Vector is automatically reduced to two dimensions for 2D (x,y) simulations.

Public Types

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

Count of macro-particles per cell per direction at initial state.

unit: none

struct OnePositionParameter

Public Members

const InCellOffset_t inCellOffset

Public Static Attributes

static constexpr uint32_t numParticlesPerCell = TYPICAL_PARTICLES_PER_CELL

Maximum number of macro-particles per cell during density profile evaluation.

Determines the weighting of a macro particle as well as the number of macro-particles which sample the evolution of the particle distribution function in phase space.

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

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

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 * particles::TYPICAL_NUM_PARTICLES_PER_MACROPARTICLE

Unit of mass.

constexpr float_64 UNIT_CHARGE = -1.0 * SI::BASE_CHARGE_SI * 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_64 TYPICAL_NUM_PARTICLES_PER_MACROPARTICLE = (SI::BASE_DENSITY_SI * SI::CELL_WIDTH_SI * SI::CELL_HEIGHT_SI * SI::CELL_DEPTH_SI) / float_64(particles::TYPICAL_PARTICLES_PER_CELL)

Typical number of particles per macro particle (= typical 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

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

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

namespace traits

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

rate calculation from given atomic data, extracted from flylite, based on FLYCHK

References:

  • Axel Huebl flylite, not yet published

    • R. Mewe. “Interpolation formulae for the electron impact excitation of ions in

      the H-, He-, Li-, and Ne-sequences.” Astronomy and Astrophysics 20, 215 (1972)

    • H.-K. Chung, R.W. Lee, M.H. Chen. “A fast method to generate collisional excitation cross-sections of

      highly charged ions in a hot dense matter” High Energy Dennsity Physics 3, 342-352 (2007)

Note

this file uses the same naming convention for updated and incident field as Solver.kernel.

Note

In this file we use camelCase “updatedField” in both code and comments to denote field E or B that is being updated (i.e. corrected) in the kernel. The other of the two fields is called “incidentField”. And for the incidentField source we explicitly use “functor” to not confuse it with the field itself. Please refer to https://picongpu.readthedocs.io/en/latest/models/total_field_scattered_field.html for theoretical background of this procedure.

namespace particles

Typedefs

using InitPipeline = pmacc::mp_list<>

InitPipeline defines in which order species are initialized.

the functors are called in order (from first to last functor). The functors must be default-constructible and take the current time iteration as the only parameter.

List of all initialization methods for particle species.