ChiantiPy.core package¶
Subpackages¶
Submodules¶
ChiantiPy.core.Continuum module¶
Continuum module
-
class
ChiantiPy.core.Continuum.
continuum
(ionStr, temperature, abundance=None, em=None, verbose=0)¶ Bases:
ChiantiPy.base._IonTrails.ionTrails
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
- Keyword Arguments
abundance (float, str, optional) – Elemental 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}\).
verbose (bool) – if True, prints additional info to the console
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()
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.
References
- 101(1,2,3)
Sutherland, R. S., 1998, MNRAS, 300, 321
- 102(1,2)
Verner & Yakovlev, 1995, A&AS, 109, 125
- 103(1,2,3,4,5,6,7)
Karzas and Latter, 1961, ApJSS, 6, 167
- 104(1,2,3,4,5)
Itoh, N. et al., 2000, ApJS, 128, 125
- 106(1,2,3,4,5,6,7,8,9,10)
Mewe, R. et al., 1986, A&AS, 65, 511
- 107(1,2,3,4)
Rybicki and Lightman, 1979, Radiative Processes in Astrophysics, (Wiley-VCH)
- 108(1,2,3,4)
Gronenschild, E.H.B.M. and Mewe, R., 1978, A&AS, 32, 283
- 109
Mao, J., Kaastra, J., Badlell, N., 2017, A&A, A&A, 599, A10 _
-
calculate_free_bound_loss
(**kwargs)¶ 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 108 and 106. Eq. 1a of 106 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 106 (see mewe_gaunt_factor for more details) and \(C_{ff}\) is the numerical constant as given in Eq. 4 of 108 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}\]
-
freeBound
(wvl, includeAbund=True, includeIoneq=True, verner=True, verbose=False)¶ Calculates the free-bound (radiative recombination) continuum emissivity of an ion. Provides emissivity in units of ergs \(\mathrm{cm}^{-2}\) \(\mathrm{s}^{-1}\) \(\mathrm{str}^{-1}\) \(\mathrm{Å}^{-1}\) for an individual ion. If includeAbund is set, the abundance is included. If includeIoneq is set, the ionization equililibrium for the given ion is included
Notes
Uses the Gaunt factors of 103 for recombination to the excited levels
Uses the photoionization cross sections of 3 to develop the free-bound cross section
Include the elemental abundance and ionization fraction by default
The specified ion is the target ion
uses the corrected version of the K-L bf gaunt factors available in CHIANTI V10
revised to calculate the bf cross, fb cross section and the maxwell energy distribution
References
-
freeBoundLoss
(includeAbund=True, includeIoneq=True, verner=True, verbose=False)¶ to calculate the free-bound (radiative recombination) energy loss rate coefficient of an ion, the ion is taken to be the target ion, including the elemental abundance and the ionization equilibrium population uses the Gaunt factors of 103 Karzas, W.J, Latter, R, 1961, ApJS, 6, 167 provides rate = erg cm^-2 s^-1 .. rubric:: Notes
Uses the Gaunt factors of 103 for recombination to the ground level
Uses the photoionization cross sections of 2 to develop the free-bound cross section
Does not include the elemental abundance or ionization fraction
The specified ion is the target ion
uses the corrected version of the K-L bf gaunt factors available in CHIANTI V10
revised to calculate the bf cross, fb cross section and the maxwell energy distribution
the Verner cross sections are not included for now
using the RESTART formulation
References
-
freeBoundLossMao
(includeAbund=False, includeIoneq=False)¶ to calculate the radiative loss rate from the parameters of 109 Mao J., Kaastra J., Badnell N.R. <Astron. Astrophys. 599, A10 (2017)>=2017A&A…599A..10M
-
freeBoundLossMewe
(**kwargs)¶ 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 108 and 106. Eq. 1a of 106 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 106 (see mewe_gaunt_factor for more details) and \(C_{ff}\) is the numerical constant as given in Eq. 4 of 108 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}\]
-
freeBoundRate
(verner=True, verbose=False)¶ Calculates the free-bound (radiative recombination) rate of an ion. Provides emissivity in units of ergs \(\mathrm{cm}^{-2}\) \(\mathrm{s}^{-1}\) \(\mathrm{str}^{-1}\) \(\mathrm{\AA}^{-1}\) for an individual ion. If includeAbund is set, the abundance is included. If includeIoneq is set, the ionization equililibrium for the given ion is included This method is only for the purposes of testing. The best rr/FB rates come from the ion class
Notes
Uses the Gaunt factors of 103 for recombination to the excited levels
Uses the photoionization cross sections of 3 to develop the free-bound cross section
Include the elemental abundance and ionization fraction by default
The specified ion is the target ion
uses the corrected version of the K-L bf gaunt factors available in CHIANTI V10
revised to calculate the bf cross, fb cross section and the maxwell energy distribution
-
freeFree
(wavelength, includeAbund=True, includeIoneq=True, **kwargs)¶ Calculates the free-free emission for a single ion. The result is returned as a dict to the FreeFree attribute. The dict has the keywords intensity, wvl, temperature, em.
The free-free emission for the given ion is calculated according Eq. 5.14a of 107, 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 104 and 101, depending on the temperature and energy regime. See itoh_gaunt_factor and sutherland_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
includeAbund (bool, optional) – If True, include the ion abundance in the final output.
includeIoneq (bool, optional) – If True, include the ionization equilibrium in the final output
-
freeFreeLoss
(includeAbund=True, includeIoneq=True, **kwargs)¶ Calculate the free-free energy loss rate of an ion. The result is returned to the FreeFreeLoss attribute.
The free-free radiative loss rate is given by Eq. 5.15a of 107. 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 ion 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 103. Note that this expression for the loss rate is just the integral over wavelength of Eq. 5.14a of 103, the free-free emission, and is expressed in units of erg \(\mathrm{cm}^3\,\mathrm{s}^{-1}\).
-
free_free_loss
(includeAbund=True, includeIoneq=True, **kwargs)¶ 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 107. 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 103. Note that this expression for the loss rate is just the integral over wavelength of Eq. 5.14a of 107, the free-free emission, and is expressed in units of erg \(\mathrm{cm}^3\,\mathrm{s}^{-1}\).
-
ioneqOne
()¶ Provide the ionization equilibrium for the selected ion as a function of temperature. Similar to but not identical to ion.ioneqOne() - the ion class needs to be able to handle the ‘dielectronic’ ions returned in self.IoneqOne
-
ioneq_one
(stage, **kwargs)¶ 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)¶ Calculates the free-free gaunt factors of 104.
An analytic fitting formulae for the relativistic Gaunt factor is given by Eq. 4 of 104,
\[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 104. 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 104
-
mewe_gaunt_factor
(**kwargs)¶ Calculate the Gaunt factor according to 106 for a single ion \(Z_z\).
Using Eq. 9 of 106, 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 106), and \(\zeta_0\) is the number of vacancies in the 0th shell and is given in Table 1 of 106. 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
-
sutherland_gaunt_factor
(wavelength)¶ Calculates the free-free gaunt factor calculations of 101.
The Gaunt factors of 101 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 using ~scipy.ndimage.map_coordinates.
-
vernerCross
(wvl)¶ Calculates the photoionization cross-section using data from 102 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 102,
\[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.
Parameters:
- wvllist, ndarray
wavelength in Angstroms
ChiantiPy.core.Ion module¶
Ion class
-
class
ChiantiPy.core.Ion.
ion
(ionStr, temperature=None, eDensity=None, pDensity='default', radTemperature=None, rStar=None, abundance=None, setup=True, em=None, verbose=0)¶ Bases:
ChiantiPy.base._IoneqOne.ioneqOne
,ChiantiPy.base._IonTrails.ionTrails
,ChiantiPy.base._SpecTrails.specTrails
The top level class for performing spectral calculations for an ion in the CHIANTI database.
- Parameters
ionStr (str) – CHIANTI notation for the given ion, e.g. ‘fe_12’ that corresponds to the Fe XII ion.
temperature (float , tuple, list, ~numpy.ndarray, optional) – Temperature array (Kelvin)
eDensity (float , tuple, list, or ~numpy.ndarray, optional) – Electron density array (\(\mathrm{cm^{-3}}\) )
pDensity (float, tuple, list or ~numpy.ndarray, optional) – Proton density (\(\mathrm{cm}^{-3}\) )
radTemperature (float or ~numpy.ndarray, optional) – Radiation black-body temperature (in Kelvin)
rStar (float or ~numpy.ndarray, optional) – Distance from the center of the star (in stellar radii)
abundance (float or str, optional) – Elemental abundance relative to Hydrogen or name of CHIANTI abundance file to use, without the ‘.abund’ suffix, e.g. ‘sun_photospheric_1998_grevesse’.
setup (bool or str, optional) – If True, run ion setup function. Otherwise, provide a limited number of attributes of the selected ion
em (float or ~numpy.ndarray, optional) – Emission Measure, for the line-of-sight emission measure (\(\mathrm{\int \, n_e \, n_H \, dl}\)) (\(\mathrm{cm}^{-5}\).), for the volumetric emission measure \(\mathrm{\int \, n_e \, n_H \, dV}\) (\(\mathrm{cm^{-3}}\)).
- Variables
IonStr (str) – Name of element plus ion, e.g. fe_12 for Fe XII
Z (int) – the nuclear charge, 26 for fe_12.
Ion (int) – the ionization stage, 12 for fe_12.
Dielectronic (bool) – true if the ion is a ‘dielectronic’ ion where the levels are populated by dielectronic recombination.
Spectroscopic (str) – the spectroscopic notation for the ion, such as Fe XII for fe_12.
Filename (str) – the complete name of the file generic filename in the CHIANTI database, such as $XUVTOP/fe/fe_12/fe_12.
Ip (float) – the ionization potential of the ion
FIP (float) – the first ionization potential of the element
Defaults (dict) – these are specified by the software unless a chiantirc file is found in ‘$HOME/.chianti’:
Notes
The keyword arguments temperature, eDensity, radTemperature, rStar, em must all be either a float or have the same dimension as the rest if specified as lists, tuples or arrays.
The Defaults dict should have the following keys:
abundfile, the elemental abundance file, unless specified in chiantirc this defaults to sun_photospheric_1998_grevesse.
ioneqfile, the ionization equilibrium file name. Unless specified in ‘chiantirc’ this is defaults to chianti. Other choices are availble in $XUVTOP/ioneq
wavelength, the units of wavelength (Angstroms, nm, or kev), unless specified in the ‘chiantirc’ this is defaults to ‘angstrom’.
flux, specified whether the line intensities are give in energy or photon fluxes, unless specified in the ‘chiantirc’ this is defaults to energy.
gui, specifies whether to use gui selection widgets (True) or to make selections on the command line (False). Unless specified in the ‘chiantirc’ this is defaults to False.
-
boundBoundLoss
(allLines=1)¶ Calculate the summed radiative loss rate for all spectral lines of the specified ion.
- Parameters
allLines (bool) – If True, include losses from both observed and unobserved lines. If False, only include losses from observed lines.
includes elemental abundance and ionization fraction.
- Returns
creates the attribute
BoundBoundLoss (dict with the keys below.) – rate : the radiative loss rate (\(\mathrm{erg \, cm^{-3}} \, \mathrm{s}^{-1}\)) per unit emission measure.
temperature : (K).
eDensity : electron density (\(\mathrm{cm^{-3}}\))
-
diCross
(energy=None, verbose=False)¶ Calculate the direct ionization cross section (cm$^2) as a function of the incident electron energy in eV, puts values into attribute DiCross
- Parameters
energy (array-like) – incident electron energy in eV
verbose (bool, int) – with verbose set to True, printing is enabled
- Variables
DiCross (dict) – keys: energy, cross
-
diRate
(ngl=20)¶ Calculate the direct ionization rate coefficient as a function of temperature (K)
-
drRate
()¶ Provide the dielectronic recombination rate coefficient as a function of temperature (K).
-
drRateLvl
(verbose=0)¶ to calculate the level resolved dielectronic rate from the higher ionization stage to the ion of interest rates are determined from autoionizing A-values the dictionary self.DrRateLvl contains rate = the dielectronic rate into an autoionizing level effRate = the dielectronic rate into an autoionizing level mutilplied by the branching ratio for a stabilizing transition totalRate = the sum of all the effRates
-
eaCross
(energy=None, verbose=False)¶ Provide the excitation-autoionization cross section.
Energy is given in eV.
-
eaDescale
()¶ Calculates the effective collision strengths (upsilon) for excitation-autoionization as a function of temperature.
-
eaRate
()¶ Calculate the excitation-autoionization rate coefficient.
-
emiss
(allLines=True)¶ Calculate the emissivities for lines of the specified ion.
units: ergs s^-1 str^-1
Does not include elemental abundance or ionization fraction
Wavelengths are sorted
set allLines = True to include unidentified lines
-
emissList
(index=None, wvlRange=None, wvlRanges=None, top=10, relative=0, outFile=0)¶ List the emissivities.
wvlRange, a 2 element tuple, list or array determines the wavelength range
Top specifies to plot only the top strongest lines, default = 10
normalize = 1 specifies whether to normalize to strongest line, default = 0
-
emissPlot
(index=None, wvlRange=None, top=10, linLog='lin', relative=0, verbose=0, plotFile=0, saveFile=0)¶ Plot the emissivities.
wvlRange, a 2 element tuple, list or array determines the wavelength range
Top specifies to plot only the top strongest lines, default = 10
linLog specifies a linear or log plot, want either lin or log, default = lin
normalize = 1 specifies whether to normalize to strongest line, default = 0
-
emissRatio
(wvlRange=None, wvlRanges=None, top=10)¶ Plot the ratio of 2 lines or sums of lines. Shown as a function of density and/or temperature. For a single wavelength range, set wvlRange = [wMin, wMax] For multiple wavelength ranges, set wvlRanges = [[wMin1,wMax1],[wMin2,wMax2], …] A plot of relative emissivities is shown and then a dialog appears for the user to choose a set of lines.
-
gofnt
(wvlRange=None, top=10, verbose=False, plot=True)¶ Calculate the ‘so-called’ G(T) function.
Given as a function of both temperature and eDensity.
Only the top( set by ‘top’) brightest lines are plotted. the G(T) function is returned in a dictionary self.Gofnt
- Note: if the default value for gui is set to False, then it is usually
necessary to invoke this method twice to get the desired result.
- Keyword Arguments
wvlRange (array-like) – the wavelength range to be considered, a two element array-type
top (int) – specifies to plot only the top strongest lines, default = 10
verbose (bool) – if True, additional information is printed to the console
plot (bool) – if True, the G(T) functionis plotted, default - True
-
intensity
(allLines=1, verbose=0)¶ Calculate the intensities for lines of the specified ion.
units: ergs cm^-3 s^-1 str^-1
includes elemental abundance and ionization fraction.
the emission measure ‘em’ is included if specified
-
intensityRatioInterpolate
(data, scale='lin', plot=False, verbose=False)¶ Take a set of data and interpolate against the IntensityRatio.
Given a set of intensities ratio data it correlates them to electron density values. This is done by interpolating against the curve computed from the IntenstyRatio (dict). For the interpolation, the B-spline representation is used from scipy.
- Parameters
data (float or array-like) – The intensities ratio data to be converted into density values.
scale ({‘lin’, ‘loglog’, ‘logx’, ‘logy’}, str, default : ‘lin’) – The scale to be used for the interpolation (default ‘lin’ indicates linear).
plot (bool, default : False) – If True the input data will be plotted on top of the IntensityRatio curve.
verbose (bool, default : False) – If True displays a message about the physical quantity that was used (i.e. Temperature or Density). Additionally, prints the interpolated data in pairs of data,value (i.e. input intensities ratio data and corresponding density values).
- Returns
Creates the attribute
IntensityRatioInterpolated (dict) – Dictionary that contains the input data and the derived density values, with keys
data : input intensities ratio (no units)
value : electron density (\(\mathrm{cm^{-3}}\))
Examples
>>> import ChiantiPy.core as ch >>> temp = 10**5.9 >>> dens = 10**(np.linspace(8,12,num=50,endpoint=True)) >>> fe12 = ch.ion('fe_12',temperature=temp,eDensity=dens) >>> fe12.intensityRatio(wvlRange=[186.,196.]) >>> # Select for example lines 186.887 and 195.119 from the GUI >>> int_ratios = np.linspace(0.1,0.8,20) >>> fe12.intensityRatioInterpolate(int_ratios) >>> dens_values = fe12.IntensityRatioInterpolated['value']
>>> # Same setup but using a 2D array of input data >>> x = 10 >>> y = 30 >>> int_ratios_map = np.random.random((x,y)) >>> fe12.intensityRatioInterpolate(int_ratios_map,verbose=True) >>> dens_values_map = fe12.IntensityRatioInterpolated['value'].reshape((x,y))
-
ionizCross
(energy=None)¶ Provides the total ionization cross section.
Notes
uses diCross and eaCross.
-
ionizRate
(ngl=20)¶ Provides the total ionization rate.
Calls diRate and eaRate.
-
p2eRatio
()¶ Calculates the proton density to electron density ratio using Eq. 7 of 1.
Notes
Uses the abundance and ionization equilibrium.
References
-
popPlot
(top=10, levels=[], scale=0, plotFile=0, outFile=0, pub=0, addTitle=None, addLegend=True)¶ Plots populations vs temperature or eDensity.
top specifies the number of the most highly populated levels to plot (the default)
or can set levels to an array such as a list to set the desired levels to plot
if scale is set, then the population, if plotted vs. density, is divided by density - only useful if plotting level populations vs density
if pub is set, the want publication plots (bw, lw=2).
if addLegend is set, a matplotlib legend is added
-
populate
(popCorrect=True, verbose=False)¶ Calculate level populations for specified ion. possible keyword arguments include temperature, eDensity, pDensity, radTemperature and rStar populate assumes that all of the population in the higher ionization stages exists only in the ground level use drPopulate() for cases where the population of various levels in the higher ionization stage figure into the calculation
-
recombRate
()¶ Provides the total recombination rate coefficient.
Calls drRate and rrRate
-
rrRate
()¶ Provide the radiative recombination rate coefficient as a function of temperature (K).
-
rrlvlDescale
(verbose=1)¶ Interpolate and extrapolate rrlvl rates. Used in level population calculations.
-
setup
(alternate_dir=None, verbose=False)¶ Setup various CHIANTI files for the ion including .wgfa, .elvlc, .scups, .psplups, .reclvl, .cilvl, and others.
- Parameters
alternate_dir (str) – directory cotaining the necessary files for a ChiantiPy ion; use to setup an ion with files not in the current CHIANTI directory
verbose (bool)
Notes
If ion is initiated with setup=False, call this method to do the setup at a later point.
-
setupIonrec
(alternate_dir=None, verbose=False)¶ Setup method for ion recombination and ionization rates.
Notes
Allows a bare-bones ion object to be setup up with just the ionization and recombination rates. For ions without a complete set of files - one that is not in the MasterList.
-
spectrum
(wavelength, filter=(<function gaussianR>, 1000.0), label=0, allLines=1)¶ Calculates the line emission spectrum for the specified ion.
Convolves the results of intensity to make them look like an observed spectrum the default filter is the gaussianR filter with a resolving power of 1000. Other choices include chianti.filters.box and chianti.filters.gaussian. When using the box filter, the width should equal the wavelength interval to keep the units of the continuum and line spectrum the same.
includes ionization equilibrium and elemental abundances
can be called multiple times to use different filters and widths uses label to keep the separate applications of spectrum sorted by the label for example, do .spectrum( …. label = ‘test1’) and do .spectrum( …. label = ‘test2’) then will get self.Spectrum.keys() = test1, test2 and self.Spectrum[‘test1’] = {‘intensity’:aspectrum, ‘wvl’:wavelength, ‘filter’:useFilter.__name__, ‘filterWidth’:useFactor}
Notes
scipy.ndimage.filters also includes a range of filters.
-
twoPhoton
(wvl, verbose=False)¶ to calculate the two-photon continuum - only for hydrogen- and helium-like ions includes the elemental abundance and the ionization equilibrium includes the emission measure if specified discussed in 105
References
-
twoPhotonEmiss
(wvl)¶ To calculate the two-photon continuum rate coefficient - only for hydrogen- and helium-like ions
-
twoPhotonLoss
()¶ to calculate the two-photon energy loss rate - only for hydrogen- and helium-like ions includes the elemental abundance and the ionization equilibrium does not include the emission measure
-
upsilonDescale
(prot=0)¶ Provides the temperatures and effective collision strengths (upsilons) set prot for proton rates otherwise, ce will be set for electron collision rates uses the new format “scups” files
ChiantiPy.core.Ioneq module¶
Ionization equilibrium class
-
class
ChiantiPy.core.Ioneq.
ioneq
(el_or_z)¶ Bases:
object
Calculation ion fractions as a function of temperature assuming ionization equilibrium.
- Parameters
el_or_z (int or str) – Atomic number or symbol
Note
When either loading or calculating a set of ion fractions, the temperature and ion fractions are returned to the Temperature and Ioneq attributes, respectively.
-
calculate
(temperature)¶ Calculate ion fractions for given temperature array using the total ionization and recombination rates.
-
load
(ioneqName=None)¶ Read temperature and ion fractions from a CHIANTI “.ioneq” file.
-
plot
(stages=0, tRange=0, yr=0, oplot=False, label=1, title=1, bw=False, semilogx=0, heightAdjust=1.3, verbose=0)¶ Plots the ionization equilibria.
self.plot(tRange=None, yr=None, oplot=False) stages = sequence of ions to be plotted, neutral == 1, fully stripped == Z+1 tRange = temperature range, yr = ion fraction range
for overplotting: oplot=”ioneqfilename” such as ‘mazzotta’ or if oplot=True or oplot=1 and a widget will come up so that a file can be selected. bw, if True, the plot is made in black and white
-
plotRatio
(stageN, stageD, tRange=0, yr=0, label=1, title=1, bw=0, semilogx=1, verbose=0)¶ Plots the ratio of the ionization equilibria of two stages of a given element
self.plotRatio(stageN, stageD) stages = sequence of ions to be plotted, neutral == 1, fully stripped == Z+1 stageN = numerator stageD = denominator tRange = temperature range, yr = ion fraction range
ChiantiPy.core.IpyMspectrum module¶
-
ChiantiPy.core.IpyMspectrum.
doAll
(inpt)¶ to process ff, fb and line inputs
-
class
ChiantiPy.core.IpyMspectrum.
ipymspectrum
(temperature, eDensity, wavelength, filter=(<function gaussianR>, 1000.0), label=None, elementList=None, ionList=None, minAbund=None, keepIons=0, doLines=1, doContinuum=1, allLines=1, em=None, abundance=None, verbose=0, timeout=0.1)¶ Bases:
ChiantiPy.base._IonTrails.ionTrails
,ChiantiPy.base._SpecTrails.specTrails
this is the multiprocessing version of spectrum for using inside an IPython Qtconsole or notebook.
be for creating an instance, it is necessary to type something like the following into a console
> ipcluster start –n=3
this is the way to invoke things under the IPython 6 notation
Calculate the emission spectrum as a function of temperature and density.
temperature and density can be arrays but, unless the size of either is one (1), the two must have the same size
the returned spectrum will be convolved with a filter of the specified width on the specified wavelength array
the default filter is gaussianR with a resolving power of 100. Other filters, such as gaussian, box and lorentz, are available in ChiantiPy.filters. When using the box filter, the width should equal the wavelength interval to keep the units of the continuum and line spectrum the same.
A selection of elements can be make with elementList a list containing the names of elements that are desired to be included, e.g., [‘fe’,’ni’]
A selection of ions can be make with ionList containing the names of the desired lines in Chianti notation, i.e. C VI = c_6
Both elementList and ionList can not be specified at the same time
a minimum abundance can be specified so that the calculation can be speeded up by excluding elements with a low abundance. With solar photospheric abundances -
setting minAbund = 1.e-4 will include H, He, C, O, Ne setting minAbund = 2.e-5 adds N, Mg, Si, S, Fe setting minAbund = 1.e-6 adds Na, Al, Ar, Ca, Ni
Setting em will multiply the spectrum at each temperature by the value of em.
em [for emission measure], can be a float or an array of the same length as the temperature/density. allLines = 1 will include lines with either theoretical or observed wavelengths. allLines=0 will include only those lines with observed wavelengths
timeout - a small but non-zero value seems to be necessary
ChiantiPy.core.Mspectrum module¶
-
class
ChiantiPy.core.Mspectrum.
mspectrum
(temperature, eDensity, wavelength, filter=(<function gaussianR>, 1000.0), label=0, elementList=None, ionList=None, minAbund=None, keepIons=0, abundance=None, doLines=1, doContinuum=1, allLines=1, em=None, proc=3, verbose=0, timeout=0.1)¶ Bases:
ChiantiPy.base._IonTrails.ionTrails
,ChiantiPy.base._SpecTrails.specTrails
this is the multiprocessing version of spectrum set proc to the desired number of processors, default=3
Calculate the emission spectrum as a function of temperature and density.
temperature and density can be arrays but, unless the size of either is one (1), the two must have the same size
the returned spectrum will be convolved with a filter of the specified width on the specified wavelength array
the default filter is gaussianR with a resolving power of 100. Other filters, such as gaussian, box and lorentz, are available in chianti.filters. When using the box filter, the width should equal the wavelength interval to keep the units of the continuum and line spectrum the same.
A selection of elements can be make with elementList a list containing the names of elements that are desired to be included, e.g., [‘fe’,’ni’]
A selection of ions can be make with ionList containing the names of the desired lines in Chianti notation, i.e. C VI = c_6
Both elementList and ionList can not be specified at the same time
a minimum abundance can be specified so that the calculation can be speeded up by excluding elements with a low abundance. With solar photospheric abundances -
setting minAbund = 1.e-4 will include H, He, C, O, Ne setting minAbund = 2.e-5 adds N, Mg, Si, S, Fe setting minAbund = 1.e-6 adds Na, Al, Ar, Ca, Ni
Setting em will multiply the spectrum at each temperature by the value of em.
em [for emission measure], can be a float or an array of the same length as the temperature/density.
allLines = 1 will include lines with either theoretical or observed wavelengths. allLines=0 will include only those lines with observed wavelengths
proc, the number of processors to use
timeout, a small but non-zero value seems to be necessary keepIons: set this to keep the ion instances that have been calculated in a dictionary self.IonInstances with the keywords being the CHIANTI-style ion names
abundance: to select a particular set of abundances, set abundance to the name of a CHIANTI abundance file, without the ‘.abund’ suffix, e.g. ‘sun_photospheric_1998_grevesse’
If set to a blank (‘’), a gui selection menu will popup and allow the selection of an set of abundances
- Parameters
temperature (float, list, ndarray) – the temperature(s) in K
eDensity (float, ndarray) – eDensity: electron density in \(\mathrm{cm^{-3}}\)
wavelength (list or ndarray) – wavelength: array of wavelengths, generally in Angstroms
elementList (list) – elementList: list of elements to include, such as ‘fe’, ‘ne’, ‘s’
ionList (list) – ionList: list of ions to include, such as ‘fe_16’, ‘ne_10’
minAbund (float) – minAbund: minimum abundance (relative to H) to include
doLines (`bool1) – doLines: if true, line intensities are calculated
doContinuum (bool) – doContinuum: if true, continuum intensities are calculated only if wavelengths are in angstroms
keepIons (bool) –
- keepIons: keep the ion instances used in the calculation
should be used with caution otherwise the bunch instance can become quite large
em (float, list, ndarray) – em: the emission measure - the integral of the electron density times the hydrogen density along the line of sight
abundance (str) –
- abuncance: the file name of the abuncance set to be used
must be one in the $XUVTOP/abund directory
allLInes (bool) – allLines: whether or not to include unobserved lines
verbose (bool) – verbose: whether to allow certain print statements
ChiantiPy.core.RadLoss module¶
-
class
ChiantiPy.core.RadLoss.
radLoss
(temperature, eDensity, elementList=None, ionList=None, minAbund=None, doContinuum=True, doLines=True, abundance=None, verbose=0, allLines=1)¶ Bases:
ChiantiPy.base._IonTrails.ionTrails
,ChiantiPy.base._SpecTrails.specTrails
Calculate the radiative emission loss rate as a function of temperature and density.
includes elemental abundances or ionization equilibria
temperature and density can be arrays but, unless the size of either is one (1), the two must have the same size
A selection of ions can be make with ionList containing the names of the desired lines in Chianti notation, i.e. C VI = c_6
a minimum abundance can be specified so that the calculation can be speeded up by excluding elements with a low abundance. With solar photospheric abundances -
setting minAbund = 1.e-4 will include H, He, C, O, Ne setting minAbund = 2.e-5 adds N, Mg, Si, S, Fe setting minAbund = 1.e-6 adds Na, Al, Ar, Ca, Ni
Setting em will multiply the spectrum at each temperature by the value of em.
em [for emission measure], can be a float or an array of the same length as the temperature/density.
- abundance: to select a particular set of abundances, set abundance to the name of a CHIANTI abundance file,
without the ‘.abund’ suffix, e.g. ‘sun_photospheric_1998_grevesse’ If set to a blank (‘’), a gui selection menu will popup and allow the selection of an set of abundances
-
radLossPlot
(doTitle=False)¶ to plot the radiative losses vs temperature
- Parameters
doTitle (bool) – if True, a title is applied to the plot. The default is for no title
ChiantiPy.core.MradLoss module¶
-
class
ChiantiPy.core.MradLoss.
mradLoss
(temperature, eDensity, elementList=None, ionList=None, minAbund=None, doContinuum=True, doLines=True, abundance=None, verbose=0, allLines=1, proc=4)¶ Bases:
ChiantiPy.base._IonTrails.ionTrails
,ChiantiPy.base._SpecTrails.specTrails
Calculate the radiative emission loss rate as a function of temperature and density.
this is the multiprocessing version of radloss
includes elemental abundances or ionization equilibria
temperature and density can be arrays but, unless the size of either is one (1), the two must have the same size
A selection of ions can be make with ionList containing the names of the desired lines in Chianti notation, i.e. C VI = c_6
a minimum abundance can be specified so that the calculation can be speeded up by excluding elements with a low abundance. With solar photospheric abundances
setting minAbund = 1.e-4 will include H, He, C, O, Ne setting minAbund = 2.e-5 adds N, Mg, Si, S, Fe setting minAbund = 1.e-6 adds Na, Al, Ar, Ca, Ni
Setting em will multiply the spectrum at each temperature by the value of em.
em (for emission measure), can be a float or an array of the same length as the temperature/density.
- abundance: to select a particular set of abundances, set abundance to the name of a
CHIANTI abundance file, without the ‘.abund’ suffix, e.g. ‘sun_photospheric_1998_grevesse’ If set to a blank (‘’), a gui selection menu will popup and allow the selection of an set of abundances
- Parameters
temperature (float, list, ndarray) – the temperature(s) in K
eDensity (float, ndarray) – eDensity: electron density in \(\mathrm{cm^{-3}}\)
elementList (list) – elementList: list of elements to include, such as ‘fe’, ‘ne’, ‘s’
ionList (list) – ionList: list of ions to include, such as ‘fe_16’, ‘ne_10’
minAbund (float) – minAbund: minimum abundance (relative to H) to include
doLines (`bool1) – doLines: if true, line intensities are calculated
doContinuum (bool) – doContinuum: if true, continuum intensities are calculated only if wavelengths are in angstroms
abundance (str) –
- abuncance: the file name of the abuncance set to be used
must be one in the $XUVTOP/abund directory
allLInes (bool) – allLines: whether or not to include unobserved lines
verbose (bool) – verbose: whether to allow certain print statements
-
radLossPlot
(title=0)¶ to plot the radiative losses vs temperature
ChiantiPy.core.Bunch module¶
-
class
ChiantiPy.core.Bunch.
bunch
(temperature, eDensity, wvlRange, elementList=None, ionList=None, minAbund=None, keepIons=False, em=None, abundance=None, verbose=False, allLines=True)¶ Bases:
ChiantiPy.base._IonTrails.ionTrails
,ChiantiPy.base._SpecTrails.specTrails
Calculate the emission line spectrum as a function of temperature and density.
‘bunch’ is very similar to ‘spectrum’ except that continuum is not calculated and the spectrum is not convolved over a filter. However, this can be done with the inherited convolve method
one of the convenient things is that all of the instantiated ion classes, determined through such keywords as ‘elementList’, ‘ionList’, and ‘minAbund’ are kept in a dictionary self.IonInstances where self.IonInstances[‘mg_7’] is the class instance of ChiantiPy.core.ion for ‘mg_7’. All its methods and attributes are available.
includes elemental abundances and ionization equilibria
the set of abundances, a file in $XUVTOP/abundance, can be set with the keyword argument ‘abundanceName’
temperature and density can be arrays but, unless the size of either is one (1), the two must have the same size
Inherited methods include ‘intensityList’, ‘intensityRatio’ (between lines of different ions), and ‘intensityRatioSave’ and ‘convolve’.
A selection of elements can be make with elementList a list containing the names of elements that are desired to be included, e.g., [‘fe’,’ni’]
A selection of ions can be make with ionList containing the names of the desired lines in Chianti notation, i.e. C VI = c_6
Both elementList and ionList can not be specified at the same time
a minimum abundance can be specified so that the calculation can be speeded up by excluding elements with a low abundance. With solar photospheric abundances -
setting minAbund = 1.e-4 will include H, He, C, O, Ne setting minAbund = 2.e-5 adds N, Mg, Si, S, Fe setting minAbund = 1.e-6 adds Na, Al, Ar, Ca, Ni
At least one of elementList, ionList, or minAbund must be set in order for ‘bunch’ to include any ions.
Setting em will multiply the spectrum at each temperature by the value of em.
em [for emission measure], can be a float or an array of the same length as the temperature/density
- Parameters
temperature (float, list, ndarray) – the temperature(s) in K
eDensity (float, ndarray) – eDensity: electron density in \(\mathrm{cm^{-3}}\)
wvlRange (2 element list or ndarray) – wvlRange: range of wavelengths to consider, generally in angstroms
elementList (list) – elementList: list of elements to include, such as ‘fe’, ‘ne’, ‘s’
ionList (list) – ionList: list of ions to include, such as ‘fe_16’, ‘ne_10’
minAbund (float) – minAbund: minimum abundance (relative to H) to include
keepIons (bool) –
- keepIons: keep the ion instances used in the calculation
should be used with caution otherwise the bunch instance can become quite large
em (float, list, ndarray) – em: the emission measure
abundance (str) –
- abuncance: the file name of the abuncance set to be used
must be one in the $XUVTOP/abund directory
allLInes (bool) – allLines: whether or not to include unobserved lines
verbose (bool) – verbose: whether to allow certain print statements
-
convolve
(wavelength=None, filter=(<function gaussianR>, 1000.0), label=None, verbose=False)¶ the first application of spectrum calculates the line intensities within the specified wavelength range and for set of ions specified
wavelength will not be used if applied to ‘spectrum’ objects
wavelength IS need for ‘bunch’ objects - in this case, the wavelength should not extend beyond the limits of the wvlRange used for the ‘bunch’ calculation
- Keyword Arguments
wavelength (‘int’, list) – if an int, the attribute ‘Wavelength’ is looked for otherwise, wavelength is used
filter (tuple) – first elements if one of the ChiantiPy.tools.filters object second element is the width appropriate to the filter
label (str) – if set, creates a Spectrum[label] attribute
verbose (bool) – if True, prints info to the terminal
ChiantiPy.core.Spectrum module¶
-
class
ChiantiPy.core.Spectrum.
spectrum
(temperature, eDensity, wavelength, filter=(<function gaussianR>, 1000.0), label=None, elementList=None, ionList=None, minAbund=None, doLines=1, doContinuum=1, em=None, keepIons=0, abundance=None, verbose=0, allLines=1)¶ Bases:
ChiantiPy.base._IonTrails.ionTrails
,ChiantiPy.base._SpecTrails.specTrails
Calculate the emission spectrum as a function of temperature and density.
one of the convenient things is that all of the instantiated ion classes, determined through such keywords as ‘elementList’, ‘ionList’, and ‘minAbund’ are kept in a dictionary self.IonInstances where self.IonInstances[‘mg_7’] is the class instance of ChiantiPy.core.ion for ‘mg_7’. All its methods and attributes are available.
includes elemental abundances and ionization equilibria
the set of abundances, a file in $XUVTOP/abundance, can be set with the keyword argument ‘abundanceName’
temperature and density can be arrays but, unless the size of either is unity (1), the two must have the same size
the returned spectrum will be convolved with a filter of the specified width on the specified wavelength array
the default filter is gaussianR with a resolving power of 1000. Other filters, such as gaussian, box and lorentz, are available in ChiantiPy.tools.filters. When using the box filter, the width should equal the wavelength interval to keep the units of the continuum and line spectrum the same.
Inherited methods include ‘intensityList’, ‘intensityRatio’ (between lines of different ions), ‘intensityRatioSave’ and ‘convolve’
A selection of elements can be make with elementList a list containing the names of elements that are desired to be included, e.g., [‘fe’,’ni’]
A selection of ions can be make with ionList containing the names of the desired lines in CHIANTI notation, i.e. C VI = c_6
Both elementList and ionList can not be specified at the same time
a minimum abundance can be specified so that the calculation can be speeded up by excluding elements with a low abundance. The default of minAbund is 1.e-6
It is necessary to specify at least an elementList, an ionList, or a minAbund to select any ions for a spectrum calculation
With solar photospheric abundances
setting minAbund = 1.e-4 will include H, He, C, O, Ne
setting minAbund = 2.e-5 adds N, Mg, Si, S, Fe
setting minAbund = 1.e-6 adds Na, Al, Ar, Ca, Ni
Setting doLines = 0 will skip the calculation of spectral lines. Setting doContinuum =0 will skip the continuum calculation.
Setting em [for emission measure] will multiply the spectrum at each temperature by the value of em.
em [for emission measure] can be a float or an array of the same length as the temperature/density
keepIons set this to keep the ion instances that have been calculated in a dictionary self.IonInstances with the keywords being the CHIANTI-style ion names
abundance - to select a particular set of abundances, set abundance to the name of a CHIANTI abundance file, without the ‘.abund’ suffix, e.g. ‘sun_photospheric_1998_grevesse’
If set to a blank (‘’), a gui selection menu will popup and allow the selection of an set of abundances
- Parameters
temperature (float, list, ndarray) – the temperature(s) in K
eDensity (float, ndarray) – eDensity: electron density in :math mathrm{cm^{-3}}
wavelength (list or ndarray) – wavelength: array of wavelengths, generally in Angstroms
elementList (list) – elementList: list of elements to include, such as ‘fe’, ‘ne’, ‘s’
ionList (list) – ionList: list of ions to include, such as ‘fe_16’, ‘ne_10’
minAbund (float) – minAbund: minimum abundance (relative to H) to include
doLines (`bool1) – doLines: if true, line intensities are calculated
doContinuum (bool) – doContinuum: if true, continuum intensities are calculated only if wavelengths are in angstroms
keepIons (bool) –
- keepIons: keep the ion instances used in the calculation
should be used with caution otherwise the bunch instance can become quite large
em (float, list, ndarray) – em: the emission measure
abundance (str) –
- abuncance: the file name of the abuncance set to be used
must be one in the $XUVTOP/abund directory
allLInes (bool) – allLines: whether or not to include unobserved lines
verbose (bool) – verbose: whether to allow certain print statements
Module contents¶
chianti.core - contains the main classes for ChiantiPy users.
This software is distributed under the terms of the ISC Software License that is found in the LICENSE file