spectrum module
This module contains codes for creating spectra from the AtomDB emissivity files
The modules is separated by the type of spectrum you wish to model. For now, the CIE, NEI and PShock models are implemented.
Roughly speaking, you load a [Modelname]Session, and you can then obtain:
- linelist :
list of lines in a certain wavelength interval)
- line_emissivity :
the emissivity of a specific line as a function of temperature, line intensity and more.
- spectrum :
the emissivity in photons cm^3 bin^-1 s^-1. This requires a response be set first
- class pyatomdb.spectrum.CIESession(linefile='$ATOMDB/apec_line.fits', cocofile='$ATOMDB/apec_coco.fits', elements=False, abundset='AG89')[source]
Bases:
object
Load and generate a collisional ionization equilibrium spectrum
- Parameters:
- linefilestring or HDUList, optional
The line emissivity data file (either name or already open)
- cocofilestring or HDUList, optional
The continuum emissivity data file (either name or already open)
- elementsiterable of int
Elements to include, listed by atomic number. if not set, include all.
- abundsetstring
The abundance set to use. Default AG89.
Examples
Create a session instance:
>>> s=CIESession()
Set up the responses, in this case a dummy response from 0.1 to 10 keV, with area 1cm^2 in each bin
>>> ebins = numpy.linspace(0.1,10,1000) >>> s.set_response(ebins, raw=True)
(Alternatively, for a real response file, s.set_response(rmffile, arf=arffile)
Turn on thermal broadening
>>> s.set_broadening(True) Will thermally broaden lines with emissivity > 1.000000e-18 ph cm3 s-1
Return spectrum at 1.0keV
>>> spec = s.return_spectrum(1.0)
spec is in photons cm^5 s^-1 bin^-1; ebins are the bin edges (so spec is 1 element shorter than ebins)
- Attributes:
- datacachedict
Any Atomdb FITS files which have to be opened are stored here
- spectraCIESpectra
Object storing the actual spectral data
- elementsiterable of int
Elements to include, listed by atomic number. if not set, include all.
- default_abundsetstring
The abundance set used for the original emissivity file calculation
- abundsetstring
The abundance set to be used for the returned spectrum
- abundsetvectorarray_like(float)
The relative abundance between default_abundset and abundset for each element
- response_setbool
Have we loaded a response (or set a dummy response)
- dolinesbool
Calculate line emission (default True)
- docontbool
Calculate Continuum emission (default True)
- dopseudobool
Calculate PseudoContinuum (weak line) emission (default True)
- broaden_limitfloat
Apply broadening to lines with epsilon > this value (ph cm3 s-1)
- thermal_broadeningbool
Apply thermal broadening to lines (default = False)
- velocity_broadeningfloat
Apply velocity broadening with this velocity (km/s). If <=0, do not apply.
- format_linelist(linelist)[source]
Print out the line list in a nice format, with all the transitions described in more detail
- Parameters:
- linelistnumpy.array(dtype=linelist)
The linelist returned from “return_linelist”
- Returns
- ——-
- strthe formatted linelist
- return_line_emissivity(Te, Z, z1, up, lo, specunit='A', teunit='keV', apply_aeff=False, apply_abund=True, log_interp=True)[source]
Get line emissivity as function of Te.
- Parameters:
- Tefloat or array(float)
Temperature(s) in keV or K
- Zint
nuclear charge of element
- z1int
ion charge +1 of ion
- upint
upper level for transition
- loint
lower level for transition
- specunit{‘Angstrom’,’keV’}
Units for wavelength or energy (a returned value)
- teunit{‘keV’ , ‘K’}
Units of Te (kev or K, default keV)
- apply_aeffbool
If true, apply the effective area to the line emissivity in the linelist to modify their intensities.
- apply_abundbool
If true, apply the abundance set in the session to the result.
- log_interpbool
Perform linear interpolation on a logT/logEpsilon grid (default), or linear. Otherwise linear
- Returns:
- retdict
Dictionary containing:
Te, tau, teunit : as input
wavelength : line wavelength (A)
energy : line energy (keV)
epsilon : emissivity in ph cm^3 s-1 (or ph cm^5 s^-1 if apply_aeff=True). First index is temperature, second is tau.
- return_linelist(Te, specrange, specunit='A', teunit='keV', apply_aeff=False, nearest=False, apply_binwidth=False)[source]
Get the list of line emissivities vs wavelengths
- Parameters:
- Tefloat
Temperature in keV or K
- specrange[float, float]
Minimum and maximum values for interval in which to search
- specunit{‘Angstrom’,’keV’}
Units for specrange
- teunit{‘keV’ , ‘K’}
Units of te (kev or K, default keV)
- apply_aeffbool
If true, apply the effective area to the lines in the linelist to modify their intensities.
- nearestbool
Return spectrum at nearest tabulated temperature, without interpolation
- apply_binwidthbool
Divide the line emissivity by the width of the bin they occupy to give emissivity per angstrom or per keV.
- Returns:
- linelistarray(dtype)
The list of lines with lambda (A), energy (keV), epsilon (ph cm3 s-1), epsilon_aeff (ph cm5 s-1) ion (string) and upper & lower levels.
- return_spectrum(te, teunit='keV', nearest=False, get_nearest_t=False, log_interp=True, dolines=True, docont=True, dopseudo=True)[source]
Get the spectrum at an exact temperature. Interpolates between 2 neighbouring spectra
Finds HDU with kT closest to desired kT in given line or coco file.
Opens the line or coco file, and looks for the header unit with temperature closest to te. Use result as index input to make_spectrum
- Parameters:
- tefloat
Temperature in keV or K
- teunit{‘keV’ , ‘K’}
Units of te (kev or K, default keV)
- nearestbool
If set, return the spectrum from the nearest tabulated temperature in the file, without interpolation
- get_nearest_tbool
If set, and nearest set, return the nearest tabulated temperature as well as the spectrum.
- log_interpbool
Perform linear interpolation on a logT/logEpsilon grid, instead of linear.
- dolinesbool
Calculate line emission (default True)
- docontbool
Calculate Continuum emission (default True)
- dopseudobool
Calculate PseudoContinuum (weak line) emission (default True)
- Returns:
- spectrumarray(float)
The spectrum in photons cm^5 s^-1 bin^-1, with the response, or photons cm^3 s^-1 bin^-1 if raw is set.
- nearest_Tfloat, optional
If get_nearest_t is set, return the actual temperature this corresponds to. Units are same as teunit
- set_abund(elements, abund)[source]
Set the elemental abundance, relative to the abundset. Defaults to 1.0 for everything
- Parameters:
- elementsint or array_like(int)
The elements to change the abundance of
- abundfloat or array_like(float)
The new abundances. If only 1 value, set all elements to this abundance Otherwise, should be of same length as elements.
- Returns:
- None
Examples
Set the abundance of iron to 0.5
>>> myspec.set_abund(26, 0.5)
Set the abundance of iron and nickel to 0.1 and 0.2 respectively
>>> myspec.set_abund([26, 28], [0.1,0.2])
Set the abundance of oxygen, neon, magnesium and iron to 0.1
>>> myspec.set_abund([8,10,12,26],0.1)
- set_abundset(abundstring=None)[source]
Set the abundance set.
- Parameters:
- abundstringstring
The abundance string (e.g. “AG89”, “uniform”. Case insensitive. See atomdb.get_abundance for list of possible abundances
- Returns:
- none
updates self.abundset and self.abundsetvector.
- set_broadening(thermal_broadening, broaden_limit=False, velocity_broadening=0.0, velocity_broadening_units='km/s', thermal_broaden_temperature=None, teunit='keV')[source]
Turn on or off thermal broadening, and the emissivity limit for thost lines
- Parameters:
- thermal_broadeningbool
If true, turn on broadening. If False, turn it off.
- broaden_limitfloat
The emissivity limit for lines to be broadened. If False, this value will not be updated.
- velocity_broadeningfloat
velocity broadening to apply. If <=0, not applied
- velocity_broadening_unitsstring
Units of velocity_broadening. ‘km/s’ is default and only value so far.
- thermal_broaden_temperaturefloat
If not None, use this temperature for all line broadening instead of plasma temperature (keV)
- Notes
- —–
- Updates attributes thermal_broadening, broaden_limit, velocity_broadening,
- velocity_broadening_units
- set_eebrems(do_eebrems)[source]
Set whether to do the electron-electron bremsstrahlung in the spectrum
- Parameters:
- do_eebremsbool
If True, include the electron-electron bremsstrahlung calculation If False, don’t (this is default, and matches XSPEC behaviour)
- set_response(rmf, arf=False, raw=False, sparse=False)[source]
Set the response. rmf, arf can either be the filenames or the opened files (latter is faster if called repeatedly)
Amends the following items:
- self.rmffilestring
The rmf file name
- self.rmfstring
The response matrix
- self.arffilestring
The arf file name
- self.arfstring
The arf data
- self.specbinsarray(float)
The spectral bins on which to calculate the spectrum (keV or A)
- self.specbin_unitsstring [‘A’,’keV’]
Units of specbins
- self.ebinsarray(float)
The spectral bins on which to calculate the spectrum (keV).
- self.ebins_outarray(float)
The spectral bins on which to return the spectrum (keV). Can be different from specbins depending on the spectrum
- self.response_setbool
A response has been loaded
- self.response_typestring
‘raw’, ‘standard’ or ‘sparse’ depending on the type implemented
- self.specbins_setbool
The spectral bins are set
- self.ebins_checksumstring
The md5checksum of the specbins
- Parameters:
- rmf: string or HDUlist or array
The response matrix file or energy bins (see raw k/w)
- arf: string or HDUlist
The ancillary response file
- rawbool
If true, the rmf variable contains the energy bin edges (keV) for all the bins, and each bin has a perfect response. This is effectively a dummy response
- sparsebool
If true, the rmf is stored as a sparse matrix and solved using sparse matrix algebra. Useful for large RMFs (e.g. XRISM). Tests show accuracy to 1 part in 10^{15}. Ignored if raw==True
- Returns:
- none
- class pyatomdb.spectrum.CIESession_RS(oscfile='$ATOMDB/apec_osc.fits', cocofile='$ATOMDB/apec_coco.fits', elements=False, abundset='AG89')[source]
Bases:
CIESession
Load and generate a collisional ionization equilibrium spectrum
- Parameters:
- linefilestring or HDUList, optional
The line emissivity data file (either name or already open)
- cocofilestring or HDUList, optional
The continuum emissivity data file (either name or already open)
- elementsiterable of int
Elements to include, listed by atomic number. if not set, include all.
- abundsetstring
The abundance set to use. Default AG89.
Examples
>> import pyatomdb,atomdb
Create a session instance:
>>> cie_rs = pyatomdb.spectrum.CIESession_RS()
Set up the responses, in this case a dummy response from 0.1 to 10 keV, with area 1cm^2 in each bin
>>> ebins = numpy.linspace(0.1,10,1000) >>> cie_rs.set_response(ebins, raw=True)
(Alternatively, for a real response file, s.set_response(rmffile, arf=arffile)
Turn on thermal broadening
>>> cie_rs.set_broadening(True) Will thermally broaden lines with emissivity > 1.000000e-18 ph cm3 s-1
Specify abundance
>>>Abund= atomdb.get_abundance(abundfile=False, abundset=’AG89’, element=[-1], datacache=False, settings=False, show=False)
Return spectrum at 1.0keV and density = 1cm^-3
>>> spec = cie_rs.return_spectrum(1.0, 1.0, Ab=Abund)
spec is in photons cm^5 s^-1 bin^-1; ebins are the bin edges (so spec is 1 element shorter than ebins)
- Attributes:
- datacachedict
Any Atomdb FITS files which have to be opened are stored here
- spectraCIESpectra
Object storing the actual spectral data
- elementsiterable of int
Elements to include, listed by atomic number. if not set, include all.
- default_abundsetstring
The abundance set used for the original emissivity file calculation
- abundsetstring
The abundance set to be used for the returned spectrum
- abundsetvectorarray_like(float)
The relative abundance between default_abundset and abundset for each element
- response_setbool
Have we loaded a response (or set a dummy response)
- dolinesbool
Calculate line emission (default True)
- docontbool
Calculate Continuum emission (default True)
- dopseudobool
Calculate PseudoContinuum (weak line) emission (default True)
- broaden_limitfloat
Apply broadening to lines with epsilon > this value (ph cm3 s-1)
- thermal_broadeningbool
Apply thermal broadening to lines (default = False)
- velocity_broadeningfloat
Apply velocity broadening with this velocity (km/s). If <=0, do not apply.
- format_linelist(linelist)
Print out the line list in a nice format, with all the transitions described in more detail
- Parameters:
- linelistnumpy.array(dtype=linelist)
The linelist returned from “return_linelist”
- Returns
- ——-
- strthe formatted linelist
- return_line_emissivity(Te, Z, z1, up, lo, specunit='A', teunit='keV', apply_aeff=False, apply_abund=True, log_interp=True)[source]
Get line emissivity as function of Te.
- Parameters:
- Tefloat or array(float)
Temperature(s) in keV or K
- Zint
nuclear charge of element
- z1int
ion charge +1 of ion
- upint
upper level for transition
- loint
lower level for transition
- specunit{‘Angstrom’,’keV’}
Units for wavelength or energy (a returned value)
- teunit{‘keV’ , ‘K’}
Units of Te (kev or K, default keV)
- apply_aeffbool
If true, apply the effective area to the line emissivity in the linelist to modify their intensities.
- apply_abundbool
If true, apply the abundance set in the session to the result.
- log_interpbool
Perform linear interpolation on a logT/logEpsilon grid (default), or linear. Otherwise linear
- Returns:
- retdict
Dictionary containing:
Te, tau, teunit : as input
wavelength : line wavelength (A)
energy : line energy (keV)
epsilon : emissivity in ph cm^3 s-1 (or ph cm^5 s^-1 if apply_aeff=True). First index is temperature, second is tau.
- return_linelist(Te, specrange, specunit='A', teunit='keV', apply_aeff=False, nearest=False, apply_binwidth=False)[source]
Get the list of line emissivities vs wavelengths
- Parameters:
- Tefloat
Temperature in keV or K
- specrange[float, float]
Minimum and maximum values for interval in which to search
- specunit{‘Angstrom’,’keV’}
Units for specrange
- teunit{‘keV’ , ‘K’}
Units of te (kev or K, default keV)
- apply_aeffbool
If true, apply the effective area to the lines in the linelist to modify their intensities.
- nearestbool
Return spectrum at nearest tabulated temperature, without interpolation
- apply_binwidthbool
Divide the line emissivity by the width of the bin they occupy to give emissivity per angstrom or per keV.
- Returns:
- linelistarray(dtype)
The list of lines with lambda (A), energy (keV), epsilon (ph cm3 s-1), epsilon_aeff (ph cm5 s-1) ion (string) and upper & lower levels.
- return_spectrum(te, N_e, Ab, teunit='keV', nearest=False, get_nearest_t=False, log_interp=True, dolines=True, docont=True, dopseudo=True)[source]
Get the spectrum at an exact temperature. Interpolates between 2 neighbouring spectra
Finds HDU with kT closest to desired kT in given line or coco file.
Opens the line or coco file, and looks for the header unit with temperature closest to te. Use result as index input to make_spectrum
- Parameters:
- tefloat
Temperature in keV or K
- teunit{‘keV’ , ‘K’}
Units of te (kev or K, default keV)
- nearestbool
If set, return the spectrum from the nearest tabulated temperature in the file, without interpolation
- get_nearest_tbool
If set, and nearest set, return the nearest tabulated temperature as well as the spectrum.
- log_interpbool
Perform linear interpolation on a logT/logEpsilon grid, instead of linear.
- dolinesbool
Calculate line emission (default True)
- docontbool
Calculate Continuum emission (default True)
- dopseudobool
Calculate PseudoContinuum (weak line) emission (default True)
- Returns:
- spectrumarray(float)
The spectrum in photons cm^5 s^-1 bin^-1, with the response, or photons cm^3 s^-1 bin^-1 if raw is set.
- nearest_Tfloat, optional
If get_nearest_t is set, return the actual temperature this corresponds to. Units are same as teunit
- set_abund(elements, abund)
Set the elemental abundance, relative to the abundset. Defaults to 1.0 for everything
- Parameters:
- elementsint or array_like(int)
The elements to change the abundance of
- abundfloat or array_like(float)
The new abundances. If only 1 value, set all elements to this abundance Otherwise, should be of same length as elements.
- Returns:
- None
Examples
Set the abundance of iron to 0.5
>>> myspec.set_abund(26, 0.5)
Set the abundance of iron and nickel to 0.1 and 0.2 respectively
>>> myspec.set_abund([26, 28], [0.1,0.2])
Set the abundance of oxygen, neon, magnesium and iron to 0.1
>>> myspec.set_abund([8,10,12,26],0.1)
- set_abundset(abundstring=None)
Set the abundance set.
- Parameters:
- abundstringstring
The abundance string (e.g. “AG89”, “uniform”. Case insensitive. See atomdb.get_abundance for list of possible abundances
- Returns:
- none
updates self.abundset and self.abundsetvector.
- set_broadening(thermal_broadening, broaden_limit=False, velocity_broadening=0.0, velocity_broadening_units='km/s', thermal_broaden_temperature=None, teunit='keV')
Turn on or off thermal broadening, and the emissivity limit for thost lines
- Parameters:
- thermal_broadeningbool
If true, turn on broadening. If False, turn it off.
- broaden_limitfloat
The emissivity limit for lines to be broadened. If False, this value will not be updated.
- velocity_broadeningfloat
velocity broadening to apply. If <=0, not applied
- velocity_broadening_unitsstring
Units of velocity_broadening. ‘km/s’ is default and only value so far.
- thermal_broaden_temperaturefloat
If not None, use this temperature for all line broadening instead of plasma temperature (keV)
- Notes
- —–
- Updates attributes thermal_broadening, broaden_limit, velocity_broadening,
- velocity_broadening_units
- set_eebrems(do_eebrems)
Set whether to do the electron-electron bremsstrahlung in the spectrum
- Parameters:
- do_eebremsbool
If True, include the electron-electron bremsstrahlung calculation If False, don’t (this is default, and matches XSPEC behaviour)
- set_response(rmf, arf=False, raw=False, sparse=False)
Set the response. rmf, arf can either be the filenames or the opened files (latter is faster if called repeatedly)
Amends the following items:
- self.rmffilestring
The rmf file name
- self.rmfstring
The response matrix
- self.arffilestring
The arf file name
- self.arfstring
The arf data
- self.specbinsarray(float)
The spectral bins on which to calculate the spectrum (keV or A)
- self.specbin_unitsstring [‘A’,’keV’]
Units of specbins
- self.ebinsarray(float)
The spectral bins on which to calculate the spectrum (keV).
- self.ebins_outarray(float)
The spectral bins on which to return the spectrum (keV). Can be different from specbins depending on the spectrum
- self.response_setbool
A response has been loaded
- self.response_typestring
‘raw’, ‘standard’ or ‘sparse’ depending on the type implemented
- self.specbins_setbool
The spectral bins are set
- self.ebins_checksumstring
The md5checksum of the specbins
- Parameters:
- rmf: string or HDUlist or array
The response matrix file or energy bins (see raw k/w)
- arf: string or HDUlist
The ancillary response file
- rawbool
If true, the rmf variable contains the energy bin edges (keV) for all the bins, and each bin has a perfect response. This is effectively a dummy response
- sparsebool
If true, the rmf is stored as a sparse matrix and solved using sparse matrix algebra. Useful for large RMFs (e.g. XRISM). Tests show accuracy to 1 part in 10^{15}. Ignored if raw==True
- Returns:
- none
- class pyatomdb.spectrum.NEISession(linefile='$ATOMDB/apec_nei_line.fits', cocofile='$ATOMDB/apec_nei_comp.fits', elements=False, abundset='AG89')[source]
Bases:
CIESession
Load and generate a Non-equilibrium ionization (NEI) spectrum
- Parameters:
- linefilestring or HDUList, optional
The line emissivity data file (either name or already open)
- cocofilestring or HDUList, optional
The continuum emissivity data file (either name or already open)
- elementsiterable of int
Elements to include, listed by atomic number. if not set, include all.
- abundsetstring
The abundance set to use. Default AG89.
Examples
Create a session instance:
>>> s = NEISession()
Set up the responses, in this case a dummy response from 0.1 to 10 keV
>>> ebins = numpy.linspace(0.1,10,1000) >>> s.set_response(ebins, raw=True)
Turn on thermal broadening
>>> s.set_broadening(True) Will thermally broaden lines with emissivity > 1.000000e-18 ph cm3 s-1
Return spectrum at 1.0keV, Ne *t = 1e11 cm^-3 s
>>> spec = s.return_spectrum(1.0, 1e11)
spec is in photons cm^3 s^-1 bin^-1; ebins are the bin edges (so spec is 1 element shorter than ebins)
- Attributes:
- datacachedict
Any Atomdb FITS files which have to be opened are stored here
- spectraNEISpectra
Object storing the actual spectral data
- elementsiterable of int
Elements to include, listed by atomic number. if not set, include all.
- default_abundsetstring
The abundance set used for the original emissivity file calculation
- abundsetstring
The abundance set to be used for the returned spectrum
- abundsetvectorarray_like(float)
The relative abundance between default_abundset and abundset for each element
- response_setbool
Have we loaded a response (or set a dummy response)
- dolinesbool
Calculate line emission (default True)
- docontbool
Calculate Continuum emission (default True)
- dopseudobool
Calculate PseudoContinuum (weak line) emission (default True)
- broaden_limitfloat
Apply broadening to lines with epsilon > this value (ph cm3 s-1)
- thermal_broadeningbool
Apply thermal broadening to lines (default = False)
- velocity_broadeningfloat
Apply velocity broadening with this velocity (km/s). If <=0, do not apply.
- format_linelist(linelist)
Print out the line list in a nice format, with all the transitions described in more detail
- Parameters:
- linelistnumpy.array(dtype=linelist)
The linelist returned from “return_linelist”
- Returns
- ——-
- strthe formatted linelist
- return_line_emissivity(Telist, taulist, Z, z1, up, lo, specunit='A', teunit='keV', apply_aeff=False, apply_abund=True, log_interp=True, init_pop='ionizing', freeze_ion_pop=False)[source]
Get line emissivity as function of Te, tau. Assumes ionization from neutral.
- Parameters:
- Telistfloat or array(float)
Temperature(s) in keV or K
- taulistfloat or array(float)
ionization timescale(s), ne * t (cm^-3 s).
- Zint
nuclear charge of element
- z1int
ion charge +1 of ion
- upint
upper level for transition
- loint
lower level for transition
- specunit{‘Angstrom’,’keV’}
Units for wavelength or energy (a returned value)
- teunit{‘keV’ , ‘K’}
Units of Telist (kev or K, default keV)
- apply_aeffbool
If true, apply the effective area to the line emissivity in the linelist to modify their intensities.
- apply_abundbool
If true, apply the abundance set in the session to the result.
- log_interpbool
Perform linear interpolation on a logT/logEpsilon grid (default), or linear.
- init_popstring or float
if ‘ionizing’ : all ionizing from neutral (so [1,0,0,0…]) if ‘recombining’: all recombining from ionized (so[…0,0,1]) if dict of arrays : the acutal fractional populations (so init_pop[6] is the array for carbon) if single float : the temperature (same units as Te)
- freeze_ion_popbool
If True, skip the ion population calculation, use init_pop as the final pop instead.
- Returns:
- retdict
Dictionary containing:
Te, tau, teunit: as input
wavelength : line wavelength (A)
energy : line energy (keV)
- epsilonemissivity in ph cm^3 s-1 (or ph cm^5 s^-1 if apply_aeff=True)
first index is temperature, second is tau. If Te or Tau was supplied as a scalar, then that index is removed
- return_linelist(Te, tau, specrange, init_pop='ionizing', specunit='A', teunit='keV', apply_aeff=False, by_ion_drv=False, nearest=False, apply_binwidth=False, log_interp=True, freeze_ion_pop=False)[source]
Get the list of line emissivities vs wavelengths
- Parameters:
- Tefloat
Temperature in keV or K
- taufloat
ionization timescale, ne * t (cm^-3 s).
- specrange[float, float]
Minimum and maximum values for interval in which to search
- specunit{‘Angstrom’,’keV’}
Units for specrange
- teunit{‘keV’ , ‘K’}
Units of te (kev or K, default keV)
- apply_aeffbool
If true, apply the effective area to the lines in the linelist to modify their intensities.
- by_ion_drvbool
If true, keep lines which are the same but have different ion_drv separate. Otherwise, merge them.
- nearestbool
calculate spectrum at nearest temperature in linelist, no interpolation. Ionization fraction calculation still exact.
- apply_binwidth
- log_interpbool
Perform linear interpolation on a logT/logEpsilon grid (default), or linear.
- freeze_ion_popbool
If True, skip the ion population calculation, use init_pop as the final pop instead.
- Returns:
- linelistarray(linelist)
The list of lines with lambda (A), energy (keV), epsilon (ph cm3 s-1), epsilon_aeff (ph cm5 s-1) ion (string) and upper & lower levels.
- return_spectrum(Te, tau, init_pop='ionizing', teunit='keV', nearest=False, log_interp=True, freeze_ion_pop=False)[source]
Get the spectrum at an exact temperature. Interpolates between 2 neighbouring spectra
Finds HDU with kT closest to desired kT in given line or coco file.
Opens the line or coco file, and looks for the header unit with temperature closest to te. Use result as index input to make_spectrum
- Parameters:
- Tefloat
Temperature in keV or K
- taufloat
ionization timescale, ne * t (cm^-3 s).
- init_popstring or float
if ‘ionizing’ : all ionizing from neutral (so [1,0,0,0…]) if ‘recombining’: all recombining from ionized (so[…0,0,1]) if dict of arrays : the acutal fractional populations (so init_pop[6] is the array for carbon) if single float : the temperature (same units as Te)
- teunit{‘keV’ , ‘K’}
Units of te (kev or K, default keV)
- nearestbool
If set, return the spectrum from the nearest tabulated temperature in the file, without interpolation
- log_interpbool
Perform linear interpolation on a logT/logEpsilon grid (default), or linear.
- freeze_ion_popbool
If True, skip the ion population calculation, use init_pop as the final pop instead.
- Returns:
- spectrumarray(float)
The spectrum in photons cm^5 s^-1 bin^-1, with the response, or photons cm^3 s^-1 bin^-1 if raw is set.
- nearest_Tfloat, optional
If nearest is set, return the actual temperature this corresponds to. Units are same as teunit
- set_abund(elements, abund)
Set the elemental abundance, relative to the abundset. Defaults to 1.0 for everything
- Parameters:
- elementsint or array_like(int)
The elements to change the abundance of
- abundfloat or array_like(float)
The new abundances. If only 1 value, set all elements to this abundance Otherwise, should be of same length as elements.
- Returns:
- None
Examples
Set the abundance of iron to 0.5
>>> myspec.set_abund(26, 0.5)
Set the abundance of iron and nickel to 0.1 and 0.2 respectively
>>> myspec.set_abund([26, 28], [0.1,0.2])
Set the abundance of oxygen, neon, magnesium and iron to 0.1
>>> myspec.set_abund([8,10,12,26],0.1)
- set_abundset(abundstring=None)
Set the abundance set.
- Parameters:
- abundstringstring
The abundance string (e.g. “AG89”, “uniform”. Case insensitive. See atomdb.get_abundance for list of possible abundances
- Returns:
- none
updates self.abundset and self.abundsetvector.
- set_broadening(thermal_broadening, broaden_limit=False, velocity_broadening=0.0, velocity_broadening_units='km/s', thermal_broaden_temperature=None, teunit='keV')
Turn on or off thermal broadening, and the emissivity limit for thost lines
- Parameters:
- thermal_broadeningbool
If true, turn on broadening. If False, turn it off.
- broaden_limitfloat
The emissivity limit for lines to be broadened. If False, this value will not be updated.
- velocity_broadeningfloat
velocity broadening to apply. If <=0, not applied
- velocity_broadening_unitsstring
Units of velocity_broadening. ‘km/s’ is default and only value so far.
- thermal_broaden_temperaturefloat
If not None, use this temperature for all line broadening instead of plasma temperature (keV)
- Notes
- —–
- Updates attributes thermal_broadening, broaden_limit, velocity_broadening,
- velocity_broadening_units
- set_eebrems(do_eebrems)
Set whether to do the electron-electron bremsstrahlung in the spectrum
- Parameters:
- do_eebremsbool
If True, include the electron-electron bremsstrahlung calculation If False, don’t (this is default, and matches XSPEC behaviour)
- set_response(rmf, arf=False, raw=False, sparse=False)
Set the response. rmf, arf can either be the filenames or the opened files (latter is faster if called repeatedly)
Amends the following items:
- self.rmffilestring
The rmf file name
- self.rmfstring
The response matrix
- self.arffilestring
The arf file name
- self.arfstring
The arf data
- self.specbinsarray(float)
The spectral bins on which to calculate the spectrum (keV or A)
- self.specbin_unitsstring [‘A’,’keV’]
Units of specbins
- self.ebinsarray(float)
The spectral bins on which to calculate the spectrum (keV).
- self.ebins_outarray(float)
The spectral bins on which to return the spectrum (keV). Can be different from specbins depending on the spectrum
- self.response_setbool
A response has been loaded
- self.response_typestring
‘raw’, ‘standard’ or ‘sparse’ depending on the type implemented
- self.specbins_setbool
The spectral bins are set
- self.ebins_checksumstring
The md5checksum of the specbins
- Parameters:
- rmf: string or HDUlist or array
The response matrix file or energy bins (see raw k/w)
- arf: string or HDUlist
The ancillary response file
- rawbool
If true, the rmf variable contains the energy bin edges (keV) for all the bins, and each bin has a perfect response. This is effectively a dummy response
- sparsebool
If true, the rmf is stored as a sparse matrix and solved using sparse matrix algebra. Useful for large RMFs (e.g. XRISM). Tests show accuracy to 1 part in 10^{15}. Ignored if raw==True
- Returns:
- none
- class pyatomdb.spectrum.PShockSession(linefile='$ATOMDB/apec_nei_line.fits', cocofile='$ATOMDB/apec_nei_comp.fits', elements=False, abundset='AG89')[source]
Bases:
NEISession
Load and generate a Parallel Shock (PShock) spectrum
- Parameters:
- linefilestring or HDUList, optional
The line emissivity data file (either name or already open)
- cocofilestring or HDUList, optional
The continuum emissivity data file (either name or already open)
- elementsiterable of int
Elements to include, listed by atomic number. if not set, include all.
- abundsetstring
The abundance set to use. Default AG89.
Examples
Create a session instance:
>>> s=PShockSession()
Set up the responses, in this case a dummy response from 0.1 to 10 keV
>>> ebins = numpy.linspace(0.1,10,1000) >>> s.set_response(ebins, raw=True)
Turn on thermal broadening
>>> s.set_broadening(True) Will thermally broaden lines with emissivity > 1.000000e-18 ph cm3 s-1
Return spectrum at 1.0keV, Tau_u = 1e11, Tau_l = 0.0 (default)
>>> spec = s.return_spectrum(1.0, 1e11)
Same with Tau_l = 1e9
>>> spec = s.return_spectrum(1.0, 1e11, tau_l = 1e9)
spec is in photons cm^3 s^-1 bin^-1; ebins are the bin edges (so spec is 1 element shorter than ebins)
- Attributes:
- datacachedict
Any Atomdb FITS files which have to be opened are stored here
- spectraPShockSpectra
Object storing the actual spectral data
- elementslist(int)
Nuclear charge of elements to include.
- default_abundsetstring
The abundance set used for the original emissivity file calculation
- abundsetstring
The abundance set to be used for the returned spectrum
- abundsetvectorarray_like(float)
The relative abundance between default_abundset and abundset for each element
- response_setbool
Have we loaded a response (or set a dummy response)
- dolinesbool
Calculate line emission (default True)
- docontbool
Calculate Continuum emission (default True)
- dopseudobool
Calculate PseudoContinuum (weak line) emission (default True)
- broaden_limitfloat
Apply broadening to lines with epsilon > this value (ph cm3 s-1)
- thermal_broadeningbool
Apply thermal broadening to lines (default = False)
- velocity_broadeningfloat
Apply velocity broadening with this velocity (km/s). If <=0, do not apply.
- format_linelist(linelist)
Print out the line list in a nice format, with all the transitions described in more detail
- Parameters:
- linelistnumpy.array(dtype=linelist)
The linelist returned from “return_linelist”
- Returns
- ——-
- strthe formatted linelist
- return_line_emissivity(Telist, tau_ulist, Z, z1, up, lo, tau_llist=0.0, specunit='A', teunit='keV', apply_aeff=False, apply_abund=True, log_interp=True, init_pop='ionizing')[source]
Get line emissivity as function of Te, tau. Assumes ionization from neutral.
- Parameters:
- Telistfloat or array(float)
Temperature(s) in keV or K
- tau_ulistfloat
ionization timescale(s), ne * t (cm^-3 s).
- Zint
nuclear charge of element
- z1int
ion charge +1 of ion
- upint
upper level for transition
- loint
lower level for transition
- tau_llistfloat
lower limit of ionization timescale(s), ne * t (cm^-3 s).
- specunit{‘Angstrom’,’keV’}
Units for wavelength or energy (a returned value)
- teunit{‘keV’ , ‘K’}
Units of Telist (kev or K, default keV)
- apply_aeffbool
If true, apply the effective area to the line emissivity in the linelist to modify their intensities.
- apply_abundbool
If true, apply the abundance set in the session to the result.
- log_interpbool
Perform linear interpolation on a logT/logEpsilon grid (default), or linear.
- Returns:
- retdict
Dictionary containing:
Te, tau, teunit: as input
wavelength : line wavelength (A)
energy : line energy (keV)
- epsilonemissivity in ph cm^3 s-1 (or ph cm^5 s^-1 if apply_aeff=True)
first index is temperature, second is tau.
- return_linelist(Te, tau_u, specrange, tau_l=0.0, init_pop='ionizing', specunit='A', teunit='keV', apply_aeff=False, by_ion_drv=False, nearest=False, apply_binwidth=False)[source]
Get the list of line emissivities vs wavelengths
- Parameters:
- Tefloat
Temperature in keV or K
- tau_ufloat
Upper limit of ionization timescale, ne * t (cm^-3 s).
- specrange[float, float]
Minimum and maximum values for interval in which to search
- tau_l
lower limit of ionization timescale, ne * t (cm^-3 s) (defaults to 0)
- init_popstring or float
if ‘ionizing’ : all ionizing from neutral (so [1,0,0,0…]) if ‘recombining’: all recombining from ionized (so[…0,0,1]) if dict of arrays : the acutal fractional populations (so init_pop[6] is the array for carbon) if single float : the temperature (same units as Te)
- specunit{‘Angstrom’,’keV’}
Units for specrange
- teunit{‘keV’ , ‘K’}
Units of te (kev or K, default keV)
- apply_aeffbool
If true, apply the effective area to the lines in the linelist to modify their intensities.
- by_ion_drvbool
If true, keep lines which are the same but have different ion_drv separate. Otherwise, merge them.
- nearest
- apply_binwidth
- Returns:
- linelistarray(dtype)
The list of lines with lambda (A), energy (keV), epsilon (ph cm3 s-1), epsilon_aeff (ph cm5 s-1) ion (string) and upper & lower levels.
- return_spectrum(Te, tau_u, tau_l=0.0, init_pop='ionizing', teunit='keV', nearest=False, log_interp=True)[source]
Get the spectrum at an exact temperature. Interpolates between 2 neighbouring spectra
Finds HDU with kT closest to desired kT in given line or coco file.
Opens the line or coco file, and looks for the header unit with temperature closest to te. Use result as index input to make_spectrum
- Parameters:
- Tefloat
Temperature in keV or K
- tau_ufloat
Upper limit of ionization timescale, ne * t (cm^-3 s).
- tau_l
lower limit of ionization timescale, ne * t (cm^-3 s) (defaults to 0)
- init_popstring or float
if ‘ionizing’ : all ionizing from neutral (so [1,0,0,0…]) if ‘recombining’: all recombining from ionized (so[…0,0,1]) if dict of arrays : the acutal fractional populations (so init_pop[6] is the array for carbon) if single float : the temperature (same units as Te)
- teunit{‘keV’ , ‘K’}
Units of te (kev or K, default keV)
- nearestbool
If set, return the spectrum from the nearest tabulated temperature in the file, without interpolation
- log_interpbool
Perform linear interpolation on a logT/logEpsilon grid (default), or linear.
- Returns:
- spectrumarray(float)
The spectrum in photons cm^5 s^-1 bin^-1, with the response, or photons cm^3 s^-1 bin^-1 if raw is set.
- set_abund(elements, abund)
Set the elemental abundance, relative to the abundset. Defaults to 1.0 for everything
- Parameters:
- elementsint or array_like(int)
The elements to change the abundance of
- abundfloat or array_like(float)
The new abundances. If only 1 value, set all elements to this abundance Otherwise, should be of same length as elements.
- Returns:
- None
Examples
Set the abundance of iron to 0.5
>>> myspec.set_abund(26, 0.5)
Set the abundance of iron and nickel to 0.1 and 0.2 respectively
>>> myspec.set_abund([26, 28], [0.1,0.2])
Set the abundance of oxygen, neon, magnesium and iron to 0.1
>>> myspec.set_abund([8,10,12,26],0.1)
- set_abundset(abundstring=None)
Set the abundance set.
- Parameters:
- abundstringstring
The abundance string (e.g. “AG89”, “uniform”. Case insensitive. See atomdb.get_abundance for list of possible abundances
- Returns:
- none
updates self.abundset and self.abundsetvector.
- set_broadening(thermal_broadening, broaden_limit=False, velocity_broadening=0.0, velocity_broadening_units='km/s', thermal_broaden_temperature=None, teunit='keV')
Turn on or off thermal broadening, and the emissivity limit for thost lines
- Parameters:
- thermal_broadeningbool
If true, turn on broadening. If False, turn it off.
- broaden_limitfloat
The emissivity limit for lines to be broadened. If False, this value will not be updated.
- velocity_broadeningfloat
velocity broadening to apply. If <=0, not applied
- velocity_broadening_unitsstring
Units of velocity_broadening. ‘km/s’ is default and only value so far.
- thermal_broaden_temperaturefloat
If not None, use this temperature for all line broadening instead of plasma temperature (keV)
- Notes
- —–
- Updates attributes thermal_broadening, broaden_limit, velocity_broadening,
- velocity_broadening_units
- set_eebrems(do_eebrems)
Set whether to do the electron-electron bremsstrahlung in the spectrum
- Parameters:
- do_eebremsbool
If True, include the electron-electron bremsstrahlung calculation If False, don’t (this is default, and matches XSPEC behaviour)
- set_response(rmf, arf=False, raw=False, sparse=False)
Set the response. rmf, arf can either be the filenames or the opened files (latter is faster if called repeatedly)
Amends the following items:
- self.rmffilestring
The rmf file name
- self.rmfstring
The response matrix
- self.arffilestring
The arf file name
- self.arfstring
The arf data
- self.specbinsarray(float)
The spectral bins on which to calculate the spectrum (keV or A)
- self.specbin_unitsstring [‘A’,’keV’]
Units of specbins
- self.ebinsarray(float)
The spectral bins on which to calculate the spectrum (keV).
- self.ebins_outarray(float)
The spectral bins on which to return the spectrum (keV). Can be different from specbins depending on the spectrum
- self.response_setbool
A response has been loaded
- self.response_typestring
‘raw’, ‘standard’ or ‘sparse’ depending on the type implemented
- self.specbins_setbool
The spectral bins are set
- self.ebins_checksumstring
The md5checksum of the specbins
- Parameters:
- rmf: string or HDUlist or array
The response matrix file or energy bins (see raw k/w)
- arf: string or HDUlist
The ancillary response file
- rawbool
If true, the rmf variable contains the energy bin edges (keV) for all the bins, and each bin has a perfect response. This is effectively a dummy response
- sparsebool
If true, the rmf is stored as a sparse matrix and solved using sparse matrix algebra. Useful for large RMFs (e.g. XRISM). Tests show accuracy to 1 part in 10^{15}. Ignored if raw==True
- Returns:
- none
- pyatomdb.spectrum.calc_ee_brems_spec(ebins, Te, dens, teunit='keV')[source]
Calculate the electron-electron bremstrahlung per-bin emissivity
- Parameters:
- ebinsarray(float)
The spectral bin edges in energy order (keV)
- Tefloat
Electron temperature (default, keV)
- densfloat
The electron density (cm^-3)
- teunitstring
Units of Te (keV by default, K also allowed)
- Returns:
- eebremsarray(float)
electron-electron bremsstrahlung in ph cm^3 s^-1 bin^-1