radiation.cu
This section contains functions defined in the radiation.cu file related to the calculation of quantities associated with radiation transport.
Functions
Functions
-
__device__ double jnu_inv(const double nu, const double Thetae, const double Ne, const double B, const double theta, cudaTextureObject_t besselTexObj)
Calculates the invariant emissivity for synchrotron radiation.
This function computes the local synchrotron emissivity \( j_\nu \) through means of the function
jnu_synchand converts it to the Lorentz-invariant form \( j_\nu / \nu^2 \).Emissivity \( j_\nu \) is determined by the local plasma properties ( \( N_e, \Theta_e, B \)) and the emission angle \( \theta \) relative to the magnetic field.
- Parameters:
nu – Photon frequency \( \nu \) in the fluid frame.
Thetae – Dimensionless electron temperature \( \Theta_e = k_B T_e / m_e c^2 \).
Ne – Electron number density \( N_e \) in cgs.
B – Magnetic field strength \( B \) in Gauss.
theta – Pitch angle \( \theta \) between the photon wavevector and the magnetic field.
besselTexObj – CUDA texture object used for fast lookup of modified Bessel functions.
- Returns:
The invariant emissivity \( j_\nu / \nu^2 \).
-
__device__ double Bnu_inv(const double nu, const double Thetae)
Calculates the invariant Planck source function for a thermal distribution.
This function computes the frequency-dependent Planck function \( B_\nu(T) \) and returns it in its Lorentz-invariant form: \( B_\nu / \nu^3 \). This quantity represents the invariant intensity of a blackbody in thermal equilibrium.
The invariant source function is defined as:
\[ \mathcal{B} = \frac{B_\nu}{\nu^3} = \frac{2h}{c^2} \frac{1}{\exp\left(\frac{h\nu}{k_B T}\right) - 1} \]where \( k_B T = m_e c^2 \Theta_e \).Note
To prevent precision loss and division by zero at low frequencies, the function uses a 4th-order Taylor expansion of \( e^x - 1 \) when \( x < 10^{-3} \):
\[ e^x - 1 \approx x + \frac{x^2}{2} + \frac{x^3}{6} + \frac{x^4}{24} \]- Parameters:
nu – Photon frequency \( \nu \) in the fluid frame.
Thetae – Dimensionless electron temperature \( \Theta_e \).
- Returns:
The invariant Planck function \( B_\nu / \nu^3 \) in CGS units.
-
__device__ double kappa_es(const double nu, const double Thetae, const double *d_table_ptr)
Calculates the electron scattering opacity in CGS units.
This function determines the scattering opacity \(\kappa_{\rm es}\) by looking up the total Compton cross-section for a given photon frequency and electron temperature. It assumes a pure hydrogen composition to convert the cross-section (area) into opacity (area per unit mass).
The opacity is calculated as: [ \kappa_{es} = \frac{\sigma_{Compton}(E_g, \Theta_{\rm e})}{m_p} ] where \( E_{\rm g} = \frac{h\nu}{m_{\rm e} c^2} \) is the dimensionless photon energy and \( m_{\rm p} \) is the proton mass.
- Parameters:
nu – Photon frequency in the fluid frame (Hz).
Thetae – Dimensionless electron temperature \( \Theta_{\rm e} \).
d_table_ptr – Pointer to the GPU memory containing the precomputed cross-section lookup table.
- Returns:
The electron scattering opacity in \( \text{cm}^2/\text{g} \).
-
__device__ double get_fluid_nu(const double X[NDIM], const double K[NDIM], const double Ucov[NDIM])
Computes the photon frequency in the local fluid frame (Hz).
This function projects the photon 4-wavevector \(K^\mu\) onto the covariant fluid 4-velocity \(U_\mu\) to determine the energy as measured by a co-moving observer.
The energy of a particle with 4-momentum \(p^\mu\) measured by an observer with 4-velocity \(u^\mu\) is given by the scalar \(E = -p^\mu u_\mu\).
Note
The energy is initially calculated in dimensionless electron rest-mass units ( \(m_e c^2\)). This is converted to physical frequency using the relation \(\nu = E \cdot (m_e c^2 / h)\).
- Parameters:
X – The 4-position of the photon (primarily used for diagnostic logging).
K – The photon 4-wavevector \(K^\mu\).
Ucov – The covariant fluid 4-velocity \(U_\mu\).
- Returns:
The physical frequency \(\nu\) in Hz.
-
__device__ double alpha_inv_scatt(const double nu, const double Thetae, const double Ne, const double *d_table_ptr)
Calculates the Lorentz invariant scattering coefficient.
This function converts the scattering opacity (cross-section per unit mass) into the invariant scattering coefficient \(\nu \alpha_{s}\).
[ \text{Result} = \nu \cdot \alpha_{s} = \nu \cdot (\rho \cdot \kappa_{es}) ] where \(\rho \approx N_e m_p\) (assuming a pure hydrogen plasma).
- Parameters:
nu – Photon frequency in the fluid frame (Hz).
Thetae – Dimensionless electron temperature \(\Theta_\rm{e} \).
Ne – Electron number density \(N_e\) (cm \(^{-3}\)).
d_table_ptr – Pointer to the Compton cross-section lookup table.
- Returns:
The invariant scattering coefficient \(\nu \alpha_{s}\) [units: Hz \( \rm{cm}^{-1} \)].
-
__device__ double alpha_inv_abs(const double nu, const double Thetae, const double Ne, const double B, double theta, cudaTextureObject_t besselTexObj)
Calculates the Lorentz invariant absorption coefficient via Kirchhoff’s Law.
This function determines the invariant absorption coefficient \(\nu \alpha_{\rm a}\) by relating the invariant emissivity to the invariant Planck function (source function).
Based on Kirchhoff’s Law of Thermal Radiation, the absorption coefficient \(\alpha_\nu\) is the ratio of emissivity to the Planck function: [ \alpha_\nu = \frac{j_\nu}{B_\nu} ] In terms of the invariant quantities provided by the internal helper functions: [ \nu \alpha_\nu = \frac{j_\nu / \nu^2}{B_\nu / \nu^3} ]
Note
A small epsilon ( \( 10^{-100} \)) is added to the denominator to prevent division by zero.
- Parameters:
nu – Photon frequency in the fluid frame (Hz).
Thetae – Dimensionless electron temperature \(\Theta_\rm{e} \).
Ne – Electron number density \(N_e\) (cm \(^{-3}\)).
B – Magnetic field strength (Gauss).
theta – Pitch angle between photon and magnetic field.
besselTexObj – CUDA texture for Bessel function lookups.
- Returns:
The invariant absorption coefficient \(\nu \alpha_{a}\) [units: Hz \( \rm{cm}^{-1} \)].
-
__device__ double get_bk_angle(const double X[NDIM], const double K[NDIM], const double Ucov[NDIM], const double Bcov[NDIM], const double B)
Calculates the pitch angle between the photon wavevector and the magnetic field.
This function determines the angle \(\theta\) between the photon’s propagation direction and the local magnetic field vector, as measured in the local fluid frame.
- Parameters:
X – Photon 4-position.
K – Photon 4-wavevector \(K^\mu\).
Ucov – Covariant fluid 4-velocity \(U_\mu\).
Bcov – Covariant magnetic field 4-vector \(B_\mu\).
B – Magnetic field strength in Gauss (CGS).
- Returns:
The pitch angle \(\theta\) in radians \([0, \pi]\).