exo_k.atm2

@author: jeremy leconte

Contains classes for atmospheric profiles and their radiative properties. That’s where forward models are computed.

Module Contents

class Atm_profile2(composition={}, psurf=None, ptop=None, logplev=None, tlev=None, Tsurf=None, Tstrat=None, grav=None, Rp=None, Mgas=None, rcp=0.28, Nlev=20, **kwargs)[source]

Bases: object

A class defining an atmospheric PT profile with some global data (gravity, etc.)

The derived class Atm handles radiative transfer calculations.

Initializes atmospheric profiles

Parameters
  • composition (dict) – Keys are molecule names and values the vmr. Vmr can be arrays of size Nlev-1 (i.e. the number of layers).

  • grav (float) – Planet surface gravity (gravity constant with altitude for now).

  • Rp (float or Astropy.unit quantity) – Planet radius. If float, Jupiter radii are assumed.

  • rcp (float) – Adiabatic lapse rate for the gas (R/cp)

  • Mgas (float, optional) – Molar mass of the gas (kg/mol). If given, overrides the molar mass computed from composition.

There are two ways to define the profile. Either define:

  • Nlev: int Number of level interfaces (Number of layers is Nlev-1)

  • psurf, Tsurf: float Surface pressure (Pa) and temperature

  • ptop: float Pressure at the top of the model (Pa)

  • Tstrat: float Stratospheric temperature

or:

  • logplev: array

  • tlev: array (same size) These will become the pressures (Pa) and temperatures at the level interfaces. This will be used to define the surface and top pressures. Nlev becomes the size of the arrays.

Warning

Layers are counted from the top down (increasing pressure order). All methods follow the same convention.

set_logPT_profile(self, log_plev, tlev)[source]

Set the logP-T profile of the atmosphere with a new one

Parameters
  • log_plev (numpy array) – Log pressure (in Pa) at the level surfaces

  • tlev (numpy array (same size)) – temperature at the level surfaces.

set_T_profile(self, tlev)[source]

Reset the temperature profile without changing the pressure levels

compute_pressure_levels(self)[source]

Computes various pressure related quantities

set_adiab_profile(self, Tsurf=None, Tstrat=None)[source]

Initializes the logP-T atmospheric profile with an adiabat with index R/cp=rcp

Parameters
  • Tsurf (float) – Surface temperature.

  • Tstrat (float, optional) – Temperature of the stratosphere. If None is given, an isothermal atmosphere with T=Tsurf is returned.

set_grav(self, grav=None)[source]

Sets the surface gravity of the planet

Parameters

grav (float) – surface gravity (m/s^2)

set_gas(self, composition_dict)[source]

Sets the composition of the atmosphere

Parameters
  • composition_dict (dictionary) – Keys are molecule names, and values are volume mixing ratios. A ‘background’ value means that the gas will be used to fill up to vmr=1 If they do not add up to 1 and there is no background gas_mix, the rest of the gas_mix is considered transparent.

  • compute_col_dens (boolean, optional) – If True, the column density per layer of the atmosphere is recomputed. This si mostly to save time when we know this will be done later on.

set_Mgas(self, Mgas=None)[source]

Sets the mean molar mass of the atmosphere.

Parameters

Mgas (float or array of size Nlay) – Mean molar mass (kg/mol). If None is give, the mmm is computed from the composition.

set_rcp(self, rcp)[source]

Sets the adiabatic index of the atmosphere

Parameters

rcp (float) – R/c_p

set_Rp(self, Rp)[source]

Sets the radius of the planet

Parameters

Rp (float) – radius of the planet (m)

set_Rstar(self, Rstar)[source]

Sets the radius of the star

Parameters

Rstar (float) – radius of the star (m)

compute_density(self)[source]

Computes the number density (m^-3) profile of the atmosphere

compute_layer_col_density(self)[source]

Computes the column number density (molecules/m^-2) per layer of the atmosphere

compute_altitudes(self)[source]

Compute altitudes of the level surfaces (zlev) and mid layers (zlay).

compute_area(self)[source]

Computes the area of the annulus covered by each layer in a transit setup.

compute_tangent_path(self)[source]

Computes a triangular array of the tangent path length (in m) spent in each layer. self.tangent_path[ilay][jlay] is the length that the ray that is tangent to the ilay layer spends in the jlay>=ilay layer (accounting for a factor of 2 due to symmetry)

class Atm2(k_database=None, cia_database=None, wn_range=None, wl_range=None, **kwargs)[source]

Bases: exo_k.atm2.Atm_profile2

Class based on Atm_profile that handles radiative trasnfer calculations.

Radiative data are accessed through the gas_mix.Gas_mix class.

Initialization method that calls Atm_Profile().__init__() and links to Kdatabase and other radiative data.

set_k_database(self, k_database=None)[source]

Change the radiative database used by the Gas_mix object handling opacities inside Atm.

See gas_mix.Gas_mix.set_k_database for details.

Parameters

