fiasco currently only supports version 8 of the CHIANTI database.

Ion#

class fiasco.Ion(ion_name, temperature: Unit('K'), abundance='sun_coronal_1992_feldman_ext', *args, **kwargs)[source]#

Bases: IonBase, ContinuumBase

Class for representing a CHIANTI ion.

The ion object is the fundamental unit of fiasco. This object contains all of the properties and methods needed to access important information about each ion from the CHIANTI database as well as compute common derived quantities.

Parameters:
  • ion_name (str or tuple) – Name of the ion. This can be either a string denoting the name or a tuple containing the atomic number and ionization stage. See parse_ion for a list of all possible input formats.

  • temperature (Quantity) – Temperature array over which to evaluate temperature dependent quantities.

  • ioneq_filename (str, optional) – Ionization equilibrium dataset

  • abundance (str or float, optional) – If a string is provided, use the appropriate abundance dataset. If a float is provided, use that value as the abundance.

  • ip_filename (str, optional) – Ionization potential dataset

Attributes Summary

abundance

Elemental abundance relative to H.

dielectronic_recombination_rate

Dielectronic recombination rate as a function of temperature.

direct_ionization_rate

Total ionization rate due to collisions as a function of temperature.

effective_collision_strength

Maxwellian-averaged collision strength, typically denoted by \(\Upsilon\), as a function of temperature.

electron_collision_deexcitation_rate

Collisional de-excitation rate coefficient for electrons.

electron_collision_excitation_rate

Collisional excitation rate coefficient for electrons.

excitation_autoionization_rate

Ionization rate due to excitation autoionization.

formation_temperature

Temperature at which ioneq is maximum.

helium_like

Is the ion helium like or not.

hydrogenic

Is the ion hydrogen-like or not.

ioneq

Ionization equilibrium data interpolated to the given temperature

ionization_rate

Total ionization rate as a function of temperature.

ip

Ionization potential.

proton_collision_deexcitation_rate

Collisional de-excitation rate coefficient for protons.

proton_collision_excitation_rate

Collisional excitation rate coefficient for protons.

radiative_recombination_rate

Radiative recombination rate as a function of temperature.

recombination_rate

Total recombination rate as a function of temperature.

transitions

Methods Summary

contribution_function(density, **kwargs)

Contribution function \(G(n_e,T)\) for all transitions.

direct_ionization_cross_section(energy)

Direct ionization cross-section as a function of energy.

emissivity(density, **kwargs)

Emissivity as a function of temperature and density for all transitions.

free_bound(wavelength[, use_verner])

Free-bound continuum emission of the recombined ion.

free_free(wavelength)

Free-free continuum emission as a function of temperature and wavelength.

gaunt_factor_free_free(wavelength)

Free-free Gaunt factor as a function of wavelength.

intensity(density, emission_measure, **kwargs)

Line-of-sight intensity computed assuming a particular column emission measure.

level_populations(density[, ...])

Energy level populations as a function of temperature and density.

next_ion()

Return an Ion instance with the next highest ionization stage.

previous_ion()

Return an Ion instance with the next lowest ionization stage.

spectrum(*args, **kwargs)

Construct the spectrum using a given filter over a specified wavelength range.

Attributes Documentation

abundance#

Elemental abundance relative to H.

dielectronic_recombination_rate#

Dielectronic recombination rate as a function of temperature.

The dielectronic recombination rate, as a function of \(T\), is computed using one of two methods. The methodology used depends on the type of dielectronic recombination rate fitting coefficients available for the particular ion in the CHIANTI atomic database.

The first method is given in Eq. 3 of Zatsarinny et al. [ZGK+03],

\[\alpha_{DR} = T^{-3/2}\sum_ic_ie^{-E_i/T}\]

where \(c_i\) and \(E_i\) are fitting coefficients stored in the CHIANTI database.

The second method is given by Eq. 5 of Shull and van Steenberg [SvanSteenberg82],

\[\alpha_{DR} = A T^{-3/2}e^{-T_0/T}(1 + B e^{-T_1/T})\]

where \(A,B,T_0,T_1\) are fitting coefficients stored in the CHIANTI database.

direct_ionization_rate#

