jnu_mixed.cu

Documentation for the file jnu_mixed.cu. This file contains functions related to the thermal synchrotron emissivity calculation.

Functions

Functions

__host__ __device__ double jnu_synch(const double nu, const double Ne, const double Thetae, const double B, const double theta)

Calculates the thermal synchrotron emissivity \(j_\nu(\nu, \theta)\) following Leung et al.(2011).

This function implements the analytical fit for the emissivity of a thermal distribution of electrons: \( j_\nu(\nu, \theta) = \frac{\sqrt{2}\pi e^2 n_e \nu_s}{3 c K_2(1/\Theta_e)} (X^{1/2} + 2^{11/12} X^{1/6})^2 \exp(-X^{1/3}) \)

Parameters:
  • nu – Photon frequency in the plasma frame.

  • Ne – Electron number density ( \(n_{\rm e}\)).

  • Thetae – Dimensionless electron temperature ( \(\Theta_{\rm e}\)).

  • B – Magnetic field strength.

  • theta – Angle between the magnetic field and the wave vector ( \(\theta\)).

  • besselTexObj – (CUDA only) Texture object for Bessel function look-up tables.

Returns:

The local emissivity \( j_\nu \) in CGS units.

__host__ __device__ double int_jnu(double Ne, double Thetae, double Bmag, double nu)

Calculates the angle-integrated thermal synchrotron emissivity \( J_\nu \) in CGS units.

This function evaluates the total energy emitted per unit time, per unit volume, and per unit frequency, integrated over all solid angles \( \int j_\nu d\Omega \).

  • Implements: \( J_\nu = \frac{\sqrt{2} e^3 n_e B \Theta_e^2}{27 m_e c^2 K_2(1/\Theta_e)} F(\Theta_e, B, \nu) \)

Parameters:
  • Ne – Electron number density \( n_e \).

  • Thetae – Dimensionless electron temperature \( \Theta_{\rm e} \).

  • Bmag – Magnetic field strength \( B \).

  • nu – Photon frequency \( \nu \) in the plasma frame.

  • besselTexObj – (CUDA only) Texture object for Bessel function look-up tables.

Returns:

The integrated emissivity in \( \text{erg} \cdot \text{s}^{-1} \cdot \text{cm}^{-3} \cdot \text{Hz}^{-1} \).

double jnu_integrand(double th, void *params)

Provides the integrand for calculating the angle-integrated thermal synchrotron emissivity.

This function is used by the GSL integrator to evaluate the integral over solid angles \( \int j_\nu \sin\theta d\theta \).

The integrand is given by:

\[ \sin^2\theta \left( (K/\sin\theta)^{1/2} + 2^{11/12}(K/\sin\theta)^{1/6} \right)^2 \exp\left[-(K/\sin\theta)^{1/3}\right] \]

Parameters:
  • th – The angle \( \theta \) between the magnetic field and the photon wave vector.

  • params – A pointer to the dimensionless frequency parameter \( K \).

Returns:

The value of the integrand at the given angle.

__host__ void init_emiss_tables(void)

Initializes lookup tables for the angle-integrated emissivity function \( F \) and the modified Bessel function \( K_2 \).

Precomputes \( F(K) \) by integrating the angular dependence of synchrotron emission over the solid angle using GSL integration and stores the second-order modified Bessel function \( K_2(1/\Theta_e) \) in log-space for rapid retrieval during the simulation.

The integration goes as:

\(F[k] = \ln \left( 4\pi \int_{0}^{\pi/2} \sin^2\theta \left[ \sqrt{\frac{K}{\sin\theta}} + 2^{11/12} \left( \frac{K}{\sin\theta} \right)^{1/6} \right]^2 \exp\left[ -\left( \frac{K}{\sin\theta} \right)^{1/3} \right] d\theta \right)\)

Returns:

void

__host__ __device__ double K2_eval(const double Thetae)

Evaluation of the modified Bessel function of the second kind, order 2, \( K_2(1/\Theta_e) \).

This function retrieves the value of \( K_2(1/\Theta_e) \) from precomputed lookup tables for efficient computation during simulations.

Note

for large \( \Theta_e \), we can use the asymptotic expansion \( K_2(1/\Theta_e) \approx 2 \Theta_e^2 \).

Note

TMAX and THETA_MIN define the valid range for \( \Theta_e \) in the lookup tables.

Parameters:
  • Thetae – Dimensionless electron temperature \( \Theta_e \).

  • besselTexObj – (CUDA only) Texture object for Bessel function look-up tables.

Returns:

The value of \( K_2(1/\Theta_e) \).

__host__ __device__ double F_eval(const double Thetae, const double Bmag, const double nu)

Evaluates the frequency-dependent component \( F(K) \) of the angle-integrated emissivity. This function calculates the dimensionless frequency parameter \( K \) and retrieves the corresponding value of the precomputed integral \( \int j_\nu d\Omega \).

  • If \( K > K_{\rm MAX} \): Returns 0 due to the exponential decay of the synchrotron spectrum.

  • If \( K < K_{\rm MIN} \): Uses a power-law approximation based on the asymptotic behavior of the emissivity at low frequencies.

  • Otherwise: Performs a linear interpolation using the precomputed F table.

Parameters:
  • Thetae – Dimensionless electron temperature \( \Theta_e \).

  • Bmag – Magnetic field strength \( B \).

  • nu – Photon frequency \( \nu \) in the plasma frame.

Returns:

The value of the integrated emissivity function \( F \).

__host__ __device__ double linear_interp_F(const double K)

Performs linear interpolation on the precomputed emissivity function \( F(K) \).

Parameters:

K – Dimensionless frequency parameter \( K \).

Returns:

The interpolated value of \( F(K) \).

__host__ __device__ double linear_interp_K2(const double Thetae)

Performs linear interpolation on the precomputed modified Bessel function \( K_2(1/\Theta_e) \).

Parameters:
  • Thetae – Dimensionless electron temperature \( \Theta_e \).

  • besselTexObj – (CUDA only) Texture object for Bessel function look-up tables.

Returns:

The interpolated value of \( K_2(1/\Theta_e) \).