Transition Radiation
The spectrally resolved far field radiation created by electrons passing through a metal foil.
Our simulation computes the transition radiation to calculate the emitted electromagnetic spectra for different observation angles.
Variable |
Meaning |
|---|---|
\(N_e\) |
Amount of real electrons |
\(\psi\) |
Azimuth angle of momentum vector from electrons to y-axis of simulation |
\(\theta\) |
Azimuth angle of observation vector |
\(\phi\) |
Polar angle between momentum vector from electrons and observation vector |
\(\omega\) |
The circular frequency of the radiation that is observed. |
\(h(\vec{r}, \vec{p})\) |
Normalized phasespace distribution of electrons |
\(g(\vec{p})\) |
Normalized momentum distribution of electrons |
\(\vec{k}\) |
Wavevector of electrons |
\(\vec{v}\) |
Velocity vector of electrons |
\(u\) |
Normalized momentum of electrons \(\beta \gamma\) |
\(\mathcal{E}\) |
Normalized electric fields generated by electrons |
\(\mathcal{N}\) |
Denominator of normalized energies |
\(F\) |
Normalized formfactor of electrons, contains phase informations |
This plugin allows to predict the emitted virtual transition radiation, which would be caused by the electrons in the simulation box passing through a virtual metal foil which is set at a specific location. The transition radiation can only be calculated for electrons at the moment.
External Dependencies
The plugin is available as soon as the openPMD API is compiled in.
.param files
In order to setup the transition radiation plugin, the transitionRadiation.param has to be configured and the radiating particles need to have the attributes weighting, momentum, position, and transitionRadiationMask (which can be added in speciesDefinition.param) as well as the flags massRatio and chargeRatio.
In transitionRadiation.param, the number of frequencies nOmega and observation directions nTheta and nPhi are defined.
Frequency range
The frequency range is set up by choosing a specific namespace that defines the frequency setup
/* choose linear frequency range */
namespace frequencies = linearFrequencies;
Currently you can choose from the following setups for the frequency range:
namespace |
Description |
|---|---|
|
linear frequency range from |
|
logarithmic frequency range from |
|
|
All three options require variable definitions in the according namespaces as described below:
For the linear frequency scale all definitions need to be in the picongpu::plugins::transitionRadiation::linearFrequencies namespace.
The number of total sample frequencies nOmega need to be defined as constexpr unsigned int.
In the sub-namespace SI, a minimal frequency omegaMin and a maximum frequency omegaMax need to be defined as constexpr float_64.
For the logarithmic frequency scale all definitions need to be in the picongpu::plugins::transitionRadiation::logFrequencies namespace.
Equivalently to the linear case, three variables need to be defined:
The number of total sample frequencies nOmega need to be defined as constexpr unsigned int.
In the sub-namespace SI, a minimal frequency omegaMin and a maximum frequency omegaMax need to be defined as constexpr float_64.
For the file-based frequency definition, all definitions need to be in the picongpu::plugins::transitionRadiation::listFrequencies namespace.
The number of total frequencies nOmega need to be defined as constexpr unsigned int and the path to the file containing the frequency values in units of \([\mathrm{s}^{-1}]\) needs to be given as constexpr char * listLocation = "/path/to/frequency_list";.
The frequency values in the file can be separated by newlines, spaces, tabs, or any other whitespace. The numbers should be given in such a way, that C++ standard std::ifstream can interpret the number e.g., as 2.5344e+16.
Note
Currently, the variable listLocation is required to be defined in the picongpu::plugins::transitionRadiation::listFrequencies namespace, even if listFrequencies is not used.
The string does not need to point to an existing file, as long as the file-based frequency definition is not used.
Observation directions
The number of observation directions nTheta and the distribution of observation directions is defined in transitionRadiation.param.
There, the function observationDirection defines the observation directions.
This function returns the x, y and z component of a unit vector pointing in the observation direction.
HDINLINE float3_X
observationDirection( const int observation_id_extern )
{
/* use the scalar index const int observation_id_extern to compute an
* observation direction (x,y,z) */
return float3_X( x , y , z );
}
Note
The transitionRadiation.param set up will be subject to further changes, since the radiationObserver.param it is based on is subject to further changes.
These might be namespaces that describe several preconfigured layouts or a functor if C++ 11 is included in the nvcc.
Foil Position
This parameter defines the position of the metal foil.
A foil position which is unequal to the current positions of the electrons adds the electrons momentum vectors onto the electrons until they reach the given y-coordinate.
This gives a rough estimate of the effect of the divergence on the electron bunch.
The position is controlled by a run time parameter which can be set with --<species>_transRad.foilPositionY, by default set to 0.
One should make sure that the foil is positioned in front of the bunch to get physically interpretable results.
The foil can also change its position together with the comoving window when a second parameter --<species>_transRad.comovingFoil is set to true.
In this case, foilPositionY defines a relative offset between the window and the foil.
If the foil position is set to zero, the foil is always positioned close to the electron bunch.
The following image shows the relation of foil and bunch position for some possible example parameters.
Macro-particle form factor
The macro-particle form factor is a method, which considers the shape of the macro particles when computing the radiation.
One can select between different macro particle shapes. Currently eight shapes are implemented. A shape can be selected by choosing one of the available namespaces:
/* choosing the 3D CIC-like macro particle shape */
namespace macroParticleFormFactor = radFormFactor_CIC_3D;
Namespace |
Description |
|---|---|
|
3D Cloud-In-Cell shape |
|
3D Triangular shaped density cloud |
|
3D Quadratic spline density shape (Piecewise Cubic Spline assignment function) |
|
Cloud-In-Cell shape in y-direction, dot like in the other directions |
|
symmetric Gauss charge distribution |
|
Gauss charge distribution according to cell size |
|
forces a completely incoherent emission by scaling the macro particle charge with the square root of the weighting |
|
forces a completely coherent emission by scaling the macro particle charge with the weighting |
Note
One should not confuse this macro-particle form factor with the form factor \(F\), which was previously mentioned. This form factor is equal to the macro-particle shape, while \(F\) contains the phase information of the whole electron bunch. Both are necessary for a physically correct transition radiation calculation.
Gamma filter
In order to consider the radiation only of particles with a gamma higher than a specific threshold.
In order to do that, the radiating particle species needs the flag transitionRadiationMask (which is initialized as false) which further needs to be manipulated, to set to true for specific (random) particles.
Using a filter functor as:
using GammaFilter = picongpu::particles::manipulators::generic::Free<
GammaFilterFunctor
>;
(see TransitionRadiation example for details)
sets the flag to true if a particle fulfills the gamma condition.
Note
More sophisticated filters might come in the near future. Therefore, this part of the code might be subject to changes.
.cfg file
For a specific (charged) species <species> e.g. e, the radiation can be computed by the following commands.
Command line option |
Description |
|---|---|
|
Gives the number of time steps between which the radiation should be calculated. |
|
Position in SI units to put a virtual foil for calculating the transition radiation. See above for more information. Default: 0.0 |
|
Flag that determines if virtual foil moves together with the comoving window. Enabled = 1. Default: 0 |
|
File name to store transition radiation in. Default: transRad |
|
Backend for openPMD output. Default: default that is used internally |
|
Optional text file output in numpy-readable format. Enabled = 1. Default: 0 |
Memory Complexity
Accelerator
two counters (float_X) and two counters (complex_X) are allocated permanently
Host
as on accelerator.
Output
Contains ASCII files in simOutput/transRad that have the total spectral intensity until the timestep specified by the filename.
Each row gives data for one observation direction (same order as specified in the observationDirection function).
The values for each frequency are separated by tabs and have the same order as specified in transitionRadiation.param.
The spectral intensity is stored in the units \(\mathrm{[Js]}\).
Analysing tools
The transition_radiation_visualizer.py in lib/python/picongpu/extra/plugins/plot_mpl can be used to analyze the radiation data after the simulation.
See transition_radiation_visualizer.py --help for more information.
It only works, if the input frequencies are divided logarithmically!
Known Issues
The output is currently only physically correct for electron passing through a metal foil.
References
Schroeder, C. B. and Esarey, E. and van Tilborg, J. and Leemans, W. P., Theory of coherent transition radiation generated at a plasma-vacuum interface, American Physical Society (2004) https://doi.org/10.1103/PhysRevE.69.016501
Downer, M. C. and Zgadzaj, R. and Debus, A. and Schramm, U. and Kaluza, M. C., Diagnostics for plasma-based electron accelerators, American Physical Society (2018), https://doi.org/10.1103/RevModPhys.90.035002
Carstens, F.-O., Synthetic characterization of ultrashort electron bunches using transition radiation, Bachelor thesis on the transition radiation plugin (2019), https://doi.org/10.5281/zenodo.3469663
Pausch, R., Quantitatively consistent computation of coherent and incoherent radiation in particle-in-cell codes — A general form factor formalism for macro-particles, Description for the effect of macro-particle shapes in particle-in-cell codes (2018), https://doi.org/10.1016/j.nima.2018.02.020