exo_k.ktable5d

@author: jeremy leconte

Module Contents

class exo_k.ktable5d.Ktable5d(*filename_filters, filename=None, path=None, p_unit='unspecified', file_p_unit='unspecified', kdata_unit='unspecified', file_kdata_unit='unspecified', remove_zeros=False, search_path=None, mol=None, **kwargs)[source]

Bases: exo_k.data_table.Data_table

A class that handles tables of k-coefficients with a variable gas. Based on the Data_table class that handles basic operations common to Ktable and Xtable.

This class is specifically designed to deal with LMDZ type ktable where there is a variable gas.

Initializes k coeff table with variable gas and supporting data from various sources (see below by order of precedence)

Parameters:
  • filename (str (optional)) – Relative or absolute name of the file to be loaded.

  • path (str) – If none of the above is specifed, path can point to a directory with a LMDZ type k coeff table. In this case, see read_LMDZ for the keywords to specify.

If there is no input, just creates an empty object to be filled later

See Ktable __init__ mehthod for documentation on p_unit, file_p_unit, kdata_unit, file_kdata_unit, remove_zeros, search_path, and mol arguments.

read_hdf5(filename=None, mol=None, wn_range=None, wl_range=None)[source]

Initializes k coeff table and supporting data from an Exomol hdf5 file

Parameters:

filename (str) – Name of the input hdf5 file

write_hdf5(filename, compression='gzip', compression_level=9, kdata_unit=None, p_unit=None)[source]

Saves data in a hdf5 format

Parameters:

filename (str) – Name of the file to be created and saved

read_LMDZ(path=None, res=None, band=None, mol=None)[source]

Initializes k coeff table and supporting data from a .dat file in a gcm friendly format.

Units are assumed to be cm^2 for kdata and mbar for pressure.

Parameters:
  • path (str) – Name of the directory with the various input files

  • res (str) – “IRxVI” where IR and VI are the numbers of bands in the infrared and visible of the k table to load.

  • band (str) – “IR” or “VI” to specify which band to load.

write_LMDZ(path, band='IR', fmt='%22.15e', write_only_metadata=False)[source]

Saves data in a LMDZ friendly format.

The gcm requires p in mbar and kdata in cm^2/molec. The conversion is done automatically.

Parameters:
  • path (str) – Name of the directory to be created and saved, the one that will contain all the necessary files

  • band (str) – The band you are computing: ‘IR’ or ‘VI’

  • fmt (str) – Fortran format for the corrk file.

  • write_only_metadata (bool, optional) – If True, only supporting files are written (T.dat, p.dat, etc.)

hires_to_ktable(path=None, filename_grid=None, logpgrid=None, tgrid=None, xgrid=None, wnedges=None, quad='legendre', order=20, g_split=0.9, weights=None, ggrid=None, mid_dw=True, write=0, mol=None, grid_p_unit='Pa', p_unit='unspecified', kdata_unit='unspecified', file_kdata_unit='unspecified', remove_zeros=False, **kwargs)[source]

Computes a k coeff table from Hires_spectrum objects.

see exo_k.ktable.Ktable.hires_to_ktable() method for details on the arguments and options.

Parameters:

xgrid (array, np.ndarray) – Input grid in vmr of the variable gas. Needed for a Ktable5d.

Warning

By default, log pressures are specified in Pa in logpgrid!!! If you want to use another unit, do not forget to specify it with the grid_p_unit keyword.

setup_interpolation(log_interp=None)[source]

Creates interpolating functions to be called later on. and loads it as attribute (inplace).

set_kdata(new_kdata)[source]

Changes kdata (inplace). this is preferred to directly accessing kdata because one could forget to run setup_interpolation().

Parameters:

new_kdata (array, np.ndarray) – New array of kdata.

interpolate_kdata(logp_array=None, t_array=None, x_array=None, log_interp=None, logp_interp=True, wngrid_limit=None, **kwargs)[source]

interpolate_kdata interpolates the kdata at on a given temperature and log pressure profile.

Parameters:
  • logp_array (Array) – log 10 pressure array to interpolate to

  • t_array (Array, same size a logp_array) – Temperature array to interpolate to

  • x_array (Array, same size a logp_array) – vmr of variable gas array to interpolate to

  • given (If floats are) –

  • 1. (they are interpreted as arrays of size) –

  • wngrid_limit (list or array, optional) – if an array is given, interpolates only within this array

  • log_interp (bool, dummy) – Dummy variable to be consistent with interpolate_kdata in data_table. Whether the interpolation is linear in kdata or in log(kdata) is actually controlled by self._settings._log_interp but only when the ktable is loaded. If you change that after the loading, you should rerun setup_interpolation().

Returns:

array of shape (logp_array.size, self.Nw , self.Ng)

The interpolated kdata.

