Contains the SpectrumFactory class, which is
the core of the RADIS Line-by-Line module.
Examples
Calculate a CO Spectrum, fetching the lines from HITRAN
# This is how you get a spectrum (see spectrum.py for front-end functions# that do just that)sf=SpectrumFactory(2125,2249.9,molecule='CO',isotope=1,cutoff=1e-30,# for faster calculations. See# `plot_linestrength_hist` for more details**kwargs)sf.fetch_databank()# auto download from HITRANs=sf.eq_spectrum(Tgas=300)s.plot('abscoeff')# opacity# Here we get some extra informations:s.plot('radiance',wunit='nm',Iunit='µW/cm2/sr/nm',# Iunit is arbitrary. Use whatever makes senseshow_points=True)# show_points to have an idea of the resolution
digraph inheritanceb7d785f131 {
bgcolor=transparent;
rankdir=LR;
size="8.0, 12.0";
"BandFactory" [URL="radis.lbl.bands.html#radis.lbl.bands.BandFactory",fillcolor=white,fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",fontsize=10,height=0.25,shape=box,style="setlinewidth(0.5),filled",target="_top",tooltip="See Also"];
"BroadenFactory" -> "BandFactory" [arrowsize=0.5,style="setlinewidth(0.5)"];
"BaseFactory" [URL="radis.lbl.base.html#radis.lbl.base.BaseFactory",fillcolor=white,fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",fontsize=10,height=0.25,shape=box,style="setlinewidth(0.5),filled",target="_top"];
"DatabankLoader" -> "BaseFactory" [arrowsize=0.5,style="setlinewidth(0.5)"];
"BroadenFactory" [URL="radis.lbl.broadening.html#radis.lbl.broadening.BroadenFactory",fillcolor=white,fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",fontsize=10,height=0.25,shape=box,style="setlinewidth(0.5),filled",target="_top",tooltip="A class that holds all broadening methods."];
"BaseFactory" -> "BroadenFactory" [arrowsize=0.5,style="setlinewidth(0.5)"];
"DatabankLoader" [URL="radis.lbl.loader.html#radis.lbl.loader.DatabankLoader",fillcolor=white,fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",fontsize=10,height=0.25,shape=box,style="setlinewidth(0.5),filled",target="_top",tooltip=".. inheritance-diagram:: radis.lbl.factory.SpectrumFactory"];
"SpectrumFactory" [URL="#radis.lbl.factory.SpectrumFactory",fillcolor=white,fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",fontsize=10,height=0.25,shape=box,style="setlinewidth(0.5),filled",target="_top",tooltip="A class to put together all functions related to loading CDSD and HITRAN"];
"BandFactory" -> "SpectrumFactory" [arrowsize=0.5,style="setlinewidth(0.5)"];
}
Notes
Debugging:
Remember one can access all locals variable in a function either by running
%debug after a crash in IPython. Or by adding the line globals().update(locals())
at any point
Performance:
Fast version: iterate over chunks of dataframe
Note that we can’t use the full dataframe because it’s too big and takes too much memory
See Performance for more details.
for Developers:
To implement new database formats, see the databases parsers in cdsd.py / hitran.py,
and the partition function interpolators / calculators methods of SpectrumFactory:
_build_partition_function_calculator() and
_build_partition_function_interpolator()
A class to put together all functions related to loading CDSD and HITRAN
databases, calculating the broadenings, and summing over all the lines.
Parameters:
wmin, wmax (float or Quantity) – a hybrid parameter which can stand for minimum (maximum) wavenumber or minimum
(maximum) wavelength depending upon the unit accompanying it. If dimensionless,
wunit is considered as the accompanying unit.
wunit ('nm', 'cm-1') – the unit accompanying wmin and wmax. Can only be passed with wmin
and wmax. Default is "cm-1".
wavenum_min, wavenum_max (float(cm^-1) or Quantity) – minimum (maximum) wavenumber to be processed in \(cm^{-1}\).
use astropy.units to specify arbitrary inverse-length units.
wavelength_min, wavelength_max (float(nm) or Quantity) – minimum (maximum) wavelength to be processed in \(nm\). This wavelength
can be in 'air' or 'vacuum' depending on the value of the parameter
medium=.
use astropy.units to specify arbitrary length units.
pressure (float(bar) or Quantity) – partial pressure of gas in bar. Default 1.01325 (1 atm).
use astropy.units to specify arbitrary pressure units.
For example, 1013.25*u.mbar.
mole_fraction (float [ 0 - 1]) – species mole fraction. Default 1. Note that the rest of the gas
is considered to be air for collisional broadening.
path_length (float(cm) or Quantity) – path length in cm. Default 1.
use astropy.units to specify arbitrary length units.
The positive or neutral atomic species (negative ions aren’t supported). It may be given in spectroscopic notation or any form that can be converted by to_conventional_name()
Default None.
isotope (int, list, str of the form '1,2', or 'all') – isotope id
For molecules, this is the isotopologue ID (sorted by relative density: (eg: 1: CO2-626, 2: CO2-636 for CO2) - see [HITRAN-2020] documentation for isotope list for all species.
For atoms, use the isotope number of the isotope (the total number of protons and neutrons in the nucleus) - use 0 to select rows where the isotope is unspecified, in which case the standard atomic weight from the periodictable module is used when mass is required.
If 'all',
all isotopes in database are used (this may result in larger computation
times!).
Default 'all'
medium ('air', 'vacuum') – propagating medium when giving inputs with 'wavenum_min', 'wavenum_max'.
Does not change anything when giving inputs in wavenumber. Default 'air'
diluent (str or dictionary) – can be a string of a single diluent or a dictionary containing diluent
name as key and its mole_fraction as value.
If left unspecified, it defaults to 'air' for molecules and atomic hydrogen ‘H’ for atoms.
For free electrons, use the symbol ‘e-’. Currently, only H, H2, H2, and e- are supported for atoms - any other diluents have no effect besides diluting the mole fractions of the other constituents.
Other Parameters:
Tref (K) – Reference temperature for calculations (linestrength temperature
correction). HITRAN database uses 296 Kelvin. Default 296 K
self_absorption (boolean) – Compute self absorption. If False, spectra are optically thin. Default True.
truncation (float (\(cm^{-1}\))) – Half-width over which to compute the lineshape, i.e. lines are truncated
on each side after truncation (\(cm^{-1}\)) from the line center.
If None, use no truncation (lineshapes spread on the full spectral range).
Default is 50\(cm^{-1}\)
Note
Large values (> 50) can induce a performance drop (computation of lineshape
typically scale as \(~truncation ^2\) ). The default 50 was
chosen to maintain a good accuracy, and still exhibit the sub-Lorentzian
behavior of most lines far (few hundreds \(cm^{-1}\)) from the line center.
neighbour_lines (float (\(cm^{-1}\))) – The calculated spectral range is increased (by neighbour_lines cm-1
on each side) to take into account overlaps from out-of-range lines.
Default is 0\(cm^{-1}\).
wstep (float (cm-1) or 'auto') – Resolution of wavenumber grid. Default 0.01 cm-1.
If 'auto', it is ensured that there
are slightly more points for each linewidth than the value of "GRIDPOINTS_PER_LINEWIDTH_WARN_THRESHOLD"
in radis.config (~/radis.json)
Note
wstep = ‘auto’ is optimized for performances while ensuring accuracy,
but is still experimental in 0.9.30. Feedback welcome!
cutoff (float (~ unit of Linestrength: cm-1/(#.cm-2))) – discard linestrengths that are lower that this, to reduce calculation
times. 1e-27 is what is generally used to generate databases such as
CDSD. If 0, no cutoff. Default 1e-27.
parsum_mode (‘full summation’, ‘tabulation’) – how to compute partition functions, at nonequilibrium or when partition
function are not already tabulated. 'fullsummation' : sums over all
(potentially millions) of rovibrational levels. 'tabulation' :
builds an on-the-fly tabulation of rovibrational levels (500 - 4000x faster
and usually accurate within 0.1%). Default fullsummation'
Note
parsum_mode= ‘tabulation’ is new in 0.9.30, and makes nonequilibrium
calculations of small spectra extremely fast. Will become the default
after 0.9.31.
pseudo_continuum_threshold (float) – if not 0, first calculate a rough approximation of the spectrum, then
moves all lines whose linestrength intensity is less than this threshold
of the maximum in a semi-continuum. Values above 0.01 can yield significant
errors, mostly in highly populated areas. 80% of the lines can typically
be moved in a continuum, resulting in 5 times faster spectra. If 0,
no semi-continuum is used. Default 0.
save_memory (boolean) – if True, removes databases calculated by intermediate functions (for
instance, delete the full database once the linestrength cutoff criteria
was applied). This saves some memory but requires to reload the database
and recalculate the linestrength for each new parameter. Default False.
export_populations ('vib', 'rovib', None) – if not None, store populations in Spectrum. Either store vibrational
populations (‘vib’) or rovibrational populations (‘rovib’). Default None
export_lines (boolean) – if True, saves details of all calculated lines in Spectrum. This is
necessary to later use line_survey(),
but can take some space. Default False.
chunksize (int, or None) – Splits the lines database in several chunks during calculation, else
the multiplication of lines over all spectral range takes too much memory
and slows the system down. Chunksize let you change the default chunk
size. If None, all lines are processed directly. Usually faster but
can create memory problems. Default None
folding_thresh (float) – Folding is a correction procedure that is applied when the lineshape is calculated with
the fft broadening method and the linewidth is comparable to wstep, that prevents
sinc(v) modulation of the lineshape. Folding continues until the lineshape intensity
is below folding_threshold. Setting to 1 or higher effectively disables folding correction.
Range: 0.0 < folding_thresh <= 1.0
Default: 1e-6
zero_padding (int) – Zero padding is used in conjunction with the fft broadening method to prevent circular
convolution at the cost of performance. When set to -1, padding is set equal to the spectrum length,
which guarantees a linear convolution.
broadening_method ("voigt_poly", "convolve", "fft") – Available methods:
"voigt_poly": polynomial approximation of a Voigt profile.
This is the fastest method in most cases (default).
"convolve": independent calculation of Doppler (Gaussian) and
collisional broadening (Lorentzian), then convolution in real
space. This is slower but is an exact definition of a Voigt and
can be used as a reference for comparison with other methods.
"fft": fast Fourier transform convolution.
Only available with optimization="simple" or optimization="min-RMS".
Because LDM convolves all lines at once,
this method is often faster than real-space convolutions on large arrays
("voigt_poly", "convolve").
This SpectrumFactory parameter can be manually adjusted a posteriori with:
By default, optimization="simple" and broadening_method="voigt_poly".
warnings (bool, or one of ['warn','error','ignore'], dict) – If one of ['warn','error','ignore'], set the default behaviour
for all warnings. Can also be a dictionary to set specific warnings only.
Example:
verbose (boolean, or int) – If False, stays quiet. If True, tells what is going on.
If >=2, gives more detailed messages (for instance, details of
calculation times). Default True.
lbfunc (callable) –
An alternative function to be used instead of the default in calculating Lorentzian broadening, which receives the following:
df: the dataframe self.df1 containing the quantities used for calculating the spectrum
pressure_atm: self.pressure in units of atmospheric pressure (1.01325 bar)
mole_fraction: self.input.mole_fraction, the mole fraction of the species for which the spectrum is being calculated
Tgas: self.input.Tgas, gas temperature in K
Tref: self.input.Tref, reference temperature for calculations in K
diluent: self._diluent, the dictionary of diluents giving the mole fraction of each
diluent_broadening_coeff: a dictionary of the broadening coefficients for each diluent
isneutral: When calculating the spectrum of an atomic species, whether or not it is neutral (always None for molecules)
Returns:
gamma_lb, shift - The total Lorentzian HWHM [\(cm^{-1}\)], and the shift [\(cm^{-1}\)] to be subtracted from the wavenumber array to account for lineshift. If setting the lineshift here is not desired, the 2nd return object can be anything for which bool(shift)==False like None. gamma_lb must be array-like but can also be a vaex expression if the dataframe type is vaex.
If unspecified, the broadening is handled by default by gamma_vald3() for atoms when using the Kurucz databank, and pressure_broadening_HWHM() for molecules.
For the NIST databank, the lbfunc parameter is compulsory as NIST doesn’t provide broadening parameters.
pfsource (string) – The source for the partition function tables for an interpolator or energy level tables for a calculator. Sources implemented so far are ‘barklem’ and ‘kurucz’ for the former, and ‘nist’ for the latter. ‘default’ is currently ‘nist’. The pfsource can be changed post-initialisation using the set_atomic_partition_functions() method. See the provided example for more details.
potential_lowering (float (cm-1/Zeff**2)) – The value of potential lowering, only relevant when pfsource is ‘kurucz’ as it depends on both temperature and potential lowering. Can be changed on the fly by setting sf.input.potential_lowering. Allowed values are typically: -500, -1000, -2000, -4000, -8000, -16000, -32000.
Again, see the provided example for more details.
fromradisimportSpectrumFactoryfromastropyimportunitsasusf=SpectrumFactory(wavelength_min=4165*u.nm,wavelength_max=4200*u.nm,isotope='1,2',truncation=10,# cm-1optimization=None,medium='vacuum',verbose=1,# more for more details)sf.load_databank('HITRAN-CO2-TEST')# predefined in ~/radis.jsons=sf.eq_spectrum(Tgas=300*u.K,path_length=1*u.cm)s.rescale_path_length(0.01)# cms.plot('radiance_noslit',Iunit='µW/cm2/sr/nm')
initial line database after loading.
If for any reason, you want to manipulate the line database manually (for instance, keeping only lines emitting
by a particular level), you need to access the df0 attribute of
SpectrumFactory.
Warning
never overwrite the df0 attribute, else some metadata may be lost in the process.
Only use inplace operations. If reducing the number of lines, add
a df0.reset_index()
For instance:
fromradisimportSpectrumFactorysf=SpectrumFactory(wavenum_min=2150.4,wavenum_max=2151.4,pressure=1,isotope=1)sf.load_databank('HITRAN-CO-TEST')sf.df0.drop(sf.df0[sf.df0.vu!=1].index,inplace=True)# keep lines emitted by v'=1 onlysf.eq_spectrum(Tgas=3000,name='vu=1').plot()
df0 contains the lines as they are loaded from the database.
df1 is generated during the spectrum calculation, after the
line database reduction steps, population calculation, and scaling of intensity and broadening parameters
with the calculated conditions.
mole_fraction (float) – database species mole fraction. If None, Factory mole fraction is used.
path_length (float or Quantity) – slab size (cm). If None, the default Factory
path_length is used.
pressure (float or Quantity) – pressure (bar). If None, the default Factory
pressure is used.
name (str) – output Spectrum name (useful in batch)
Other Parameters:
device_id (int, str) – Select the GPU device. If int, specifies the device index, which is printed for convenience during GPU initialization with backend=’vulkan’ (default).
If str, return the first device that includes the specified string (case-insensitive). If not found, return the device at index 0.
default = 0
exit_gpu (bool) – Specifies whether the GPU app should be exited after producing the spectrum. Usually this is undesirable, because the GPU
computations start to benefit after the first spectrum is produced by calling s.recalc_gpu(). See also recalc_gpu()
default = False
backend (str) – Since version 0.16, only 'vulkan' backend is supported.
In previous versions, 'gpu-cuda' and 'cpu-cuda' were available to switch to a CUDA backend,
but this has been deprecated in favor of the Vulkan backend.
.. warning:: – The backend parameter is deprecated. Only the Vulkan backend is supported.
Fit an experimental spectrum with an arbitrary model and an arbitrary
number of fit parameters. This method calls fit_legacy()
which is still functional. However, we recommend using fit_spectrum().
Parameters:
s_exp (Spectrum) – experimental spectrum. Should have only spectral array only. Use
take(), e.g:
sf.fit_legacy(s_exp.take('transmittance'))
model (func -> Spectrum) – a line-of-sight model returning a Spectrum. Example :
Tvib12Tvib3Trot_NonLTEModel()
Calculate emission spectrum in non-equilibrium case. Calculates
absorption with broadened linestrength and emission with broadened
Einstein coefficient.
Parameters:
Tvib (float) – vibrational temperature [K]
can be a tuple of float for the special case of more-than-diatomic
molecules (e.g: CO2); only applicable for molecules, not atoms
Trot (float) – rotational temperature [K]; only applicable for molecules, not atoms
Ttrans (float) – translational temperature [K]. If None, translational temperature is
taken as rotational temperature (valid at 1 atm for times above ~ 2ns
which is the RT characteristic time)
Telec (float) – electronic temperature [K]; only implemented for atoms, not molecules
mole_fraction (float) – database species mole fraction. If None, Factory mole fraction is used.
diluent (str or dictionary) – can be a string of a single diluent or a dictionary containing diluent
name as key and its mole_fraction as value
path_length (float or Quantity) – slab size (cm). If None, the default Factory
path_length is used.
pressure (float or Quantity) – pressure (bar). If None, the default Factory
pressure is used.
Other Parameters:
vib_distribution ('boltzmann', 'treanor') – vibrational distribution
rot_distribution ('boltzmann') – rotational distribution
overpopulation (dict, or None) –
add overpopulation factors for given levels:
{level:overpopulation_factor}
name (str) – output Spectrum name (useful in batch)
Multi-vibrational temperature. Below we compare non-LTE spectra of CO2 where all
vibrational temperatures are equal, or where the bending and symmetric modes are in
equilibrium with rotation
fromradisimportSpectrumFactorysf=SpectrumFactory(wavenum_min=2000,wavenum_max=3000,molecule="CO2",isotope="1,2,3",)sf.fetch_databank("hitemp",load_columns='noneq')# nonequilibrium between bending+symmetric and asymmetric modes :s1=sf.non_eq_spectrum(Tvib=(600,600,2000),Trot=600,path_length=1,pressure=1)# all vibrational temperatures are equal :s2=sf.non_eq_spectrum(Tvib=(2000,2000,2000),Trot=600,path_length=1,pressure=1)
Calculate total power emitted in equilibrium or non-equilibrium case
in the optically thin approximation: it sums all emission integral over
the total spectral range.
Warning
this is a fast implementation that doesnt take into account
the contribution of lines outside the given spectral range. It is valid for spectral ranges
surrounded by no lines, and spectral ranges much broadened than the typical
line broadening (~ 1-10 cm-1 in the infrared)
If what you’re looking for is an accurate simulation on a narrow spectral range
you better calculate the spectrum (that does take all of that into account)
and integrate it with get_power()
Parameters:
Tgas (float) – equilibrium temperature [K]
If doing a non equilibrium case it should be None. Use Ttrans for
translational temperature
Tvib (float) – vibrational temperature [K]
Trot (float) – rotational temperature [K]
Ttrans (float) – translational temperature [K]. If None, translational temperature is
taken as rotational temperature (valid at 1 atm for times above ~ 2ns
which is the RT characteristic time)
mole_fraction (float) – database species mole fraction. If None, Factory mole fraction is used.
path_length (float) – slab size (cm). If None, Factory mole fraction is used.
unit (str) – output unit. Default 'mW/cm2/sr'
Returns:
float – see unit=.
Return type:
Returns total power density in mW/cm2/sr (unless different unit is chosen),
predict_time(self) uses the input parameters like Spectral Range, Number of lines, wstep,
truncation to predict the estimated calculation time for the Spectrum
broadening step(bottleneck step) for the current optimization and broadening_method. The formula
for predicting time is based on benchmarks performed on various parameters for different optimization,
broadening_method and deriving its time complexity.