continuum¶
-
class
ChiantiPy.core.
continuum
(ion_string, temperature, abundance=None, em=None)[source] [edit on github]¶ Bases:
object
The top level class for continuum calculations. Includes methods for the calculation of the free-free and free-bound continua.
Parameters: ionStr :
str
CHIANTI notation for the given ion, e.g. ‘fe_12’ that corresponds to the Fe XII ion.
temperature : array-like
In units of Kelvin
abundance :
float
orstr
, optionalElemental abundance relative to Hydrogen or name of CHIANTI abundance file, without the ‘.abund’ suffix, e.g. ‘sun_photospheric_1998_grevesse’.
em : array-like, optional
Line-of-sight emission measure (\(\int\mathrm{d}l\,n_en_H\)), in units of \(\mathrm{cm}^{-5}\), or the volumetric emission measure (\(\int\mathrm{d}V\,n_en_H\)) in units of \(\mathrm{cm}^{-3}\).
Notes
The methods for calculating the free-free and free-bound emission and losses return their result to an attribute. See the respective docstrings for more information.
Examples
>>> import ChiantiPy.core as ch >>> import numpy as np >>> temperature = np.logspace(4,9,20) >>> cont = ch.continuum('fe_15',temperature) >>> wavelength = np.arange(1,1000,10) >>> cont.freeFree(wavelength) >>> cont.freeBound(wavelength, include_abundance=True, include_ioneq=False) >>> cont.calculate_free_free_loss() >>> cont.calculate_free_bound_loss()
Methods Summary
calculate_free_bound_loss
(**kwargs)Calculate the free-bound energy loss rate of an ion. calculate_free_free_loss
(**kwargs)Calculate the free-free energy loss rate of an ion. freeBound
(wavelength[, include_abundance, …])Calculate the free-bound emission of an ion. freeBoundLoss
(**kwargs)Calculate the free-bound energy loss rate of an ion. freeFree
(wavelength[, include_abundance, …])Calculates the free-free emission for a single ion. freeFreeLoss
(**kwargs)Calculate the free-free energy loss rate of an ion. ioneqOne
()Provide the ionization equilibrium for the selected ion as a function of temperature. ioneq_one
(stage, **kwargs)Calculate the equilibrium fractional ionization of the ion as a function of temperature. itoh_gaunt_factor
(wavelength)Calculates the free-free gaunt factors of [R15]. karzas_cross_section
(photon_energy, …)Calculate the photoionization cross-sections using the Gaunt factors of [R16]. klgfbInterp
(wvl, n, l)A Python version of the CHIANTI IDL procedure karzas_xs. mewe_gaunt_factor
(**kwargs)Calculate the Gaunt factor according to [R18] for a single ion \(Z_z\). sutherland_gaunt_factor
(wavelength)Calculates the free-free gaunt factor calculations of [R19]. verner_cross_section
(photon_energy)Calculates the photoionization cross-section using data from [R20] for transitions to the ground state. Methods Documentation
-
calculate_free_bound_loss
(**kwargs)[source] [edit on github]¶ Calculate the free-bound energy loss rate of an ion. The result is returned to the
free_bound_loss
attribute.The free-bound loss rate can be calculated by integrating the free-bound emission over the wavelength. This is difficult using the expression in
calculate_free_bound_emission
so we instead use the approach of [R21] and [R22]. Eq. 1a of [R22] can be integrated over wavelength to get the free-bound loss rate,\[\frac{dW}{dtdV} = C_{ff}\frac{k}{hc}T^{1/2}G_{fb},\]in units of erg \(\mathrm{cm}^3\,\mathrm{s}^{-1}\) where \(G_{fb}\) is the free-bound Gaunt factor as given by Eq. 15 of [R22] (see
mewe_gaunt_factor
for more details) and \(C_{ff}\) is the numerical constant as given in Eq. 4 of [R21] and can be written in terms of the fine structure constant \(\alpha\),\[C_{ff}\frac{k}{hc} = \frac{8}{3}\left(\frac{\pi}{6}\right)^{1/2}\frac{h^2\alpha^3}{\pi^2}\frac{k_B}{m_e^{3/2}} \approx 1.43\times10^{-27}\]References
[R21] (1, 2, 3) Gronenschild, E.H.B.M. and Mewe, R., 1978, A&AS, 32, 283 [R22] (1, 2, 3, 4) Mewe, R. et al., 1986, A&AS, 65, 511
-
calculate_free_free_loss
(**kwargs)[source] [edit on github]¶ Calculate the free-free energy loss rate of an ion. The result is returned to the
free_free_loss
attribute.The free-free radiative loss rate is given by Eq. 5.15a of [R23]. Writing the numerical constant in terms of the fine structure constant \(\alpha\),
\[\frac{dW}{dtdV} = \frac{4\alpha^3h^2}{3\pi^2m_e}\left(\frac{2\pi k_B}{3m_e}\right)^{1/2}Z^2T^{1/2}\bar{g}_B\]where where \(Z\) is the nuclear charge, \(T\) is the electron temperature, and \(\bar{g}_{B}\) is the wavelength-averaged and velocity-averaged Gaunt factor. The Gaunt factor is calculated using the methods of [R24]. Note that this expression for the loss rate is just the integral over wavelength of Eq. 5.14a of [R23], the free-free emission, and is expressed in units of erg \(\mathrm{cm}^3\,\mathrm{s}^{-1}\).
References
[R23] (1, 2, 3) Rybicki and Lightman, 1979, Radiative Processes in Astrophysics, (Wiley-VCH) [R24] (1, 2) Karzas and Latter, 1961, ApJSS, 6, 167
-
freeBound
(wavelength, include_abundance=True, include_ioneq=True, use_verner=True, **kwargs)[source] [edit on github]¶ Calculate the free-bound emission of an ion. The result is returned as a 2D array to the
free_bound_emission
attribute.The total free-bound continuum emissivity is given by,
\[\frac{dW}{dtdVd\lambda} = \frac{1}{4\pi}\frac{2}{hk_Bc^3m_e\sqrt{2\pi k_Bm_e}}\frac{E^5}{T^{3/2}}\sum_i\frac{\omega_i}{\omega_0}\sigma_i^{bf}\exp\left(-\frac{E - I_i}{k_BT}\right)\]where \(E=hc/\lambda\) is the photon energy, \(\omega_i\) and \(\omega_0\) are the statistical weights of the \(i^{\mathrm{th}}\) level of the recombined ion and the ground level of the recombing ion, respectively, \(\sigma_i^{bf}\) is the photoionization cross-section, and \(I_i\) is the ionization potential of level \(i\). This expression comes from Eq. 12 of [R27]. For more information about the free-bound continuum calculation, see Peter Young’s notes on free-bound continuum.
The photoionization cross-sections are calculated using the methods of [R26] for the transitions to the ground state and [R25] for all other transitions. See
verner_cross_section
andkarzas_cross_section
for more details.The free-bound emission is in units of erg \(\mathrm{cm}^3\mathrm{s}^{-1}\mathrm{\mathring{A}}^{-1}\mathrm{str}^{-1}\). If the emission measure has been set, the units will be multiplied by \(\mathrm{cm}^{-5}\) or \(\mathrm{cm}^{-3}\), depending on whether it is the line-of-sight or volumetric emission measure, respectively.
Parameters: wavelength : array-like
In units of angstroms
include_abundance :
bool
, optionalIf True, include the ion abundance in the final output.
include_ioneq :
bool
, optionalIf True, include the ionization equilibrium in the final output
use_verner :
bool
, optionalIf True, cross-sections of ground-state transitions using [R26], i.e.
verner_cross_section
Raises: ValueError
If no .fblvl file is available for this ion
References
[R25] (1, 2) Karzas and Latter, 1961, ApJSS, 6, 167 [R26] (1, 2, 3) Verner & Yakovlev, 1995, A&AS, 109, 125 [R27] (1, 2) Young et al., 2003, ApJSS, 144, 135
-
freeBoundLoss
(**kwargs)[source] [edit on github]¶ Calculate the free-bound energy loss rate of an ion. The result is returned to the
free_bound_loss
attribute.The free-bound loss rate can be calculated by integrating the free-bound emission over the wavelength. This is difficult using the expression in
calculate_free_bound_emission
so we instead use the approach of [R28] and [R29]. Eq. 1a of [R29] can be integrated over wavelength to get the free-bound loss rate,\[\frac{dW}{dtdV} = C_{ff}\frac{k}{hc}T^{1/2}G_{fb},\]in units of erg \(\mathrm{cm}^3\,\mathrm{s}^{-1}\) where \(G_{fb}\) is the free-bound Gaunt factor as given by Eq. 15 of [R29] (see
mewe_gaunt_factor
for more details) and \(C_{ff}\) is the numerical constant as given in Eq. 4 of [R28] and can be written in terms of the fine structure constant \(\alpha\),\[C_{ff}\frac{k}{hc} = \frac{8}{3}\left(\frac{\pi}{6}\right)^{1/2}\frac{h^2\alpha^3}{\pi^2}\frac{k_B}{m_e^{3/2}} \approx 1.43\times10^{-27}\]References
[R28] (1, 2, 3) Gronenschild, E.H.B.M. and Mewe, R., 1978, A&AS, 32, 283 [R29] (1, 2, 3, 4) Mewe, R. et al., 1986, A&AS, 65, 511
-
freeFree
(wavelength, include_abundance=True, include_ioneq=True, **kwargs)[source] [edit on github]¶ Calculates the free-free emission for a single ion. The result is returned as a dict to the
FreeFree
attribute. The dict has the keywordsintensity
,wvl
,temperature
,em
.The free-free emission for the given ion is calculated according Eq. 5.14a of [R30], substituting \(\nu=c/\lambda\), dividing by the solid angle, and writing the numerical constant in terms of the fine structure constant \(\alpha\),
\[\frac{dW}{dtdVd\lambda} = \frac{c}{3m_e}\left(\frac{\alpha h}{\pi}\right)^3\left(\frac{2\pi}{3m_ek_B}\right)^{1/2}\frac{Z^2}{\lambda^2T^{1/2}}\exp{\left(-\frac{hc}{\lambda k_BT}\right)}\bar{g}_{ff},\]where \(Z\) is the nuclear charge, \(T\) is the electron temperature in K, and \(\bar{g}_{ff}\) is the velocity-averaged Gaunt factor. The Gaunt factor is estimated using the methods of [R31] and [R32], depending on the temperature and energy regime. See
itoh_gaunt_factor
andsutherland_gaunt_factor
for more details.The free-free emission is in units of erg \(\mathrm{cm}^3\mathrm{s}^{-1}\mathrm{\mathring{A}}^{-1}\mathrm{str}^{-1}\). If the emission measure has been set, the units will be multiplied by \(\mathrm{cm}^{-5}\) or \(\mathrm{cm}^{-3}\), depending on whether it is the line-of-sight or volumetric emission measure, respectively.
Parameters: wavelength : array-like
In units of angstroms
include_abundance :
bool
, optionalIf True, include the ion abundance in the final output.
include_ioneq :
bool
, optionalIf True, include the ionization equilibrium in the final output
References
[R30] (1, 2) Rybicki and Lightman, 1979, Radiative Processes in Astrophysics, (Wiley-VCH) [R31] (1, 2) Itoh, N. et al., 2000, ApJS, 128, 125 [R32] (1, 2) Sutherland, R. S., 1998, MNRAS, 300, 321
-
freeFreeLoss
(**kwargs)[source] [edit on github]¶ Calculate the free-free energy loss rate of an ion. The result is returned to the
free_free_loss
attribute.The free-free radiative loss rate is given by Eq. 5.15a of [R33]. Writing the numerical constant in terms of the fine structure constant \(\alpha\),
\[\frac{dW}{dtdV} = \frac{4\alpha^3h^2}{3\pi^2m_e}\left(\frac{2\pi k_B}{3m_e}\right)^{1/2}Z^2T^{1/2}\bar{g}_B\]where where \(Z\) is the nuclear charge, \(T\) is the electron temperature, and \(\bar{g}_{B}\) is the wavelength-averaged and velocity-averaged Gaunt factor. The Gaunt factor is calculated using the methods of [R34]. Note that this expression for the loss rate is just the integral over wavelength of Eq. 5.14a of [R33], the free-free emission, and is expressed in units of erg \(\mathrm{cm}^3\,\mathrm{s}^{-1}\).
References
[R33] (1, 2, 3) Rybicki and Lightman, 1979, Radiative Processes in Astrophysics, (Wiley-VCH) [R34] (1, 2) Karzas and Latter, 1961, ApJSS, 6, 167
-
ioneqOne
()[source] [edit on github]¶ Provide the ionization equilibrium for the selected ion as a function of temperature. Similar to but not identical to ion.ioneqOne() returned in self.IoneqOne
-
ioneq_one
(stage, **kwargs)[source] [edit on github]¶ Calculate the equilibrium fractional ionization of the ion as a function of temperature.
Uses the
ChiantiPy.core.ioneq
module and does a first-order spline interpolation to the data. An ionization equilibrium file can be passed as a keyword argument,ioneqfile
. This can be passed through as a keyword argument to any of the functions that uses the ionization equilibrium.Parameters: stage : int
Ionization stage, e.g. 25 for Fe XXV
-
itoh_gaunt_factor
(wavelength)[source] [edit on github]¶ Calculates the free-free gaunt factors of [R35].
An analytic fitting formulae for the relativistic Gaunt factor is given by Eq. 4 of [R35],
\[g_{Z} = \sum^{10}_{i,j=0}a_{ij}t^iU^j\]where,
\[t = \frac{1}{1.25}(\log_{10}{T} - 7.25),\ U = \frac{1}{2.5}(\log_{10}{u} + 1.5),\]\(u=hc/\lambda k_BT\), and \(a_{ij}\) are the fitting coefficients and are read in using
ChiantiPy.tools.io.itohRead
and are given in Table 4 of [R35]. These values are valid for \(6<\log_{10}(T)< 8.5\) and \(-4<\log_{10}(u)<1\).See also
ChiantiPy.tools.io.itohRead
- Read in Gaunt factor coefficients from [R35]
References
[R35] (1, 2, 3, 4, 5) Itoh, N. et al., 2000, ApJS, 128, 125
-
karzas_cross_section
(photon_energy, ionization_potential, n, l)[source] [edit on github]¶ Calculate the photoionization cross-sections using the Gaunt factors of [R36].
The free-bound photoionization cross-section is given by,
\[\sigma_i^{bf} = 1.077294\times8065.54\times10^{16}\left(\frac{I_i}{hc}\right)^2\left(\frac{hc}{E}\right)^3\frac{g_{bf}}{n_i},\]where \(I_i\) is the ionization potential of the \(i^{\mathrm{th}}\) level, \(E\) is the photon energy, \(g_{bf}\) is the Gaunt factor calculated according to [R36], and \(n_i\) is the principal quantum number of the \(i^{\mathrm{th}}\) level. \(\sigma_i^{bf}\) is units of \(\mathrm{cm}^{2}\). This expression is given by Eq. 13 of [R37]. For more information on the photoionization cross-sections, see Peter Young’s notes on free-bound continuum.
Parameters: photon_energy : array-like
ionization_potential :
float
n :
int
l :
int
References
[R36] (1, 2, 3) Karzas and Latter, 1961, ApJSS, 6, 167 [R37] (1, 2) Young et al., 2003, ApJSS, 144, 135
-
klgfbInterp
(wvl, n, l)[source] [edit on github]¶ A Python version of the CHIANTI IDL procedure karzas_xs.
Interpolates free-bound gaunt factor of Karzas and Latter, (1961, Astrophysical Journal Supplement Series, 6, 167) as a function of wavelength (wvl).
-
mewe_gaunt_factor
(**kwargs)[source] [edit on github]¶ Calculate the Gaunt factor according to [R38] for a single ion \(Z_z\).
Using Eq. 9 of [R38], the free-bound Gaunt factor for a single ion can be written as,
\[G_{fb}^{Z,z} = \frac{E_H}{k_BT}\mathrm{Ab}(Z)\frac{N(Z,z)}{N(Z)}f(Z,z,n)\]where \(E_H\) is the ground-state potential of H, \(\mathrm{Ab}(Z)\) is the elemental abundance, \(\frac{N(Z,z)}{N(Z)}\) is the fractional ionization, and \(f(Z,z,n)\) is given by Eq. 10 and is approximated by Eq 16 as,
\[f(Z,z,n) \approx f_2(Z,z,n_0) = 0.9\frac{\zeta_0z_0^4}{n_0^5}\exp{\left(\frac{E_Hz_0^2}{n_0^2k_BT}\right)} + 0.42\frac{z^4}{n_0^{3/2}}\exp{\left(\frac{E_Hz^2}{(n_0 + 1)^2k_BT}\right)}\]where \(n_0\) is the principal quantum number, \(z_0\) is the effective charge (see Eq. 7 of [R38]), and \(\zeta_0\) is the number of vacancies in the 0th shell and is given in Table 1 of [R38]. Here it is calculated in the same manner as in fb_rad_loss.pro of the CHIANTI IDL library. Note that in the expression for \(G_{fb}\), we have not included the \(N_H/n_e\) factor.
Raises: ValueError
If no .fblvl file is available for this ion
References
[R38] (1, 2, 3, 4, 5) Mewe, R. et al., 1986, A&AS, 65, 511
-
sutherland_gaunt_factor
(wavelength)[source] [edit on github]¶ Calculates the free-free gaunt factor calculations of [R39].
The Gaunt factors of [R39] are read in using
ChiantiPy.tools.io.gffRead
as a function of \(u\) and \(\gamma^2\). The data are interpolated to the appropriate wavelength and temperature values usingmap_coordinates
.References
[R39] (1, 2, 3) Sutherland, R. S., 1998, MNRAS, 300, 321
-
verner_cross_section
(photon_energy)[source] [edit on github]¶ Calculates the photoionization cross-section using data from [R40] for transitions to the ground state.
The photoionization cross-section can be expressed as \(\sigma_i^{fb}=F(E/E_0)\) where \(F\) is an analytic fitting formula given by Eq. 1 of [R40],
\[F(y) = ((y-1)^2 + y_w^2)y^{-Q}(1 + \sqrt{y/y_a})^{-P},\]where \(E\) is the photon energy, \(n\) is the principal quantum number, \(l\) is the orbital quantum number, \(Q = 5.5 + l - 0.5P\), and \(\sigma_0,E_0,y_w,y_a,P\) are fitting paramters. These can be read in using
ChiantiPy.tools.io.vernerRead
.References
[R40] (1, 2, 3) Verner & Yakovlev, 1995, A&AS, 109, 125
-