Plugins

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

Typedefs

using ChargeDensity_Seq = deriveField::CreateEligible_t<VectorAllSpecies, deriveField::derivedAttributes::ChargeDensity>

FieldTmp output (calculated at runtime) *******************************.

Those operations derive scalar field quantities from particle species at runtime. Each value is mapped per cell. Some operations are identical up to a constant, so avoid writing those twice to save storage.

you can choose any of these particle to grid projections:

  • Density: particle position + shape on the grid

  • BoundElectronDensity: density of bound electrons note: only makes sense for partially ionized ions

  • ChargeDensity: density * charge note: for species that do not change their charge state, this is the same as the density times a constant for the charge

  • Energy: sum of kinetic particle energy per cell with respect to shape

  • EnergyDensity: average kinetic particle energy per cell times the particle density note: this is the same as the sum of kinetic particle energy divided by a constant for the cell volume

  • Momentum: sum of chosen component of momentum per cell with respect to shape

  • MomentumDensity: average chosen component of momentum per cell times the particle density note: this is the same as the sum of the chosen momentum component divided by a constant for the cell volume

  • LarmorPower: radiated Larmor power (species must contain the attribute momentumPrev1)

for debugging:

  • MidCurrentDensityComponent: density * charge * velocity_component

  • Counter: counts point like particles per cell

  • MacroCounter: counts point like macro particles per cell

combined attributes: These attributes are defined as a function of two primary derived attributes.

  • AverageAttribute: template to compute a per particle average of any derived value

  • RelativisticDensity: equals to 1/gamma^2 * n. Where gamma is the Lorentz factor for the mean kinetic energy in the cell and n ist the usual number density. Useful for Faraday Rotation calculation.

    Example use:

    using AverageEnergy_Seq = deriveField::CreateEligible_t<
        VectorAllSpecies,
        deriveField::combinedAttributes::AverageAttribute<deriveField::derivedAttributes::Energy>>;
    using RelativisticDensity_Seq
        = deriveField::CreateEligible_t<PIC_Electrons, deriveField::combinedAttributes::RelativisticDensity>;
    

Filtering: You can create derived fields from filtered particles. Only particles passing the filter will contribute to the field quantity. For that you need to define your filters in particleFilters.param and pass a filter as the 3rd (optional) template argument in CreateEligible_t.

Example: This will create charge density field for all species that are eligible for the this attribute and the chosen filter.

using ChargeDensity_Seq
= deriveField::CreateEligible_t< VectorAllSpecies,
deriveField::derivedAttributes::ChargeDensity, filter::FilterOfYourChoice>;

using EnergyDensity_Seq = deriveField::CreateEligible_t<VectorAllSpecies, deriveField::derivedAttributes::EnergyDensity>
using FieldTmpSolvers = MakeSeq_t<ChargeDensity_Seq, EnergyDensity_Seq>

FieldTmpSolvers groups all solvers that create data for FieldTmp ******.

FieldTmpSolvers is used in

See also

FieldTmp to calculate the exchange size

using NativeFileOutputFields = MakeSeq_t<FieldE, FieldB>

FileOutputFields: Groups all Fields that shall be dumped.

Possible native fields: FieldE, FieldB, FieldJ

using FileOutputFields = MakeSeq_t<NativeFileOutputFields, FieldTmpSolvers>
using FileOutputParticles = VectorAllSpecies

FileOutputParticles: Groups all Species that shall be dumped **********.

hint: to disable particle output set to using FileOutputParticles = MakeSeq_t< >;

isaac.param

Definition which native fields and density fields of (filtered) particles will be visualizable with ISAAC.

ISAAC is an in-situ visualization library with which the PIC simulation can be observed while it is running avoiding the time consuming writing and reading of simulation data for the classical post processing of data.

ISAAC can directly visualize natives fields like the E or B field, but density fields of particles need to be calculated from PIConGPU on the fly which slightly increases the runtime and the memory consumption. Every particle density field will reduce the amount of memory left for PIConGPUs particles and fields.

