exo_k.atable

@author: jeremy leconte

A class to handle continuum absorption (CIA)

Module Contents

class exo_k.atable.Atable(*filename_filters, filename=None, aerosol_name=None, search_path=None, mks=False, remove_zeros=False, N_per_line=5, wn_range=None, wl_range=None)[source]

Bases: exo_k.util.spectral_object.Spectral_object

A class to handle aerosol optical properties in table form.

Initialization for Atables.

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 ‘*’.

  • 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.

_init_empty()[source]

Initializes attributes to none.

read_LMDZ(filename, aerosol_name=None, N_per_line=5)[source]

Reads LMDZ like optical properties files.

Parameters:

filename (str) – Name of the file to be read.

_read_arrays(file, Nvalue, Narray, N_per_line=5)[source]

Reads an array in an optical property LMDZ file. Assumes that the arrays are arranged 5 values per line.

Parameters:
  • file (file stream) – File to be read.

  • Nvalue (int) – Number of values to be read in each array.

  • Narray (int) – Number of arrays to be read.

  • N_per_line (int) – Number of values per lines

Returns:

Array

A numpy array with the values.

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

Reads hdf5 cia files and load temperature, wavenumber, and absorption coefficient grid.

Parameters:

filename (str) – Name of the file to be read.

write_hdf5(filename)[source]

Writes hdf5 cia files.

Parameters:

filename (str) – Name of the file to be written.

sample(wngrid, remove_zeros=False, use_grid_filter=False, sample_all_vars=True, **kwargs)[source]

Method to re sample an Atable to a new grid of wavenumbers (in place)

Parameters:
  • wngrid (array, np.ndarray) – new wavenumber grid (cm-1)

  • use_grid_filter (boolean, optional) – If true, the table is sampled only within the boundaries of its current wavenumber grid. The coefficients are set to zero elswere (except if remove_zeros is set to True). If false, the values at the boundaries are used when sampling outside the grid.

  • sample_all_vars (boolean, optional) – Whether to sample the single_scattering albedo and asymmetry_factor as well.

sample_cp(wngrid, **kwargs)[source]

Creates a copy of the object before resampling it.

Parameters:

details. (See sample method for) –

Returns:

Atable object

the re-sampled Atable

interpolate_optical_properties(r_array=None, var_type=0, log_interp=None, wngrid_limit=None)[source]

interpolate_cia interpolates the kdata at on a given temperature profile.

Parameters:
  • r_array (float or array) – Effective radius array to interpolate to. If a float is given, it is interpreted as an array of size 1.

  • var_type (int) –

    type of data to interpolate:
    • 0 is extinction coefficient

    • 1 is single scattering albedo

    • 2 is asymmetry factor

  • wngrid_limit (array, np.ndarray, optional) – If an array is given, interpolates only within this array.

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

cross_section(r_array, wngrid_limit=None, log_interp=None)[source]

Computes the cross section due to the aerosol in area per particles.

Parameters:
  • r_array (float or array) – Effective radius array to interpolate to. If a float is given, it is interpreted as an array of size 1.

  • wngrid_limit (array, np.ndarray, optional) – If an array is given, interpolates only within this array.

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

absorption_coefficient(r_array, n_density, wngrid_limit=None, log_interp=None)[source]

Computes the opacity due to the aerosols.

Warning

If n_density is the number density of aerosol particles (in m^-3), then the result is the absorption coefficient in m^-1.

If n_density is the ratio of the number density of aerosol particles normalized to the number density of gas molecules (dimless), then the result is a cross section (in m^2) for aerosols normalized “per gas molecule”. This latter choice allows us to directly add it to gaseous cross sections later-on without further normalizations. This is what needs to be used when coupling with the atmospheric model.

Parameters:
  • r_array (float or 1d array) – Effective radius array to interpolate to. If a float is given, it is interpreted as an array of size 1.

  • n_density (float or 1d array (same dim as r_array)) – Number density of aerosol or ratio of aerosol to gas particle number density (see above).

  • wngrid_limit (array, np.ndarray, optional) – If an array is given, interpolates only within this array.

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

optical_properties(r_array, n_density, wngrid_limit=None, log_interp=None, compute_all_opt_prop=True)[source]

Computes all the optical properties for the aerosols.

Warning

If n_density is the number density of aerosol particles (in m^-3), then kdata is the absorption coefficient in m^-1.

If n_density is the ratio of the number density of aerosol particles normalized to the number density of gas molecules (dimless), then kdata is a cross section (in m^2) for aerosols normalized “per gas molecule”. This latter choice allows us to directly add it to gaseous cross sections later-on without further normalizations. This is what needs to be used when coupling with the atmospheric model.

Parameters:
  • r_array (float or 1d array) – Effective radius array to interpolate to. If a float is given, it is interpreted as an array of size 1.

  • n_density (float or 1d array (same dim as r_array)) – Number density of aerosol or ratio of aerosol to gas particle number density (see above).

  • wngrid_limit (array, np.ndarray, optional) – If an array is given, interpolates only within this array.

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

plot_spectrum(ax, r=1e-06, x_axis='wls', xscale=None, yscale=None, var_type=0, **kwarg)[source]

Plot the spectrum for a given point

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

  • r (float) – Effective radius (m)

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

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

convert_to_mks()[source]

Converts units to MKS

rindex(r)[source]

Finds the index corresponding to the given radius r (units must be the same as the ktable)

remove_zeros(deltalog_min_value=0.0)[source]

Finds zeros in the ext_coeff and set them to (10.^-deltalog_min_value) times the minimum positive value in the table. This is to be able to work in logspace.

copy()[source]

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

classmethod concatenate_spectral_range(atab1, atab2, verbose=False)[source]

Combine the instance of atable with another one for which there is some overlap of the spectral range.

exo_k.atable.combine_tables(aerosol_name, *atables)[source]

Combine several tables representing the same aerosol type for different wavelength range.

Deprecated, use concatenate_spectral_range

exo_k.atable.interp_data(data_to_interp, rind, rweight, Nw, wngrid_filter, log_interp)[source]