Total ionization rate due to collisions as a function of temperature.

The ionization rate due to the collisions with free electrons can be written as the integral of the velocity-weighted collisional cross-section over the Maxwell-Boltzmann distribution. Following Section 3.5.1 of Del Zanna and Mason [DZM18], this can be written as,

\[C^I = \sqrt{\frac{8}{\pi m_e}}(k_BT)^{-3/2}\int_I^{\infty}\mathrm{d}E\,E\sigma_I(E)\exp{\left(-\frac{E}{k_BT}\right)}\]

where \(E\) is the energy of the incident electron, \(I\) is the ionization energy of the initially bound electron, and \(\sigma_I\) is the ionization cross-section.

Making the substitution \(x=(E-I)/k_BT\), the above integral can be rewritten as,

\[\begin{split}\begin{aligned} C^I = \sqrt{\frac{8k_BT}{\pi m_e}}\exp{\left(-\frac{I}{k_BT}\right)}&\left(\int_0^{\infty}\mathrm{d}x\,x\sigma_{I}(k_BTx+I)e^{-x} \right. \\ &\left. + \frac{I}{k_BT}\int_0^{\infty}\mathrm{d}x\,\sigma_{I}(k_BTx+I)e^{-x}\right). \end{aligned}\end{split}\]

Each of these integrals is of the form such that they can be evaluated using Gauss-Laguerre quadrature. Note that there is a typo in the expression for the ionization rate integral in Eq. 32 of Del Zanna and Mason [DZM18]. The details of the ionization cross-section calculation can be found in direct_ionization_cross_section.

See also

direct_ionization_cross_section

Calculation of \(\sigma_I\) as a function of \(E\).

effective_collision_strength#

Maxwellian-averaged collision strength, typically denoted by \(\Upsilon\), as a function of temperature.

According to Eq. 4.11 of Phillips et al. [PFL08], \(\Upsilon\) is given by,

\[\Upsilon = \int_0^\infty\mathrm{d}\left(\frac{E}{k_BT}\right)\,\Omega_{ji}\exp{\left(-\frac{E}{k_BT}\right)}\]

where \(\Omega_{ji}\) is the collision strength. These Maxwellian-averaged collision strengths are stored in dimensionless form in CHIANTI and are rescaled to the appropriate temperature.

See also

fiasco.util.burgess_tully_descale

Descale and interpolate \(\Upsilon\).

electron_collision_deexcitation_rate#

Collisional de-excitation rate coefficient for electrons.

According to Eq. 4.12 of Phillips et al. [PFL08], the rate coefficient for collisional de-excitation is given by,

\[C^d_{ji} = I_Ha_0^2\sqrt{\frac{8\pi}{mk_B}}\frac{\Upsilon}{\omega_jT^{1/2}},\]

where \(j,i\) are the upper and lower level indices, respectively, \(I_H\) is the ionization potential for H, \(a_0\) is the Bohr radius, \(\Upsilon\) is the effective collision strength, and \(\omega_j\) is the statistical weight of the level \(j\).

See also

electron_collision_excitation_rate

Excitation rate due to collisions

effective_collision_strength

Maxwellian-averaged collision strength, \(\Upsilon\)

electron_collision_excitation_rate#

Collisional excitation rate coefficient for electrons.

The rate coefficient for collisional excitation is given by,

\[C^e_{ij} = \frac{\omega_j}{\omega_i}C^d_{ji}\exp{\left(-\frac{k_BT_e}{\Delta E_{ij}}\right)}\]

where \(j,i\) are the upper and lower level indices, respectively, \(\omega_j,\omega_i\) are the statistical weights of the upper and lower levels, respectively, and \(\Delta E_{ij}\) is the energy of the transition [PFL08].

Parameters:

deexcitation_rate (Quantity, optional) – Optionally specify deexcitation rate to speedup calculation

See also

electron_collision_deexcitation_rate

De-excitation rate due to collisions

excitation_autoionization_rate#

Ionization rate due to excitation autoionization.

Following Eq. 4.74 of Phillips et al. [PFL08], the excitation autoionization rate is given by,

