radis.api.exomolapi module

Molecular database (MDB) class

  • MdbExomol is the MDB for ExoMol

Initial code borrowed from the Exojax code (which you should also have a look at !), by @HajimeKawahara, under MIT License.

class MdbExomol(path, molecule, database=None, local_databases=None, name='EXOMOL-{molecule}', nurange=[0.0, inf], margin=0.0, crit=-inf, bkgdatm='Air', broadf=True, engine='vaex', verbose=True, cache=True, skip_optional_data=True)[source]

Bases: DatabaseManager

molecular database of ExoMol

MdbExomol is a class for ExoMol.

Parameters:
  • path (str) – path for Exomol data directory/tag. For instance, “/home/CO/12C-16O/Li2015”

  • nurange (array) – wavenumber range list (cm-1) or wavenumber array

  • margin (float) – margin for nurange (cm-1)

  • crit (float) – line strength lower limit for extraction

  • bkgdatm (str) – background atmosphere for broadening. e.g. H2, He,

  • broadf (bool) – if False, the default broadening parameters in .def file is used

Other Parameters:
  • engine (str) – which memory mapping engine to use : ‘vaex’, ‘pytables’ (HDF5), ‘feather’

  • skip_optional_data (bool) – If False, fetch all fields which are marked as available in the ExoMol definition file. If True, load only the first 4 columns of the states file (“i”, “E”, “g”, “J”). The structure of the columns above 5 depend on the the definitions file (*.def) and the Exomol version. If skip_optional_data=False, two errors may occur:

Notes

The trans/states files can be very large. For the first time to read it, we convert it to the feather or hdf5-format. After the second-time, we use the feather/hdf5 format instead.

Examples

# Init database, download files if needed.
mdb = MdbExomol(
    local_path,
    molecule=molecule,
    name=databank_name,
    local_databases=local_databases,
    # nurange=[load_wavenum_min, load_wavenum_max],
    engine="vaex",
)

# Get cache files to load :
mgr = mdb.get_datafile_manager()
local_files = [mgr.cache_file(f) for f in mdb.trans_file]

# Load files
df = mdb.load(
    local_files,
    columns=columns_exomol,
    lower_bound=([('nu_lines', load_wavenum_min)] if load_wavenum_min else []) + ([("Sij0", mdb.crit)] if not np.isneginf(mdb.crit) else []),
    upper_bound=([('nu_lines', load_wavenum_max)] if load_wavenum_max else []),
    output="jax", # or "pytables", "vaex"
)

Compare CO cross-sections from HITRAN, HITEMP, GEISA, and ExoMol

Compare CO cross-sections from HITRAN, HITEMP, GEISA, and ExoMol

Calculate a spectrum from ExoMol

Calculate a spectrum from ExoMol

DataFrame columns

nu_lines (nd array): line center (cm-1) Sij0 (nd array): line strength at T=Tref (cm) dev_nu_lines (np array): line center in device (cm-1) logsij0 (np array): log line strength at T=Tref A (np array): Einstein A coefficient elower (np array): the lower state energy (cm-1) gpp (np array): statistical weight jlower (np array): J_lower jupper (np array): J_upper n_Tref (np array): temperature exponent alpha_ref (np array): alpha_ref (gamma0) n_Tref_def: default temperature exponent in .def file, used for jlower not given in .broad alpha_ref_def: default alpha_ref (gamma0) in .def file, used for jlower not given in .broad

References

See also

MdbExomol

QT_interp(T)[source]

interpolated partition function

Parameters:

T (temperature)

Return type:

Q(T) interpolated in jnp.array

download(molec, extension, numtag=None)[source]

Downloading Exomol files

Parameters:
  • molec (like “12C-16O__Li2015”)

  • extension (extension list e.g. [“.pf”,”.def”,”.trans.bz2”,”.states.bz2”,”.broad”])

  • numtag (number tag of transition file if exists. e.g. “11100-11200”)

