PIC Core

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 picongpuABSORBER_CELLS[3][2] = { {32, 32}, {32, 32}, {32, 32} }

Defines the size of the absorbing zone (in cells)

unit: none

constexpr float_X picongpuABSORBER_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 picongpumovePoint = 0.9

When to start moving the co-moving window.

Slide point model: A virtual photon starts at t=0 at the lower end of the global simulation box in y-direction of the simulation. When it reaches movePoint % of the global simulation box, the co-moving window starts to move with the speed of light.

Note
global simulation area: there is one additional “hidden” row of gpus at the y-front, when you use the co-moving window. 1.0 would correspond to: start moving exactly when the above described “virtual photon” from the lower end of the box’ Y-axis reaches the beginning of this “hidden” row of GPUs.

namespace picongpuSI

Variables

constexpr float_64 picongpu::SIDELTA_T_SI = 0.8e-16

Duration of one timestep unit: seconds.

constexpr float_64 picongpu::SICELL_WIDTH_SI = 0.1772e-6

equals X unit: meter

constexpr float_64 picongpu::SICELL_HEIGHT_SI = 0.4430e-7

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

constexpr float_64 picongpu::SICELL_DEPTH_SI = CELL_WIDTH_SI

equals Z unit: meter

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 picongpusimDim = SIMDIM

components.param

Select the laser profile and the field solver here.

Defines

ENABLE_CURRENT

enable (1) or disable (0) current calculation (deprecated)

namespace picongpu

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

namespace laserProfile

Laser Profile Selection:

  • laserNone : no laser init
  • laserGaussianBeam : Gaussian beam (focusing)
  • laserPulseFrontTilt : Gaussian beam with a tilted pulse envelope in ‘x’ direction
  • laserWavepacket : wavepacket (Gaussian in time and space, not focusing)
  • laserPlaneWave : a plane wave (Gaussian in time)
  • laserPolynom : a polynomial laser envelope

Adjust the settings of the selected profile in laser.param

namespace fieldSolver

Field Solver Selection:

  • fieldSolverYee : standard Yee solver
  • fieldSolverLehe: Num. Cherenkov free field solver in a chosen direction
  • fieldSolverDirSplitting: Sentoku’s Directional Splitting Method
  • fieldSolverNone: disable the vacuum update of E and B

For development purposes:

  • fieldSolverYeeNative : generic version of fieldSolverYee (need more shared memory per GPU and is slow)

Adjust the settings of the selected field solver in fieldSolver.param

fieldSolver.param

Configure the selected field solver method.

You can set/modify Maxwell solver specific options in the section of each “FieldSolver”.

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< simDim >:
    • default for staggered grids/Yee-scheme
    • updates E
  • Binomial< simDim >: 2nd order Binomial filter
    • smooths the current before assignment in staggered grid
    • updates E & breaks local charge conservation slightly
  • NoneDS< simDim >:
    • experimental assignment for all-centered/directional splitting
    • updates E & B at the same time

namespace picongpu

namespace picongpufieldSolverDirSplitting

Typedefs

using picongpu::fieldSolverDirSplitting::CurrentInterpolation = typedef currentInterpolation::NoneDS< simDim >
namespace picongpufieldSolverLehe

Lehe Solver The solver proposed by R.

Lehe et al in Phys. Rev. ST Accel. Beams 16, 021301 (2013)

Typedefs

using picongpu::fieldSolverLehe::CherenkovFreeDir = typedef CherenkovFreeDirection_Y

Distinguish the direction where numerical Cherenkov Radiation by moving particles shall be suppressed.

using picongpu::fieldSolverLehe::CurrentInterpolation = typedef currentInterpolation::None< simDim >
namespace picongpufieldSolverNone

Typedefs

using picongpu::fieldSolverNone::CurrentInterpolation = typedef currentInterpolation::None< simDim >
namespace picongpufieldSolverYee

Typedefs