\[\alpha_{EA} = \frac{h^2}{(2\pi m_e)^{3/2}}(k_BT)^{-1/2}\sum_{lj}\Upsilon^{EA}_{lj}\exp{\left(-\frac{\Delta E_{lj}}{k_BT}\right)}\]

where \(\Upsilon^{EA}\) is the thermally-averaged excitation autoionization cross-section as stored in CHIANTI and includes the additional \(\omega_j\) multiplicity factor compared to the expression in Phillips et al. [PFL08]. The sum is taken over inelastic collisions to level \(j\) from a level \(l\) below the ionization threshold. Additionally, note that the constant has been rewritten in terms of \(h\) rather than \(I_H\) and \(a_0\).

formation_temperature#

Temperature at which ioneq is maximum. This is a useful proxy for the temperature at which lines for this ion are formed.

helium_like#

Is the ion helium like or not.

Notes

This is True if \(Z - z = 2\).

hydrogenic#

Is the ion hydrogen-like or not.

Notes

This is True if \(Z - z = 1\).

ioneq#

Ionization equilibrium data interpolated to the given temperature

Interpolated the pre-computed ionization fractions stored in CHIANTI to the temperature of the ion. Returns NaN where interpolation is out of range of the data. For computing ionization equilibrium outside of this temperature range, it is better to use the ionization and recombination rates.

Note

The cubic interpolation is performed in log-log space using a Piecewise Cubic Hermite Interpolating Polynomial with PchipInterpolator. This helps to ensure smoothness while reducing oscillations in the interpolated ionization fractions.

ionization_rate#

Total ionization rate as a function of temperature.

The total ionization rate, as a function of temperature, for a given ion is the sum of the direct ionization and excitation autoionization rates such that,

\[\alpha_{I} = \alpha_{DI} + \alpha_{EA}\]
ip#

Ionization potential.

proton_collision_deexcitation_rate#

Collisional de-excitation rate coefficient for protons.

As in the electron case, the proton collision de-excitation rate is given by,

\[C^{d,p}_{ji} = \frac{\omega_i}{\omega_j}\exp{\left(\frac{E}{k_BT}\right)}C^{e,p}_{ij}\]

where \(C^{e,p}_{ji}\) is the excitation rate due to collisions with protons.

Note that \(T\) is technically the proton temperature. In the case of a thermal plasma, the electron and proton temperatures are equal, \(T_e=T_p\). See Section 4.9.4 of Phillips et al. [PFL08] for additional information on proton collision rates.

See also

proton_collision_excitation_rate

Excitation rate due to collisions with protons

proton_collision_excitation_rate#

Collisional excitation rate coefficient for protons.

These excitation rates are stored in CHIANTI and then rescaled to the appropriate temperatures using the method of Burgess and Tully [BT92].

radiative_recombination_rate#

Radiative recombination rate as a function of temperature.

The recombination rate due to interaction with the ambient radiation field is calculated using a set of fit parameters using one of two methods. The methodology used depends on the type of radiative recombination rate fitting coefficients available for the particular ion in the CHIANTI atomic database.

The first method is given in Eq. 4 of Verner and Ferland [VF96] and Eq. 1 of Badnell [Bad06],

\[\alpha_{RR} = A(\sqrt{T/T_0}(1 + \sqrt{T/T_0})^{1-B}(1 + \sqrt{T/T_1})^{1+B})^{-1}\]

where \(A,B,T_0,T_1\) are fitting coefficients provided for each ion in the CHIANTI atomic database. In some cases, the fitting coefficient \(B\) is also modified as,

\[B \to B + Ce^{-T_2/T}\]

where \(C\) and \(T_2\) are additional fitting coefficients (see Eq. 2 of Badnell [Bad06]).

The second method is given by Eq. 4 of Shull and van Steenberg [SvanSteenberg82] and Eq. 1 of Verner and Ferland [VF96],

\[\alpha_{RR} = A(T/T_0)^{-\eta}\]

where \(A\) and \(\eta\) are fitting parameters provided in the CHIANTI atomic database and \(T_0=10^4\) K.

recombination_rate#

Total recombination rate as a function of temperature.

The total recombination rate, as a function of temperature, for a given ion is the sum of the radiative and dielectronic recombination rates such that,

