radis.spectrum package

Submodules

Module contents

Spectrum-class module for post-processing.

PerfectAbsorber(s: radis.spectrum.spectrum.Spectrum) radis.spectrum.spectrum.Spectrum[source]

Makes a new Spectrum with the same transmittance/absorbance as Spectrum s, but with radiance set to 0. Useful to get contribution of different slabs in line-of-sight calculations (see example).

Parameters

s (Spectrum) – Spectrum object

Returns

s_trSpectrum object, with only the transmittance, absorbance and/or abscoeff part of s, where radiance_noslit , emisscoeff and emissivity_noslit (if they exist) have been set to 0

Return type

Spectrum

Examples

Let’s say you have a total line of sight:

s_los = s1 > s2 > s3

If you now want to get the contribution of s2 to the line-of-sight emission, you can do:

(s2 > PerfectAbsorber(s3)).plot('radiance_noslit')

And the contribution of s1 would be:

(s1 > PerfectAbsorber(s2>s3)).plot('radiance_noslit')

See more examples in Line-of-Sight module

Radiance(s: radis.spectrum.spectrum.Spectrum) radis.spectrum.spectrum.Spectrum[source]

Returns a new Spectrum with only the radiance component of s

Parameters

s (Spectrum) – Spectrum object

Returns

s_trSpectrum object, with only radiance defined

Return type

Spectrum

Examples

This function is useful to use Spectrum algebra operations:

s = calc_spectrum(...)   # contains emission & absorption arrays
rad = Radiance(s)        # contains radiance array only
rad -= 0.1    # arithmetic operation is applied to Radiance only

Equivalent to:

rad = s.take('radiance')
Radiance_noslit(s: radis.spectrum.spectrum.Spectrum) radis.spectrum.spectrum.Spectrum[source]

Returns a new Spectrum with only the radiance_noslit component of s

Parameters

s (Spectrum) – Spectrum object

Returns

s_trSpectrum object, with only radiance_noslit defined

Return type

Spectrum

Examples

This function is useful to use Spectrum algebra operations:

s = calc_spectrum(...)   # contains emission & absorption arrays
rad = Radiance_noslit(s) # contains 'radiance_noslit' array only
rad -= 0.1    # arithmetic operation is applied to Radiance_noslit only

Equivalent to:

rad = s.take('radiance_noslit')
class Spectrum(quantities, units=None, conditions=None, cond_units=None, populations=None, lines=None, wunit=None, name=None, references={}, check_wavespace=True, **kwargs)[source]

Bases: object

This class holds results calculated with the SpectrumFactory calculation, with other radiative codes, or experimental data. It can be used to plot different quantities a posteriori, or manipulate output units (for instance convert a spectral radiance per wavelength units to a spectral radiance per wavenumber).

See more information on how to generate, edit or combine Spectrum objects on the Spectrum object guide.