using picongpu::fieldSolverYee::CurrentInterpolation = typedef currentInterpolation::None< simDim >
namespace picongpufieldSolverYeeNative

Typedefs

using picongpu::fieldSolverYeeNative::CurrentInterpolation = typedef currentInterpolation::None< simDim >

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 picongpudensityProfiles

Typedefs

using picongpu::densityProfiles::Gaussian = typedef GaussianImpl< GaussianParam >
using picongpu::densityProfiles::Homogenous = typedef HomogenousImpl
using picongpu::densityProfiles::LinearExponential = typedef LinearExponentialImpl< LinearExponentialParam >
using picongpu::densityProfiles::GaussianCloud = typedef GaussianCloudImpl< GaussianCloudParam >
using picongpu::densityProfiles::SphereFlanks = typedef SphereFlanksImpl<SphereFlanksParam>
using picongpu::densityProfiles::FromHDF5 = typedef FromHDF5Impl< FromHDF5Param >
using picongpu::densityProfiles::FreeFormula = typedef 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 = float_X(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 picongpu::densityProfilesFreeFormulaFunctor

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 picongpuSI

Variables

constexpr float_64 picongpu::SIBASE_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

pusher.param

Configure particle pushers.

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

namespace picongpu

namespace picongpuparticlePusherAxel

Enums

enum picongpu::particlePusherAxelTrajectoryInterpolationType

Values:

picongpu::particlePusherAxelLINEAR = 1u
picongpu::particlePusherAxelNONLINEAR = 2u

Variables

constexpr TrajectoryInterpolationType picongpu::particlePusherAxelTrajectoryInterpolation = LINEAR

laser.param

Configure laser profiles.

namespace picongpu

namespace picongpulaser

Variables

constexpr uint32_t picongpu::laserinitPlaneY = 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

namespace picongpulaserGaussianBeam

Enums

enum picongpu::laserGaussianBeamPolarisationType

Available polarisation types.

Values:

picongpu::laserGaussianBeamLINEAR_X = 1u
picongpu::laserGaussianBeamLINEAR_Z = 2u
picongpu::laserGaussianBeamCIRCULAR = 4u

Functions

picongpu::laserGaussianBeam::PMACC_CONST_VECTOR(float_X, MODENUMBER+ 1, LAGUERREMODES, 1. 0)

Variables

constexpr float_64 picongpu::laserGaussianBeamPULSE_INIT = 20.0

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

unit: none

constexpr float_X picongpu::laserGaussianBeamLASER_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 picongpu::laserGaussianBeamMODENUMBER = 0

Use only the 0th Laguerremode for a standard Gaussian.

constexpr PolarisationType picongpu::laserGaussianBeamPolarisation = CIRCULAR

Polarization selection.

namespace picongpu::laserGaussianBeamSI

Variables

constexpr float_64 picongpu::laserGaussianBeam::SIWAVE_LENGTH_SI = 0.8e-6

unit: meter

constexpr float_64 picongpu::laserGaussianBeam::SIUNITCONV_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 picongpu::laserGaussianBeam::SIAMPLITUDE_SI = 1.738e13

unit: W / m^2

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

constexpr float_64 picongpu::laserGaussianBeam::SIPULSE_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 picongpu::laserGaussianBeam::SIW0_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 picongpu::laserGaussianBeam::SIFOCUS_POS_SI = 4.62e-5

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

namespace picongpulaserNone

No Laser initialization.

namespace picongpu::laserNoneSI

Variables

constexpr float_64 picongpu::laserNone::SIWAVE_LENGTH_SI = 0.0

unit: meter

constexpr float_64 picongpu::laserNone::SIAMPLITUDE_SI = 0.0

unit: Volt /meter

constexpr float_64 picongpu::laserNone::SIPULSE_LENGTH_SI = 0.0
namespace picongpulaserPlaneWave

plane wave (use periodic boundaries!)

no transverse spacial envelope based on the electric potential Phi = Phi_0 * exp(0.5 * (x-x_0)^2 / sigma^2) * cos(k*(x - x_0) - phi) by applying -grad Phi = -d/dx Phi = E(x) we get: E = -Phi_0 * exp(0.5 * (x-x_0)^2 / sigma^2) * [k*sin(k*(x - x_0) - phi) + x/sigma^2 * cos(k*(x - x_0) - phi)]

This approach ensures that int_{-infinity}^{+infinity} E(x) = 0 for any phase if we have no transverse profile as we have with this plane wave train

Since PIConGPU requires a temporally defined electric field, we use: t = x/c and (x-x_0)/sigma = (t-t_0)/tau and k*(x-x_0) = omega*(t-t_0) with omega/k = c and tau * c = sigma and get: E = -Phi_0*omega/c * exp(0.5 * (t-t_0)^2 / tau^2) * [sin(omega*(t - t_0) - phi) + t/(omega*tau^2) * cos(omega*(t - t_0) - phi)] and define: E_0 = -Phi_0*omega/c integrationCorrectionFactor = t/(omega*tau^2)

Please consider: 1) The above formulae does only apply to a Gaussian envelope. If the plateau length is not zero, the integral over the volume will only vanish if the plateau length is a multiple of the wavelength. 2) Since we define our envelope by a sigma of the laser intensity, tau = PULSE_LENGTH * sqrt(2)