\[\alpha_{R} = \alpha_{RR} + \alpha_{DR}\]

Important

For most ions, this total recombination rate is computed by summing the outputs of the radiative_recombination_rate and dielectronic_recombination_rate methods. However, for some ions, total recombination rate data is available in the so-called .trparams files. For these ions, the output of this method will not be equal to the sum of the dielectronic_recombination_rate and radiative_recombination_rate method. As such, when computing the total recombination rate, this method should always be used.

transitions#

Methods Documentation

contribution_function(density: Unit('1 / cm3'), **kwargs)[source]#

Contribution function \(G(n_e,T)\) for all transitions.

The contribution function for ion \(k\) of element \(X\) for a particular transition \(ij\) is given by,

\[G_{ij} = mathrm{Ab}(X)f_{X,k}N_jA_{ij}\Delta E_{ij}\frac{1}{n_e},\]

Note that the contribution function is often defined in differing ways by different authors. The contribution function is defined as above in Young et al. [YDL+16].

The corresponding wavelengths can be retrieved with,

ion.transitions.wavelength[~ion.transitions.is_twophoton]

Important

The ratio \(n_H/n_e\), which is often approximated as \(n_H/n_e\approx0.83\), is explicitly not included here. This means that when computing an intensity with the result of this function, the accompanying emission measure is \(\mathrm{EM}=\mathrm{d}hn_Hn_e\) rather than \(n_e^2\).

Parameters:
  • density (Quantity) – Electron number density

  • couple_density_to_temperature (bool, optional) – If True, the density will vary along the same axis as temperature in the computed level populations. The number of densities must be the same as the number of temperatures. This is useful, for example, when computing the level populations at constant pressure and is also much faster than computing the level populations along an independent density axis. By default, this is set to False.

Returns:

g (Quantity) – A (l, m, k) shaped quantity, where l is the number of temperatures, m is the number of densities, and k is the number of transitions corresponding to the transition wavelengths described above. If couple_density_to_temperature=True, then m=1 and l represents the number of temperatures and densities.

direct_ionization_cross_section(energy: Unit('erg'))[source]#

Direct ionization cross-section as a function of energy.

The direction ionization cross-section is calculated one of two ways. See Dere [Der07], Sections 3.1 and 3.2 for details. For H-like ions with \(Z\ge6\) and He-like ions with \(Z\ge10\), the cross-section is computed according to the method of Fontes et al. [FSZ99],

\[\sigma_I = B\frac{\pi a_0^2}{I^2}Q_R\]

where \(B=1\) for H-like ions and \(B=2\) for He-like ions, \(I\) is the ionization energy (expressed in Ry), \(a_0\) is the Bohr radius, and \(Q_R\) is a reduced cross-section which can be approximated by the fitting formula given in Eqs. 2.10, 2.11, and 2.12 of Fontes et al. [FSZ99].

For all other ions, the cross-section is computed according to the method of Dere [Der07] which employs a scaling similar to that used by Burgess and Tully [BT92]. Rearranging Eq. 3 of Dere [Der07],

\[\sigma_I = \frac{\Sigma (\log{u} + 1)}{uI^2}\]

where \(u=E/I\) is the energy of the incident electron scaled by ionization potential and \(\Sigma\) is the scaled cross-section which is defined over,

\[U = 1 - \frac{\log{f}}{\log{u - 1 + f}}\]

where \(f\) is a fitting parameter. \(U,f,\Sigma\) are all stored in the CHIANTI database such that \(\sigma_I\) can be computed for a given \(E\).

emissivity(density: Unit('1 / cm3'), **kwargs)[source]#

Emissivity as a function of temperature and density for all transitions.

The emissivity is given by the expression,

\[\epsilon(n_e,T) = G(n_e,T)n_Hn_e\]

where \(G\) is the contribution function, \(n_H\) is the H (or proton) density, \(n_e\) is the electron density, and \(T\) is the temperature. Note that, like the contribution function, emissivity is often defined in in differing ways by different authors. Here, we use the definition of the emissivity as given by Eq. 3 of Young et al. [YDL+16].

Note

