Berger, M.J., Seltzer, S.M., Maeda, K.,
Some new results on electron transport in the atmosphere,
Journal of Atmospheric and Terrestrial Physics, v36, i4, pp. 591–617,
April 1974,
doi: 10.1016/0021-9169(74)90085-3
By default, the log(coefficients) are interpolated wrt. log(energy)
and log(zm) using scipy.interpolated.Rbf.
The default “multiquadric” should work fine, if not consider
using “thin-plate” splines.
Parameters:
energy (array_like (M, ...)) – Energy E_0 of the mono-energetic electron beam [keV].
A scalar (0-D) value is promoted to 1-D with one element.
flux (array_like (M,...)) – Energy flux Q_0 of the mono-energetic electron beam [keV / cm² / s¹].
scale_height (array_like (N, ...)) – The atmospheric scale heights [cm].
rho (array_like (N, ...)) – The atmospheric mass density [g / cm³].
ens (array_like (I,), optional) – The energies (one axis) of the coefficient array,
used to interpolate the coefficients to energy.
Defaults to the Berger et al., 1974 coefficients.
zm_p_en (array_like (J,), optional) – The atmospheric depth (the other axis) of the coefficient array,
used to interpolate the coefficients to z = scale_height * rhos.
Defaults to the Berger et al., 1974 coefficients.
coeffs (array_like, (J, I), optional) – The bremsstrahlung energy dissipation coefficients.
Defaults to the Berger et al., 1974 coefficients.
fillna (float or None, optional (default None)) – Value to use for nan values in coeffs, None skips them.
log3 (bool, optional (default True)) – Interpolate the coefficients as log(ens)-log(zm)-log(coeff)
instead of a linear variant.
rbf (str or callable, optional (default "multiquadric")) – Radial basis functions to use for scipy.interpolate.Rbf.
Returns:
a_br – A scalar (0-D) energy is promoted to 1-D, and the result will
have shape (1, N), not (N,).
Energy dissipation rate, units: [keV cm⁻³ s⁻¹]
A. Brekke, J. R. Doupnik, P. M. Banks,
Incoherent scatter measurements of E region conductivities and currents in the auroral zone,
J. Geophys. Res., 79(25), 3773–3790, Sept. 1974,
doi: 10.1029/JA079i025p03773
James F. Vickrey, Richard R. Vondrak, Stephen J. Matthews,
The diurnal and latitudinal variation of auroral zone ionospheric conductivity,
J. Geophys. Res., 86(A1), 65–75, Jan. 1981,
doi: 10.1029/JA086iA01p00065
R. M. Robinson, R. R. Vondrak, K. Miller, T. Dabbs, D. Hardy,
On calculating ionospheric conductances from the flux and energy of precipitating electrons,
J. Geophys. Res. Space Phys., 92(A3), 2565–2569, Mar. 1987,
doi: 10.1029/JA092iA03p02565
scale_height (array_like (N,...)) – The atmospheric scale height(s) [cm].
rho (array_like (N,...)) – The atmospheric densities [g / cm³], corresponding to the scale heights.
pij (array_like (8, 4), optional) – Polynomial coefficents for the electron energy dissipation
per atmospheric depth. Default: None (as given in the reference).
Integrate Fang et al., 2010 over a Maxwellian spectrum
Integrates the mono-energetic parametrization from Fang et al., 2010 [10]
over a Maxwellian spectrum with characteristic energy energy and
total energy flux flux.
Parameters:
energy (float or array_like (M,...)) – Characteristic energy E_0 [keV] of the Maxwellian distribution.
flux (float or array_like (M,...)) – Integrated energy flux Q_0 [keV / cm² / s¹]
scale_height (float or array_like (N,...)) – The atmospheric scale heights [cm].
rho (float or array_like (N,...)) – The atmospheric mass density [g / cm³]
bounds (tuple, optional) – (min, max) [keV] of the integration range to integrate the Maxwellian.
Make sure that this is appropriate to encompass the spectrum.
Default: (0.1, 300.)
nsteps (int, optional) – Number of integration steps, default: 128.
pij (array_like (8, 4), optional) – Polynomial coefficents for the electron energy dissipation
per atmospheric depth. Default: None (as given in the reference).
Atmospheric electron energy dissipation from Fang et al., 2010
Parametrization for mono-energetic electrons [11].
Parameters:
energy (array_like (M,...)) – Energy E_0 of the mono-energetic electron beam [keV].
flux (array_like (M,...)) – Energy flux Q_0 of the mono-energetic electron beam [keV / cm² / s¹].
scale_height (array_like (N,...)) – The atmospheric scale heights [cm].
rho (array_like (N,...)) – The atmospheric mass densities [g / cm³], corresponding to the scale heights.
pij (array_like (8, 4), optional) – Polynomial coefficents for the electron energy dissipation
per atmospheric depth. Default: None (as given in the reference).
Integrates over the mono-energetic parametrization q from [12]
using the given differential particle spectrum phi:
\(\int_\text{spec} \phi(E) q(E, Q) E \text{d}E\)
Parameters:
ens (array_like (M,...)) – Central (bin) energies of the spectrum
dfluxes (array_like (M,...)) – Differential particle fluxes in the given bins
scale_height (array_like (N,...)) – The atmospheric scale heights
rho (array_like (N,...)) – The atmospheric densities, corresponding to the
scale heights.
pij (array_like (8, 4), optional) – Polynomial coefficents for the electron energy dissipation
per atmospheric depth. Default: None (as given in the reference).
axis (int, optional) – The axis to use for integration, default: -1 (last axis).
Fang, X., Lummerzheim, D., and Jackman, C. H. (2013),
Proton impact ionization and a fast calculation method,
J. Geophys. Res. Space Physics, 118, 5369–5378, doi:10.1002/jgra.50484.
Proton ionization parametrization by Fang et al., 2013
Parametrization for mono-energetic protons as described in [15].
Parameters:
energy (array_like (M,...)) – Energy E_0 of the mono-energetic proton beam [keV].
flux (array_like (M,...)) – Energy flux Q_0 of the mono-energetic proton beam [keV / cm² / s¹].
scale_height (array_like (N,...)) – The atmospheric scale heights [cm].
rho (array_like (N,...)) – The atmospheric mass densities [g / cm³], corresponding to the scale heights.
pij (array_like (12, 4), optional) – Polynomial coefficents for the proton energy dissipation
per atmospheric depth. Default: None (as given in the reference).
Integrate mono-energetic parametrization over a spectrum
Integrates the mono-energetic parametrization over a spectrum given by a
functional dependence with characteristic energy energy and total energy
flux flux.
Parameters:
energy (float or array_like (M,...)) – Characteristic energy E_0 [keV] of the spectral distribution.
flux (float or array_like (M,...)) – Integrated energy flux Q_0 [keV / cm² / s¹]
scale_height (float or array_like (N,...)) – The atmospheric scale heights [cm].
rho (float or array_like (N,...)) – The atmospheric mass density [g / cm³]
ediss_func (callable) – Mono-energetic energy dissipation function to integrate.
ediss_kws (dict-like, optional) – Optional keyword arguments to pass to the mono-energetic
energy dissipation function. Default: None
bounds (tuple, optional) – (min, max) [keV] of the integration range to integrate the Maxwellian.
Make sure that this is appropriate to encompass the spectrum.
Default: (0.1, 300.)
nsteps (int, optional) – Number of integration steps, default: 128.
Standard exponential distribution with
\(\lambda\) = 1 / en_0 or \(\beta\) = en_0.
normalized to unit number flux, i.e.
\(\int_0^\infty \phi(E) \text{d}E = 1\).
\[\phi(E, E_0) = E / E_0^2 \cdot \exp\{-E / E_0\}\]
Equal to a standard Gamma distribution with
\(\alpha\) = 2 and \(\beta\) = 1 / en_0,
or
\(k\) = 2 and \(\theta\) = en_0.
Normalized to \(\int_0^\infty \phi(E) \text{d}E = 1\).
\[\phi(E, E_0) = E / 2E_0^3 \cdot \exp\{-E / E_0\}\]
Equal to a standard Gamma distribution with
\(\alpha\) = 3 and \(\beta\) = 1 / en_0,
or
\(k\) = 3 and \(\theta\) = en_0.
The total precipitating energy flux is fixed to 1 keV cm-2 s-1,
multiply by Q_0 [keV cm-2 s-1] to scale the particle flux.
Normalized to \(\int_0^\infty \phi(E) E \text{d}E = 1\).
Scales to arbitrary energy flux \(Q\) via multiplication:
\(\tilde\phi = Q \cdot \phi\).
The minus-sign (-) and is used for the high-energy tail variant,
and the plus-sign (+) for the low-energy tail variant.
The exponent gamma needs to be set appropriately,
< -1 for het, and > 1 for let.
Normalized to \(\int_{E_0}^\infty \phi(E) E \text{d}E = 1\)
for the high-energy tail version, and to
\(\int_0^{E_0} \phi(E) E \text{d}E = 1\)
for the low-energy tail version.
Scales to arbitrary energy flux \(Q\) via multiplication:
\(\tilde\phi = Q \cdot \phi\).
en_0 (float, optional) – Characteristic energy in [keV], i.e. mode of the distribution.
Default: 10 keV
gamma (float, optional) – Exponent of the power-law distribution, in [keV].
het (bool, optional (default True)) – Return a high-energy tail (true) for en > en_0,
or low-energy tail (false) for en < en_0.
Adjusts the normalization accordingly.
Returns:
phi – Hemispherical differential particle flux at en in [keV-1 cm-2 s-1]
([keV-2] scaled by unit energy flux).
The minus-sign (-) and is used for the high-energy tail variant,
and the plus-sign (+) for the low-energy tail variant.
The exponent gamma needs to be set appropriately,
< -1 for het, and > 1 for let.
The “high-energy tail” version (het = True)
resembles a Pareto distribution with
scale parameter \(x_m\) = en_0
and shape parameter \(\alpha\) = -(gamma + 1).
Adapted from Strickland et al., 1993 [26] and
normalized to unit particle flux:
\(\int_{E_0}^\infty \phi(E) \text{d}E = 1\)
for the high-energy tail version, and
\(\int_0^{E_0} \phi(E) \text{d}E = 1\)
for the low-energy tail version.
en_0 (float, optional) – Characteristic energy in [keV], i.e. mode of the distribution.
Default: 10 keV
gamma (float, optional) – Exponent of the power-law distribution, in [keV].
het (bool, optional) – Return a high-energy tail (het, default: true) for en > en_0,
or low-energy tail (false) for en < en_0.
Adjusts the normalization accordingly.
Returns:
phi – Normalized differential hemispherical number flux at en in [keV-1 cm-2 s-1]
([keV] or scaled by 1 keV-2 cm-2 s-1, e.g.).
Atmospheric ionization from auroral particle precipitation
Bundles some of the parametrizations for middle and upper atmospheric
ionization and recombination rates for precipitating
auroral (100 eV–30 keV) and radiation-belt (30 keV–1 MeV) electrons.