To get best performance, ISAAC defines an exponential amount of different visualization kernels for every combination of (at runtime) activated fields. So furthermore a lot of fields will increase the compilation time.

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 isaacP

Typedefs

using Particle_Seq = VectorAllSpecies

Intermediate list of native particle species of PIConGPU which shall be visualized.

using Native_Seq = MakeSeq_t<FieldE, FieldB, FieldJ>

Intermediate list of native fields of PIConGPU which shall be visualized.

using Density_Seq = deriveField::CreateEligible_t<Particle_Seq, deriveField::derivedAttributes::Density>

Intermediate list of particle species, from which density fields shall be created at runtime to visualize them.

You can create such densities from filtered particles by passing a particle filter as the third template argument (filter::All by default). Don’t forget to add your filtered densities to the Fields_Seq below.

using Density_Seq_Filtered = deriveField::CreateEligible_t<Particle_Seq,
    deriveField::derivedAttributes::Density, filter::All>;

using Fields_Seq = MakeSeq_t<Native_Seq, Density_Seq>

Compile time sequence of all fields which shall be visualized.

Basically joining Native_Seq and Density_Seq.

using VectorFields_Seq = Native_Seq

Compile time sequence of all fields which shall be visualized.

Basically joining Native_Seq and Density_Seq.

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

Functions

HDINLINE float2_X mapYawPitchToNormedRange (const float_X yaw, const float_X pitch, const float_X maxYaw, const float_X maxPitch)

Map yaw and pitch into [0,1] respectively.

These ranges correspond to the normalized histogram range of the calorimeter (0: first bin, 1: last bin). Out-of-range values are mapped to the first or the last bin.

Useful for fine tuning the spatial calorimeter resolution.

Parameters:
  • yaw – -maxYaw…maxYaw

  • pitch – -maxPitch…maxPitch

  • maxYaw – maximum value of angle yaw

  • maxPitch – maximum value of angle pitch

Returns:

Two values within [-1,1]

radiation.param

Definition of frequency space, number of observers, filters, form factors and window functions of the radiation plugin.

All values set here determine what the radiation plugin will compute. The observation direction is defined in a seperate file radiationObserver.param. On the comand line the plugin still needs to be called for each species the radiation should be computed for.

Defines

PIC_VERBOSE_RADIATION

radiation verbose level: 0=nothing, 1=physics, 2=simulation_state, 4=memory, 8=critical

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 plugins
namespace radiation

Typedefs

using RadiationParticleFilter = picongpu::particles::manipulators::generic::Free<GammaFilterFunctor>

filter to (de)select particles for the radiation calculation

to activate the filter:

  • goto file speciesDefinition.param

  • add the attribute radiationMask to the particle species

struct GammaFilterFunctor

select particles for radiation example of a filter for the relativistic Lorentz factor gamma

Public Functions

template<typename T_Particle> inline HDINLINE void operator() (T_Particle &particle)

Public Static Attributes

static constexpr float_X radiationGamma = 5.0

Gamma value above which the radiation is calculated.

namespace frequencies_from_list

Variables

constexpr const char *listLocation = "/path/to/frequency_list"

path to text file with frequencies

constexpr unsigned int N_omega = 2048

number of frequency values to compute if frequencies are given in a file [unitless]

namespace linear_frequencies

Variables

constexpr unsigned int N_omega = 2048

number of frequency values to compute in the linear frequency [unitless]

namespace SI

Variables

constexpr float_64 omega_min = 0.0

mimimum frequency of the linear frequency scale in units of [1/s]

constexpr float_64 omega_max = 1.06e16

maximum frequency of the linear frequency scale in units of [1/s]

namespace log_frequencies

Variables

constexpr unsigned int N_omega = 2048

number of frequency values to compute in the logarithmic frequency [unitless]

namespace SI

Variables

constexpr float_64 omega_min = 1.0e14

mimimum frequency of the logarithmic frequency scale in units of [1/s]

constexpr float_64 omega_max = 1.0e17

maximum frequency of the logarithmic frequency scale in units of [1/s]