k_database (Kdatabase object) – New Kdatabase to use.

set_cia_database(self, cia_database=None)[source]

Change the CIA database used by the Gas_mix object handling opacities inside Atm.

See gas_mix.Gas_mix.set_cia_database for details.

Parameters

cia_database (CIAdatabase object) – New CIAdatabase to use.

set_spectral_range(self, wn_range=None, wl_range=None)[source]

Sets the spectral range in which computations will be done by specifying either the wavenumber (in cm^-1) or the wavelength (in micron) range.

See gas_mix.Gas_mix.set_spectral_range for details.

spectral_integration(self, spectral_var)[source]

Spectrally integrate an array, taking care of whether we are dealing with corr-k or xsec data.

Parameters

spectral_var (array) – array to integrate

Returns

var – array integrated over wavenumber (and g-space if relevant)

Return type

array

opacity(self, rayleigh=False, **kwargs)[source]

Computes the opacity of each of the layers.

See gas_mix.Gas_mix.cross_section for details.

source_function(self, integral=True, source=True)[source]

Compute the blackbody source function (Pi*Bnu) for each layer of the atmosphere.

Parameters
  • integral (boolean, optional) –

    • If true, the black body is integrated within each wavenumber bin.

    • If not, only the central value is used. False is faster and should be ok for small bins, but True is the correct version.

  • source (boolean, optional) – If False, the source function is put to 0 (for solar absorption calculations)

setup_emission_caculation(self, mu_eff=0.5, rayleigh=False, integral=True, source=True, **kwargs)[source]

Computes all necessary quantities for emission calculations (opacity, source, etc.)

emission_spectrum(self, integral=True, mu0=0.5, mu_quad_order=None, **kwargs)[source]

Returns the emission flux at the top of the atmosphere (in W/m^2/cm^-1)

Parameters

integral (boolean, optional) –

  • If true, the black body is integrated within each wavenumber bin.

  • If not, only the central value is used. False is faster and should be ok for small bins, but True is the correct version.

Other Parameters
  • mu0 (float) – Cosine of the quadrature angle use to compute output flux

  • mu_quad_order (int) – If an integer is given, the emission intensity is computed for a number of angles and integrated following a gauss legendre quadrature rule of order mu_quad_order.

Returns

A spectrum with the Spectral flux at the top of the atmosphere (in W/m^2/cm^-1)

Return type

Spectrum object

emission_spectrum_quad(self, integral=True, mu_quad_order=3, **kwargs)[source]

Returns the emission flux at the top of the atmosphere (in W/m^2/cm^-1) using gauss legendre qudrature of order mu_quad_order

Parameters

integral (boolean, optional) –

  • If true, the black body is integrated within each wavenumber bin.

  • If not, only the central value is used. False is faster and should be ok for small bins, but True is the correct version.

Returns

A spectrum with the Spectral flux at the top of the atmosphere (in W/m^2/cm^-1)

Return type

Spectrum object

emission_spectrum_2stream(self, integral=True, mu0=0.5, method='toon', dtau_min=1e-10, mid_layer=False, rayleigh=False, flux_top_dw=None, source=True, compute_kernel=False, **kwargs)[source]

Returns the emission flux at the top of the atmosphere (in W/m^2/cm^-1)

Parameters

integral (boolean, optional) –

  • If true, the black body is integrated within each wavenumber bin.

  • If not, only the central value is used. False is faster and should be ok for small bins, but True is the correct version.

Other Parameters
  • mu0 (float) – Cosine of the quadrature angle use to compute output flux

  • dtau_min (float) – If the optical depth in a layer is smaller than dtau_min, dtau_min is used in that layer instead. Important as too transparent layers can cause important numerical rounding errors.

Returns

A spectrum with the Spectral flux at the top of the atmosphere (in W/m^2/cm^-1)

Return type

Spectrum object

compute_kernel(self, solve_2stream_nu, epsilon=0.01, mid_layer=False, mu0=0.5, per_unit_mass=True, integral=True, **kwargs)[source]

Compute the Jacobian matrix d Heating[lev=i] / d T[lev=j]

flux_divergence(self, internal_flux=0.0, per_unit_mass=True)[source]

Computes the divergence of the net flux in the layers (used to compute heating rates).

emission_spectrum_2stream() needs to be ran first.

Parameters
  • internal_flux (float) – flux coming from below in W/m^2

  • per_unit_mass (bool) – If True, the heating rates are normalized by the mass of each layer (result in W/kg).

Returns

  • net (array) – Net fluxes at level surfaces

  • H (array) – Heating rate in each layer (Difference of the net fluxes). Positive means heating. The last value is the net flux impinging on the surface + the internal flux.

bolometric_fluxes(self, internal_flux=0.0, per_unit_mass=True)[source]

Computes the bolometric fluxes at levels and the divergence of the net flux in the layers (used to compute heating rates).

emission_spectrum_2stream() needs to be ran first.

