Betatron radiation

This plugin allows to calculate on-fly the synchrotron radiation assuming a simplified description valid for the relativistic particles and the strong-wiggler regime. These assumptions can be formulated as the conditions on the Lorentz factor of the emitting particle \(\gamma_p\gg 1\), and the angular excursion of it’s oscillations is \(\sigma_{\theta}\gg 1/\gamma_p\). This is a typical case for the laser-plasma betatron or non-linear Compton source.

For the assumptions above, radiation emitted by an electron during time \(\Delta t\) can be described by the dimensionless spectral density distribution as:

\[\frac{\partial E_{rad} }{ \partial \hbar \omega_{ph}} = \frac{P_{rad}\; \Delta t }{ \hbar\omega_c } \; S\left(\frac{\omega}{\omega_c}\right) \,,\]

where the instantaneous power follows from the relativistic Larmor formula:

\[P_{rad} = \frac{e^2 \gamma^6}{6 \pi\epsilon_0 c} \; [|\dot{\boldsymbol \beta}|^2 - (\boldsymbol \beta\times\dot{\boldsymbol \beta})^2] \,.\]

The spectral profile is approximated with the standard function of \(\xi=\omega/\omega_c\):

\[S(x) = \frac{9\sqrt{3}}{8\pi} \; x\; \int_x^{\infty} K_{5/3}(\xi)\; d\xi'\,,\]

where critical frequency is calculated from with curvature radius along trajectory:

\[\omega_c = \frac{3}{2}\;\gamma^3\; \cfrac{c }{ \rho} = \frac{3}{2}\;\gamma^3 \frac{\sqrt{ \left|\dot{\boldsymbol \beta}\right|^2 - \left(\boldsymbol \beta \cdot \dot{\boldsymbol\beta} \right)^2 / |\boldsymbol \beta|^2} } {|\boldsymbol \beta|^2} \,,\]

and motion quantities are calculated from the motion equations with fields:

\[ \begin{align}\begin{aligned}\boldsymbol\beta = \mathbf{u} / \gamma \,,\;\; \partial_t \mathbf{u} = - \frac{e}{m c} ( \boldsymbol E + \boldsymbol\beta \times c \mathbf{B} )\,,\\\partial_t \gamma = - \frac{e}{m c} \; \boldsymbol\beta \cdot \mathbf{E} \,,\;\; \partial_t \boldsymbol\beta = \cfrac{1}{\gamma} \;(\partial_t \mathbf{u} - \boldsymbol\beta \; \partial_t \gamma)\,.\end{aligned}\end{align} \]

Radiation is emitted in the direction along electron propagation \(\boldsymbol u/|\boldsymbol u|\) into the angular profile approximated as \(\propto (1-\beta\cos\theta)^{-5}\). The field is calculated in the discrete 3D spectral-angular domain \((\theta_x, \theta_y, \hbar\omega)\), and the particles are deposited using fast bilinear interpolations, so to account for the changing angular shapes and emission directions we add to emission direction a random angular spread \(\delta\theta \sim \mathcal N(2^{-3/2} \gamma^{-1})\), and with time averages to the proper profile.

Optionally the plugin can deduce the photon momentum from the emitting electron in order to simulate the classical radiation reaction.

Warning

Currently this plugin only works in the lab-frame simulations.

Plugin is activated for the chosen particle specie (electron) with the following method:

Particles.activate_synchrotron(photon_energy_axis, theta_x_axis, theta_y_axis, gamma_cutoff=10.0, radiation_reaction=False, x_max=20, nSamples=2048, boost=None)[source]

Activate synchrotron radiation.

Parameters:
  • photon_energy_axis (tuple) – Parameters for the photon energy axis provided as (photon_energy_min, photon_energy_max, N_photon_energy), where photon_energy_min and photon_energy_max are floats in Joules and N_photon_energy is integer

  • theta_x_axis (tuple) – Parameters for the x-elevation angle axis provided as (theta_x_min, theta_x_max, N_theta_x), where theta_x_min and theta_x_max are floats in (rad) and N_theta_x is integer

  • theta_y_axis (tuple) – Parameters for the y-elevation angle axis provided as (theta_y_min, theta_y_max, N_theta_y), where theta_y_min and theta_y_max are floats in radians and N_theta_y is integer

  • gamma_cutoff (float, optional) – Minimal particle gamma factor for which radiation is calculated

  • radiation_reaction (bool) – Whether to consider radiation reaction on the electrons

  • x_max (float, optional) – Extent of the sampling used for the spectral profile function

  • nSamples (integer, optional) – number of sampling points for the spectral profile function

The calculated radiation data in the dimensionless units (\(\partial^2 E_{rad} / \partial \hbar \omega_{ph} \partial \theta\)) can be added as a standard openPMD diagnostics:

class fbpic.openpmd_diag.SynchrotronRadiationDiagnostic(period=None, dt_period=None, species={}, comm=None, write_dir=None, iteration_min=0, iteration_max=inf)[source]

Initialize the synchrotron radiation diagnostic

Parameters:
  • period (int, optional) – The period of the diagnostics, in number of timesteps. (i.e. the diagnostics are written whenever the number of iterations is divisible by period). Specify either this or dt_period

  • dt_period (float (in seconds), optional) – The period of the diagnostics, in physical time of the simulation. Specify either this or period

  • species (a dictionary of Particles objects) – The object that is written (e.g. elec) is assigned to the particle name of this species. (e.g. {“electrons”: elec }). All species must have synchrotron radiation activated with the same specral-angular grids.

  • comm (an fbpic BoundaryCommunicator object or None) – If this is not None, the data is gathered on the first proc, and the guard cells are removed from the output. Otherwise, each proc writes its own data, including guard cells (Make sure to use different write_dir in this case)

  • write_dir (string, optional) – The POSIX path to the directory where the results are to be written. If none is provided, this will be the path of the current working directory

  • iteration_min (ints, optional) – The iterations between which data should be written (iteration_min is inclusive, iteration_max is exclusive)

  • iteration_max (ints, optional) – The iterations between which data should be written (iteration_min is inclusive, iteration_max is exclusive)