namespace parameters

Variables

constexpr unsigned int N_observer = 256

number of observation directions

namespace radFormFactor_CIC_1Dy
namespace radFormFactor_CIC_3D

correct treatment of coherent and incoherent radiation from macro particles

Choose different form factors in order to consider different particle shapes for radiation

  • radFormFactor_CIC_3D … CIC charge distribution

  • radFormFactor_TSC_3D … TSC charge distribution

  • radFormFactor_PCS_3D … PCS charge distribution

  • radFormFactor_CIC_1Dy … only CIC charge distribution in y

  • radFormFactor_Gauss_spherical … symmetric Gauss charge distribution

  • radFormFactor_Gauss_cell … Gauss charge distribution according to cell size

  • radFormFactor_incoherent … only incoherent radiation

  • radFormFactor_coherent … only coherent radiation

namespace radFormFactor_coherent
namespace radFormFactor_Gauss_cell
namespace radFormFactor_Gauss_spherical
namespace radFormFactor_incoherent
namespace radFormFactor_PCS_3D
namespace radFormFactor_TSC_3D
namespace radiationNyquist

selected mode of frequency scaling:

options:

  • linear_frequencies

  • log_frequencies

  • frequencies_from_list

Variables

constexpr float_32 NyquistFactor = 0.5

Nyquist factor: fraction of the local Nyquist frequency above which the spectra is set to zero should be in (0, 1).

namespace radWindowFunctionGauss
namespace radWindowFunctionHamming
namespace radWindowFunctionNone
namespace radWindowFunctionTriangle

add a window function weighting to the radiation in order to avoid ringing effects from sharpe boundaries default: no window function via radWindowFunctionNone

Choose different window function in order to get better ringing reduction radWindowFunctionTriangle radWindowFunctionHamming radWindowFunctionTriplett radWindowFunctionGauss radWindowFunctionNone

namespace radWindowFunctionTriplett

radiationObserver.param

This file defines a function describing the observation directions.

It takes an integer index from [ 0, picongpu::parameters::N_observer ) and maps it to a 3D unit vector in R^3 (norm=1) space that describes the observation direction in the PIConGPU cartesian coordinate system.

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 plugins
namespace radiation
namespace radiation_observer

Functions

HDINLINE vector_64 observationDirection (const int observation_id_extern)

Compute observation angles.

This function is used in the Radiation plug-in kernel to compute the observation directions given as a unit vector pointing towards a ‘virtual’ detector

This default setup is an example of a 2D detector array. It computes observation directions for 2D virtual detector field with its center pointing toward the +y direction (for theta=0, phi=0) with observation angles ranging from theta = [angle_theta_start : angle_theta_end] phi = [angle_phi_start : angle_phi_end ] Every observation_id_extern index moves the phi angle from its start value toward its end value until the observation_id_extern reaches N_split. After that the theta angle moves further from its start value towards its end value while phi is reset to its start value.

The unit vector pointing towards the observing virtual detector can be described using theta and phi by: x_value = sin(theta) * cos(phi) y_value = cos(theta) z_value = sin(theta) * sin(phi) These are the standard spherical coordinates.

The example setup describes an detector array of 16x16 detectors ranging from -pi/8= -22.5 degrees to +pi/8= +22.5 degrees for both angles with the center pointing toward the y-axis (laser propagation direction).

Parameters:

observation_id_extern – int index that identifies each block on the GPU to compute the observation direction

Returns:

unit vector pointing in observation direction type: vector_64

png.param

Defines

EM_FIELD_SCALE_CHANNEL1
EM_FIELD_SCALE_CHANNEL2
EM_FIELD_SCALE_CHANNEL3
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 scale_image = 1.0
constexpr bool scale_to_cellsize = true
constexpr bool white_box_per_GPU = false
namespace visPreview

Unnamed Group

DINLINE float_X preChannel1 (const float3_X &field_B, const float3_X &field_E, const float3_X &field_Current)

Calculate values for png channels for given field values.

