radis.lbl.factory moduleΒΆ
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() # autodownload from HITRAN
s = 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 sense
show_points=True) # show_points to have an idea of the resolution
See radis.test.lbl.test_factory
for more examples
Routine ListingΒΆ
(the following sections are marked on the slider bar in the Spyder IDE with the # XXX symbol)
PUBLIC METHODS
eq_spectrum()
>>> calc equilibrium spectrumnon_eq_spectrum()
>>> calc non equilibrium spectrumoptically_thin_power()
>>> get total power (equilibrium or non eq)
Most methods are written in inherited class with the following inheritance scheme:
DatabankLoader
> BaseFactory
>
BroadenFactory
> BandFactory
>
SpectrumFactory
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()
- class SpectrumFactory(wmin=None, wmax=None, wunit=cm - 1, wavenum_min=None, wavenum_max=None, wavelength_min=None, wavelength_max=None, Tref=296, pressure=1.01325, mole_fraction=1, path_length=1, wstep=0.01, molecule=None, isotope='all', medium='air', truncation=50, neighbour_lines=0, pseudo_continuum_threshold=0, self_absorption=True, chunksize=None, optimization='simple', folding_thresh=1e-06, zero_padding=-1, broadening_method='voigt', cutoff=0, parsum_mode='full summation', verbose=True, warnings=True, save_memory=False, export_populations=None, export_lines=False, emulate_gpu=False, **kwargs)[source]ΒΆ
Bases:
BandFactory
A class to put together all functions related to loading CDSD / HITRAN databases, calculating the broadenings, and summing over all the lines.
- Parameters
wmin, wmax (
float
orQuantity
) β 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)
orQuantity
) β minimum (maximum) wavenumber to be processed in \(cm^{-1}\). use astropy.units to specify arbitrary inverse-length units.wavelength_min, wavelength_max (
float(nm)
orQuantity
) β minimum (maximum) wavelength to be processed in \(nm\). This wavelength can be in'air'
or'vacuum'
depending on the value of the parametermedium=
. use astropy.units to specify arbitrary length units.pressure (
float(bar)
orQuantity
) β partial pressure of gas in bar. Default1.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. Default1
. Note that the rest of the gas is considered to be air for collisional broadening.path_length (
float(cm)
orQuantity
) β path length in cm. Default1
. use astropy.units to specify arbitrary length units.molecule (
int
,str
, orNone
) β molecule id (HITRAN format) or name. IfNone
, the molecule can be infered from the database files being loaded. See the list of supported molecules inMOLECULES_LIST_EQUILIBRIUM
andMOLECULES_LIST_NONEQUILIBRIUM
. DefaultNone
.isotope (
int
,list
,str
of the form'1,2'
, or'all'
) β isotope id (sorted by relative density: (eg: 1: CO2-626, 2: CO2-636 for CO2). See HITRAN documentation for isotope list for all species. 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'
- 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. DefaultTrue
.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. IfNone
, use no truncation (lineshapes spread on the full spectral range). Default is300
\(cm^{-1}\)Note
Large values (>
50
) can induce a performance drop (computation of lineshape typically scale as \(~truncation ^2\) ). The default300
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 is0
\(cm^{-1}\).βwstep (float (cm-1) or
'auto'
) β Resolution of wavenumber grid. Default0.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"
inradis.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. If0
, no cutoff. Default1e-27
.parsum_mode (βfull summationβ, βtabulationβ) β how to compute partition functions, at nonequilibrium or when partition function are not already tabulated.
'full summation'
: 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%). Defaultfull summation'
Note
parsum_mode= βtabulationβ is new in 0.9.30, and makes nonequilibrium calculations of small spectra extremelly 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. If0
, no semi-continuum is used. Default0
.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 & recalculate the linestrength for each new parameter. DefaultFalse
.export_populations (
'vib'
,'rovib'
,None
) β if not None, store populations in Spectrum. Either store vibrational populations (βvibβ) or rovibrational populations (βrovibβ). DefaultNone
export_lines (boolean) β if
True
, saves details of all calculated lines in Spectrum. This is necessary to later useline_survey()
, but can take some space. DefaultFalse
.chunksize (int, or
None
) β Splits the lines database in several chuncks 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 chunck size. IfNone
, all lines are processed directly. Usually faster but can create memory problems. DefaultNone
optimization (
"simple"
,"min-RMS"
,None
) β If either"simple"
or"min-RMS"
LDM optimization for lineshape calculation is used: -"min-RMS"
: weights optimized by analytical minimization of the RMS-error (See: [Spectral-Synthesis-Algorithm]) -"simple"
: weights equal to their relative position in the gridIf using the LDM optimization, broadening method is automatically set to
'fft'
. IfNone
, no lineshape interpolation is performed and the lineshape of all lines is calculated.Refer to [Spectral-Synthesis-Algorithm] for more explanation on the LDM method for lineshape interpolation.
Default
"min-RMS"
folding_thresh (float) β Folding is a correction procedure thet is applied when the lineshape is calculated with the
fft
broadening method and the linewidth is comparable towstep
, that prevents sinc(v) modulation of the lineshape. Folding continues until the lineshape intensity is belowfolding_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.Range: 0 <= zero_padding <= len(w), or zero_padding = -1 Default: -1
broadening_method (
"voigt"
,"convolve"
,"fft"
) β Calculates broadening with a direct voigt approximation (βvoigtβ) or by convoluting independantly calculated Doppler and collisional broadening (βconvolveβ). First is much faster, 2nd can be used to compare results. This SpectrumFactory parameter can be manually adjusted a posteriori with:sf = SpectrumFactory(...) sf.params.broadening_method = 'voigt'
Fast fourier transform
'fft'
is only available if using the LDM lineshape calculationoptimization
. Because the LDM convolves all lines at the same time, and thus operates on large arrays,'fft'
becomes more appropriate than convolutions in real space ('voit'
,'convolve'
)By default, use
"fft"
for anyoptimization
, and"voigt"
if optimization isNone
.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:warnings = {'MissingSelfBroadeningWarning':'ignore', 'NegativeEnergiesWarning':'ignore', 'HighTemperatureWarning':'ignore'}
See
default_warning_status
for more information.verbose (boolean, or int) β If
False
, stays quiet. IfTrue
, tells what is going on. If>=2
, gives more detailed messages (for instance, details of calculation times). DefaultTrue
.
Examples
An example using
SpectrumFactory
,load_databank()
, theSpectrum
methods, andunits
from radis import SpectrumFactory from astropy import units as u sf = SpectrumFactory(wavelength_min=4165 * u.nm, wavelength_max=4200 * u.nm, isotope='1,2', truncation=10, # cm-1 optimization=None, medium='vacuum', verbose=1, # more for more details ) sf.load_databank('HITRAN-CO2-TEST') # predefined in ~/radis.json s = sf.eq_spectrum(Tgas=300 * u.K, path_length=1 * u.cm) s.rescale_path_length(0.01) # cm s.plot('radiance_noslit', Iunit='Β΅W/cm2/sr/nm')
Refer to the online Examples for more cases.
GPU Accelerated SpectraUse Custom AbundancesReal-time GPU Accelerated Spectra (Interactive)
Real-time GPU Accelerated Spectra (Interactive)Spectrum DatabaseMulti-temperature Fit1 temperature fitNotes
High-level wrapper to SpectrumFactory:
Main Methods:
For advanced use:
Inputs and parameters can be accessed a posteriori with :
input
: physical inputparams
: computational parametersmisc
: miscallenous parameters (donβt change output)
See also
- eq_spectrum(Tgas, mole_fraction=None, path_length=None, pressure=None, name=None) Spectrum [source]ΒΆ
Generate a spectrum at equilibrium.
- Parameters
Tgas (float or
Quantity
) β Gas temperature (K)mole_fraction (float) β database species mole fraction. If None, Factory mole fraction is used.
path_length (float or
Quantity
) β slab size (cm). IfNone
, the default Factorypath_length
is used.pressure (float or
Quantity
) β pressure (bar). IfNone
, the default Factorypressure
is used.name (str) β output Spectrum name (useful in batch)
- Returns
s β
Returns a
Spectrum
object- Return type
Examples
- ::
from radis import SpectrumFactory sf = SpectrumFactory( wavenum_min=2900, wavenum_max=3200, molecule=βOHβ, wstep=0.1, ) sf.fetch_databank(βhitempβ)
s1 = sf.eq_spectrum(Tgas=300, path_length=1, pressure=0.1) s2 = sf.eq_spectrum(Tgas=500, path_length=1, pressure=0.1)
References
- 1
RADIS doc: Spectrum how to?
See also
- eq_spectrum_gpu(Tgas, mole_fraction=None, path_length=None, pressure=None, name=None, emulate=False) Spectrum [source]ΒΆ
Generate a spectrum at equilibrium with calculation of lineshapes and broadening done on the GPU.
Note
This method requires CUDA compatible hardware to execute. For more information on how to setup your system to run GPU-accelerated methods using CUDA and Cython, check GPU Spectrum Calculation on RADIS
- Parameters
Tgas (float or
Quantity
) β Gas temperature (K)mole_fraction (float) β database species mole fraction. If None, Factory mole fraction is used.
path_length (float or
Quantity
) β slab size (cm). IfNone
, the default Factorypath_length
is used.pressure (float or
Quantity
) β pressure (bar). IfNone
, the default Factorypressure
is used.name (str) β output Spectrum name (useful in batch)
- Other Parameters
emulate (bool) β if
True
, execute the GPU code on the CPU (useful for development)- Returns
s β
Returns a
Spectrum
object- Return type
Examples
GPU Accelerated SpectraReal-time GPU Accelerated Spectra (Interactive)
Real-time GPU Accelerated Spectra (Interactive)See also
- eq_spectrum_gpu_interactive(var='transmittance', slit_FWHM=0.0, plotkwargs={}, *vargs, **kwargs) Spectrum [source]ΒΆ
- Parameters
var (TYPE, optional) β DESCRIPTION. The default is βtransmittanceβ.
slit_FWHM (TYPE, optional) β DESCRIPTION. The default is 0.0.
*vargs (TYPE) β arguments forwarded to
eq_spectrum_gpu()
**kwargs (dict) β arguments forwarded to
eq_spectrum_gpu()
In particular, seeemulate=
.plotkwargs (dict) β arguments forwarded to
plot()
- Returns
s β DESCRIPTION.
- Return type
TYPE
Examples
from radis import SpectrumFactory from radis.tools.plot_tools import ParamRange sf = SpectrumFactory(2200, 2400, # cm-1 molecule='CO2', isotope='1,2,3', wstep=0.002, ) sf.fetch_databank('hitemp') s = sf.eq_spectrum_gpu_interactive(Tgas=ParamRange(300.0,2000.0,1200.0), #K pressure=0.2, #bar mole_fraction=0.1, path_length=ParamRange(0,10,2.0), #cm slit_FWHM=ParamRange(0,1,0), #cm emulate=False, # if True, execute the GPU code on the CPU )
- fit_spectrum(s_exp, model, fit_parameters, bounds={}, plot=False, solver_options={'maxiter': 300}, **kwargs) Union[Spectrum, OptimizeResult] [source]ΒΆ
Fit an experimental spectrum with an arbitrary model and an arbitrary number of fit parameters.
- Parameters
s_exp (Spectrum) β experimental spectrum. Should have only spectral array only. Use
take()
, e.g:sf.fit_spectrum(s_exp.take('transmittance'))
model (func -> Spectrum) β a line-of-sight model returning a Spectrum. Example :
Tvib12Tvib3Trot_NonLTEModel()
fit_parameters (dict) β
example:
{fit_parameter:initial_value}
bounds (dict, optional) β
example:
{fit_parameter:[min, max]}
fixed_parameters (dict) β
fixed parameters given to the model. Example:
fit_spectrum(fixed_parameters={"vib_distribution":"treanor"})
- Other Parameters
plot (bool) β if
True
, plot spectra as they are computed; and plot the convergence of the residual. DefaultFalse
solver_options (dict) β parameters forwarded to the solver. More info in
minimize
Example:{"maxiter": (int) max number of iteration default ``300``, }
kwargs (dict) β forwarded to
fit_spectrum()
- Returns
s_best (Spectrum) β best spectrum
res (OptimizeResults) β output of
minimize
Examples
See a one-temperature fit example and a non-LTE fit example
Multi-temperature Fit1 temperature fitMore advanced tools for interactive fitting of multi-dimensional, multi-slabs spectra can be found in
fitroom
.See also
- generate_perf_profile()[source]ΒΆ
Generate a visual/interactive performance profile diagram for the last calculated spectrum, with
tuna
.Requires a
profiler
key with in Spectrum.conditions, andtuna
installed.Note
print_perf_profile()
is an ascii-version which does not requiretuna
.Examples
sf = SpectrumFactory(...) sf.eq_spectrum(...) sf.generate_perf_profile()
See typical output in https://github.com/radis/radis/pull/325
Note
You can also profile with
tuna
directly:python -m cProfile -o program.prof your_radis_script.py tuna your_radis_script.py
See also
- non_eq_spectrum(Tvib, Trot, Ttrans=None, mole_fraction=None, path_length=None, pressure=None, vib_distribution='boltzmann', rot_distribution='boltzmann', overpopulation=None, name=None) Spectrum [source]ΒΆ
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)
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 or
Quantity
) β slab size (cm). IfNone
, the default Factorypath_length
is used.pressure (float or
Quantity
) β pressure (bar). IfNone
, the default Factorypressure
is used.
- Other Parameters
vib_distribution (
'boltzmann'
,'treanor'
) β vibrational distributionrot_distribution (
'boltzmann'
) β rotational distributionoverpopulation (dict, or
None
) βadd overpopulation factors for given levels:
{level:overpopulation_factor}
name (str) β output Spectrum name (useful in batch)
- Returns
s β
Returns a
Spectrum
object- Return type
Examples
from radis import SpectrumFactory sf = SpectrumFactory( wavenum_min=2000, wavenum_max=3000, molecule="CO", isotope="1,2,3", ) sf.fetch_databank("hitemp", load_columns='noneq') s1 = sf.non_eq_spectrum(Tvib=2000, Trot=600, path_length=1, pressure=0.1) s2 = sf.non_eq_spectrum(Tvib=2000, Trot=600, path_length=1, pressure=0.1)
Multi-vibrational temperature. Below we compare non-LTE spectra of CO2 where all vibrational temperatures are equal, or where the bending & symmetric modes are in equilibrium with rotation
from radis import SpectrumFactory sf = 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)
References
- 1
RADIS doc: Spectrum how to?
See also
- optically_thin_power(Tgas=None, Tvib=None, Trot=None, Ttrans=None, vib_distribution='boltzmann', rot_distribution='boltzmann', mole_fraction=None, path_length=None, unit='mW/cm2/sr')[source]ΒΆ
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 broaded 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),
See also
- predict_time()[source]ΒΆ
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.
Benchmarks: https://anandxkumar.github.io/Benchmark_Visualization_GSoC_2021/
Complexity vs Calculation Time Visualizations LBL>Voigt: LINK, DIT>Voigt: LINK, DIT>FFT: LINK
- Returns
float
- Return type
Predicted time in seconds.
- print_perf_profile(number_format='{:.3f}', precision=16)[source]ΒΆ
Prints Profiler output dictionary in a structured manner for the last calculated spectrum
Examples
sf.print_perf_profile() # output >> spectrum_calculation 0.189s ββββββββββββββββ check_line_databank 0.000s check_non_eq_param 0.042s βββ fetch_energy_5 0.015s β calc_weight_trans 0.008s reinitialize 0.002s copy_database 0.000s memory_usage_warning 0.002s reset_population 0.000s calc_noneq_population 0.041s βββ part_function 0.035s ββ population 0.006s scaled_non_eq_linestrength 0.005s map_part_func 0.001s corrected_population_se 0.003s calc_emission_integral 0.006s applied_linestrength_cutoff 0.002s calc_lineshift 0.001s calc_hwhm 0.007s generate_wavenumber_arrays 0.001s calc_line_broadening 0.074s ββββββ precompute_LDM_lineshapes 0.012s LDM_Initialized_vectors 0.000s LDM_closest_matching_line 0.001s LDM_Distribute_lines 0.001s LDM_convolve 0.060s βββββ others 0.001s calc_other_spectral_quan 0.003s generate_spectrum_obj 0.000s others -0.016s
- Other Parameters
precision (int, optional) β total number of blocks. Default 16.
See also