exo_k.data_table

@author: jeremy leconte

Module Contents

class exo_k.data_table.Data_table[source]

Bases: exo_k.util.spectral_object.Spectral_object

An abstract class that will serve as a basis for Ktable and Xtable. This class includes all the interpolation and remapping methods.

Initializes all attributes to None

filename = None
mol = None
isotopolog_id = 0
pgrid = None
logpgrid = None
tgrid = None
kdata = None
logk = None
Np = None
Nt = None
Ng = None
DOI = 'unknown'
sampling_method = 'unknown'
Date_ID
Nx = None
p_unit = 'unspecified'
kdata_unit = 'unspecified'
_settings
finalize_init(p_unit='unspecified', file_p_unit='unspecified', kdata_unit='unspecified', file_kdata_unit='unspecified', remove_zeros=False)[source]

Common code at the end of the initialization of inheriting classes put here to avoid duplicates

copy_attr(other, cp_kdata=False)[source]

Copy attributes from other

Parameters:
  • other (Data_table) – Data_table object that will be copied

  • cp_kdata (bool, optional) – If False, only metadata are copied

remove_zeros(deltalog_min_value=10.0)[source]

Finds zeros in the kdata and set them to (10.**-deltalog_min_value) times the minimum positive value in the table (inplace).

This is to be able to work in logspace.

Parameters:

deltalog_min_value (float)

convert_p_unit(p_unit='unspecified', file_p_unit='unspecified')[source]

Converts pressure to a new unit (inplace).

Parameters:
  • p_unit (str) – String identifying the pressure units to convert to (e.g. ‘bar’, ‘Pa’, ‘mbar’, or any pressure unit recognized by the astropy.units library). If =’unspecified’, no conversion is done.

  • file_p_unit (str, optional) – String to specify the current pressure unit if it is unspecified or if you have reasons to believe it is wrong (e.g. you just read a file where you know that the pressure grid and the pressure unit do not correspond)

convert_kdata_unit(kdata_unit='unspecified', file_kdata_unit='unspecified')[source]

Converts kdata to a new unit (inplace).

Parameters:
  • kdata_unit (str) – String to identify the units to convert to. Accepts ‘cm^2’, ‘m^2’ or any surface unit recognized by the astropy.units library. If =’unspecified’, no conversion is done. In general, kdata should be kept in ‘per number’ or ‘per volume’ units (as opposed to ‘per mass’ units) as composition will always be assumed to be a number or volume mixing ratio. Opacities per unit mass are not supported yet. Note that you do not need to specify the ‘/molec’ or ‘/molecule’ in the unit.

  • file_kdata_unit (str) – String to specify the current kdata unit if it is unspecified or if you have reasons to believe it is wrong (e.g. you just read a file where you know that the kdata grid and the kdata unit do not correspond)

convert_to_mks()[source]

Converts units to MKS (inplace).

interpolate_kdata(logp_array=None, t_array=None, x_array=1.0, 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. If a volume mixing ratio profile (x_array) is given, the cross section computed for the species is multiplied by x_array to account for the ‘dilution’ of the opacity.

Parameters:
  • logp_array (array, np.ndarray) – log 10 pressure array to interpolate to

  • t_array (array, np.ndarray, same size as logp_array) – Temperature array to interpolate to

  • x_array (None) – Volume mixing ratio array used to renormalize the cross section.

If floats are given, they are interpreted as arrays of size 1.

Parameters:
  • wngrid_limit (list or array) – if an array of to values is given, interpolates only within this array

  • log_interp (bool, optional) – Whether the interpolation is linear in kdata or in log(kdata)

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 temperature and log pressure grid (inplace).

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

  • t_array (Array) – temperature array to interpolate to

  • x_array (dummy argument to be consistent with interpolate_kdata in Ktable5d)

Whether the interpolation is linear in kdata or in log10(kdata) is controlled by self._settings._log_interp

pindex(p)[source]

Finds and returns the index corresponding to the given pressure p (units must be the same as the ktable)

tindex(t)[source]

Finds and returns the index corresponding to the given temperature t (in K)

plot_spectrum(ax, p=1e-05, t=200.0, x=1.0, g=None, x_axis='wls', xscale=None, yscale=None, **kwarg)[source]

Plot the spectrum for a given point

Parameters:
  • ax (pyplot.Axes) – A pyplot axes instance where to put the plot.

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

  • t (float) – Temperature (K)

  • g (float) – Gauss point

  • x (float) – Mixing ratio of the species

  • x_axis (str, optional) – If ‘wls’, x axis is wavelength. Wavenumber otherwise.

  • x/yscale (str, optional) – If ‘log’ log axes are used.

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

Dummy function to be defined in inheriting classes

vmr_normalize(x_self)[source]

Rescales kdata to account for the fact that the gas is not a pure species

Parameters:

x_self (float or array of shape (`self.Np,self.Nt)`) – The volume mixing ratio of the species.

Returns:

array

The vmr normalized kdata (x_self*self.kdata).

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

Method to create a new Data_table where the kdata of ‘self’ are:

  • randomly mixed with ‘other’ for a Ktable.

  • added to ‘other’ for an Xtable

Parameters:
  • other (same class as self) – A Data_table object to be mixed with. Dimensions should be the same as self.

  • x_self (float or array, 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:

Data_table

A new table for the mix

set_kdata(new_kdata)[source]

Changes kdata (inplace). this is preferred to directly accessing kdata because this method checks that the array has the right dimensions

Parameters:

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

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(wngrid_left=None, wngrid_right=None, wnedges_left=None, wnedges_right=None, remove_zeros=False)[source]

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

Parameters:
  • wngrid_left (array, np.ndarray) – Array of wavenumbers to add to the small wn end of the table.

  • wngrid_right (array, np.ndarray) – Array of wavenumbers to add to the high wn end of the table.

Warning

There should not be any overlap between wngrid_left, wngrid_right, and the current wavenumber grid of the table.

Parameters:

remove_zeros (bool (optional)) – Whether zeros in the resulting table should be removed using remove_zeros().

bin_down_cp(wnedges=None, **kwargs)[source]

Creates a copy of the instance and bins it down using the methods in Ktable or Xtable.

See exo_k.ktable.Ktable.bin_down() or exo_k.xtable.Xtable.bin_down() for details on parameters.

Returns:

Ktable or Xtable object

The binned down table.

property shape
Returns the shape that self.kdata should have to be compatible with
the parameter grid.
change_molecule_name(mol_name)[source]

Changes name of the molecule (self.mol attribute).

Parameters:

mol_name (str) – New molecule name

property molar_mass
Computes molar mass from molecule name
blackbody(Temperature, **kwargs)[source]

See Spectral_object.blackbody

write_hdf5_common(f, compression='gzip', compression_level=9, p_unit=None)[source]

Method that writes datasets and attributes that are common to X and Ktables.

Parameters:

f (h5py file instance) – The file to write that is created in daughter classes

toLogK()[source]

Changes kdata to log 10.

Changes kdata back from log to linear scale.