Parameters
  • internal_flux (float) – flux coming from below in W/m^2

  • per_unit_mass (bool) – If True, the heating rates are normalized by the mass of each layer (result in W/kg).

Returns

  • up (array) – Upward fluxes at level surfaces

  • dw (array) – Downward fluxes at level surfaces

  • net (array) – Net fluxes at level surfaces

  • H (array) – Heating rate in each layer (Difference of the net fluxes). Positive means heating. The last value is the net flux impinging on the surface + the internal flux.

exp_minus_tau(self)[source]

Sums Exp(-tau) over gauss points

exp_minus_tau_g(self, g_index)[source]

Sums Exp(-tau) over gauss points

surf_bb(self, integral=True)[source]

Computes the surface black body flux (in W/m^2/cm^-1)

Parameters

integral (boolean, optional) –

  • If true, the black body is integrated within each wavenumber bin.

  • If not, only the central value is used. False is faster and should be ok for small bins, but True is the correct version.

Returns

Spectral flux at the surface (in W/m^2/cm^-1)

Return type

Spectrum object

top_bb(self, integral=True)[source]

Computes the top of atmosphere black body flux (in W/m^2/cm^-1)

Parameters

integral (boolean, optional) –

If true, the black body is integrated within each wavenumber bin. If not, only the central value is used.

False is faster and should be ok for small bins, but True is the correct version.

Returns

Spectral flux of a bb at the temperature at the top of atmosphere (in W/m^2/cm^-1)

Return type

Spectrum object

transmittance_profile(self, **kwargs)[source]

Computes the transmittance profile of an atmosphere, i.e. Exp(-tau) for each layer of the model. Real work done in the numbafied function path_integral_corrk/xsec depending on the type of data.

transmission_spectrum(self, normalized=False, Rstar=None, **kwargs)[source]

Computes the transmission spectrum of the atmosphere. In general (see options below), the code returns the transit depth:

\[\delta_\nu=(\pi R_p^2+\alpha_\nu)/(\pi R_{star}^2),\]

where

\[\alpha_\nu=2 \pi \int_0^{z_{max}} (R_p+z)*(1-e^{-\tau_\nu(z)) d z.\]
Parameters
  • Rstar (float, optional) –

    Radius of the host star. Does not need to be given here if as already been specified as an attribute of the self.Atm object. If specified, the result is the transit depth:

    \[\delta_\nu=(\pi R_p^2+\alpha_\nu)/(\pi R_{star}^2).\]

  • normalized (boolean, optional) –

    Used only if self.Rstar and Rstar are None:

    • If True, the result is normalized to the planetary radius:

      \[\delta_\nu=1+\frac{\alpha_\nu}{\pi R_p^2}.\]
    • If False,

      \[\delta_\nu=\pi R_p^2+\alpha_\nu.\]

Returns

The transit spectrum (see above for normalization options).

Return type

array

heating_rate(self, Fin=1.0, Tstar=5570.0, szangle=60.0, **kwargs)[source]

Computes the heating rate in the atmosphere

Parameters
  • Fin (float) – Bolometric stellar flux at the top of atmopshere (W/m^2).

  • Tstar (float) – Stellar temperature

  • szangle (float) – Solar zenith angle

Returns

Heating rate in each atmospheric layer (K/s).

Return type

array

emission_spectrum_exp_integral(self, integral=True, **kwargs)[source]

Computes the emission flux at the top of the atmosphere (in W/m^2/cm^-1)

Warning

This method uses a formulation with exponential integrals that is not yet perceived as accurate enough by the author. It is left for testing only.

Parameters

integral (boolean, optional) –

  • If true, the black body is integrated within each wavenumber bin.

  • If not, only the central value is used. False is faster and should be ok for small bins, but True is the correct version.

Returns

A spectrum with the Spectral flux at the top of the atmosphere (in W/m^2/cm^-1)

Return type

Spectrum object

convadj(timestep, Nlev, tlev, exner, dp_lay, verbose=False)[source]

Computes the heating rates needed to adjust unstable regions of a given atmosphere to a convectively neutral T profile on a given timestep.

Important

The layers here are not the same as the layers in the radiative transfer!

In the radiative transfer, a layer is the volume between two pressure level boundaries (plev).

Here, a layer is centered on a pressure level (plev) with its temperature tlev. Therefore, dp_lay[0]=play[0]-plev[0], dp_lay[i]=play[i]-play[i-1], and dp_lay[-1]=psurf-play[-1].

Parameters
  • timestep (float) – Duration of the adjustment in seconds.

  • Nlev (int) – Number of atmospheric layers

  • tlev (array) – Temperatures of the atmospheric layers

  • exner (array) –

    Exner function computed at the layer centers ((p/psurf)**rcp)

    \[\Pi=(p / p_{s})^{R/c_p}\]

  • dp_lay (array) – Pressure difference between the bottom and top of each layer

Returns

Heating rate in each atmospheric layer (K/s).

Return type

array