Parameters:
  • field_B – normalized magnetic field value

  • field_E – normalized electric field value

  • field_Current – normalized electric current value (note - not current density)

DINLINE float_X preChannel2 (const float3_X &field_B, const float3_X &field_E, const float3_X &field_Current)
DINLINE float_X preChannel3 (const float3_X &field_B, const float3_X &field_E, const float3_X &field_Current)

Variables

constexpr float_64 customNormalizationSI[3] = {5.0e12 / SI::SPEED_OF_LIGHT_SI, 5.0e12, 15.0}

SI values to be used for Custom normalization.

The order of normalization values is: B, E, current (note - current, not current density). This variable must always be defined, but has no effect for other normalization types.

constexpr float_X preParticleDens_opacity = 0.25_X
constexpr float_X preChannel1_opacity = 1.0_X
constexpr float_X preChannel2_opacity = 1.0_X
constexpr float_X preChannel3_opacity = 1.0_X

pngColorScales.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 colorScales
namespace blue

Functions

HDINLINE void addRGB (float3_X &img, const float_X value, const float_X opacity)
namespace gray

Functions

HDINLINE void addRGB (float3_X &img, const float_X value, const float_X opacity)
namespace grayInv

Functions

HDINLINE void addRGB (float3_X &img, const float_X value, const float_X opacity)
namespace green

Functions

HDINLINE void addRGB (float3_X &img, const float_X value, const float_X opacity)
namespace none

Functions

HDINLINE void addRGB (const float3_X &, const float_X, const float_X)
namespace red

Functions

HDINLINE void addRGB (float3_X &img, const float_X value, const float_X opacity)

transitionRadiation.param

Definition of frequency space, number of observers, filters and form factors of the transition radiation plugin.

All values set here determine what the radiation plugin will compute. On the comand line the plugin still needs to be called for each species the transition radiation should be computed for.

Defines

PIC_VERBOSE_RADIATION

Uses the same verbose level schemes as the radiation plugin.

radiation verbose level: 0=nothing, 1=physics, 2=simulation_state, 4=memory, 8=critical

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 plugins
namespace radiation
namespace radFormFactor_CIC_1Dy
namespace radFormFactor_CIC_3D

correct treatment of coherent and incoherent radiation from macro particles

Choose different form factors in order to consider different particle shapes for radiation

  • radFormFactor_CIC_3D … CIC charge distribution

  • radFormFactor_TSC_3D … TSC charge distribution

  • radFormFactor_PCS_3D … PCS charge distribution

  • radFormFactor_CIC_1Dy … only CIC charge distribution in y

  • radFormFactor_Gauss_spherical … symmetric Gauss charge distribution

  • radFormFactor_Gauss_cell … Gauss charge distribution according to cell size

  • radFormFactor_incoherent … only incoherent radiation

  • radFormFactor_coherent … only coherent radiation

namespace radFormFactor_coherent
namespace radFormFactor_Gauss_cell
namespace radFormFactor_Gauss_spherical
namespace radFormFactor_incoherent
namespace radFormFactor_PCS_3D
namespace radFormFactor_TSC_3D
namespace transitionRadiation

Typedefs

using GammaFilter = picongpu::particles::manipulators::generic::Free<GammaFilterFunctor>

filter to (de)select particles for the radiation calculation

to activate the filter:

  • goto file speciesDefinition.param

  • add the attribute transitionRadiationMask to the particle species

Warning

Do not remove this filter. Otherwise still standing electrons would generate NaNs in the output of the plugin and transition radiation is scaling proportionally to gamma^2.

Functions

HDINLINE float3_X observationDirection (const int observation_id_extern)

Compute observation angles.

This function is used in the transition radiation plugin kernel to compute the observation directions given as a unit vector pointing towards a ‘virtual’ detector

This default setup is an example of a 2D detector array. It computes observation directions for 2D virtual detector field with its center pointing toward the +y direction (for theta=0, phi=0) with observation angles ranging from theta = [angle_theta_start : angle_theta_end] phi = [angle_phi_start : angle_phi_end ] Every observation_id_extern index moves the phi angle from its start value toward its end value until the observation_id_extern reaches N_split. After that the theta angle moves further from its start value towards its end value while phi is reset to its start value.