Enums

enum picongpu::laserPlaneWavePolarisationType

Available polarisation types.

Values:

picongpu::laserPlaneWaveLINEAR_X = 1u
picongpu::laserPlaneWaveLINEAR_Z = 2u
picongpu::laserPlaneWaveCIRCULAR = 4u

Variables

constexpr float_64 picongpu::laserPlaneWaveRAMP_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 picongpu::laserPlaneWaveLASER_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 picongpu::laserPlaneWavePolarisation = LINEAR_X

Polarization selection.

namespace picongpu::laserPlaneWaveSI

Variables

constexpr float_64 picongpu::laserPlaneWave::SIWAVE_LENGTH_SI = 0.8e-6

unit: meter

constexpr float_64 picongpu::laserPlaneWave::SIUNITCONV_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 picongpu::laserPlaneWave::SI_A0 = 1.5

unit: W / m^2

unit: none

constexpr float_64 picongpu::laserPlaneWave::SIAMPLITUDE_SI = _A0 * UNITCONV_A0_to_Amplitude_SI

unit: Volt /meter

constexpr float_64 picongpu::laserPlaneWave::SILASER_NOFOCUS_CONSTANT_SI = 13.34e-15

unit: Volt /meter

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

constexpr float_64 picongpu::laserPlaneWave::SIPULSE_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)

namespace picongpulaserPolynom

not focusing polynomial laser pulse

no phase shifts, just spacial envelope

Variables

constexpr float_X picongpu::laserPolynomLASER_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

namespace picongpu::laserPolynomSI

Variables

constexpr float_64 picongpu::laserPolynom::SIWAVE_LENGTH_SI = 0.8e-6

unit: meter

constexpr float_64 picongpu::laserPolynom::SIUNITCONV_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 picongpu::laserPolynom::SIAMPLITUDE_SI = 1.738e13

unit: W / m^2

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

constexpr float_64 picongpu::laserPolynom::SIPULSE_LENGTH_SI = 4.0e-15

Pulse length: PULSE_LENGTH_SI = total length of polynamial laser pulse Rise time = 0.5 * PULSE_LENGTH_SI Fall time = 0.5 * PULSE_LENGTH_SI in order to compare to a gaussian pulse: rise time = sqrt{2} * T_{FWHM} unit: seconds.

constexpr float_64 picongpu::laserPolynom::SIW0x_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 picongpu::laserPolynom::SIW0z_SI = W0x_SI
namespace picongpulaserPulseFrontTilt

Enums

enum picongpu::laserPulseFrontTiltPolarisationType

Available polarisation types.

Values:

picongpu::laserPulseFrontTiltLINEAR_X = 1u
picongpu::laserPulseFrontTiltLINEAR_Z = 2u
picongpu::laserPulseFrontTiltCIRCULAR = 4u

Variables

constexpr float_64 picongpu::laserPulseFrontTiltPULSE_INIT = 20.0

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

constexpr float_X picongpu::laserPulseFrontTiltLASER_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 picongpu::laserPulseFrontTiltPolarisation = LINEAR_X

Polarization selection.

namespace picongpu::laserPulseFrontTiltSI

Variables

constexpr float_64 picongpu::laserPulseFrontTilt::SIWAVE_LENGTH_SI = 0.8e-6

unit: meter

constexpr float_64 picongpu::laserPulseFrontTilt::SIUNITCONV_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 picongpu::laserPulseFrontTilt::SIAMPLITUDE_SI = 1.738e13

unit: Volt / meter

constexpr float_64 picongpu::laserPulseFrontTilt::SIPULSE_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 picongpu::laserPulseFrontTilt::SIW0_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 picongpu::laserPulseFrontTilt::SIFOCUS_POS_SI = 4.62e-5

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

constexpr float_64 picongpu::laserPulseFrontTilt::SITILT_X_SI = 0

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

namespace picongpulaserWavepacket

not focusing wavepaket with spacial gaussian envelope

no phase shifts, just spacial envelope including correction to laser formular derived from vector potential, so the integration along propagation direction gives 0 this is important for few-cycle laser pulses

Enums

enum picongpu::laserWavepacketPolarisationType

Available polarisation types.

Values:

picongpu::laserWavepacketLINEAR_X = 1u
picongpu::laserWavepacketLINEAR_Z = 2u
picongpu::laserWavepacketCIRCULAR = 4u

Variables

constexpr float_64 picongpu::laserWavepacketRAMP_INIT = 20.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 float_X picongpu::laserWavepacketLASER_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 picongpu::laserWavepacketPolarisation = LINEAR_X

Polarization selection.

namespace picongpu::laserWavepacketSI

Variables

constexpr float_64 picongpu::laserWavepacket::SIWAVE_LENGTH_SI = 0.8e-6

unit: meter

constexpr float_64 picongpu::laserWavepacket::SIUNITCONV_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 picongpu::laserWavepacket::SIAMPLITUDE_SI = 1.738e13

unit: W / m^2

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

constexpr float_64 picongpu::laserWavepacket::SILASER_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 constexprant area between the up and downramp unit: seconds.

constexpr float_64 picongpu::laserWavepacket::SIPULSE_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 picongpu::laserWavepacket::SIW0_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, 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 picongpu::laserWavepacket::SIW0_Z_SI = W0_X_SI

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 picongpuparticles

Variables

constexpr float_X picongpu::particlesMIN_WEIGHTING = 10.0

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

unit: none

constexpr uint32_t picongpu::particlesTYPICAL_PARTICLES_PER_CELL = 2

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 picongpu::particlesmanipulators

Typedefs

using picongpu::particles::manipulators::AssignXDrift = typedef DriftImpl< DriftParam, nvidia::functors::Assign >

definition of manipulator that assigns a dirft in X

using picongpu::particles::manipulators::AddTemperature = typedef TemperatureImpl< TemperatureParam, nvidia::functors::Add >
using picongpu::particles::manipulators::AssignXDriftToLowerHalfXPosition = typedef IfRelativeGlobalPositionImpl< IfRelativeGlobalPositionParam, AssignXDrift >

definition of a relative position selection that assigns a drift in X

using picongpu::particles::manipulators::DoubleWeighting = typedef FreeImpl< DoubleWeightingFunctor >

definition of a free particle manipulator: double weighting

typedef FreeRngImpl<RandomEnabledRadiationFunctor, nvidia::rng::distributions::Uniform_float> picongpu::particles::manipulatorsRandomEnabledRadiation
using picongpu::particles::manipulators::RandomPosition = typedef RandomPositionImpl<>

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 picongpu::particles::manipulatorsDoubleWeightingFunctor