remap_logPT(logp_array=None, t_array=None, x_array=None)[source]

remap_logPT re-interpolates the kdata on a new temprature and log pressure grid (inplace). If an x_array is specified, the interpolation occurs along the vmr axis as well.

Parameters:
  • logp_array (Array) – log 10 pressure array to interpolate to

  • t_array (Array) – temperature array to interpolate to

  • x_array (Array) – vmr of variable gas array to interpolate to

Warning

Whether the interpolation is linear in kdata or in log(kdata) is controlled by self._settings._log_interp but only when the ktable is loaded. If you change that after the loading, you should rerun setup_interpolation().

copy(cp_kdata=True)[source]

Creates a new instance of Ktable5d object and (deep) copies data into it

Parameters:

cp_kdata (bool, optional) – If false, the kdata table is not copied and only the structure and metadata are.

Returns:

Ktable

A new Ktable5d instance with the same structure as self.

gindex(g)[source]

Finds the index corresponding to the given g

xindex(x)[source]

Finds the index corresponding to the given x

spectrum_to_plot(p=1e-05, t=200.0, x=1.0, g=None)[source]

provide the spectrum for a given point to be plotted

Parameters:
  • p (float) – Pressure (Ktable pressure unit)

  • t (float) – Temperature(K)

  • g (float) – Gauss point

  • x (float) – Mixing ratio of the species

plot_distrib(ax, p=1e-05, t=200.0, wl=1.0, x=1.0, xscale=None, yscale='log', **kwarg)[source]

Plot the distribution for a given point

Parameters:
  • p (float) – Pressure (Ktable pressure unit)

  • t (float) – Temperature(K)

  • wl (float) – Wavelength (micron)

combine_with(other, x_self=None, x_other=None, **kwargs)[source]

Method to create a new Ktable5d where the kdata of ‘self’ are randomly mixed with ‘other’ (that must be a Ktable).

The main purpose is to add the opacity of a trace species to the background gas of the Ktable5d instance.

Warning

Because:

  • the opacity from the background and variable gases cannot be isolated,

  • The values of the array for the vmr of the variable gas (self.xgrid) are not modified (diluted),

the treatment here is valid only if x_other << 1.

For this reason, x_self should be either left to None, or 1-x_other depending exactly what you want to do. But if you put `x_self`<<1, you are on your own.

Parameters:
  • other (Ktable) – A Ktable object to be mixed with. Dimensions should be the same as self (except for xgrid).

  • x_self (float only, optional) – Volume mixing ratio of self.

  • x_other (float or array, optional) – Volume mixing ratio of the species to be mixed with (other).

If either x_self or x_other are set to None (default), the cross section of the species in question are considered to be already normalized with respect to the mixing ratio.

Returns:

Ktable5d

A new Ktable5d with the opacity of the new species added.

bin_down(wnedges=None, weights=None, ggrid=None, remove_zeros=False, num=300, use_rebin=False, write=0)[source]

Method to bin down a kcoeff table to a new grid of wavenumbers (inplace).

Parameters:
  • wnedges (array, np.ndarray) – Edges of the new bins of wavenumbers (cm-1) onto which the kcoeff should be binned down. if you want Nwnew bin in the end, wnedges.size must be Nwnew+1 wnedges[0] should be greater than self.wnedges[0] (JL20 not sure anymore) wnedges[-1] should be lower than self.wnedges[-1]

  • weights (array, np.ndarray, optional) – Desired weights for the resulting Ktable.

  • ggrid (array, np.ndarray, optional) – Desired g-points for the resulting Ktable. Must be consistent with provided weights. If not given, they are taken at the midpoints of the array given by the cumulative sum of the weights

clip_spectral_range(wn_range=None, wl_range=None)[source]

Limits the data to the provided spectral range (inplace):

  • Wavenumber in cm^-1 if using wn_range argument

  • Wavelength in micron if using wl_range

extend_spectral_range(**kwargs)[source]

Extends the spectral range of an existing table (inplace). The new bins are filled with zeros (except if remove_zeros=True).

See exo_k.data_table.Data_table.extend_spectral_range() method for details on the arguments and options.

remap_g(ggrid=None, weights=None)[source]

Method to resample a kcoeff table to a new g grid (inplace).

Parameters:
  • ggrid (array, np.ndarray) – New grid of abcissas for quadrature

  • weights (array, np.ndarray) – New grid of weights

exo_k.ktable5d.read_Qdat(filename)[source]

Reads Q.dat files LMDZ style and extract the vmr grid.

Parameters:

filename (str) – Path to file to read.

Returns:

background_mol_names: list

list of names of molecules in background gas

var_mol: str

Name of variable molecule

Nx: int

Size of xgrid

xgrid: array, np.ndarray

grid of vmr for variable gas