Plugins

fileOutput.param

namespace picongpu

Typedefs

using picongpu::ChargeDensity_Seq = typedef 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
  • MomentumComponent: ratio between a selected momentum component and the absolute momentum with respect to shape
  • 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

using picongpu::EnergyDensity_Seq = typedef deriveField::CreateEligible_t< VectorAllSpecies, deriveField::derivedAttributes::EnergyDensity >
using picongpu::MomentumComponent_Seq = typedef deriveField::CreateEligible_t< VectorAllSpecies, deriveField::derivedAttributes::MomentumComponent< 0 > >
using picongpu::FieldTmpSolvers = typedef MakeSeq_t< ChargeDensity_Seq, EnergyDensity_Seq, MomentumComponent_Seq >

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

FieldTmpSolvers is used in

See
FieldTmp to calculate the exchange size

using picongpu::NativeFileOutputFields = typedef MakeSeq_t< FieldE, FieldB >

FileOutputFields: Groups all Fields that shall be dumped.

Possible native fields: FieldE, FieldB, FieldJ

using picongpu::FileOutputFields = typedef MakeSeq_t< NativeFileOutputFields, FieldTmpSolvers >
using picongpu::FileOutputParticles = typedef 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 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
namespace picongpuisaacP

Typedefs

using picongpu::isaacP::Native_Seq = typedef MakeSeq_t< FieldE, FieldB, FieldJ >

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

using picongpu::isaacP::Density_Seq = typedef deriveField::CreateEligible_t< VectorAllSpecies, deriveField::derivedAttributes::Density >

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

using picongpu::isaacP::Fields_Seq = typedef MakeSeq_t< Native_Seq, Density_Seq >

Compile time sequence of all fields which shall be visualized.

Basically the join of Native_Seq and Density_Seq.

particleCalorimeter.param

namespace picongpu
namespace picongpuparticleCalorimeter

Functions

HDINLINE float2_X picongpu::particleCalorimeter::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.

Return
Two values within [-1,1]
Parameters
  • yaw: -maxYaw…maxYaw
  • pitch: -maxPitch…maxPitch
  • maxYaw: maximum value of angle yaw
  • maxPitch: maximum value of angle pitch

particleMerger.param

namespace picongpu
namespace picongpuplugins
namespace picongpu::pluginsparticleMerging

Variables

constexpr size_t picongpu::plugins::particleMergingMAX_VORONOI_CELLS = 128

maximum number of active Voronoi cells per supercell.

If the number of active Voronoi cells reaches this limit merging events are dropped.

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
namespace picongpuparameters

Variables

constexpr unsigned int picongpu::parametersN_observer = 256

number of observation directions

namespace picongpurad_frequencies_from_list

Variables

constexpr char picongpu::rad_frequencies_from_listlistLocation[] = "/path/to/frequency_list"

path to text file with frequencies

constexpr unsigned int picongpu::rad_frequencies_from_listN_omega = 2048

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

namespace picongpurad_linear_frequencies

Variables

constexpr unsigned int picongpu::rad_linear_frequenciesN_omega = 2048

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

namespace picongpu::rad_linear_frequenciesSI

Variables

constexpr float_64 picongpu::rad_linear_frequencies::SIomega_min = 0.0

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

constexpr float_64 picongpu::rad_linear_frequencies::SIomega_max = 1.06e16

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

namespace picongpurad_log_frequencies

Variables

constexpr unsigned int picongpu::rad_log_frequenciesN_omega = 2048

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

namespace picongpu::rad_log_frequenciesSI

Variables

constexpr float_64 picongpu::rad_log_frequencies::SIomega_min = 1.0e14

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

constexpr float_64 picongpu::rad_log_frequencies::SIomega_max = 1.0e17

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

namespace picongpuradFormFactor_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 picongpuradiation

Typedefs

using picongpu::radiation::RadiationParticleFilter = typedef 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 picongpu::radiationGammaFilterFunctor

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

Public Functions

template <typename T_Particle>
HDINLINE void picongpu::radiation::GammaFilterFunctor::operator()(T_Particle & particle)

Public Static Attributes

constexpr float_X picongpu::radiation::GammaFilterFunctorradiationGamma = 5.0

Gamma value above which the radiation is calculated.

namespace picongpuradiationNyquist

selected mode of frequency scaling:

options:

  • rad_linear_frequencies
  • rad_log_frequencies
  • rad_frequencies_from_list

Variables

constexpr float_32 picongpu::radiationNyquistNyquistFactor = 0.5

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

namespace picongpuradWindowFunctionTriangle

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

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
namespace picongpuradiation_observer

Functions

HDINLINE vector_64 picongpu::radiation_observer::observation_direction(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).

Return
unit vector pointing in observation direction type: vector_64
Parameters
  • observation_id_extern: int index that identifies each block on the GPU to compute the observation direction

png.param

Defines

EM_FIELD_SCALE_CHANNEL1
EM_FIELD_SCALE_CHANNEL2
EM_FIELD_SCALE_CHANNEL3
namespace picongpu

Variables

constexpr float_64 picongpuscale_image = 1.0
constexpr bool picongpuscale_to_cellsize = true
constexpr bool picongpuwhite_box_per_GPU = false
namespace picongpuvisPreview

Functions

DINLINE float_X picongpu::visPreview::preChannel1(const float3_X & field_B, const float3_X & field_E, const float3_X & field_J)
DINLINE float_X picongpu::visPreview::preChannel2(const float3_X & field_B, const float3_X & field_E, const float3_X & field_J)
DINLINE float_X picongpu::visPreview::preChannel3(const float3_X & field_B, const float3_X & field_E, const float3_X & field_J)

Variables

constexpr float_X picongpu::visPreviewpreParticleDens_opacity = 0.25_X
constexpr float_X picongpu::visPreviewpreChannel1_opacity = 1.0_X
constexpr float_X picongpu::visPreviewpreChannel2_opacity = 1.0_X
constexpr float_X picongpu::visPreviewpreChannel3_opacity = 1.0_X

pngColorScales.param

namespace picongpu
namespace picongpucolorScales
namespace picongpu::colorScalesblue

Functions

HDINLINE void picongpu::colorScales::blue::addRGB(float3_X & img, const float_X value, const float_X opacity)
namespace picongpu::colorScalesgray

Functions

HDINLINE void picongpu::colorScales::gray::addRGB(float3_X & img, const float_X value, const float_X opacity)
namespace picongpu::colorScalesgrayInv

Functions

HDINLINE void picongpu::colorScales::grayInv::addRGB(float3_X & img, const float_X value, const float_X opacity)
namespace picongpu::colorScalesgreen

Functions

HDINLINE void picongpu::colorScales::green::addRGB(float3_X & img, const float_X value, const float_X opacity)
namespace picongpu::colorScalesnone

Functions

HDINLINE void picongpu::colorScales::none::addRGB(const float3_X &, const float_X, const float_X)
namespace picongpu::colorScalesred

Functions

HDINLINE void picongpu::colorScales::red::addRGB(float3_X & img, const float_X value, const float_X opacity)