exo_k.util
¶
@author: jeremy leconte __init__ module to load util subpackage
Submodules¶
Package Contents¶
-
PI
= 3.14159265359¶
-
PLANCK
= 6.62606957e-34¶
-
C_LUM
= 299792458¶
-
KBOLTZ
= 1.380648813e-23¶
-
SIG_SB
= 5.670367e-08¶
-
PLANCK_CST1MKS
¶
-
PLANCK_CST1
¶
-
PLANCK_CST2
¶
-
PLANCK_CST1_lamb
¶
-
Bnu
(nu, T)[source]¶ Computes the Planck law in wavenumber domain.
- Parameters
nu (float) – The wavenumber in cm^-1.
T (float) – Temperature (K).
- Returns
Bnu in W/m^2/str/cm^-1 (not SI units).
- Return type
float
sigma*T^4=Int(Bnu dnu) where nu is expressed in cm^-1
-
Bnu_integral
(nu_edges, T)[source]¶ Computes the integral of the Planck function in wavenumber bins.
- Parameters
nu_edges (array) – Edges of the wavenumber bins in cm^-1.
T (float) – Temperature (K).
- Returns
Bnu_integral in W/m^2/str.
- Return type
array
sigma*T^4/PI=sum of Bnu_integral
-
Bnu_integral_num
(nu_edges, T, n=30)[source]¶ Computes the integral of the Planck function in wavenumber bins.
- Parameters
nu_edges (array) – Edges of the wavenumber bins in cm^-1.
T (float) – Temperature (K).
n (int) – number of terms to take into account in the black body calculation.
- Returns
Bnu_integral_num: Integral of the source function at temperatures T_array inside the bins in W/m^2/str
- Return type
Array of shape (T_array.size,nu_edges.size-1)
sigma*T^4=sum of Bnu_integral_num
-
Bnu_integral_array
(nu_edges, T_array, Nw, Nt, n=30)[source]¶ Computes the integral of the Planck function in wavenumber bins.
- Parameters
nu_edges (array) – Edges of the wavenumber bins in cm^-1.
T (array) – Array of temperature (K).
Nw (int) – Number of spectral bins.
Nt (int) – Number of temperatures.
n (int, optional) – number of terms to take into account in the black body calculation.
- Returns
Bnu_integral_num: Integral of the source function at temperatures T_array inside the bins in W/m^2/str
- Return type
Array of shape (T_array.size,nu_edges.size-1)
sigma*T^4=sum of Bnu_integral_num
-
Bnu_integral_old
(nu_edges, T)[source]¶ Deprecated.
Computes the integral of the Planck function in wavenumber bins. nu_edges is an array of the edges of the wavenumber bins in cm^-1.
sigma*T^4=sum of Bnu_integral
Bnu_integral is in W/m^2/str
-
Bmicron
(lamb, T)[source]¶ Computes the Planck law in wavelength domain.
- Parameters
lamb (float) – The wavelength in microns.
T (float) – Temperature (K).
- Returns
- float
Bmicron in W/m^2/str/micron (not SI units).
sigma*T^4=Int(Bmicron d lamb) where lamb is expressed in microns
-
BBnu_Stellar_Spectrum
(nu, T, Flux)[source]¶ Computes the outgoing spectral flux of a black body at temperature T so that its bolometric flux is equal to the input ‘Flux’
- Parameters
nu (float) – The wavenumber in cm^-1.
T (float) – Temperature (K).
Flux (float) – bolometric flux for rescaling.
- Returns
BBnu_Stellar_Spectrum is in W/m^2/cm^-1 (not SI units).
- Return type
float
sigma*T^4=Int(Bnu dnu) where nu is expressed in cm^-1
-
BBmicron_Stellar_Spectrum
(lamb, T, Flux)[source]¶ Computes the outgoing spectral flux of a black body at temperature T so that its bolometric flux is equal to the input ‘Flux’
- Parameters
lamb (float) – The wavelength in microns.
T (float) – Temperature (K).
Flux (float) – bolometric flux for rescaling.
- Returns
BBmicron_Stellar_Spectrum is in W/m^2/micron (not SI units).
- Return type
float
sigma*T^4=Int(Bmicron d lamb) where lamb is expressed in microns.
-
wavelength_grid_R
(lambda_min, lambda_max, R)[source]¶ Creates a wavelength grid starting from lambda_min and ending at lambda_max with a resolution of R (roughly).
- Parameters
lambda_min/lambda_max (int or float) – Min/max wavelength (in micron).
R (int or float) – Resolution.
- Returns
Grid of wavelength.
- Return type
array
-
wavenumber_grid_R
(nu_min, nu_max, R)[source]¶ Creates a wavenumber grid starting from nu_min and ending at nu_max with a resolution of R (roughly).
- Parameters
nu_min/nu_max (int or float) – Min/max wavenumber (in cm^-1).
R (int or float) – Resolution.
- Returns
Grid of wavenumber.
- Return type
array
-
rad_prop_corrk
(dcol_density, opacity_prof, mu0)[source]¶ Computes the optical depth of each of the layers for the opacity given for every wavelength and g point.
- Parameters
dcol_density (array) – Column density per layer (molecules/m^2/layer).
opacity_prof (array) – effective cross section per molecule (m^2/molecule).
mu0 (cosine of the equivalent zenith angle for the rays.) – e.g. Cos(60deg)=0.5 is chosen (See, for example, chap 4 of Pierrehumbert 2010)
- Returns
tau (Array) – cumulative optical depth from the top (with zeros at the top of the first layer)
dtau (Array) – optical depth of each layer for each wavenumber.
-
rad_prop_xsec
(dcol_density, opacity_prof, mu0)[source]¶ Computes the optical depth of each of the layers for the opacity given for every wavelength.
- Parameters
dcol_density (array) – Column density per layer (molecules/m^2/layer).
opacity_prof (array) – effective cross section per molecule (m^2/molecule).
mu0 (cosine of the equivalent zenith angle for the rays.) – e.g. Cos(60deg)=0.5 is chosen (See, for example, chap 4 of Pierrehumbert 2010)
- Returns
tau (Array) – cumulative optical depth from the top (with zeros at the top of the first layer)
dtau (Array) – optical depth of each layer for each wavenumber.
-
path_integral_corrk
(Nlay, Nw, Ng, tangent_path, density_prof, opacity_prof, weights)[source]¶ Computes the transmittance for each layer of the atmosphere for k-coefficients.
- Parameters
Nlay (int) – Number of layers.
Nw (int) – Number of wavenumber bins.
Ng (int) – Number of gauss points.
tangent_path (List of arrays) – Triangular array of tangent path for each layer.
density_prof (array) – Array with the number density of molecules in the atmosphere.
opacity_prof (array) – Effective cross section of the atmsophere for each (layer,wavenumber,g-point).
weights (array) – Weights for the quadrature in g-space.
- Returns
transmittance – Transmittance for each layer and wavenumber bin (Exp(-tau_sigma)).
- Return type
array
-
path_integral_xsec
(Nlay, Nw, tangent_path, density_prof, opacity_prof)[source]¶ Computes the transmittance for each layer of the atmosphere from cross sections.
- Parameters
Nlay (int) – Number of layers.
Nw (int) – Number of wavenumber bins.
tangent_path (List of arrays) – Triangular array of tangent path for each layer.
density_prof (array) – Array with the number density of molecules in the atmosphere.
opacity_prof (array) – Effective cross section of the atmsophere for each (layer,wavenumber).
- Returns
transmittance – Transmittance for each layer and wavenumber bin (Exp(-tau_sigma)).
- Return type
array
-
G
= 6.67384e-11¶
-
RSOL
= 695500000.0¶
-
RJUP
= 71492000.0¶
-
RJUP_POL
= 69911000.0¶
-
PI
= 3.14159265359¶
-
MSOL
= 1.98847542e+30¶
-
MJUP
= 1.898e+27¶
-
AU
= 149597871000.0¶
-
KBOLTZ
= 1.380648813e-23¶
-
RGP
= 8.31446¶
-
PLANCK
= 6.62606957e-34¶
-
C_LUM
= 299792458¶
-
SIG_SB
= 5.670367e-08¶
-
N_A
= 6.022140857e+23¶
-
tiny
= 1e-13¶
-
OneMintiny
¶
-
ktable_long_name_attributes
¶
-
nemesis_hitran_id_numbers
¶
-
class
Spectrum
(value=None, wns=None, wnedges=None, filename=None, from_taurex=False, dataset='native_spectrum')[source]¶ Bases:
exo_k.util.spectral_object.Spectral_object
A class defining a Spectrum object to plot and manipulate.
Instanciate with a value, bin centers, and bin edges. Can also load a Taurex spectrum if filename is provided.
-
plot_spectrum
(self, ax, per_wavenumber=True, x_axis='wls', xscale=None, yscale=None, **kwarg)[source]¶ Plot the spectrum
- Parameters
ax (
pyplot.Axes
) – A pyplot axes instance where to put the plot.per_wavenumber (bool, optional) – Defines the units of spectral flux density. False converts to per wavelength units.
x_axis (str, optional) – If ‘wls’, x axis is wavelength. Wavenumber otherwise.
x/yscale (str, optional) – If ‘log’ log axes are used.
-
bin_down
(self, wnedges)[source]¶ Bins down the spectrum to a new grid of wnedges by conserving area.
- Parameters
wnedges (array) – Wavenumbers of the bin edges to be used
-
bin_down_cp
(self, wnedges)[source]¶ Returns a new binned down spectrum to a grid of wnedges by conserving area.
- Parameters
wnedges (array) – Wavenumbers of the bin edges to be used
- Returns
Binned down spectrum
- Return type
-
property
total
(self)¶ Defines the weighted sum over the spectrum
-
read_hdf5
(self, filename=None)[source]¶ Reads data in a hdf5 file
- Parameters
filename (str) – Name of the file to be created and saved
-
-
exception
EndOfFile
[source]¶ Bases:
Exception
Error for an end of file
Initialize self. See help(type(self)) for accurate signature.
-
select_kwargs
(kwargs, filter_keys_list=[])[source]¶ Function to select only some keyword arguments from a kwargs dict to pass to a function afterward.
- Parameters
kwargs (dict) – Dictionary of keyword arguments and values.
filter_keys_list (list of str) – Names of the keys to select.
- Returns
filtered_kwargs – A dictionary with only the selected keys (if they were present in kwargs).
- Return type
dict
Examples
>>> def func(allo=None): >>> print(allo) >>> >>> def bad_global_func(**kwargs): >>> print(kwargs) >>> func(**kwargs) >>> >>> def good_global_func(**kwargs): >>> print(kwargs) >>> func(**select_kwargs(kwargs,['allo']))
>>> bad_global_func(allo=3.,yeah=None) {'allo': 3.0, 'yeah': None} TypeError: func() got an unexpected keyword argument 'yeah'
>>> good_global_func(allo=3.,yeah=None) {'allo': 3.0, 'yeah': None} 3.0
-
create_fname_grid_Kspectrum_LMDZ
(Np, Nt, Nx=None, suffix='', nb_digit=3)[source]¶ Creates a grid of filenames consistent with Kspectrum/LMDZ from the number of pressure, temperatures (, and vmr) points (respectively) in the grid.
- Parameters
Nt (Np,) – Number of Pressure and temperature points.
Nx (int, optional) – Number of vmr points if there is a variable gas
suffix (str, optional) – String to add behind the ‘k001’
nb_digit (int, optional) – Total number of digit including the zeros (‘k001’ is 3)
- Returns
List of the files in the right order and formating to be given to
exo_k.ktable.Ktable.hires_to_ktable()
orexo_k.xtable.Xtable.hires_to_xtable()
.- Return type
list of str
Examples
>>> exo_k.create_fname_grid_Kspectrum_LMDZ(2,3) array([['k001', 'k002', 'k003'], ['k004', 'k005', 'k006']], dtype='<U4')
>>> exo_k.create_fname_grid_Kspectrum_LMDZ(2,3,suffix='.h5') array([['k001.h5', 'k002.h5', 'k003.h5'], ['k004.h5', 'k005.h5', 'k006.h5']], dtype='<U7')
-
convert_exo_transmit_to_hdf5
(file_in, file_out, mol='unspecified')[source]¶ Converts exo_transmit like spectra to hdf5 format for speed and space.
- Parameters
file_in (str) – Initial exo_transmit filename.
file_out (str) – Name of the final hdf5 file to be created.
-
convert_old_kspectrum_to_hdf5
(file_in, file_out, skiprows=0)[source]¶ Deprecated.
Converts kspectrum like spectra to hdf5 format for speed and space
- Parameters
file_in (str) – Initial kspectrum filename.
file_out (str) – Name of the final hdf5 file to be created.
-
class
Molar_mass
[source]¶ Bases:
exo_k.util.singleton.Singleton
A class to compute the molar mass (in kg/mol) of regular molecules with a name written in a regular way (e.g. CO2, H2O, etc.). This class can also store the molar mass of custom gases with arbitrary names (for example: My_gas, earth_background).
-
class
Ktable
(*filename_filters, filename=None, xtable=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.ktable_io.Ktable_io
A class that handles 4D tables of k-coefficients.
Based on the
Data_table
class that handles basic operations common to Xtable.Based on the
Ktable_io
class that incorporates all io routines (read/write_xxx for all the xxx supported formats).Initializes k coeff table 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.
filename_filters (sequence of string) – As many strings as necessary to uniquely define a file in the global search path defined in
Settings
. This path will be searched for a file with all the filename_filters in the name. The filename_filters can contain ‘*’.xtable (Xtable object) – If no filename nor filename_filters are provided, this xtable object will be used to create a ktable. In this case, wavenumber bins must be given with the wnedges keyword.
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 none of the parameters above is specified, just creates an empty object to be filled later.
- Other Parameters
p_unit (str, optional) – 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)
kdata_unit (str, optional) – String to identify the unit 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, optional) – 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)
remove_zeros (boolean, optional) – If True, the zeros in the kdata table are replaced by a value 10 orders of magnitude smaller than the smallest positive value
mol (str, optional) – The name of the gas or molecule described by the
Ktable
search_path (str, optional) – If search_path is provided, it locally overrides the global _search_path in
Settings
and only files in search_path are returned.
-
property
shape
(self)¶ Returns the shape of self.kdata
-
xtable_to_ktable
(self, xtable=None, wnedges=None, weights=None, ggrid=None, quad='legendre', order=20, g_split=0.9, mid_dw=True, write=0, remove_zeros=False)[source]¶ Fills the
Ktable
object with a k-coeff table computed from aXtable
object (inplace).The p and kcorr units are inherited from the
Xtable
object.- Parameters
xtable (
Xtable
) – input Xtable object instancewnedges (Array) – edges of the wavenumber bins to be used to compute the corrk
- Other Parameters
weights (array, optional) – If weights are provided, they are used instead of the legendre quadrature.
quad (string, optional) – Type of quadrature used. Default is ‘legendre’. Also available: ‘split-legendre’ which uses half quadrature points between 0. and g_split and half between g_split and 1.
order (Integer, optional) – Order of the Gauss legendre quadrature used. Default is 20.
g_split (float, optional) – Used only if quad=’split-legendre’. See above.
mid_dw (boolean, optional) –
If True, the Xsec values in the high resolution xtable data are assumed to cover a spectral interval that is centered around the corresponding wavenumber value. The first and last Xsec values are discarded.
If False, each interval runs from the wavenumber value to the next one. The last Xsec value is dicarded.
-
hires_to_ktable
(self, path=None, filename_grid=None, logpgrid=None, tgrid=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 (inplace).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.
- Parameters
path (String) – directory with the input files
filename_grid (array of str with shape (logpgrid.size,tgrid.size)) – Names of the input high-res spectra. If None, the files are assumed to follow Kspectrum/LMDZ convention, i.e. be of the type ‘k001’, ‘k002’, etc. See
create_fname_grid_Kspectrum_LMDZ()
for possible additional keyword arguments.logpgrid (array) – Grid in log(pressure) of the input. Default unit is Pa, but can be changed with the grid_p_unit keyword.
tgrid (array) – Grid in temperature of the input.
wnedges (array) – Edges of the wavenumber bins to be used to compute the corr-k
- Other Parameters
weights (array, optional) – If weights are provided, they are used instead of the legendre quadrature.
quad (string, optional) – Type of quadrature used. Default is ‘legendre’. Also available: ‘split-legendre’ which uses half quadrature points between 0. and g_split and half between g_split and 1.
order (Integer, optional) – Order of the Gauss legendre quadrature used. Default is 20.
g_split (float, optional) – Used only if quad=’split-legendre’. See above.
mid_dw (boolean, optional) –
If True, the Xsec values in the high resolution xsec data are assumed to cover a spectral interval that is centered around the corresponding wavenumber value. The first and last Xsec values are discarded.
If False, each interval runs from the wavenumber value to the next one. The last Xsec value is dicarded.
mol (string, optional) – Give a name to the molecule. Useful when used later in a Kdatabase to track molecules.
p_unit (str, optional) – Pressure unit to convert to.
grid_p_unit (str, optional) – Unit of the specified logpgrid.
kdata_unit (str, optional) – Kdata unit to convert to.
file_kdata_unit (str, optional) – Kdata unit in input files.
See
create_fname_grid_Kspectrum_LMDZ()
orHires_spectrum
for a list of additional arguments that can be provided to those funtion through **kwargs.
-
copy
(self, cp_kdata=True, ktab5d=False)[source]¶ Creates a new instance of
Ktable
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.
ktab5d (bool, optional) – If true, creates a Ktable5d object with the same structure. Data are not copied.
- Returns
A new
Ktable
orKtable5d
instance with the same structure as self.- Return type
-
spectrum_to_plot
(self, 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
(self, 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)
-
RandOverlap
(self, other, x_self, x_other, write=0, use_rebin=False)[source]¶ Method to randomly mix the opacities of 2 species (self and other).
- Parameters
If one of these is None, the kcoeffs of the species in question are considered to be already normalized with respect to the mixing ratio.
- Returns
array of k-coefficients for the mix.
- Return type
array
-
bin_down
(self, 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) – 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, optional) – Desired weights for the resulting Ktable.
ggrid (array, 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
-
class
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.-
property
shape
(self)¶ Returns the shape of self.kdata
-
read_hdf5
(self, filename=None, mol=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
(self, 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
(self, 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
(self, 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
(self, 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.- Other Parameters
xgrid (array) – 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
(self, log_interp=None)[source]¶ Creates interpolating functions to be called later on. and loads it as attribute (inplace).
-
set_kdata
(self, 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) – New array of kdata.
-
interpolate_kdata
(self, logp_array=None, t_array=None, x_array=None, log_interp=None, wngrid_limit=None)[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
floats are given, they are interpreted as arrays of size 1. (If) –
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
The interpolated kdata.
- Return type
array of shape (logp_array.size, self.Nw , self.Ng)
-
remap_logPT
(self, 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).
- 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
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
(self, cp_kdata=True)[source]¶ Creates a new instance of
Ktable5d
object and (deep) copies data into it
-
spectrum_to_plot
(self, 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
(self, 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
(self, 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 aKtable
).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
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
A new Ktable5d with the opacity of the new species added.
- Return type
-
bin_down
(self, 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) – 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, optional) – Desired weights for the resulting Ktable.
ggrid (array, 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
(self, 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
(self, **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.
-
class
Xtable
(*filename_filters, filename=None, p_unit='unspecified', file_p_unit='unspecified', kdata_unit='unspecified', file_kdata_unit='unspecified', remove_zeros=False, search_path=None, mol=None)[source]¶ Bases:
exo_k.data_table.Data_table
A class that handles tables of cross sections.
Initializes cross section table and supporting data from a file based on its extension.
- Parameters
filename (str) – Relative or absolute path to the input file.
filename_filters (sequence of string) – As many strings as necessary to uniquely define a file in the global search path defined in
Settings
. This path will be searched for a file with all the filename_filters in the name. The filename_filters can contain ‘*’.
If there is no filename or filename_filters provided, 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 keywords.-
read_hdf5
(self, filename=None, mol=None)[source]¶ Initializes k coeff table and supporting data from an Exomol hdf5 file
- Parameters
filename (str) – Name of the input hdf5 file
mol (str, optional) – Overrides the name of the molecule to be put in the
Xtable
object.
-
write_hdf5
(self, filename, compression='gzip', compression_level=9, kdata_unit=None, p_unit=None, exomol_units=False)[source]¶ Saves data in a hdf5 format
- Parameters
filename (str) – Name of the file to be created and saved
exomol_units (bool (optional)) – If True, data are converted back to cm^2 and bar units before being written.
-
read_exo_transmit
(self, filename, mol=None)[source]¶ Creates an xsec object from an exo_transmit like spectra. See https://github.com/elizakempton/Exo_Transmit or Kempton et al. (2016) for details. Pressures are expected to be in Pa and cross sections in m^2/molecule
- Parameters
filename (str) – Name of the input file.
mol (str) – Overrides the name of the molecule to be put in the
Xtable
object.
-
hires_to_xtable
(self, path=None, filename_grid=None, logpgrid=None, tgrid=None, write=0, mol=None, grid_p_unit='Pa', p_unit='unspecified', kdata_unit='unspecified', file_kdata_unit='unspecified', **kwargs)[source]¶ Loads an Xtable from high-resolution spectra (inplace).
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.
see
exo_k.ktable.Ktable.hires_to_ktable()
method for details on the arguments and options.
-
bin_down
(self, wnedges=None, remove_zeros=False, write=0)[source]¶ Method to bin down a xsec table to a new grid of wavenumbers (inplace).
- Parameters
wnedges (array) – Edges of the new bins of wavenumbers (cm-1) onto which the xsec should be binned down. if you want Nwnew bin in the end, wngrid.size must be Nwnew+1 wnedges[0] should be greater than self.wnedges[0] wnedges[-1] should be lower than self.wnedges[-1]
remove_zeros (bool, optional) – If True, remove zeros in kdata.
-
sample
(self, wngrid, remove_zeros=False, log_interp=None)[source]¶ Method to re sample a xsec table to a new grid of wavenumbers (inplace).
- Parameters
wngrid (array) – Location of the new wavenumbers points (cm-1)
-
spectrum_to_plot
(self, 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)
x (float) – Volume mixing ratio of the species
g (is unused but here to be consistent with the method in data_table) –
-
copy
(self, cp_kdata=True)[source]¶ Creates a new instance of
Xtable
object and (deep) copies data into it
-
property
shape
(self)¶ Returns the shape of self.kdata
-
class
Hires_spectrum
(filename, file_kdata_unit='unspecified', kdata_unit='unspecified', mult_factor=None, binary=False, **kwargs)[source]¶ Bases:
exo_k.util.spectral_object.Spectral_object
A class defining a Hires_spectrum object.
Reads a high-resolution spectrum from a file (either hdf5 or ascii).
- Parameters
filename (str) – Full pathname to the file. Extension defines the format used.
file_kdata_unit (str) – Specifies the unit for the opacity data in the file. This is needed for ascii formats as the units are not known. The type of quantity may differ whether we are handling cross sections (surface) or absorption coefficients (inverse length)
kdata_unit (str) – Unit to convert to.
mult_factor (float) – A multiplicative factor that can be applied to kdata (for example to correct for any dilution effect, or specific conversion).
see
read_ascii()
for additional arguments to use with ascii files-
read_ascii
(self, filename, data_type=None, skiprows=0, wn_column=None, kdata_column=None)[source]¶ Read native kspectrum format
- Parameters
filename (str) – Initial hires-spectrum filename.
data_type ('xsec' or 'abs_coeff') – Whether the data read are cross-sections or absorption coefficients.
skiprows (int, optional) – Number of header lines to skip. For the latest Kspectrum format, the header is skipped automatically.
wn_column/kdata_column (int, optional) – Number of column to be read for wavenumber and kdata in python convention (0 is first, 1 is second, etc.)
-
read_binary
(self, filename, mass_amu=None)[source]¶ Reads spectra file in binary format (petitRADTRANS style)
Assumed to be in cm^2/g with wavelength in cm.
Will be automatically converted to cm^2/molecule and wns in cm^-1 (unless conversion to mks is requested).
-
convert_kdata_unit
(self, 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)
-
plot_spectrum
(self, ax, x_axis='wls', xscale=None, yscale=None, x=1.0, **kwarg)[source]¶ Plot the spectrum
- Parameters
ax (
pyplot.Axes
) – A pyplot axes instance where to put the plot.x_axis (str, optional) – If ‘wls’, x axis is wavelength. Wavenumber otherwise.
x/yscale (str, optional) – If ‘log’ log axes are used.
-
convert_data_type
(self, pressure, temperature, kdata_unit=None, convert_to=None)[source]¶ Converts from one data_type (cross sections or absorption coefficents) to the other (inplace).
Conversion to mks is done by default if a conversion takes place and no kdata_unit is specified.
- Parameters
pressure (float) – Pressure used for the conversion (in Pa)
temperature (float) – Temperature used for the conversion (in K)
kdata_unit (str (optional)) – Unit to use for the output
convert_to (str ('xsec' or 'abs_coeff', optional)) – Data type to convert to. Nothing is done if convert_to is equal to self.data_type. If None, converts to the other type.
-
PI
= 3.14159265359¶
-
create_table
(logpgrid=None, tgrid=None, xgrid=None, wngrid=None, wnedges=None, ggrid=None, weights=None, p_unit='Pa', kdata_unit='m^2/molecule', value=None, mol=None)[source]¶ Create a new table with the required dimensions padded with zeros.
The class of the object created will depend on the dimensions required.
- Parameters
logpgrid (array) – Grid in log p
tgrid (array) – Grid in T
xgrid (array) – Grid in volume mixing ratio
wngrid (array) – Grid wavenumber
ggrid (array) – Grid g-space
weights (array) – Grid of quadrature weights
p_unit (str (optional)) – Unit of pressure (default Pa)
kdata_unit (str (optional)) – Unit of cross section (default is m^2/molecule)
- Returns
new_table – A new table
- Return type
Xtable, Ktable, or Ktable5d
-
hires_to_ktable
(filename_grid=None, xgrid=None, **kwargs)[source]¶ Emulates
exo_k.ktable.Ktable.hires_to_ktable()
andexo_k.ktable5d.Ktable5d.hires_to_ktable()
as functions and not methods so that the user can call them without creating a Ktable first.See those methods for details on the available arguments and options.
-
hires_to_xtable
(**kwargs)[source]¶ Emulates
exo_k.xtable.Xtable.hires_to_xtable()
as a function and not a method so that the user can call it without creating an Xtable first.See
hires_to_xtable()
for details on the available arguments and options.
-
convert_kspectrum_to_hdf5
(file_in, file_out=None, **kwargs)[source]¶ Converts kspectrum like spectra to hdf5 format for speed and space. Helper function. Real work done in
Hires_spectrum
__init__ funtion.Go there to see all the available arguments and options.
- Parameters
file_in (str) – Initial kspectrum filename.
file_out (str) – Name of the final hdf5 file to be created. If not provided, ‘file_in.h5’ will be used.
-
create_fname_grid
(base_string, logpgrid=None, tgrid=None, xgrid=None, p_kw=None, t_kw=None, x_kw=None)[source]¶ Creates a grid of filenames from an array of pressures, temperatures ( and vmr if there is a variable gas).
- Parameters
base_string (str) – Generic name of the spectra files with specific keywords to be replaced by the relevant numerical values
logpgrid (Array) – Grid in log(pressure/Pa) of the input
tgrid (Array) – Grid in temperature of the input
xgrid (Array) – Input grid in vmr of the variable gas
Warning
The result of this function is much more predictable if the values in the above arrays are given as integers or directly strings. If you want to use floats anyway, good luck.
- Parameters
p_kw (str) –
t_kw (str) –
x_kw (str) – The pattern string that will be recognized as keywords between {} in base_string (See examples).
Examples
>>> logpgrid=[1,2] >>> tgrid=np.array([100.,150.,200.]) >>> file_grid=exo_k.create_fname_grid('spectrum_CO2_1e{logp}Pa_{t}K.hdf5', logpgrid=logpgrid,tgrid=tgrid,p_kw='logp',t_kw='t') array([['spectrum_CO2_1e1Pa_100K.hdf5', 'spectrum_CO2_1e1Pa_150K.hdf5', 'spectrum_CO2_1e1Pa_200K.hdf5'], ['spectrum_CO2_1e2Pa_100K.hdf5', 'spectrum_CO2_1e2Pa_150K.hdf5', 'spectrum_CO2_1e2Pa_200K.hdf5']], dtype='<U28')
-
finalize_LMDZ_dir
(corrkname, IRsize, VIsize)[source]¶ Creates the right links for a LMDZ type directory to be read by the LMDZ generic GCM.
You will need to create a proper Q.dat before using with the LMDZ GCM.
Important
You must have run
exo_k.ktable.Ktable.write_LMDZ()
orexo_k.ktable5d.Ktable5d.write_LMDZ()
for both of your IR and VI channels beforehand.- Parameters
corrkname (str) – Path to the directory with the LMDZ ktable to finalize
IRsize (int) – Number of IR spectral bins
VIsize (int) – Number of VI spectral bins
-
mmr_to_number_density
(mmr, gas_density, r_eff, condensate_density)[source]¶ Converts a mass mixing ratio (mmr or q) in a number density of particles (in number per unit volume)
- Parameters
mmr (float or array) – Mass mixing ratio (in kg per kg of air)
gas_density (float or array) – Density of the gas (in kg/m^3)
r_eff (float or array) – Effective radius of the particles
condensate_density (float or array) – Density of the constituent of the condensed particles (in kg/m^3)