thermal_bremsstrahlung
- plasmapy.formulary.radiation.thermal_bremsstrahlung(
- frequencies: Annotated[Quantity, Unit('Hz')],
- n_e: Annotated[Quantity, Unit('1 / m3')],
- T_e: Annotated[Quantity, Unit('K')],
- n_i: Annotated[Quantity, Unit('1 / m3')] = None,
- ion: str | int | integer | Particle | CustomParticle | Quantity = 'p+',
- kmax: Annotated[Quantity, Unit('m')] = None,
Calculate the bremsstrahlung emission spectrum for a Maxwellian plasma in the Rayleigh-Jeans limit \(ℏ ω ≪ k_B T_e\).
\[\frac{dP}{dω} = \frac{8 \sqrt{2}}{3\sqrt{π}} \bigg ( \frac{e^2}{4 π ε_0} \bigg )^3 \bigg ( m_e c^2 \bigg )^{-\frac{3}{2}} \bigg ( 1 - \frac{ω_{pe}^2}{ω^2} \bigg )^\frac{1}{2} \frac{Z_i^2 n_i n_e}{\sqrt(k_B T_e)} E_1(y)\]where \(E_1\) is the exponential integral
\[E_1 (y) = - \int_{-y}^∞ \frac{e^{-t}}{t}dt\]and \(y\) is the dimensionless argument
\[y = \frac{1}{2} \frac{ω^2 m_e}{k_{max}^2 k_B T_e}\]where \(k_{max}\) is a maximum wavenumber approximated here as \(k_{max} = 1/λ_B\) where \(λ_B\) is the electron de Broglie wavelength.
- Parameters:
frequencies (
Quantity) – Array of frequencies over which the bremsstrahlung spectrum will be calculated (convertible to Hz).n_e (
Quantity) – Electron number density in the plasma (convertible to m-3).T_e (
Quantity) – Temperature of the electrons (in K or convertible to eV).n_i (
Quantity, optional) – Ion number density in the plasma (convertible to m-3). Defaults to the quasineutral condition \(n_i = n_e / Z\).ion (particle-like, default:
"p+") – An instance ofParticle, or a string convertible toParticle.kmax (
Quantity) – Cutoff wavenumber (convertible to radians per meter). Defaults to the inverse of the electron de Broglie wavelength.
- Returns:
spectrum – Computed bremsstrahlung spectrum over the frequencies provided.
- Return type:
Notes
For details, see Bekefi [1966].
Examples
>>> import astropy.units as u >>> import numpy as np >>> thermal_bremsstrahlung(10**15 * u.Hz, 1e10 * u.cm**-3, 2e7 * u.K) # solar flare <Quantity 8.17560238e-23 kg / (m s2)> >>> thermal_bremsstrahlung( ... 10 ** np.arange(15, 16, 0.1) * u.Hz, 1e22 * u.cm**-3, 1e2 * u.eV ... ) <Quantity [ 79.59052452, 117.73282254, 127.85119908, 127.12505588, 121.01549498, 112.02367743, 101.45553309, 90.04503155, 78.23475796, 66.32227273] kg / (m s2)> >>> thermal_bremsstrahlung( ... 1e17 * u.Hz, 1e16 * u.cm**-3, 1e4 * u.eV, ion="Fe-56 12+" ... ) <Quantity 2.16932808e-10 kg / (m s2)>