The H number density, \(n_H\), is computed using density combined with the output of proton_electron_ratio.

Parameters:
  • density (Quantity) – Electron number density.

  • couple_density_to_temperature (bool, optional) – If True, the density will vary along the same axis as temperature in the computed level populations. The number of densities must be the same as the number of temperatures. This is useful, for example, when computing the level populations at constant pressure and is also much faster than computing the level populations along an independent density axis. By default, this is set to False.

Returns:

Quantity – A (l, m, k) shaped quantity, where l is the number of temperatures, m is the number of densities, and k is the number of transitions corresponding to the transition wavelengths described in contribution_function. If couple_density_to_temperature=True, then m=1 and l represents the number of temperatures and densities.

See also

contribution_function

Calculate contribution function, \(G(n,T)\)

free_bound(wavelength: Unit('Angstrom'), use_verner=True)[source]#

Free-bound continuum emission of the recombined ion.

Note

This does not include the ionization fraction or abundance factors. Unlike the equivalent IDL routine, the output here is not expressed per steradian and as such the factor of \(1/4\pi\) is not included.

When an electron is captured by an ion of charge \(z+1\) (the recombining ion), it creates a an ion of charge \(z\) (the recombined ion) and produces a continuum of emission called the free-bound continuum. The emission of the recombined ion is given by,

\[C_{fb}(\lambda, T) = \frac{2}{hc^3(k_B m_e)^{3/2}\sqrt{2\pi}}\frac{E^5}{T^{3/2}}\sum_i\frac{\omega_i}{\omega_0}\sigma_i^{\mathrm{bf}}\exp{\left(-\frac{E-I_i}{k_BT}\right)}\]

where \(E\) is the energy of the outgoing photon, \(\omega_i,\omega_0\) are the statastical weights of the \(i\)-th level of the recombined ion and the ground level of the recombining ion, respectively, \(\sigma_i^{\mathrm{bf}}\) is the free-bound cross-section, and \(I_i\) is the energy required to ionize the recombined ion from level \(i\). A detailed derivation of this formula can be found in CHIANTI Technical Report No. 12.

For ground state transitions, the photoionization cross-section \(\sigma_i^{\mathrm{bf}}\) is evaluated using Eq. 1 of Verner and Yakovlev [VY95] if use_verner is set to True. For all other transitions, and in all cases if use_verner is set to False, \(\sigma_i^{\mathrm{bf}}\) is evaluated using the method of Karzas and Latter [KL61].

Parameters:
  • wavelength (Quantity) –

  • use_verner (bool, optional) – If True, evaluate ground-state cross-sections using method of Verner and Yakovlev [VY95].

free_free(wavelength: Unit('Angstrom'))[source]#

Free-free continuum emission as a function of temperature and wavelength.

Note

This does not include ionization fraction or abundance factors.

Free-free emission, also known as bremsstrahlung (or “braking radiation”), is produced when an ion interacts with a free electron, reduces the momentum of the free electron, and, by conservation of energy and momentum, produces a photon. According to Eq. 4.114 of Phillips et al. [PFL08] the free-free emission produced by a thermal distribution of electrons as a function of temperature and wavelength is given by,

\[C_{ff}(\lambda,T_e) = \frac{c}{3m_e}\left(\frac{\alpha h}{\pi}\right)^3\sqrt{\frac{2\pi}{3m_ek_B}}\frac{z^2}{\lambda^2T_e^{1/2}}\exp{\left(-\frac{hc}{\lambda k_BT_e}\right)}\langle g_{ff}\rangle,\]

where \(\alpha\) is the fine-structure constant, \(z\) is the charge of the ion, and \(\langle g_{ff}\rangle\) is the velocity-averaged free-free Gaunt factor.

Parameters:

wavelength (Quantity) –

Notes

The result does not include ionization equilibrium or abundance factors.

See also

gaunt_factor_free_free

fiasco.IonCollection.free_free

Includes abundance and ionization equilibrium

gaunt_factor_free_free(wavelength: Unit('Angstrom'))[source]#

Free-free Gaunt factor as a function of wavelength.