Unary particle manipulator: double each weighting.

Public Functions

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

Parameter for a particle drift assignment.

Public Members

const DriftParam_direction_t picongpu::particles::manipulators::DriftParamdirection

Public Static Attributes

constexpr float_64 picongpu::particles::manipulators::DriftParamgamma = 1.0
struct picongpu::particles::manipulatorsIfRelativeGlobalPositionParam

Parameters for an assignment in a relative position selection.

Public Static Attributes

constexpr float_X picongpu::particles::manipulators::IfRelativeGlobalPositionParamlowerBound = 0.0
constexpr float_X picongpu::particles::manipulators::IfRelativeGlobalPositionParamupperBound = 0.5
constexpr uint32_t picongpu::particles::manipulators::IfRelativeGlobalPositionParamdimension = 0
struct picongpu::particles::manipulatorsRandomEnabledRadiationFunctor

Public Functions

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

Parameter for a temperature assignment.

Public Static Attributes

constexpr float_64 picongpu::particles::manipulators::TemperatureParamtemperature = 0.0
namespace picongpu::particlesstartPosition

Typedefs

using picongpu::particles::startPosition::Random = typedef RandomImpl< RandomParameter >

definition of random particle start

using picongpu::particles::startPosition::Quiet = typedef QuietImpl< QuietParam >

definition of quiet particle start

using picongpu::particles::startPosition::OnePosition = typedef 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 picongpu::particles::startPositionOnePositionParameter

Public Members

const InCellOffset_t picongpu::particles::startPosition::OnePositionParameterinCellOffset

Public Static Attributes

constexpr uint32_t picongpu::particles::startPosition::OnePositionParameternumParticlesPerCell = TYPICAL_PARTICLES_PER_CELL

Count of particles per cell at initial state.

unit: none

struct picongpu::particles::startPositionQuietParam

Public Types

using picongpu::particles::startPosition::QuietParamnumParticlesPerDimension = 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 picongpu::particles::startPositionRandomParameter

Public Static Attributes

constexpr uint32_t picongpu::particles::startPosition::RandomParameternumParticlesPerCell = TYPICAL_PARTICLES_PER_CELL

Count of particles per cell at initial state.

unit: none

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 picongpu::UsedParticleShape = typedef 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 CICShape = particles::shapes::CIC;

using picongpu::UsedField2Particle = typedef FieldToParticleInterpolation< UsedParticleShape, AssignedTrilinearInterpolation >

define which interpolation method is used to interpolate fields to particles

using picongpu::UsedParticleCurrentSolver = typedef 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)
  • currentSolver::ZigZag< SHAPE > : particle shapes - CIC, TSC, PCS, P4S (1st to 4th order)

using picongpu::UsedParticlePusher = typedef particles::pusher::Boris

particle pusher configuration

Define 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 development purposes:

  • particles::pusher::Axel : a pusher developed at HZDR during 2011 (testing)
  • particles::pusher::Free : free propagation, ignore fields (= free stream model)
  • particles::pusher::Photon : propagate with c in direction of normalized mom.

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

picongpualias(position)

relative (to cell origin) in-cell position of a particle With this definition we not define any type like float3,double3,…

This is only a name without a specialization

picongpuvalue_identifier(uint64_t, particleId, IdProvider<simDim>::getNewId())
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.))

momentum at (previous) timestep t-1

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

weighting of the macro particle

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(float_X, boundElectrons, float_X(0.0))

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

picongpuvalue_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.

picongpualias(shape)

alias for particle shape

See
species.param

picongpualias(particlePusher)

alias for particle pusher

See
species.param

picongpualias(ionizers)

alias for particle ionizers

See
ionizer.param

picongpualias(ionizationEnergies)

alias for ionization energy container

See
ionizationEnergies.param

picongpualias(synchrotronPhotons)

alias for synchrotronPhotons alias for ion species used for bremsstrahlung

See
speciesDefinition.param