The unit vector pointing towards the observing virtual detector can be described using theta and phi by: x_value = sin(theta) * cos(phi) y_value = cos(theta) z_value = sin(theta) * sin(phi) These are the standard spherical coordinates.

The example setup describes an detector array of 128X128 detectors ranging from 0 to pi for the azimuth angle theta and from 0 to 2 pi for the polar angle phi.

If the calculation is only supposed to be done for a single azimuth or polar angle, it will use the respective minimal angle.

Parameters:

observation_id_extern – int index that identifies each block on the GPU to compute the observation direction

Returns:

unit vector pointing in observation direction type: float3_X

struct GammaFilterFunctor

example of a filter for the relativistic Lorentz factor gamma

Public Functions

template<typename T_Particle> inline HDINLINE void operator() (T_Particle &particle)

Public Static Attributes

static constexpr float_X filterGamma = 5.0

Gamma value above which the radiation is calculated, must be positive.

namespace linearFrequencies

units for linear frequencies distribution for transition radiation plugin

Variables

constexpr unsigned int nOmega = 512

number of frequency values to compute in the linear frequency [unitless]

namespace SI

Variables

constexpr float_64 omegaMin = 0.0

mimimum frequency of the linear frequency scale in units of [1/s]

constexpr float_64 omegaMax = 1.06e16

maximum frequency of the linear frequency scale in units of [1/s]

namespace listFrequencies

units for frequencies from list for transition radiation calculation

Variables

constexpr char listLocation[] = "/path/to/frequency_list"

path to text file with frequencies

constexpr unsigned int nOmega = 512

number of frequency values to compute if frequencies are given in a file [unitless]

namespace logFrequencies

units for logarithmic frequencies distribution for transition radiation plugin

Variables

constexpr unsigned int nOmega = 256

number of frequency values to compute in the logarithmic frequency [unitless]

namespace SI

Variables

constexpr float_64 omegaMin = 1.0e13

mimimum frequency of the logarithmic frequency scale in units of [1/s]

constexpr float_64 omegaMax = 1.0e17

maximum frequency of the logarithmic frequency scale in units of [1/s]

namespace parameters

selected mode of frequency scaling:

unit for foil position

options:

  • linearFrequencies

  • logFrequencies

  • listFrequencies correct treatment of coherent radiation from macro particles

These formfactors are the same as in the radiation plugin! Choose different form factors in order to consider different particle shapes for radiation

  • picongpu::plugins::radiation::radFormFactor_CIC_3D … CIC charge distribution

  • ::picongpu::plugins::radiation::radFormFactor_TSC_3D … TSC charge distribution

  • ::picongpu::plugins::radiation::radFormFactor_PCS_3D … PCS charge distribution

  • ::picongpu::plugins::radiation::radFormFactor_CIC_1Dy … only CIC charge distribution in y

  • ::picongpu::plugins::radiation::radFormFactor_Gauss_spherical … symmetric Gauss charge distribution

  • ::picongpu::plugins::radiation::radFormFactor_Gauss_cell … Gauss charge distribution according to cell size

  • ::picongpu::plugins::radiation::radFormFactor_incoherent … only incoherent radiation

  • ::picongpu::plugins::radiation::radFormFactor_coherent … only coherent radiation

Variables

constexpr unsigned int nPhi = 128

Number of observation directions.

If nPhi or nTheta is equal to 1, the transition radiation will be calculated for phiMin or thetaMin respectively.

constexpr unsigned int nTheta = 128
constexpr unsigned int nObserver = nPhi * nTheta
constexpr float_64 thetaMin = 0.0
constexpr float_64 thetaMax = picongpu::PI
constexpr float_64 phiMin = 0.0
constexpr float_64 phiMax = 2 * picongpu::PI
namespace SI

Variables

constexpr float_64 foilPosition = 0.0