The free-free Gaunt factor is calculated from a lookup table of temperature averaged free-free Gaunt factors from Table 2 of Sutherland [Sut98] as a function of \(\log{\gamma^2},\log{u}\), where \(\gamma^2=Z^2\mathrm{Ry}/k_BT\) and \(u=hc/\lambda k_BT\).

For the regime, \(6<\log_{10}(T)< 8.5\) and \(-4<\log_{10}(u)<1\), the above prescription is replaced with the fitting formula of Itoh et al. [ISK+00] for the relativistic free-free Gaunt factor. This is given by Eq. 4 of Itoh et al. [ISK+00],

\[g_{ff} = \sum_{i,j=0}^{10} a_{ij}t^iU^j,\]

where \(t=(\log{T} - 7.25)/1.25\) and \(U=(\log{u} + 1.5)/2.5\).

intensity(density: Unit('1 / cm3'), emission_measure: Unit('1 / cm5'), **kwargs)[source]#

Line-of-sight intensity computed assuming a particular column emission measure.

The intensity along the line-of-sight can be written as,

\[I = \frac{1}{4\pi}\int\mathrm{d}T,G(n,T)n_Hn_e\frac{dh}{dT}\]

which, in the isothermal approximation, can be simplified to,

\[I(T_0) \approx \frac{1}{4\pi}G(n,T_0)\mathrm{EM}(T_0)\]

where,

\[\mathrm{EM}(T) = \int\mathrm{d}h\,n_Hn_e\]

is the column emission measure.

Parameters:
  • density (Quantity) – Electron number density

  • emission_measure (Quantity) – Column emission measure. Must be either a scalar, an array of length 1, or an array with the same length as temperature. Note that it is assumed that the emission measure is the product of the H and electron density.

  • couple_density_to_temperature (bool, optional) – If True, the density will vary along the same axis as temperature. The number of densities must be the same as the number of temperatures. This is useful, for example, when computing the intensities at constant pressure and is also much faster than computing the intensity along an independent density axis. By default, this is set to False.

Returns:

Quantity – A (l, m, k) shaped quantity, where l is the number of temperatures, m is the number of densities, and k is the number of transitions corresponding to the transition wavelengths described in contribution_function. If couple_density_to_temperature=True, then m=1 and l represents the number of temperatures and densities.

level_populations(density: Unit('1 / cm3'), include_protons=True, include_level_resolved_rate_correction=True, couple_density_to_temperature=False)[source]#

Energy level populations as a function of temperature and density.

Compute the level populations of the given ion as a function of temperature and density. This is done by solving the homogeneous linear system of equations describing the processes that populate and depopulate each energy level of each ion. Section 3 of Young et al. [YDL+16] provides a brief description of this set of equations.

Parameters:
  • density (Quantity) –

  • include_protons (bool, optional) – If True (default), include proton excitation and de-excitation rates.

  • include_level_resolved_rate_correction (bool, optional) – If True (default), include the level-resolved ionization and recombination rate correction in the resulting level populations as described in Section 2.3 of Landi et al. [LDZY+06].

  • couple_density_to_temperature (bool, optional) – If True, the density will vary along the same axis as temperature in the computed level populations and the number of densities must be the same as the number of temperatures. This is useful, for example, when computing the level populations at constant pressure and is also much faster than computing the level populations along an independent density axis. By default, this is set to False.

Returns:

Quantity – A (l, m, n) shaped quantity, where l is the number of temperatures, m is the number of densities, and n is the number of energy levels. If couple_density_to_temperature=True, then m=1 and l represents the number of temperatures and densities.

next_ion()[source]#

Return an Ion instance with the next highest ionization stage.

For example, if the current instance is Fe XII (+11), this method returns an instance of Fe XIII (+12). All other input arguments remain the same.

previous_ion()[source]#

Return an Ion instance with the next lowest ionization stage.

For example, if the current instance is Fe XII (+11), this method returns an instance of Fe XI (+10). All other input arguments remain the same.

spectrum(*args, **kwargs)[source]#

Construct the spectrum using a given filter over a specified wavelength range.

All arguments are passed directly to fiasco.IonCollection.spectrum.

See also

fiasco.IonCollection.spectrum

Compute spectrum for multiple ions

intensity

Compute LOS intensity for all transitions