picongpualias(bremsstrahlungPhotons)

alias for photon species used for bremsstrahlung

picongpualias(interpolation)

alias for particle to field interpolation

See
species.param

picongpualias(current)

alias for particle current solver

See
species.param

picongpualias(atomicNumbers)

alias for particle flag: atomic numbers

See
ionizer.param
  • only reasonable for atoms / ions / nuclei

picongpualias(effectiveNuclearCharge)

alias for particle flag: effective nuclear charge

See
ionizer.param
  • only reasonable for atoms / ions / nuclei

picongpualias(massRatio)

alias for particle mass ratio

mass ratio between base particle default value: 1.0 if unset

See
speciesConstants.param SI::BASE_MASS_SI and a user defined species

picongpualias(chargeRatio)

alias for particle charge ratio

charge ratio between base particle default value: 1.0 if unset

See
speciesConstants.param SI::BASE_CHARGE_SI and a user defined species

picongpualias(densityRatio)

alias for particle density ratio

density ratio between default density default value: 1.0 if unset

See
density.param SI::BASE_DENSITY_SI and a user defined species

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 picongpuGAMMA_THRESH = float_X(1.005)

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 picongpuGAMMA_INV_SQUARE_RAD_THRESH = float_X(0.18)

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 picongpuSI

Variables

constexpr float_64 picongpu::SIBASE_MASS_SI = ELECTRON_MASS_SI

base particle mass

reference for massRatio in speciesDefinition.param

unit: kg

constexpr float_64 picongpu::SIBASE_CHARGE_SI = ELECTRON_CHARGE_SI

base particle charge

reference for chargeRatio in speciesDefinition.param

unit: C

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 picongpu::DefaultParticleAttributes = typedef MakeSeq_t< position< position_pic >, momentum, weighting >

describe attributes of a particle

using picongpu::ParticleFlagsPhotons = typedef bmpl::vector< particlePusher< particles::pusher::Photon >, shape< UsedParticleShape >, interpolation< UsedField2Particle >, massRatio< MassRatioPhotons >, chargeRatio< ChargeRatioPhotons > >
using picongpu::PIC_Photons = typedef Particles< bmpl::string< 'p', 'h' >, ParticleFlagsPhotons, DefaultParticleAttributes >
using picongpu::ParticleFlagsElectrons = typedef bmpl::vector< particlePusher< UsedParticlePusher >, shape< UsedParticleShape >, interpolation< UsedField2Particle >, current< UsedParticleCurrentSolver >, massRatio< MassRatioElectrons >, chargeRatio< ChargeRatioElectrons > >
using picongpu::PIC_Electrons = typedef Particles< bmpl::string< 'e' >, ParticleFlagsElectrons, DefaultParticleAttributes >
using picongpu::ParticleFlagsIons = typedef bmpl::vector< particlePusher< UsedParticlePusher >, shape< UsedParticleShape >, interpolation< UsedField2Particle >, current< UsedParticleCurrentSolver >, massRatio< MassRatioIons >, chargeRatio< ChargeRatioIons >, densityRatio< DensityRatioIons >, atomicNumbers< ionization::atomicNumbers::Hydrogen_t > >
using picongpu::PIC_Ions = typedef Particles< bmpl::string< 'i' >, ParticleFlagsIons, DefaultParticleAttributes >
using picongpu::VectorAllSpecies = typedef 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)

speciesInitialization.param

Available species functors in src/picongpu/include/particles/InitFunctors.hpp.

  • CreateDensity<T_DensityFunctor, T_PositionFunctor, T_SpeciesType> Create particle distribution based on a density profile and an in-cell positioning. Fills a particle species (fillAllGaps() is called).
    See
    density.param
    Template Parameters
    • T_DensityFunctor: unary lambda functor with density description,

namespace picongpu

namespace picongpuparticles

Typedefs

using picongpu::particles::InitPipeline = typedef mpl::vector<>

InitPipeline defines in which order species are initialized.

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