Parameters
  • quantities (dict of tuples {'quantity':(w, a)} or dict {'wavelength/wavenumber': w, quantity': a}) – where quantities are spectral quantities (absorbance, radiance, etc.) and wavenum is in \(cm^{-1}\) or \(nm\) (see waveunit) Example:

    # w, k, I are numpy arrays for wavenumbers, absorption coefficient, and radiance.
    from radis import Spectrum
    s = Spectrum({"wavenumber":w, "abscoeff":k, "radiance_noslit":I}, wunit='cm-1',
                 units={"radiance_noslit":"mW/cm2/sr/nm", "abscoeff":"cm-1"})
    

    Or:

    s = Spectrum({"abscoeff":(w,k), "radiance_noslit":(w,I)},
                 wunit="cm-1"
                 units={"radiance_noslit":"mW/cm2/sr/nm", "abscoeff":"cm-1"})
    

    See also: from_array() and from_txt()

  • units (dict) – units for quantities

Other Parameters
  • conditions (dict) – physical conditions and calculation parameters

  • wunit ('nm', 'cm-1', 'nm_vac' or None) – wavelength in air ('nm'), wavenumber ('cm-1'), or wavelength in vacuum ('nm_vac'). If None, 'wavespace' must be defined in conditions. Quantities should be evenly distributed along this space for fast convolution with the slit function

  • cond_units (dict) – units for conditions

Other Parameters
  • name (str, or None) – Give a name to this Spectrum object (automatically used in plots; useful for multislab configurations). Default None

  • populations (dict) – a dictionary of all species, and levels. Should be compatible with other radiative codes such as Specair output. Suggested format: {molecules: {isotopes: {elec state: rovib levels}}} e.g:

    {'CO2':{1: 'X': df}}   # with df a Pandas Dataframe
    
  • lines (pandas Dataframe) – all lines in databank (necessary for using line_survey()). Warning if you want to play with the lines content: The signification of columns in lines may be specific to a database format. Plus, some additional columns may have been added by the calculation (e.g: Ei and S for emission integral and linestrength in SpectrumFactory). Refer to the code to know what they mean (and their units)

  • references (dict) – a dict of doi of references used to compute this object. Automatically returned with the full bibtex entry by cite() It can also be set a posteriori. Example

    s = Spectrum()
    s.references = {"10.1016/j.jqsrt.2010.05.001": "HITEMP-2010 database",
                  "10.1016/j.jqsrt.2018.09.027":["calculation", "post-processing"],  # RADIS main paper. Automatically added
                  "10.1016/j.jqsrt.2020.107476":"DIT algorithm"}
                  )
    s.cite()
    
    Returns :
    Used for DIT algorithm
    ----------------------
    @article{van_den_Bekerom_2021,
        doi = {10.1016/j.jqsrt.2020.107476},
        url = {https://doi.org/10.1016%2Fj.jqsrt.2020.107476},
        year = 2021,
        month = {mar},
        publisher = {Elsevier {BV}},
        volume = {261},
        pages = {107476},
        author = {D.C.M. van den Bekerom and E. Pannier},
        title = {A discrete integral transform for rapid spectral synthesis},
        journal = {Journal of Quantitative Spectroscopy and Radiative Transfer}
    }
    
    Used for HITEMP-2010 database
    -----------------------------
    @article{Rothman_2010,
        doi = {10.1016/j.jqsrt.2010.05.001},
        url = {https://doi.org/10.1016%2Fj.jqsrt.2010.05.001},
        year = 2010,
        month = {oct},
        publisher = {Elsevier {BV}},
        volume = {111},
        number = {15},
        pages = {2139--2150},
        author = {L.S. Rothman and I.E. Gordon and R.J. Barber and H. Dothe and R.R. Gamache and A. Goldman and V.I. Perevalov and S.A. Tashkun and J. Tennyson},
        title = {{HITEMP}, the high-temperature molecular spectroscopic database},
        journal = {Journal of Quantitative Spectroscopy and Radiative Transfer}
    }
    
    Used for calculation, post-processing
    -------------------------------------
    @article{Pannier_2019,
        doi = {10.1016/j.jqsrt.2018.09.027},
        url = {https://doi.org/10.1016%2Fj.jqsrt.2018.09.027},
        year = 2019,
        month = {jan},
        publisher = {Elsevier {BV}},
        volume = {222-223},
        pages = {12--25},
        author = {Erwan Pannier and Christophe O. Laux},
        title = {{RADIS}: A nonequilibrium line-by-line radiative code for {CO}2 and {HITRAN}-like database species},
        journal = {Journal of Quantitative Spectroscopy and Radiative Transfer}
    }
    

    </details>

  • warnings (boolean) – if True, test if inputs are valid, e.g, spectra are evenly distributed in wavelength, and raise a warning if not. Note that this take ~ 3.5 ms for a 20k points spectrum, when the rest of the creation process is only ~ 1.8ms (makes it 3 times longer, and can be a problem if hundreds of spectra are created in a row). Default True

Examples

Manipulate a Spectrum calculated by RADIS:

s = calc_spectrum(2125, 2300, Tgas=2000, databank='CDSD')
s.print_conditions()
s.plot('absorbance')
s.line_survey(overlay='absorbance')
s.plot('radiance_noslit', wunits='cm-1', Iunits='W/m2/sr/cm-1')
s.apply_slit(5)
s.plot('radiance')
w, t = s.get('transmittance_noslit')  # for use in multi-slabs configs

Any tuple of numpy arrays (w, I) can also be converted into a Spectrum object from the Spectrum class directly, or using the calculated_spectrum() function. All the following methods are equivalent:

from radis import Spectrum, calculated_spectrum
s1 = calculated_spectrum(w, I, wunit='nm', Iunit='mW/cm2/sr/nm')
s2 = Spectrum.from_array(w, I, 'radiance_noslit',
                       wunit='nm', unit='mW/cm2/sr/nm')
s3 = Spectrum({'radiance_noslit': (w, I)},
              units={'radiance_noslit':'mW/cm2/sr/nm'},
              wunit='nm')

See more examples in the [Spectrum] page.

Spectrum objects can be stored, retrieved, rescaled, resampled:

from radis import load_spec
s = load_spec('co_calculation.spec')
s.rescale_path_length(0.5)                  # calculate for new path_length
s.rescale_mole_fraction(0.02)   # calculate for new mole fraction
s.resample(w_new)               # resample on new wavespace
s.store('co_calculation2.spec')

Notes

Implementation:

quantities are stored in the self._q dictionary. They are better accessed with the get() method that deals with units and wavespace

Wavebase:

quantites are stored either in wavenum or wavelength base, but this doesnt matter as they are retrieved / plotted with the get() and plot() methods which have units as input arguments

conditions[source]

Stores computation / measurement conditions

Type

dict

c[source]

convenience wrapper to conditions:

s.c["calculation_time"] is s.conditions["calculation_time"]
>> True
Type

dict

populations[source]

Stores molecules, isotopes, electronic states and vibrational or rovibrational populations

Type

dict

References

Spectrum

See the Spectrum object page

apply_slit()[source]

Apply an instrumental slit function to all quantities in Spectrum. Slit function can be generated with usual shapes (see shape=) or imported from an experimental slit function (path to a text file or numpy array of shape n*2). Convoluted spectra are cut on the edge compared to non-convoluted spectra, to remove side effects. See mode= to change this behaviour.

Warning with units: read about 'unit' and 'return_unit' parameters.

Parameters
  • slit_function (float or str or array) –

    If float:

    generate slit function with FWHM of slit function (in nm or cm-1 depending on unit=). A (top, base) tuple of (float, float) is required when asking for a trapezoidal slit function.

    If .txt:

    import experimental slit function from .txt file: format must be 2-columns with wavelengths and intensity (doesn’t have to be normalized)

    If array:

    format must be 2-columns with wavelengths and intensity (doesn’t have to be normalized)

  • unit ('nm' or 'cm-1') – unit of slit_function (FWHM, or imported file)

  • shape ('triangular', 'trapezoidal', 'gaussian', or any of SLIT_SHAPES) –

    which shape to use when generating a slit. Will call,

    respectively, triangular_slit(), trapezoidal_slit(), gaussian_slit(). Default ‘triangular’

  • center_wavespace (float, or None) – center of slit when generated (in unit). Not used if slit is imported.

  • norm_by ('area', 'max') – normalisation type: - 'area' normalizes the slit function to an area

    of 1. It conserves energy, and keeps the same units.

    • 'max' normalizes the slit function to a maximum of 1. The convoluted spectrum units change (they are multiplied by the spectrum waveunit, e.g: a radiance non convoluted in mW/cm2/sr/nm on a wavelength (nm). range will yield a convoluted radiance in mW/cm2/sr. Note that the slit is set to 1 in the Spectrum wavespace (i.e: a Spectrum calculated in cm-1 will have a slit set to 1 in cm-1).

    Default 'area'

  • mode ('valid', 'same') –

    'same' returns output of same length as initial spectra,

    but boundary effects are still visible. 'valid' returns output of length len(spectra) - len(slit) + 1, for which lines outside of the calculated range have no impact. Default 'valid'.

Other Parameters
  • auto_recenter_crop (bool) – if True, recenter slit and crop zeros on the side when importing an experimental slit. Default True. See recenter_slit(), crop_slit()

  • plot_slit (boolean) – if True, plot slit

  • store (boolean) – if True, store slit in the Spectrum object so it can be retrieved with get_slit() and plot with plot_slit(). Default True

  • slit_dispersion (func of (lambda, in 'nm'), or None) – spectrometer reciprocal function : dλ/dx(λ) (in nm) If not None, then the slit_dispersion function is used to correct the slit function for the whole range. Can be important if slit function was measured far from the measured spectrum (e.g: a slit function measured at 632.8 nm will look broader at 350 nm because the spectrometer dispersion is higher at 350 nm. Therefore it should be corrected) Default None

    Warning

    slit dispersion function is assumed to be given in nm if your spectrum is stored in cm-1 the wavenumbers are converted to wavelengths before being forwarded to the dispersion function

    See test_auto_correct_dispersion() for an example of the slit dispersion effect.

    A Python implementation of the slit dispersion:

    >>> def f(lbd):
    >>>    return  w/(2*f)*(tan(Φ)+sqrt((2*d/m/(w*1e-9)*cos(Φ))^2-1))
    

    Theoretical / References:

    >>> dλ/dx ~ d/mf    # at first order
    >>> dλ/dx = w/(2*f)*(tan(Φ)+sqrt((2*d/m/(w)*cos(Φ))^2-1))  # cf
    

    with:

    • Φ: spectrometer angle (°)

    • f: focal length (mm)

    • m: order of dispersion

    • d: grooves spacing (mm) = 1/gr with gr in (gr/mm)

    See Laux 1999 “Experimental study and modeling of infrared air plasma radiation” for more information

  • slit_dispersion_warning_threshold (float) – if not None, check that slit dispersion is about constant (< threshold change) on the calculated range. Default 0.01 (1%). See offset_dilate_slit_function()

  • inplace (bool) – if True, adds convolved arrays directly in the Spectrum. If False, returns a new Spectrum with only the convolved arrays. Note: if you want a new Spectrum with both the convolved and non convolved quantities, use

    s.copy().apply_slit()
    
  • *args, **kwargs – are forwarded to slit generation or import function

  • verbose (bool) – print stuff

  • energy_threshold (float) – tolerance fraction when resampling. Default 1e-3 (0.1%) If areas before and after resampling differ by more than that an error is raised.

Returns

Spectrum – Allows chaining. If inplace=False, return a new Spectrum with the new spectral arrays only.

Return type

same Spectrum, with new spectral arrays.

Notes

Units:

the slit function is first converted to the wavespace (wavelength/wavenumber) that the Spectrum is stored in, and applied to the spectral quantities in their native wavespace.

Implementation:

convolve_with_slit() is applied to all quantities in get_vars() that ends with _noslit. Generate a triangular instrumental slit function (or any other shape depending of shape=) with base slit_function_base (Uses the central wavelength of the spectrum for the slit function generation)

We deal with several special cases (which makes the code a little heavy, but the method very versatile):

  • when slit unit and spectrum unit arent the same

  • when spectrum is not evenly spaced

Examples

s.apply_slit(1.2, 'nm')

To manually apply the slit to a particular quantity use:

wavenum, quantity = s['quantity']
s['convolved_quantity'] = convolve_slit(wavenum, quantity,
    slit_function_base)

See convolve_with_slit() for more details on Units and Normalization

The slit is made considering the “center wavelength” which is the mean wavelength of the full spectrum you are applying it to.

Examples using apply_slit :

c[source]
cite(format='bibentry')[source]

Prints bibliographic references used to compute this spectrum, as stored in the references dictionary.

Parameters

format (default 'bibentry'. See more in habanero.content_negotiation())

Examples

from radis import calc_spectrum
s = calc_spectrum(
    1900,
    2300,  # cm-1
    molecule="CO",
    isotope="1,2,3",
    pressure=1.01325,  # bar
    Tvib=2000,  #
    Trot=300,
    mole_fraction=0.1,
    path_length=1,  # cm
    databank="hitran",
)
s.cite()
Returns :
Used for algorithm
------------------
@article{van_den_Bekerom_2021,
    doi = {10.1016/j.jqsrt.2020.107476},
    url = {https://doi.org/10.1016%2Fj.jqsrt.2020.107476},
    year = 2021,
    month = {mar},
    publisher = {Elsevier {BV}},
    volume = {261},
    pages = {107476},
    author = {D.C.M. van den Bekerom and E. Pannier},
    title = {A discrete integral transform for rapid spectral synthesis},
    journal = {Journal of Quantitative Spectroscopy and Radiative Transfer}
}

Used for calculation, rovibrational energies
--------------------------------------------
@article{Pannier_2019,
    doi = {10.1016/j.jqsrt.2018.09.027},
    url = {https://doi.org/10.1016%2Fj.jqsrt.2018.09.027},
    year = 2019,
    month = {jan},
    publisher = {Elsevier {BV}},
    volume = {222-223},
    pages = {12--25},
    author = {Erwan Pannier and Christophe O. Laux},
    title = {{RADIS}: A nonequilibrium line-by-line radiative code for {CO}2 and {HITRAN}-like database species},
    journal = {Journal of Quantitative Spectroscopy and Radiative Transfer}
}

Used for data retrieval
-----------------------
@article{Ginsburg_2019,
    doi = {10.3847/1538-3881/aafc33},
    url = {https://doi.org/10.3847%2F1538-3881%2Faafc33},
    year = 2019,
    month = {feb},
    publisher = {American Astronomical Society},
    volume = {157},
    number = {3},
    pages = {98},
    author = {Adam Ginsburg and Brigitta M. Sip{\H{o}}cz and C. E. Brasseur and Philip S. Cowperthwaite and Matthew W. Craig and Christoph Deil and James Guillochon and Giannina Guzman and Simon Liedtke and Pey Lian Lim and Kelly E. Lockhart and Michael Mommert and Brett M. Morris and Henrik Norman and Madhura Parikh and Magnus V. Persson and Thomas P. Robitaille and Juan-Carlos Segovia and Leo P. Singer and Erik J. Tollerud and Miguel de Val-Borro and Ivan Valtchanov and Julien Woillez and},
    title = {astroquery: An Astronomical Web-querying Package in Python},
    journal = {The Astronomical Journal}
}

Used for line database
----------------------
@article{Gordon_2017,
    doi = {10.1016/j.jqsrt.2017.06.038},
    url = {https://doi.org/10.1016%2Fj.jqsrt.2017.06.038},
    year = 2017,
    month = {dec},
    publisher = {Elsevier {BV}},
    volume = {203},
    pages = {3--69},
    author = {I.E. Gordon and L.S. Rothman and C. Hill and R.V. Kochanov and Y. Tan and P.F. Bernath and M. Birk and V. Boudon and A. Campargue and K.V. Chance and B.J. Drouin and J.-M. Flaud and R.R. Gamache and J.T. Hodges and D. Jacquemart and V.I. Perevalov and A. Perrin and K.P. Shine and M.-A.H. Smith and J. Tennyson and G.C. Toon and H. Tran and V.G. Tyuterev and A. Barbe and A.G. Cs{'{a}}sz{'{a}}r and V.M. Devi and T. Furtenbacher and J.J. Harrison and J.-M. Hartmann and A. Jolly and T.J. Johnson and T. Karman and I. Kleiner and A.A. Kyuberis and J. Loos and O.M. Lyulin and S.T. Massie and S.N. Mikhailenko and N. Moazzen-Ahmadi and H.S.P. Müller and O.V. Naumenko and A.V. Nikitin and O.L. Polyansky and M. Rey and M. Rotger and S.W. Sharpe and K. Sung and E. Starikova and S.A. Tashkun and J. Vander Auwera and G. Wagner and J. Wilzewski and P. Wcis{\l}o and S. Yu and E.J. Zak},
    title = {The {HITRAN}2016 molecular spectroscopic database},
    journal = {Journal of Quantitative Spectroscopy and Radiative Transfer}
}

Used for partition function
---------------------------
@article{Gamache_2021,
    doi = {10.1016/j.jqsrt.2021.107713},
    url = {https://doi.org/10.1016%2Fj.jqsrt.2021.107713},
    year = 2021,
    month = {sep},
    publisher = {Elsevier {BV}},
    volume = {271},
    pages = {107713},
    author = {Robert R. Gamache and Bastien Vispoel and Michaël Rey and Andrei Nikitin and Vladimir Tyuterev and Oleg Egorov and Iouli E. Gordon and Vincent Boudon},
    title = {Total internal partition sums for the {HITRAN}2020 database},
    journal = {Journal of Quantitative Spectroscopy and Radiative Transfer}
}
@article{Kochanov_2016,
    doi = {10.1016/j.jqsrt.2016.03.005},
    url = {https://doi.org/10.1016%2Fj.jqsrt.2016.03.005},
    year = 2016,
    month = {jul},
    publisher = {Elsevier {BV}},
    volume = {177},
    pages = {15--30},
    author = {R.V. Kochanov and I.E. Gordon and L.S. Rothman and P. Wcis{\l}o and C. Hill and J.S. Wilzewski},
    title = {{HITRAN} Application Programming Interface ({HAPI}): A comprehensive approach to working with spectroscopic data},
    journal = {Journal of Quantitative Spectroscopy and Radiative Transfer}
}

Used for spectroscopic constants
--------------------------------
@article{Guelachvili_1983,
    doi = {10.1016/0022-2852(83)90203-5},
    url = {https://doi.org/10.1016%2F0022-2852%2883%2990203-5},
    year = 1983,
    month = {mar},
    publisher = {Elsevier {BV}},
    volume = {98},
    number = {1},
    pages = {64--79},
    author = {G. Guelachvili and D. de Villeneuve and R. Farrenq and W. Urban and J. Verges},
    title = {Dunham coefficients for seven isotopic species of {CO}},
    journal = {Journal of Molecular Spectroscopy}
}

</details>

See also

RefTracker

compare_with(other, spectra_only=False, plot=True, wunit='default', verbose=True, rtol=1e-05, ignore_nan=False, ignore_outliers=False, normalize=False, **kwargs)[source]

Compare Spectrum with another Spectrum object.

Parameters
  • other (type Spectrum) – another Spectrum to compare with

  • spectra_only (boolean, or str) – if True, only compares spectral quantities (in the same waveunit) and not lines or conditions. If str, compare a particular quantity name. If False, compare everything (including lines and conditions and populations). Default False

  • plot (boolean) – if True, use plot_diff to plot all quantities for the 2 spectra and the difference between them. Default True.

  • wunit ("nm", "cm-1", "default") – in which wavespace to compare (and plot). If "default", natural wavespace of first Spectrum is taken.

  • rtol (float) – relative difference to use for spectral quantities comparison

  • ignore_nan (boolean) – if True, nans are ignored when comparing spectral quantities

  • ignore_outliers (boolean, or float) –

    if not False, outliers are discarded. i.e, output is determined by:

    out = (~np.isclose(I, Ie, rtol=rtol, atol=0)).sum()/len(I) < ignore_outliers
    
  • normalize (bool) – Normalize the spectra to be ploted

Other Parameters

kwargs (dict) – arguments are forwarded to plot_diff()

Returns

equals – return True if spectra are equal (respective to tolerance defined by rtol and other input conditions)

Return type

boolean

Examples

Compare two Spectrum objects, or specifically the transmittance:

s1.compare_with(s2)
s1.compare_with(s2, 'transmittance')

Note that you can also simply use s1 == s2, that uses compare_with() internally:

s1 == s2       # will return True or False
cond_units[source]
conditions[source]

computation conditions, or experimetnal parameters, or any metadata you need to store with the Spectrum object.

Type

dict

copy(copy_lines=True, quantity='all')[source]

Returns a copy of this Spectrum object (performs a smart deepcopy)

Parameters
  • copy_lines (bool) – default True

  • quantity (‘all’, or one of ‘radiance_noslit’, ‘absorbance’, etc.) – if not ‘all’, copy only one quantity. Default 'all'

  • .. minigallery:: radis.spectrum.spectrum.Spectrum.copy – :add-heading:

crop(wmin=None, wmax=None, wunit='default', inplace=True)[source]

Crop spectrum to wmin-wmax range in wunit (inplace)

Parameters
  • wmin, wmax (float, or None) – boundaries of spectral range (in wunit)

  • wunit ('nm', 'cm-1', 'nm_vac') – which waveunit to use for wmin, wmax. If default: use the default Spectrum wavespace defined with get_waveunit().

Other Parameters

inplace (bool) – if True, modifies the Spectrum object directly. Else, returns a copy. Default True.

Returns

s – Cropped Spectrum. If inplace=True, Spectrum has been updated directly anyway. Allows chaining

Return type

Spectrum

Examples

Crop to experimental Spectrum, and compare:

from radis import calc_spectrum, load_spec, plot_diff
s = calc_spectrum(...)
s_exp = load_spec('typical_result.spec')
s.crop(s_exp.get_wavelength.min(), s_exp.get_wavelength.max(), 'nm')
plot_diff(s_exp, s)

See also

radis.spectrum.operations.crop()

MergeSlabs()

if used with ``resample=’full’,

out, this, to

file[source]
classmethod from_array(w, I, quantity, wunit, unit, waveunit=None, *args, **kwargs)[source]

Construct Spectrum from 2 arrays.

Parameters
  • w, I (array) – waverange and vector

  • quantity (str) – spectral quantity name

  • wunit ('nm', 'cm-1', 'nm_vac') – unit of waverange: wavelength in air ('nm'), wavenumber ('cm-1'), or wavelength in vacuum ('nm_vac').

  • unit (str) – spectral quantity unit (arbitrary). Ex: ‘mW/cm2/sr/nm’ for radiance_noslit

  • *args, **kwargs – see Spectrum doc

Other Parameters
  • conditions (dict) – physical conditions and calculation parameters

  • cond_units (dict) – units for conditions

  • populations (dict) – a dictionary of all species, and levels. Should be compatible with other radiative codes such as Specair output. Suggested format: {molecules: {isotopes: {elec state: rovib levels}}} e.g:

    {'CO2':{1: 'X': df}}   # with df a Pandas Dataframe
    
  • lines (pandas Dataframe) – all lines in databank (necessary for using line_survey()). Warning if you want to play with the lines content: The signification of columns in lines may be specific to a database format. Plus, some additional columns may have been added by the calculation (e.g: Ei and S for emission integral and linestrength in SpectrumFactory). Refer to the code to know what they mean (and their units)

Returns

s – creates a Spectrum object

Return type

Spectrum

Examples

Create a spectrum:

from radis import Spectrum
s = Spectrum.from_array(w, I, 'radiance_noslit',
                       wunit='nm', unit='mW/cm2/sr/nm')

To create a spectrum with absorption and emission components (e.g: radiance_noslit and transmittance_noslit, or emisscoeff and abscoeff) call the Spectrum class directly. Ex:

from radis import Spectrum
s = Spectrum({'abscoeff': (w, A), 'emisscoeff': (w, E)},
             units={'abscoeff': 'cm-1', 'emisscoeff':'W/cm2/sr/nm'},
             wunit='nm')
classmethod from_txt(file, quantity, wunit, unit, waveunit=None, *args, **kwargs)[source]

Construct Spectrum from txt file.

Parameters
  • file (str) – file name

  • quantity (str) – spectral quantity name

  • wunit ('nm', 'cm-1', 'nm_vac') – unit of waverange: wavelength in air ('nm'), wavenumber ('cm-1'), or wavelength in vacuum ('nm_vac').

  • unit (str) – spectral quantity unit

  • *args, **kwargs – the following inputs are forwarded to loadtxt: 'delimiter', 'skiprows' The rest if forwarded to Spectrum. see Spectrum doc

Other Parameters
  • delimiter (',', etc.) – see numpy.loadtxt()

  • skiprows (int) – see numpy.loadtxt()

  • argsort (bool) – sorts the arrays in file by wavespace. Convenient way to load a file where points have been manually added at the end. Default False.

  • *Optional Spectrum parameters*

  • conditions (dict) – physical conditions and calculation parameters

  • cond_units (dict) – units for conditions

  • populations (dict) – a dictionary of all species, and levels. Should be compatible with other radiative codes such as Specair output. Suggested format: {molecules: {isotopes: {elec state: rovib levels}}} e.g:

    {'CO2':{1: 'X': df}}   # with df a Pandas Dataframe
    
  • lines (pandas Dataframe) – all lines in databank (necessary for using line_survey()). Warning if you want to play with the lines content: The signification of columns in lines may be specific to a database format. Plus, some additional columns may have been added by the calculation (e.g: Ei and S for emission integral and linestrength in SpectrumFactory). Refer to the code to know what they mean (and their units)

Returns

s – creates a Spectrum object

Return type

Spectrum

Examples

Generate an experimental spectrum from txt. In that example the delimiter key is forwarded to loadtxt():

from radis import Spectrum
s = Spectrum.from_txt('spectrum.csv', 'radiance', wunit='nm',
                          unit='W/cm2/sr/nm', delimiter=',')

To create a spectrum with absorption and emission components (e.g: radiance_noslit and transmittance_noslit, or emisscoeff and abscoeff) call the Spectrum class directly. Ex:

from radis import Spectrum
s = Spectrum({'abscoeff': (w, A), 'emisscoeff': (w, E)},
             units={'abscoeff': 'cm-1', 'emisscoeff':'W/cm2/sr/nm'},
             wunit='nm')

Notes

Internally, the numpy loadtxt() function is used and transposed:

w, I = np.loadtxt(file).T

You can use 'delimiter' and ‘skiprows' as arguments.

generate_perf_profile()[source]

Generate a visual/interactive performance profile diagram using tuna

Examples

s = calc_spectrum(...)
s.generate_perf_profile()

See typical output in https://github.com/radis/radis/pull/325

https://user-images.githubusercontent.com/16088743/128018032-6049be72-1881-46ac-9d7c-1ed89f9c4f42.png

Note

You can also profile with tuna directly:

python -m cProfile -o program.prof your_radis_script.py
tuna your_radis_script.py
get()[source]

Retrieve a spectral quantity from a Spectrum object. You can select wavespace unit, intensity unit, or propagation medium.

Parameters
  • var (variable (‘absorbance’, ‘transmittance’, etc.)) – Should be a defined quantity among CONVOLUTED_QUANTITIES or NON_CONVOLUTED_QUANTITIES. To get the full list of quantities defined in this Spectrum object use the get_vars() method.

  • wunit ('nm', 'cm', 'nm_vac'.) – wavespace unit: wavelength in air ('nm'), wavenumber ('cm-1'), or wavelength in vacuum ('nm_vac'). if "default", default unit for waveunit is used. See get_waveunit().

  • Iunit (unit for variable var) – if "default", default unit for quantity var is used. See the units attribute. For var="radiance", one can use per wavelength (~ ‘W/m2/sr/nm’) or per wavenumber (~ ‘W/m2/sr/cm-1’) units

Other Parameters
  • copy (bool) – if True, returns a copy of the stored quantity (modifying it wont change the Spectrum object). Default True.

  • trim_nan (bool) – if True, removes nan on the sides of the spectral array (and corresponding wavespace). Default False.

Returns

w, I – wavespace, quantity (ex: wavelength, radiance_noslit). For numpy users, note that these are copies (values) of the Spectrum quantity and not a view (reference): if you modify them the Spectrum is not changed

Return type

array-like

Examples

Get transmittance in cm-1:

w, I = s.get('transmittance_noslit', wunit='cm-1')

Get radiance (in wavelength in air):

_, R = s.get('radiance_noslit', wunit='nm', Iunit='W/cm2/sr/nm')
get_conditions()[source]

Get all physical / computational parameters.

get_integral(var, wunit='default', Iunit='default', **kwargs)[source]

Returns integral of variable ‘var’ over waverange.

Parameters
  • var (str) – spectral quantity to integate

  • wunit (str) – over which waverange to integrated. If default, use the default Spectrum wavespace defined with get_waveunit().

  • Iunit (str) – default 'default'

    Warning

    this is the unit of the quantity, not the unit of the integral. Don’t forget to multiply by wunit

Other Parameters

kwargs (**dict) – forwarded to get()

Returns

integral – integral in [Iunit]*[wunit]

Return type

float

get_name()[source]

Return Spectrum name.

If not defined, returns either the file name if Spectrum was loaded from a file, or the 'spectrum{id}' with the Python id object

get_populations(molecule=None, isotope=None, electronic_state=None)[source]

Return populations that are featured in the spectrum, either as upper or lower levels.

Parameters
  • molecule (str, or None) – if None, only one molecule must be defined. Else, an error is raised

  • isotope (int, or None) – isotope number. if None, only one isotope must be defined. Else, an error is raised

  • electronic_state (str) – if None, only one electronic state must be defined. Else, an error is raised

Returns

  • pandas dataframe of levels, where levels are the index,

  • and ‘Evib’ and ‘nvib’ are featured

Notes

Structure:

{molecule: {isotope: {electronic_state: {'vib': pandas Dataframe,    # (copy of) vib levels
                                         'rovib': pandas Dataframe,  # (copy of) rovib levels
                                         'Ia': float    # isotopic abundance
                                         }}}}

(If Spectrum generated with RADIS, structure should match that of SpectrumFactory.get_populations())

get_power(unit='mW/cm2/sr')[source]

Returns integrated radiance (no slit) power density.

Parameters

Iunit (str) – power unit.

Returns

P – radiated power in unit

Return type

float

Examples

::

s.get_power(‘W/cm2/sr’)

get_quantities(which=None)[source]

Returns all spectral quantities stored in this object (convoluted or non convoluted). Wrapper to get_vars()

get_radiance(Iunit='mW/cm2/sr/nm', copy=True)[source]

Return radiance in whatever unit, and can even convert from ~1/nm to ~1/cm-1 (and the other way round)

Other Parameters

copy (boolean) – if True, returns a copy of the stored waverange (modifying it wont change the Spectrum object). Default True.

get_radiance_noslit(Iunit='mW/cm2/sr/nm', copy=True)[source]

Return radiance (non convoluted) in whatever unit, and can even convert from ~1/nm to ~1/cm-1 (and the other way round)

Other Parameters

copy (boolean) – if True, returns a copy of the stored waverange (modifying it wont change the Spectrum object). Default True.

get_rovib_levels(molecule=None, isotope=None, electronic_state=None, first=None)[source]

Return rovibrational levels calculated in the spectrum (energies, populations)

Parameters
  • molecule (str, or None) – if None, only one molecule must be defined. Else, an error is raised

  • isotope (int, or None) – isotope number. if None, only one isotope must be defined. Else, an error is raised

  • electronic_state (str) – if None, only one electronic state must be defined. Else, an error is raised

  • first (int, or ‘all’ or None) – only show the first N levels. If None or ‘all’, all levels are shown

Returns

out – pandas dataframe of levels, where levels are the index, and ‘Evib’ and ‘nvib’ are featured

Return type

pandas DataFrame

get_slit(wunit='same')[source]

Get slit function that was applied to the Spectrum.

Returns

wslit, Islit – slit function with wslit in Spectrum waveunit. See get_waveunit()

Return type

array

get_vars(which=None)[source]

Returns all spectral quantities stored in this object (convoluted or non convoluted)

get_vib_levels(molecule=None, isotope=None, electronic_state=None, first=None)[source]

Return vibrational levels in the spectrum (energies, populations)

Parameters
  • molecule (str, or None) – if None, only one molecule must be defined. Else, an error is raised

  • isotope (int, or None) – isotope number. if None, only one isotope must be defined. Else, an error is raised

  • electronic_state (str) – if None, only one electronic state must be defined. Else, an error is raised

  • first (int, or ‘all’ or None) – only show the first N levels. If None or ‘all’, all levels are shown

Returns

out – pandas dataframe of levels, where levels are the index, and ‘Evib’ and ‘nvib’ are featured

Return type

pandas DataFrame

get_wavelength(medium='air', which=None, copy=True)[source]

Return wavelength in defined medium.

Parameters

medium ('air', 'vacuum') – returns wavelength as seen in air, or vacuum. Default 'air'. See vacuum2air(), air2vacuum()

Other Parameters

copy (boolean) – if True, returns a copy of the stored waverange (modifying it wont change the Spectrum object). Default True.

Returns

w – (a copy of) spectrum wavelength for convoluted or non convoluted quantities

Return type

array_like

get_wavenumber(which=None, copy=True)[source]

Return wavenumber (if the same for all quantities)

Other Parameters

copy (boolean) – if True, returns a copy of the stored waverange (modifying it wont change the Spectrum object). Default True.

Returns

w – (a copy of) spectrum wavenumber for convoluted or non convoluted quantities

Return type

array_like

get_waveunit()[source]

Returns whether this spectrum is defined in wavelength (nm) or wavenumber (cm-1)

has_nan(ignore_wavespace=True) bool[source]
Parameters

s (Spectrum) – radis Spectrum.

Returns

b – returns whether Spectrum has nan

Return type

bool

Note

print(s) will also show which spectral quantities have ````nan.

is_at_equilibrium(check='warn', verbose=False)[source]

Returns whether this spectrum is at (thermal) equilibrium. Reads the thermal_equilibrium key in Spectrum conditions. It does not imply chemical equilibrium (mole fractions are still arbitrary)

If they are defined, also check that the following assertions are True:

Tvib = Trot = Tgas self_absorption = True overpopulation = None

If they are not, still trust the value in Spectrum conditions, but raise a warning.

Other Parameters
  • check ('warn', 'error', 'ignore') – what to do if Spectrum conditions dont match the given equilibrium state: raise a warning, raise an error, or just ignore and dont even check. Default 'warn'.

  • verbose (bool) – if True, print why is the spectrum is not at equilibrium, if applicable.

is_optically_thin()[source]

Returns whether the spectrum is optically thin, based on the value on the self_absorption key in conditions.

If not given, raises an error

line_survey(overlay=None, wunit='default', writefile=None, cutoff=None, *args, **kwargs)[source]

Plot Line Survey (all linestrengths used for calculation) Output in Plotly (html)

Parameters
  • spec (Spectrum) – result from SpectrumFactory calculation (see spectrum.py)

  • overlay (‘absorbance’, ‘transmittance’, ‘radiance’, etc… or list of the above, or None) – overlay Linestrength with specified variable calculated in spec. Get the full list with the get_vars() method. Default None.

  • wunit ('default', 'nm', 'cm-1', 'nm_vac',) – wavelength air, wavenumber, or wavelength vacuum. If 'default', Spectrum get_waveunit() is used.

  • medium ({‘air’, ‘vacuum’, ‘default’}) – Choose whether wavelength are shown in air or vacuum. If 'default' lines are shown as stored in the spectrum.

Other Parameters
  • writefile (str) – if not None, a valid filename to save the plot under .html format. If None, use the fig object returned to show the plot.

  • kwargs:: dict – Other inputs are passed to LineSurvey(). Example below (see LineSurvey() documentation for more details):

  • Iunit (hitran, splot) – Linestrength output units:

    • hitran: (cm-1/(molecule/cm-2))

    • splot: (cm-1/atm) (Spectraplot units 2)

    Note: if not None, cutoff criteria is applied in this unit. Not used if plot is not ‘S’

  • barwidth (float) – With of bars in LineSurvey. Default 0.07

Returns

  • fig (a Plotly figure.) – If using a Jupyter Notebook, the plot will appear. Else, use writefile to export to an html file.

  • Plot in Plotly. See Output in [1]_

Examples

An example using the SpectrumFactory to generate a spectrum:

from radis import SpectrumFactory
sf = SpectrumFactory(
                     wavenum_min=2380,
                     wavenum_max=2400,
                     mole_fraction=400e-6,
                     path_length=100,  # cm
                     isotope=[1],
                     export_lines=True,    # required for LineSurvey!
                     db_use_cached=True)
sf.load_databank('HITRAN-CO2-TEST')
s = sf.eq_spectrum(Tgas=1500)
s.apply_slit(0.5)
s.line_survey(overlay='radiance_noslit', barwidth=0.01)

See the output in Examples

References

1

RADIS Online Documentation (LineSurvey)

2

SpectraPlot

lines[source]
max()[source]

Maximum of the Spectrum, if only one spectral quantity is available:

s.max()

Else, use Radiance(), Radiance_noslit(), Transmittance() or Transmittance_noslit()

Radiance(s).max()
min()[source]

Minimum of the Spectrum, if only one spectral quantity is available

s.min()

Else, use Radiance(), Radiance_noslit(), Transmittance() or Transmittance_noslit()

Radiance(s).min()
name[source]
normalize(normalize_how='max', wrange=(), wunit=None, inplace=False, force=False, verbose=False)[source]

Normalise the Spectrum, if only one spectral quantity is available.

Parameters
  • normalize_how ('max', 'area', 'mean') – how to normalize. 'max' is the default but may not be suited for very noisy experimental spectra. 'area' will normalize the integral to 1. 'mean' will normalize by the mean amplitude value

  • wrange (tuple) – if not empty, normalize on this range

  • wunit ("nm", "cm-1", "nm_vac") – unit of the normalisation range above. If None, use the spectrum default waveunit.

  • inplace (bool) – if True, changes the Spectrum.

Other Parameters

force (boolean) – By default, normalizing some parametres such as transmittance is forbidden because considered non-physical. Use force=True if you really want to.

Examples

::

s.normalize(“max”, (4200, 4800), inplace=True).plot()

offset(offset: radis.spectrum.spectrum.Spectrum, unit: float, inplace: str = True) radis.spectrum.spectrum.Spectrum[source]

Offset the spectrum by a wavelength or wavenumber (inplace)

Parameters
  • offset (float) – Constant to add to all quantities in the Spectrum.

  • unit (‘nm’ or ‘cm-1’) – unit for offset

Other Parameters

inplace (bool) – if True, modifies the Spectrum object directly. Else, returns a copy. Default True.

Returns

s – Offset Spectrum. If inplace=True, Spectrum has been updated directly anyway. Allows chaining.

Return type

Spectrum

Examples

::

s.offset(5, ‘nm’)

plot()[source]

Plot a Spectrum object.

Parameters
  • var (variable (absorbance, transmittance, transmittance_noslit, etc.)) – For full list see get_vars(). If None, plot the first thing in the Spectrum. Default None.

  • wunit ('default', 'nm', 'cm-1', 'nm_vac',) – wavelength air, wavenumber, or wavelength vacuum. If 'default', Spectrum get_waveunit() is used.

  • Iunit (unit for variable) – if default, default unit for quantity var is used. for radiance, one can use per wavelength (~ W/m2/sr/nm) or per wavenumber (~ W/m2/sr/cm-1) units

Other Parameters
  • show_points (boolean) – show calculated points. Default True.

  • nfig (int, None, or ‘same’) – plot on a particular figure. ‘same’ plots on current figure. For instance:

    s1.plot()
    s2.plot(nfig='same')
    
  • show_medium (bool, 'vacuum_only') – if True and wunit are wavelengths, plot the propagation medium in the xaxis label ([air] or [vacuum]). If 'vacuum_only', plot only if wunit=='nm_vac'. Default 'vacuum_only' (prevents from inadvertently plotting spectra with different propagation medium on the same graph).

  • yscale (‘linear’, ‘log’) – plot yscale

  • normalize (boolean, or tuple.) – option to normalize quantity to 1 (ex: for radiance). Default False

  • plot_by_parts (bool) – if True, look for discontinuities in the wavespace and plot the different parts without connecting lines. Useful for experimental spectra produced by overlapping step-and-glue. Additional parameters read from kwargs : split_threshold and cutwings. See more in split_and_plot_by_parts().

  • force (bool) – plotting on an existing figure is forbidden if labels are not the same. Use force=True to ignore that.

  • show (bool) – show figure. Default False. Will still show the figure in interactive mode, e.g, %matplotlib inline in a Notebook.

  • show_ruler (bool) – if True, add a ruler tool to the Matplotlib toolbar.

    Warning

    still experimental in 0.9.30 ! Try it, feedback welcome !

  • **kwargs (**dict) – kwargs forwarded as argument to plot (e.g: lineshape attributes: lw=3, color='r')

Returns

line plot

Return type

line

Examples

Plot an experimental_spectrum() in arbitrary units:

s = experimental_spectrum(..., Iunit='mW/cm2/sr/nm')
s.plot(Iunit='W/cm2/sr/cm-1')

See more examples in the plot Spectral quantities page.

plot_populations(what=None, nunit='', correct_for_abundance=False, **kwargs)[source]

Plots vib populations if given and format is valid.

Parameters
  • what (‘vib’, ‘rovib’, None) – if None plot everything

  • nunit (‘’, ‘cm-3’) – plot either in a fraction of vibrational levels, or a molecule number in in cm-3

  • correct_for_abundance (boolean) – if True, multiplies each population by the isotopic abundance (as it is done during the calculation of emission integral)

  • kwargs (**dict) – are forwarded to the plot

plot_slit(wunit=None, waveunit=None)[source]

Plot slit function that was applied to the Spectrum.

If dispersion was used (see apply_slit()) the different slits are built again and plotted too (dotted).

Parameters

wunit ('nm', 'cm-1', or None) – plot slit in wavelength or wavenumber. If None, use the unit the slit in which the slit function was given. Default None

Returns

fix, ax – figure and ax

Return type

matplotlib objects

populations[source]
print_conditions(**kwargs)[source]

Prints all physical / computational parameters.

Parameters

kwargs (dict) – refer to print_conditions()

Examples

::

s.print_conditions()

You can also simply print the Spectrum object directly:

print(s)
print_perf_profile(number_format='{:.3f}', precision=16)[source]

Prints Profiler output dictionary in a structured manner.

Example

Spectrum.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_DLM_lineshapes        0.012s
            DLM_Initialized_vectors          0.000s
            DLM_closest_matching_line        0.001s
            DLM_Distribute_lines             0.001s
            DLM_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.

references[source]
resample(w_new, unit='same', out_of_bounds='nan', energy_threshold=0.005, print_conservation=False, inplace=True, if_conflict_drop=None, **kwargs)[source]

Resample spectrum over a new wavelength/wavenumber range.

Warning

This may result in information loss. Resampling is done with oversampling and spline interpolation. These parameters can be adjusted, and energy conservation ensured with the appropriate parameters.

To minimize information loss, always resample the high-resolution spectrum over the low-resolution spectrum, i.e.

s_highres.resample(s_lowres)

Fills with 'nan' or transparent medium (transmittance 1, radiance 0) when out of bound (see out_of_bounds)

Parameters
  • w_new (array, or Spectrum) – new wavespace to resample the spectrum on. Must be inclosed in the current wavespace (we won’t extrapolate) One can also give a Spectrum directly:

    s1.resample(s2.get_wavenumber())
    s1.resample(s2)            # also valid
    
  • unit ('same', 'nm', 'cm-1', 'nm_vac') – unit of new wavespace. It 'same' it is assumed to be the current waveunit. Default 'same'. The spectrum waveunit is changed to this unit after resampling (i.e: a spectrum calculated and stored in cm-1 but resampled in nm will be stored in nm from now on). If 'nm', wavelength in air. If 'nm_vac', wavelength in vacuum.

  • out_of_bounds ('transparent', 'nan', 'error') – what to do if resampling is out of bounds. 'transparent': fills with transparent medium. ‘nan’: fill with nan. 'error': raises an error. Default 'nan'

  • medium ('air', 'vacuum', or 'default') – in which medium is the new waverange is calculated if it is given in ‘nm’. Ignored if unit=’cm-1’

Other Parameters
  • energy_threshold (float or None) – if energy conservation (integrals on the intersecting range) is above this threshold, raise an error. If None, dont check for energy conservation Default 5e-3 (0.5%)

  • print_conservation (boolean) – if True, prints energy conservation. Default False.

  • inplace (boolean) – if True, modifies the Spectrum object directly. Else, returns a copy. Default True.

  • **kwargs (**dict) – all other arguments are sent to radis.misc.signal.resample()

Returns

Spectrum – object has been modified anyway.

Return type

resampled Spectrum object. If using inplace=True, the Spectrum

Examples

rescale_mole_fraction(new_mole_fraction, old_mole_fraction=None, inplace=True, ignore_warnings=False, force=False, verbose=True)[source]

Update spectrum with new molar fraction Convoluted values (with slit) are dropped in the process.

Parameters
  • new_mole_fraction (float) – new mole fraction

  • old_mole_fraction (float, or None) – if None, current mole fraction (conditions[‘mole_fraction’]) is used

Other Parameters
  • inplace (boolean) – if True, modifies the Spectrum object directly. Else, returns a copy. Default True.

  • force (boolean) – if False, won’t allow rescaling to 0 (not to loose information). Default False

Returns

s – Cropped Spectrum. If inplace=True, Spectrum has been updated directly anyway. Allows chaining

Return type

Spectrum

Examples

::

s.rescale_mole_fraction(0.2)

Notes

Implementation:

similar to rescale_path_length() but we have to scale abscoeff & emisscoeff Note that this is valid only for small changes in mole fractions. Then, the change in line broadening becomes significant

rescale_path_length(new_path_length, old_path_length=None, inplace=True, force=False)[source]

Rescale spectrum to new path length. Starts from absorption coefficient and emission coefficient, and solves the RTE again for the new path length Convoluted values (with slit) are dropped in the process.

Parameters
  • new_path_length (float) – new path length

  • old_path_length (float, or None) – if None, current path length (conditions[‘path_length’]) is used

Other Parameters
  • inplace (boolean) – if True, modifies the Spectrum object directly. Else, returns a copy. Default True.

  • force (boolean) – if False, won’t allow rescaling to 0 (not to loose information). Default False

Returns

s – Cropped Spectrum. If inplace=True, Spectrum has been updated directly anyway. Allows chaining

Return type

Spectrum

Examples

::
for path in [0.1, 10, 100]:

s.rescale_path_length(10, inplace=False).plot(nfig=’same’)

Notes

Implementation:

To deal with all the input cases, we first make a list of what has to be recomputed, and what has to be recalculated

save(*args, **kwargs)[source]

Alias to Spectrum.store.

See Spectrum.store for documentation

savetxt(filename, var, wunit='default', Iunit='default')[source]

Export spectral quantity var to filename.

(note that by doing this you will loose additional information, such

as the calculation conditions or the units. You better save a Spectrum object under a .spec file with store() and load it afterwards with load_spec())

Parameters
  • filename (str) – file name

  • var (str) – which spectral variable ot export

Other Parameters

wunit, Iunit, medium (str) – see get() for more information

Notes

Export variable as:

np.savetxt(filename, np.vstack(self.get(var, wunit=wunit, Iunit=Iunit,
                                    medium=medium)).T, header=header)
sort(inplace=True)[source]

Sort the Spectrum by wavelength / wavenumber.

Parameters

inplace (bool, optional) – if True, modifies the Spectrum object directly. The default is True.

Returns

s – sorted Spectrum. If inplace=True, Spectrum has been updated directly anyway. Allows chaining.

Return type

Spectrum

store(path, discard=['lines', 'populations'], compress=True, add_info=None, add_date=None, if_exists_then='error', verbose=True)[source]

Save a Spectrum object in JSON format. Object can be recovered with load_spec(). If many Spectrum are saved in a same folder you can view their properties with the SpecDatabase structure.

Parameters
  • path (path to folder (database) or file) – if a folder, file is saved to database and name is generated automatically. if a file name, then Spectrum is saved to this file and the later formatting options dont apply

  • file (str) – explicitely give a filename to save

  • compress (boolean) – if False, save under text format, readable with any editor. if True, saves under binary format. Faster and takes less space. If 2, removes all quantities that can be regenerated with s.update(), e.g, transmittance if abscoeff and path length are given, radiance if emisscoeff and abscoeff are given in non-optically thin case, etc. Default True.

  • add_info (list) – append these parameters and their values if they are in conditions example:

    add_info = ['Tvib', 'Trot']
    
  • discard (list of str) – parameters to exclude, for instance to save some memory for instance Default [lines, populations]: retrieved Spectrum looses the line_survey() ability, and plot_populations() (but it saves tons of memory!)

  • if_exists_then ('increment', 'replace', 'error', 'ignore') – what to do if file already exists. If increment an incremental digit is added. If replace file is replaced (yeah). If 'ignore' no file is created. If error (or anything else) an error is raised. Default error

Returns

Return type

Returns filename used

Notes

If many spectra are stored in a folder, it may be time to set up a SpecDatabase structure to easily see all Spectrum conditions and get Spectrum that suits specific parameters.

Implementation:

Shouldnt rely on a Database. One may just want to store/load a Spectrum once.

Examples

Store a spectrum in compressed mode, regenerate quantities after loading:

from radis import load_spec
s.store('test.spec', compress=True)   # s is a Spectrum
s2 = load_spec('test.spec')
s2.update()                           # regenerate missing quantities
take(var, copy_lines=False)[source]
Parameters

var (str) – spectral quantity

Returns

s – same Spectrum with only the var spectral quantity

Return type

Spectrum

Examples

Use it to chain other commands

s.take('radiance').normalize().plot()
trim(inplace=True)[source]

Remove nan common to all arrays on each side of the Spectrum.

Returns a smaller Spectrum (inplace or not).

Returns

s – directly anyway. Allows chaining.

Return type

Spectrum : trimmed Spectrum. If inplace=True, Spectrum has been updated

units[source]

units for spectral quantities.

Type

dict

update(quantity='all', optically_thin='default', verbose=True)[source]

Calculate missing quantities: ex: if path_length and emisscoeff are given, recalculate radiance_noslit.

Parameters
  • spec (Spectrum)

  • quantity (str) – name of the spectral quantity to recompute. If ‘same’, only the quantities in the Spectrum are recomputed. If ‘all’, then all quantities that can be derived are recomputed. Default ‘all’.

  • optically_thin (True, False, or ‘default’) – determines whether to calculate radiance with or without self absorption. If ‘default’, the value is determined from the self_absorption key in Spectrum.conditions. If not given, False is taken. Default ‘default’ Also updates the self_absorption value in conditions (creates it if doesnt exist)

Transmittance(s: radis.spectrum.spectrum.Spectrum) radis.spectrum.spectrum.Spectrum[source]

Returns a new Spectrum with only the transmittance component of s

Parameters

s (Spectrum) – Spectrum object

Returns

s_trSpectrum object, with only the transmittance, absorbance and/or abscoeff part of s, where radiance_noslit , emisscoeff and emissivity_noslit (if they exist) have been set to 0

Return type

Spectrum

Examples

This function is useful to use Spectrum algebra operations:

s = calc_spectrum(...)   # contains emission & absorption arrays
tr = Transmittance(s)    # contains 'radiance_noslit' array only
tr -= 0.1    # arithmetic operation is applied to Transmittance only

Equivalent to:

rad = s.take('transmittance')
Transmittance_noslit(s: radis.spectrum.spectrum.Spectrum) radis.spectrum.spectrum.Spectrum[source]

Returns a new Spectrum with only the transmittance_noslit component of s

Parameters

s (Spectrum) – Spectrum object

Returns

s_trSpectrum object, with only transmittance_noslit defined

Return type

Spectrum

Examples

This function is useful to use Spectrum algebra operations:

s = calc_spectrum(...)   # contains emission & absorption arrays
tr = Transmittance_noslit(s) # contains 'radiance_noslit' array only
tr -= 0.1    # arithmetic operation is applied to Transmittance_noslit only

Equivalent to:

rad = s.take('transmittance_noslit')
calculated_spectrum(w, I, wunit='nm', Iunit='mW/cm2/sr/nm', conditions=None, cond_units=None, populations=None, name=None) radis.spectrum.spectrum.Spectrum[source]

Convert (w, I) into a Spectrum object that has unit conversion, plotting and slit convolution capabilities.

Parameters
  • w (np.array) – wavelength, or wavenumber

  • I (np.array) – intensity (no slit)

  • wunit ('nm', 'cm-1', 'nm_vac') – wavespace unit: wavelength in air ('nm'), wavenumber ('cm-1'), or wavelength in vacuum ('nm_vac'). Default 'nm'.

  • Iunit (str) – intensity unit (can be ‘counts’, ‘mW/cm2/sr/nm’, etc…). Default ‘mW/cm2/sr/nm’ (note that non-convoluted Specair spectra are in ‘mW/cm2/sr/µm’)

Other Parameters
  • conditions (dict) – (optional) calculation conditions to be stored with Spectrum. Default None

  • cond_units (dict) – (optional) calculation conditions units. Default None

  • populations (dict) – populations to be stored in Spectrum. Default None

  • name (str) – (optional) give a name

Examples

# w, I are numpy arrays for wavelength and radiance
from radis import calculated_spectrum
s = calculated_spectrum(w, I, wunit='nm', Iunit='W/cm2/sr/nm')     # creates 'radiance_noslit'
experimental_spectrum(w, I, wunit='nm', Iunit='counts', conditions={}, cond_units=None, name=None) radis.spectrum.spectrum.Spectrum[source]

Convert (w, I) into a Spectrum object that has unit conversion and plotting capabilities. Convolution is not available as the spectrum is assumed to have be measured experimentally (hence it is already convolved with the slit function)

Parameters
  • w (np.array) – wavelength, or wavenumber

  • I (np.array) – intensity

  • wunit ('nm', 'cm-1', 'nm_vac') – wavespace unit: wavelength in air ('nm'), wavenumber ('cm-1'), or wavelength in vacuum ('nm_vac'). Default 'nm'.

  • Iunit (str) – intensity unit (can be ‘counts’, ‘mW/cm2/sr/nm’, etc…). Default ‘counts’ (default Winspec output)

Other Parameters
  • conditions (dict) – (optional) calculation conditions to be stored with Spectrum

  • cond_units (dict) – (optional) calculation conditions units

  • name (str) – (optional) give a name

Examples

Load and plot an experimental spectrum:

from numpy import loadtxt
from radis import experimental_spectrum
w, I = loadtxt('my_file.txt').T    # assuming 2 columns
s = experimental_spectrum(w, I, Iunit='mW/cm2/sr/nm')             # creates 'radiance'
s.plot()
get_diff(s1: radis.spectrum.spectrum.Spectrum, s2: radis.spectrum.spectrum.Spectrum, var, wunit='default', Iunit='default', resample=True, diff_window=0)[source]

Get the difference between 2 spectra.

Basically returns w1, I1 - I2 where (w1, I1) and (w2, I2) are the values of s1 and s2 for variable var. (w2, I2) is linearly interpolated if needed.

\[dI = I_1 - I_2\]
Parameters
  • s1, s2 (Spectrum objects) – 2 spectra to compare.

  • var (str) – spectral quantity (ex: 'radiance', 'transmittance'…)

  • wunit ('nm', 'cm-1', 'nm_vac') – waveunit to compare in: wavelength air, wavenumber, wavelength vacuum

  • Iunit (str) – if 'default' use s1 unit for variable var

  • medium (‘air’, ‘vacuum’, default’) – propagating medium to compare in (if in wavelength)

Other Parameters
  • resample (bool) – if not True, wavelength must be equals. Else, resample s2 on s1 if needed.

  • diff_window (int) – If non 0, calculates diff by offsetting s1 by diff_window number of units on either side, and returns the minimum. Compensates for experimental errors on the w axis. Default 0. (look up code for more details…)

Returns

w1, Idiff – difference interpolated on the wavespace range of the first Spectrum

Return type

array

Notes

Uses curve_substract() internally

get_distance(s1: radis.spectrum.spectrum.Spectrum, s2: radis.spectrum.spectrum.Spectrum, var, wunit='default', Iunit='default', resample=True)[source]

Get a regularized Euclidian distance between two spectra s1 and s2

This regularized Euclidian distance minimizes the effect of a small shift in wavelength between the two spectra

\[D(w_1)[i] = \sqrt{ \sum_j (\hat{I_1}[i] - \hat{I_2}[j] )^2 + (\hat{w_1}[i] - \hat{w_2}[j])^2}\]

Where values are normalized as:

\[\hat{A} = \frac{A}{max(A) - min(A)}\]

If waveranges dont match, s2 is interpolated over s1.

Warning

This is a distance on both the waverange and the intensity axis. It may be used to compensate for a small offset in your experimental spectrum (due to wavelength calibration, for instance) but can lead to wrong fits easily. Plus, it is very cost-intensive! Better use get_residual() for an automatized procedure.

Parameters
  • s1, s2 (Spectrum objects) – Spectrum

  • var (str) – spectral quantity

  • wunit ('nm', 'cm-1', 'nm_vac') – waveunit to compare in: wavelength air, wavenumber, wavelength vacuum

  • Iunit (str) – if 'default' use s1 unit for variable var

  • medium (‘air’, ‘vacuum’, default’) – propagating medium to compare in (if in wavelength)

Notes

Uses curve_distance() internally

get_ratio(s1: radis.spectrum.spectrum.Spectrum, s2: radis.spectrum.spectrum.Spectrum, var, wunit='default', Iunit='default', resample=True)[source]

Get the ratio between two spectra Basically returns w1, I1 / I2 where (w1, I1) and (w2, I2) are the values of s1 and s2 for variable var. (w2, I2) is linearly interpolated if needed.

\[R = I_1 / I_2\]
Parameters
  • s1, s2 (Spectrum objects) – Spectrum

  • var (str) – spectral quantity

  • wunit ('nm', 'cm-1', 'nm_vac') – waveunit to compare in: wavelength air, wavenumber, wavelength vacuum

  • Iunit (str) – if 'default' use s1 unit for variable var

Notes

Uses curve_divide() internally

get_residual(s1: radis.spectrum.spectrum.Spectrum, s2: radis.spectrum.spectrum.Spectrum, var, norm='L2', ignore_nan=False, diff_window=0, normalize=False, normalize_how='max', wunit='default') float[source]

Returns L2 norm of s1 and s2

For I1, I2, the values of variable var in s1 and s2, respectively, residual is calculated as:

For L2 norm:

\[res = \frac{\sqrt{\sum_i {(s_1[i]-s_2[i])^2}}}{N}.\]

For L1 norm:

\[res = \frac{\sqrt{\sum_i {|s_1[i]-s_2[i]|}}}{N}.\]
Parameters
  • s1, s2 (Spectrum objects) – if not on the same range, s2 is resampled on s1.

  • var (str) – spectral quantity

  • norm (‘L2’, ‘L1’) – which norm to use

Other Parameters
  • ignore_nan (boolean) – if True, ignore nan in the difference between s1 and s2 (ex: out of bound) when calculating residual. Default False. Note: get_residual will still fail if there are nan in initial Spectrum.

  • normalize (bool, or tuple) – if True, normalize the two spectra before calculating the residual. If a tuple (ex: (4168, 4180)), normalize on this range only. The unit is that of the first Spectrum by default (use wunit to change). Ex:

    s_exp   # in 'nm'
    s_calc  # in 'cm-1'
    get_residual(s_exp, s_calc, normalize=(4178, 4180))  # interpreted as 'nm'
    
  • normalize_how ('max', 'area', 'mean') – how to normalize. 'max' is the default but may not be suited for very noisy experimental spectra. 'area' will normalize the integral to 1. 'mean' will normalize by the mean amplitude value

  • wunit (str) – used if normalized is a range

Notes

0 values for I1 yield nans except if I2 = I1 = 0

when s1 and s2 dont have the size wavespace range, they are automatically resampled through get_diff on ‘s1’ range

Implementation of L2 norm:

np.sqrt((dI**2).sum())/len(dI)

Implementation of L1 norm:

np.abs(dI).sum()/len(dI)
get_residual_integral(s1: radis.spectrum.spectrum.Spectrum, s2: radis.spectrum.spectrum.Spectrum, var: str, ignore_nan: bool = False) float[source]

Returns integral of the difference between two spectra s1 and s2, relatively to the integral of spectrum s1.

Compared to get_residual(), this tends to cancel the effect of the gaussian noise of an experimental spectrum.

res = trapz(I2-I1, w1) / trapz(I1, w1)

Note: when the considered variable is transmittance or transmittance_noslit, the upper integral is used (up to 1) to normalize the integral difference

res = trapz(I2-I1, w1) / trapz(1-I1, w1)
Parameters
  • s1, s2 (Spectrum objects) – Spectrum

  • var (str) – spectral quantity

Other Parameters

ignore_nan (boolean) – if True, ignore nan in the difference between s1 and s2 (ex: out of bound) when calculating residual. Default False. Note: get_residual_integral will still fail if there are nan in initial Spectrum.

Notes

For I1, I2, the values of ‘var’ in s1 and s2, respectively, residual is calculated as:

res = trapz(I2-I1, w1) / trapz(I1, w1)

0 values for I1 yield nans except if I2 = I1 = 0

when s1 and s2 dont have the size wavespace range, they are automatically resampled through get_diff on ‘s1’ range

plot_diff(s1: radis.spectrum.spectrum.Spectrum, s2: radis.spectrum.spectrum.Spectrum, var=None, wunit='default', Iunit='default', resample=True, method='diff', diff_window=0, show_points=False, label1=None, label2=None, figsize=None, title=None, nfig=None, normalize=False, verbose=True, save=False, show=True, show_residual=False, lw_multiplier=1, diff_scale_multiplier=1, discard_centile=0, plot_medium='vacuum_only', legendargs={'loc': 'best'}, show_ruler=False)[source]

Plot two spectra, and the difference between them. method= allows you to plot the absolute difference, ratio, or both.

If waveranges dont match, s2 is interpolated over s1.

Parameters
  • s1, s2 (Spectrum objects)

  • var (str, or None) – spectral quantity to plot (ex: 'abscoeff'). If None, plot the first one in the Spectrum from 'radiance', 'radiance_noslit', 'transmittance', etc.

  • wunit ('default', 'nm', 'cm-1', 'nm_vac') – wavespace unit: wavelength air, wavenumber, wavelength vacuum. If 'default', use first spectrum wunit

  • Iunit (str) – if 'default', use first spectrum unit

  • method ('distance', 'diff', 'ratio', or list of them.) – If 'diff', plot difference of the two spectra. If 'distance', plot Euclidian distance (note that units are meaningless then) If 'ratio', plot ratio of two spectra Default 'diff'.

    Warning

    with 'distance', calculation scales as ~N^2 with N the number of points in a spectrum (against ~N with 'diff'). This can quickly override all memory.

    Can also be a list:

    method=['diff', 'ratio']
    
  • normalize (bool) – Normalize the spectra to be ploted

Other Parameters
  • diff_window (int) – If non 0, calculates diff by offsetting s1 by diff_window number of units on either side, and returns the minimum. Kinda compensates for experimental errors on the w axis. Default 0. (look up code to understand…)

  • show_points (boolean) – if True, make all points appear with ‘o’

  • label1, label2 (str) – curve names

  • figsize – figure size

  • nfig (int, str) – figure number of name

  • title (str) – title

  • verbose (boolean) – if True, plot stuff such as rescale ratio in normalize mode. Default True

  • save (str) – Default is False. By default won’t save anything, type the path of the destination if you want to save it (format in the name).

  • show (Bool) – Default is True. Will show the plots : bad if there are more than 20.

  • show_residual (bool) – if True, calculates and shows on the graph the residual in L2 norm. See get_residual(). diff_window is used in the residual calculation too. normalize has no effect.

  • diff_scale_multiplier (float) – dilate the diff plot scale. Default 1

  • discard_centile (int) – if not 0, discard the firsts and lasts centile when setting the limits of the diff window. Example:

    discard_centile=1     #  --> discards the smallest 1% and largest 1%
    discard_centile=10    #  --> discards the smallest 10% and largest 10%
    

    Useful to remove spikes in a ratio, for instance. Note that this does not change the values of the residual. It’s just a plot feature. Default 0

  • plot_medium (bool, 'vacuum_only') – if True and wunit are wavelengths, plot the propagation medium in the xaxis label ([air] or [vacuum]). If 'vacuum_only', plot only if wunit=='nm_vac'. Default 'vacuum_only' (prevents from inadvertently plotting spectra with different propagation medium on the same graph).

  • legendargs (dict) – format arguments forwarded to the legend

  • show_ruler (bool) – if True, add a ruler tool to the Matplotlib toolbar.

    Warning

    still experimental in 0.9.30 ! Try it, feedback welcome !

Returns

  • fig (figure) – fig

  • [ax0, ax1] (axes) – spectra and difference axis

Examples

Simple use:

from radis import plot_diff
plot_diff(s10, s50)                # s10, s50 are two spectra

Advanced use, plotting the total power in the label, and getting the figure and axes handle to edit them afterwards:

Punit = 'mW/cm2/sr'
fig, axes = plot_diff(s10, s50, 'radiance_noslit', figsize=(18,6),
      label1='brd 10 cm-1, P={0:.2f} {1}'.format(s10.get_power(unit=Punit),Punit),
      label2='brd 50 cm-1, P={0:.2f} {1}'.format(s50.get_power(unit=Punit),Punit)
      )
# modify fig, axes..

See an example in Compare two Spectra, which produces the output below:

https://radis.readthedocs.io/en/latest/_images/cdsd4000_vs_hitemp_3409K.svg

If you wish to plot in a logscale, it can be done in the following way:

fig, [ax0, ax1] = plot_diff(s0, s1, normalize=False, verbose=False)
ylim0 = ax0.get_ybound()
ax0.set_yscale("log")
ax0.set_ybound(ylim0)

See also

get_diff(), get_ratio(), get_distance(), get_residual(), compare_with()

transmittance_spectrum(w, T, wunit='nm', Tunit='', conditions=None, cond_units=None, name=None) radis.spectrum.spectrum.Spectrum[source]

Convert (w, I) into a Spectrum object that has unit conversion, plotting and slit convolution capabilities.

Parameters
  • w (np.array) – wavelength, or wavenumber

  • T (np.array) – transmittance (no slit)

  • wunit ('nm', 'cm-1', 'nm_vac') – wavespace unit: wavelength in air ('nm'), wavenumber ('cm-1'), or wavelength in vacuum ('nm_vac'). Default 'nm'.

  • Iunit (str) – intensity unit. Default "" (adimensionned)

Other Parameters
  • conditions (dict) – (optional) calculation conditions to be stored with Spectrum

  • cond_units (dict) – (optional) calculation conditions units

  • name (str) – (optional) give a name

Examples

# w, T are numpy arrays for wavelength and transmittance
from radis import transmittance_spectrum
s2 = transmittance_spectrum(w, T, wunit='nm')                       # creates 'transmittance_noslit'