Notes

The download URL is written in exojax.utils.url.

qr_interp(T)[source]

interpolated partition function ratio

Parameters:

T (temperature)

Return type:

qr(T)=Q(T)/Q(Tref) interpolated in jnp.array

set_broadening_coef(df, alpha_ref_def=None, n_Texp_def=None, output=None, add_columns=True)[source]

setting broadening parameters

Parameters:
  • df (Data Frame)

  • alpha_ref (set default alpha_ref and apply it. None=use self.alpha_ref_def)

  • n_Texp_def (set default n_Texp and apply it. None=use self.n_Texp_def)

  • add_columns (adds alpha_ref and n_Texp columns to df)

Return type:

None. Store values in Data Frame.

to_partition_function_tabulator()[source]

Generate a PartFuncExoMol object

check_bdat(bdat)[source]

checking codes in .broad Args:

bdat: exomol .broad data given by exomolapi.read_broad

Returns:

code level: None, a0, a1, other codes unavailable currently,

e2s(molname_exact)[source]

convert the exact molname (used in ExoMol) to the simple molname

Args:

molname_exact: the exact molname

Returns:

simple molname

Examples:

>>> print(e2s("12C-1H4"))
>>> CH4
>>> print(e2s("23Na-16O-1H"))
>>> NaOH
>>> print(e2s("HeH_p"))
>>> HeH_p
>>> print(e2s("trans-31P2-1H-2H")) #not working
>>> Warning: Exact molname  trans-31P2-1H-2H cannot be converted to simple molname
>>> trans-31P2-1H-2H
get_exomol_database_list(molecule, isotope_full_name=None)[source]

Parse ExoMol website and return list of available databases, and recommended database

Parameters:
Return type:

list of databases, database recommended by ExoMol

Examples

Get CH4 from ExoMol :

databases, recommended = get_exomol_database_list("CH4", "12C-1H4")
>>> ['xsec-YT10to10', 'YT10to10', 'YT34to10'], 'YT34to10'

Or combine with get_exomol_full_isotope_name() to get the isopologue (sorted by terrestrial abundance)

from radis.api.exomolapi import get_exomol_database_list, get_exomol_full_isotope_name
databases, recommended = get_exomol_database_list("CH4", get_exomol_full_isotope_name("CH4", 1))
>>> ['xsec-YT10to10', 'YT10to10', 'YT34to10'], 'YT34to10'
get_exomol_full_isotope_name(molecule, isotope)[source]

Get full isotope name for molecule and isotope number

Parameters:
  • molecule (str)

  • isotope (int) – terrestrial abundance

Examples

get_exomol_full_isotope_name("CH4", 1)
>>> '12C-1H4'
get_list_of_known_isotopes(molecule)[source]

find all isotopes until error ensues

make_j2b(bdat, alpha_ref_default=0.07, n_Texp_default=0.5, jlower_max=None)[source]

compute j2b (code a0, map from jlower to alpha_ref)

Args:

bdat: exomol .broad data given by exomolapi.read_broad alpha_ref_default: default value n_Texp_default: default value jlower_max: maximum number of jlower

Returns:

j2alpha_ref[jlower] provides alpha_ref for jlower j2n_Texp[jlower] provides nT_exp for jlower

make_jj2b(bdat, j2alpha_ref_def, j2n_Texp_def, jupper_max=None)[source]

compute jj2b (code a1, map from (jlower, jupper) to alpha_ref and n_Texp)

Args:

bdat: exomol .broad data given by exomolapi.read_broad j2alpha_ref_def: default value from a0 j2n_Texp_def: default value from a0 jupper_max: maximum number of jupper

Returns:

jj2alpha_ref[jlower,jupper] provides alpha_ref for (jlower, jupper) jj2n_Texp[jlower,jupper] provides nT_exp for (jlower, jupper)

Note:

The pair of (jlower, jupper) for which broadening parameters are not given, jj2XXX contains None.

pickup_gE(states, trans, dic_def, skip_optional_data=True)[source]

extract g_upper (gup), E_lower (elower), and J_lower and J_upper from states DataFrame and insert them into the transition DataFrame.

Parameters:
  • states (states DataFrame - the i, E, g, J are in the 4 first columns)

  • trans (transition numpy array)

  • dic_def (Informations about additional quantum labels)

  • skip_optional_data (bool . If True fetch all quantum labels in dic_def[‘quantum_labels’] from states into transitions (_l for lower, _u for upper states))

Return type:

A, nu_lines, elower, gup, jlower, jupper, mask, **quantum_labels

Notes

We first convert pandas DataFrame to ndarray. The state counting numbers in states DataFrame is used as indices of the new array for states (newstates). We remove the state count numbers as the column of newstate, i.e. newstates[:,k] k=0: E, 1: g, 2: J. Then, we can directly use the state counting numbers as mask.

States by default has columns

#       i      E             g    J    v
0       1      0.0           1    0    0
1       2      1.448467      3    1    0
2       3      4.345384      5    2    0
3       4      8.690712      7    3    0
read_broad(broadf)[source]

Reading broadening file (.broad) :Parameters: broadf (.broad file)

Return type:

broadening info in bdat form (pandas), defined by this instance.

Notes

See Table 16 in https://arxiv.org/pdf/1603.05890.pdf

read_def(deff)[source]

Exomol IO for a definition file

Parameters:

deff (definition file)

Returns:

  • n_Texp (float) – temperature exponent

  • alpha_ref (float) – broadening parameter

  • molmass (float) – molecular mass

  • numinf (List[float]) – limit points ([w(0), w(1), ..., w(n)], n+1 elements) defining the spectral ranges appearing in the name of *.trans.bz2 files (["w(0)-w(1)", "w(1)-w(2)", ..., w(n-1)-w(n)], n elements)

  • numtag (List[str]) – tag for wavelength ranges.

  • Note – For some molecules, ExoMol provides multiple trans files. numinf and numtag are the ranges and identifiers for the multiple trans files.

read_pf(pff)[source]

Exomol IO for partition file

Note:

T=temperature QT=partition function

Args:

pff: partition file

Returns:

partition data in pandas DataFrame

read_states(statesf, dic_def, engine='vaex', skip_optional_data=True, print_states=False)[source]

Exomol IO for a state file

Notes

States f format

i=state counting number
E=state energy
g=state degeneracy
J=total angular momentum

See Table 11 in https://arxiv.org/pdf/1603.05890.pdf

Parameters:
  • statesf (str) – state file

  • dic_def (dict) – Info from def file to read extra quantum numbers

  • engine (str) – parsing engine to use (‘vaex’, ‘csv’)

  • skip_optional_data (bool) – If False, fetch all fields which are marked as available in the ExoMol definition file. If True, load only the first 4 columns of the states file (“i”, “E”, “g”, “J”). The structure of the columns above 5 depend on the the definitions file (*.def) and the Exomol version. If skip_optional_data=False, two errors may occur:

  • print_states (bool) – print some info. When .def file looks inconsistent with .states file, turn ON and report them in Issue.

Returns:

  • states data in pandas DataFrame

  • If 'vaex', also writes a local hdf5 file statesf.with_suffix('.hdf5')

References

read_trans(transf, engine='vaex')[source]

Exomol IO for a transition file

Notes

Transf format

i_upper=Upper state counting number
i_lower=Lower state counting number
A=Einstein coefficient in s-1
nu_lines=transition wavenumber in cm-1

See Table 12 in https://arxiv.org/pdf/1603.05890.pdf [Exomol-2016]

Parameters:
  • transf (transition file)

  • engine (parsing engine to use (‘vaex’, ‘csv’))

Return type:

transition data in vaex/pandas DataFrame