POSEIDON.atmosphere

Functions for calculating atmospheric temperature, mixing ratio, and other profiles.

Module Contents

Functions

compute_T_Madhu(P, a1, a2, log_P1, log_P2, log_P3, ...)

Computes the temperature profile for an atmosphere using a re-arranged

compute_T_slope(P, T_phot, Delta_T_arr[, log_P_phot, ...])

Computes the temperature profile for an atmosphere using the 'slope' P-T

compute_T_field_gradient(P, T_bar_term, Delta_T_term, ...)

Creates the 3D temperature profile defined in MacDonald & Lewis (2022).

compute_T_field_two_gradients(P, T_bar_term_high, ...)

Extension of 'compute_T_field_gradient' with a second gradient in the

compute_X_field_gradient(P, log_X_state, N_sectors, ...)

Creates the 4D abundance profile array storing X(species, layer, sector, zone).

compute_X_field_two_gradients(P, log_X_state, ...[, ...])

Extension of 'compute_X_field_gradient' with a second gradient in the

add_bulk_component(P, X_param, N_species, N_sectors, ...)

Concatenates mixing ratios of the bulk species to the parametrised mixing

radial_profiles(P, T, g_0, R_p, P_ref, R_p_ref, mu, ...)

Solves the equation of hydrostatic equilibrium [ dP/dr = -G*M*rho/r^2 ]

radial_profiles_constant_g(P, T, g_0, P_ref, R_p_ref, ...)

Solves the equation of hydrostatic equilibrium [ dP/dr = -G*M*rho/r^2 ]

mixing_ratio_categories(P, X, N_sectors, N_zones, ...)

Sort mixing ratios into those of active species, collision-induced

compute_mean_mol_mass(P, X, N_species, N_sectors, ...)

Computes the mean molecular mass in each atmospheric column.

count_atoms(molecule)

Count how many atoms of each element are contained in a molecule.

elemental_ratio(included_species, X, element_1, element_2)

Calculate the abundance ratio between any two elements in the atmosphere.

profiles(P, R_p, g_0, PT_profile, X_profile, PT_state, ...)

Main function to calculate the vertical profiles in each atmospheric

POSEIDON.atmosphere.compute_T_Madhu(P, a1, a2, log_P1, log_P2, log_P3, T_set, P_set)

Computes the temperature profile for an atmosphere using a re-arranged form of the P-T profile parametrisation in Madhusudhan & Seager (2009).

Parameters:
  • P (np.array of float) – Atmosphere pressure array (bar).

  • a1 (float) – Alpha_1 parameter (encodes slope in layer 1).

  • a2 (float) – Alpha_2 parameter encodes slope in layer 2).

  • log_P1 (float) – Pressure of layer 1-2 boundary.

  • log_P2 (float) – Pressure of inversion.

  • log_P3 (float) – Pressure of layer 2-3 boundary.

  • T_set (float) – Atmosphere temperature reference value at P = P_set (K).

  • P_set (float) – Pressure whether the temperature parameter T_set is defined (bar).

Returns:

Temperature of each layer as a function of pressure (K). Only the first axis is used for this 1D profile.

Return type:

T (3D np.array of float)

POSEIDON.atmosphere.compute_T_slope(P, T_phot, Delta_T_arr, log_P_phot=0.5, log_P_arr=[-3.0, -2.0, -1.0, 0.0, 1.0, 1.5, 2.0])

Computes the temperature profile for an atmosphere using the ‘slope’ P-T profile parametrisation defined in Piette & Madhusudhan (2021).

Note: The number of temperature difference parameters is the same as the

number of pressure points (including the photosphere). For the default values, we have: [Delta_T (10-1mb), Delta_T (100-10mb), Delta_T (1-0.1b),

Delta_T (3.2-1b), Delta_T (10-3.2b), Delta_T (32-10b), Delta_T (100-32b)], where ‘b’ = bar.

Parameters:
  • P (np.array of float) – Atmosphere pressure array (bar).

  • T_phot (float) – Temperature at the photosphere (located at log_P_phot) (K).

  • Delta_T_arr (np.array of float) – Temperature differences between each pressure point (K).

  • log_P_phot (float) – Photosphere pressure (default: 3.16 bar).

  • log_P_arr (list) – Pressures where the temperature difference parameters are defined (log bar).

Returns:

Temperature of each layer as a function of pressure (K). Only the first axis is used for this 1D profile.

Return type:

T (3D np.array of float)

POSEIDON.atmosphere.compute_T_field_gradient(P, T_bar_term, Delta_T_term, Delta_T_DN, T_deep, N_sectors, N_zones, alpha, beta, phi, theta, P_deep=10.0, P_high=1e-05)

Creates the 3D temperature profile defined in MacDonald & Lewis (2022).

Note: This profile assumes a linear in log-pressure temperature gradient

between P_deep and P_high (these anchor points are fixed such that the photosphere generally lies within them). For pressures above and below the considered range, the temperature is isothermal.

Parameters:
  • P (np.array of float) – Atmosphere pressure array (bar).

  • T_bar_term (float) – Average terminator plane temperature at P_high.

  • Delta_T_term (float) – Temperature difference between the evening and morning terminators at P_high(K).

  • Delta_T_DN (float) – Temperature difference between the dayside and nightside at P_high (K).

  • T_deep (float) – Global deep temperature at P_deep.

  • N_sectors (int) – Number of azimuthal sectors.

  • N_zones (int) – Number of zenith zones.

  • alpha (float) – Terminator opening angle (degrees).

  • beta (float) – Day-night opening angle (degrees).

  • phi (np.array of float) – Mid-sector angles (radians).

  • theta (np.array of float) – Mid-zone angles (radians).

  • P_deep (float) – Anchor point in deep atmosphere below which the atmosphere is homogenous.

  • P_high (float) – Anchor point high in the atmosphere above which all columns are isothermal.

Returns:

Temperature of each layer as a function of pressure, sector, and zone (K).

Return type:

T (3D np.array of float)

POSEIDON.atmosphere.compute_T_field_two_gradients(P, T_bar_term_high, T_bar_term_mid, Delta_T_term_high, Delta_T_term_mid, Delta_T_DN_high, Delta_T_DN_mid, log_P_mid, T_deep, N_sectors, N_zones, alpha, beta, phi, theta, P_deep=10.0, P_high=1e-05)

Extension of ‘compute_T_field_gradient’ with a second gradient in the photosphere region. The two gradients connect at a new parameter: P_mid.

Parameters:
  • P (np.array of float) – Atmosphere pressure array (bar).

  • T_bar_term_high (float) – Average terminator plane temperature at P_high.

  • T_bar_term_mid (float) – Average terminator plane temperature at P_mid.

  • Delta_T_term_high (float) – Temperature difference between the evening and morning terminators at P_high (K).

  • Delta_T_term_mid (float) – Temperature difference between the evening and morning terminators at P_mid (K).

  • Delta_T_DN_high (float) – Temperature difference between the dayside and nightside at P_high (K).

  • Delta_T_DN_mid (float) – Temperature difference between the dayside and nightside at P_mid (K).

  • log_P_mid (float) – log10 of pressure where the two gradients switch (bar).

  • T_deep (float) – Global deep temperature at P_deep.

  • N_sectors (int) – Number of azimuthal sectors.

  • N_zones (int) – Number of zenith zones.

  • alpha (float) – Terminator opening angle (degrees).

  • beta (float) – Day-night opening angle (degrees).

  • phi (np.array of float) – Mid-sector angles (radians).

  • theta (np.array of float) – Mid-zone angles (radians).

  • P_deep (float) – Anchor point in deep atmosphere below which the atmosphere is homogenous.

  • P_high (float) – Anchor point high in the atmosphere above which all columns are isothermal.

Returns:

Temperature of each layer as a function of pressure, sector, and zone (K).

Return type:

T (3D np.array of float)

POSEIDON.atmosphere.compute_X_field_gradient(P, log_X_state, N_sectors, N_zones, param_species, species_has_profile, alpha, beta, phi, theta, P_deep=10.0, P_high=1e-05)

Creates the 4D abundance profile array storing X(species, layer, sector, zone).

The functional dependence is the same as defined above in the function ‘compute_T_field_gradient’ (see also MacDonald & Lewis (2022)), with the exception that not all chemical species have a vertical profile. Any chemical species not in the array ‘species_has_profile’ have constant mixing ratios with altitude.

Parameters:
  • P (np.array of float) – Atmosphere pressure array (bar).

  • log_X_state (2D np.array of float) – Mixing ratio state array.

  • N_sectors (int) – Number of azimuthal sectors.

  • N_zones (int) – Number of zenith zones.

  • param_species (np.array of str) – Chemical species with parametrised mixing ratios.

  • species_has_profile (np.array of int) – Array with an integer ‘1’ if a species in ‘param_species’ has a gradient profile, or ‘0’ for an constant mixing ratio with altitude.

  • alpha (float) – Terminator opening angle (degrees).

  • beta (float) – Day-night opening angle (degrees).

  • phi (np.array of float) – Mid-sector angles (radians).

  • theta (np.array of float) – Mid-zone angles (radians).

  • P_deep (float) – Anchor point in deep atmosphere below which all columns are isochemical.

  • P_high (float) – Anchor point high in the atmosphere above which all columns are isochemical.

Returns:

Mixing ratios in each layer as a function of pressure, sector, and zone.

Return type:

X_profiles (4D np.array of float)

POSEIDON.atmosphere.compute_X_field_two_gradients(P, log_X_state, N_sectors, N_zones, param_species, species_has_profile, alpha, beta, phi, theta, P_deep=10.0, P_high=1e-05)

Extension of ‘compute_X_field_gradient’ with a second gradient in the photosphere region. The two gradients connect at P_mid.

Parameters:
  • P (np.array of float) – Atmosphere pressure array (bar).

  • log_X_state (2D np.array of float) – Mixing ratio state array.

  • N_sectors (int) – Number of azimuthal sectors.

  • N_zones (int) – Number of zenith zones.

  • param_species (np.array of str) – Chemical species with parametrised mixing ratios.

  • species_has_profile (np.array of int) – Array with an integer ‘1’ if a species in ‘param_species’ has a gradient profile, or ‘0’ for an constant mixing ratio with altitude.

  • alpha (float) – Terminator opening angle (degrees).

  • beta (float) – Day-night opening angle (degrees).

  • phi (np.array of float) – Mid-sector angles (radians).

  • theta (np.array of float) – Mid-zone angles (radians).

  • P_deep (float) – Anchor point in deep atmosphere below which all columns are isochemical.

  • P_high (float) – Anchor point high in the atmosphere above which all columns are isochemical.

Returns:

Mixing ratios in each layer as a function of pressure, sector, and zone.

Return type:

X_profiles (4D np.array of float)

POSEIDON.atmosphere.add_bulk_component(P, X_param, N_species, N_sectors, N_zones, bulk_species, He_fraction)

Concatenates mixing ratios of the bulk species to the parametrised mixing ratios, forming the full mixing ratio array (i.e. sums to 1).

Note: For H2 and He as bulk species, the output array has the mixing ratio

profile of H2 as the first element and He second. For other bulk species (e.g. N2), that species occupies the first element.

Parameters:
  • P (np.array of float) – Atmosphere pressure array (bar).

  • X_param (4D np.array of float) – Mixing ratios of the parametrised chemical species in each layer as a function of pressure, sector, and zone.

  • N_sectors (int) – Number of azimuthal sectors.

  • N_zones (int) – Number of zenith zones.

  • bulk_species (np.array of str) – Bulk species dominating atmosphere (e.g. [‘H2’, ‘He’]).

  • He_fraction (float) – Assumed H2/He ratio (0.17 default corresponds to the solar ratio).

Returns:

Same as X_param, but with the bulk species appended.

Return type:

X (4D np.array of float)

POSEIDON.atmosphere.radial_profiles(P, T, g_0, R_p, P_ref, R_p_ref, mu, N_sectors, N_zones)

Solves the equation of hydrostatic equilibrium [ dP/dr = -G*M*rho/r^2 ] to compute the radius in each atmospheric layer.

Note: g is taken as an inverse square law with radius by assuming the

enclosed planet mass at a given radius is M_p. This assumes most mass is in the interior (negligible atmosphere mass).

Parameters:
  • P (np.array of float) – Atmosphere pressure array (bar).

  • T (3D np.array of float) – Temperature profile (K).

  • g_0 (float) – Gravitational field strength at white light radius (m/s^2).

  • R_p (float) – Observed white light planet radius (m).

  • P_ref (float) – Reference pressure (bar).

  • R_p_ref (float) – Planet radius corresponding to reference pressure (m).

  • mu (3D np.array of float) – Mean molecular mass (kg).

  • N_sectors (int) – Number of azimuthal sectors comprising the background atmosphere.

  • N_zones (int) – Number of zenith zones comprising the background atmosphere.

Returns:

Number density profile (m^-3). r (3D np.array of float):

Radial distance profile (m).

r_up (3D np.array of float):

Upper layer boundaries (m).

r_low (3D np.array of float):

Lower layer boundaries (m).

dr (3D np.array of float):

Layer thicknesses (m).

Return type:

n (3D np.array of float)

POSEIDON.atmosphere.radial_profiles_constant_g(P, T, g_0, P_ref, R_p_ref, mu, N_sectors, N_zones)

Solves the equation of hydrostatic equilibrium [ dP/dr = -G*M*rho/r^2 ] to compute the radius in each atmospheric layer.

Note: This version of the solver assumes the gravitational field strength

is constant with altitude (for testing purposes). The standard ‘radial_profiles’ function should be used for any real calculation.

Parameters:
  • P (np.array of float) – Atmosphere pressure array (bar).

  • T (3D np.array of float) – Temperature profile (K).

  • g_0 (float) – Gravitational field strength at white light radius (m/s^2).

  • P_ref (float) – Reference pressure (bar).

  • R_p_ref (float) – Planet radius corresponding to reference pressure (m).

  • mu (3D np.array of float) – Mean molecular mass (kg).

  • N_sectors (int) – Number of azimuthal sectors comprising the background atmosphere.

  • N_zones (int) – Number of zenith zones comprising the background atmosphere.

Returns:

Number density profile (m^-3). r (3D np.array of float):

Radial distant profile (m).

r_up (3D np.array of float):

Upper layer boundaries (m).

r_low (3D np.array of float):

Lower layer boundaries (m).

dr (3D np.array of float):

Layer thicknesses (m).

Return type:

n (3D np.array of float)

POSEIDON.atmosphere.mixing_ratio_categories(P, X, N_sectors, N_zones, included_species, active_species, CIA_pairs, ff_pairs, bf_species)

Sort mixing ratios into those of active species, collision-induced absorption (CIA), free-free opacity, and bound-free opacity.

Parameters:
  • P (np.array of float) – Atmosphere pressure array (bar).

  • X (4D np.array of float) – Mixing ratio profile

  • N_sectors (int) – Number of azimuthal sectors comprising the background atmosphere.

  • N_zones (int) – Number of zenith zones comprising the background atmosphere.

  • included_species (np.array of str) – List of chemical species included in the model (including bulk species).

  • active_species (np.array of str) – Spectroscopically active chemical species (see supported_opac.py).

  • CIA_pairs (np.array of str) – Collisionally-induced absorption (CIA) pairs.

  • ff_pairs (np.array of str) – Free-free absorption pairs.

  • bf_species (np.array of str) – Bound-free absorption species.

Returns:

Mixing ratios of active species. X_CIA (5D np.array of float):

Mixing ratios of CIA pairs

X_ff (5D np.array of float):

Mixing ratios of free-free pairs.

X_bf (4D np.array of float):

Mixing ratios of bound-free species.

Return type:

X_active (4D np.array of float)

POSEIDON.atmosphere.compute_mean_mol_mass(P, X, N_species, N_sectors, N_zones, masses_all)

Computes the mean molecular mass in each atmospheric column.

Parameters:
  • P (np.array of float) – Atmosphere pressure array (bar).

  • X (4D np.array of float) – Mixing ratio profile.

  • N_species (int) – Number of chemical species.

  • N_sectors (int) – Number of azimuthal sectors comprising the background atmosphere.

  • N_zones (int) – Number of zenith zones comprising the background atmosphere.

  • masses_all (np.array of float) – Masses of each chemical species (atomic mass units, u)

Returns:

Mean molecular mass (kg).

Return type:

mu (3D np.array of float)

POSEIDON.atmosphere.count_atoms(molecule)

Count how many atoms of each element are contained in a molecule.

Parameters:

molecule (str) – Name of the molecule (e.g. ‘H2O’, ‘CaOH’, ‘HC(CH3)3’).

Returns:

Dictionary containing element counts (e.g. for H2O: {‘H’: 2, ‘O’: 1}).

Return type:

counts (dict)

POSEIDON.atmosphere.elemental_ratio(included_species, X, element_1, element_2)

Calculate the abundance ratio between any two elements in the atmosphere.

Example: to compute the C/O ratio, use element_1 = ‘C’ and element_2 = ‘O’.

Parameters:
  • included_species (np.array of str) – List of all chemical species included in the model.

  • X (4D np.array of float) – Mixing ratio profiles.

  • element_1 (str) – First element in ratio.

  • element_2 (str) – First element in ratio.

Returns:

Abundance ratio in each layer, sector, and zone.

Return type:

element_ratio (3D np.array of float)

POSEIDON.atmosphere.profiles(P, R_p, g_0, PT_profile, X_profile, PT_state, P_ref, R_p_ref, log_X_state, included_species, bulk_species, param_species, active_species, CIA_pairs, ff_pairs, bf_species, N_sectors, N_zones, alpha, beta, phi, theta, species_vert_gradient, He_fraction, T_input, X_input, P_param_set, constant_gravity=False, chemistry_grid=None)

Main function to calculate the vertical profiles in each atmospheric column. The profiles cover the temperature, number density, mean molecular mass, layer radial extents, and mixing ratio arrays.

Notes: Most profiles are 3D arrays with format (N_layers, N_sectors, N_zones).

The mixing ratio profiles are 4D (separate profiles for each species) or 5D (CIA or free-free pairs). The layer index starts from the base of the atmosphere. The sector index starts from the evening terminator. The zone index starts from the dayside.

Parameters:
  • P (np.array of float) – Atmosphere pressure array (bar).

  • R_p (float) – Observed white light planet radius (m).

  • g_0 (float) – Gravitational field strength at white light radius (m/s^2).

  • PT_profile (str) – Chosen P-T profile parametrisation (Options: isotherm / gradient / two-gradients / Madhu / slope / file_read).

  • X_profile (str) – Chosen mixing ratio profile parametrisation (Options: isochem / gradient / two-gradients / file_read).

  • PT_state (np.array of float) – P-T profile state array.

  • P_ref (float) – Reference pressure (bar).

  • R_p_ref (float) – Planet radius corresponding to reference pressure (m).

  • log_X_state (2D np.array of float) – Mixing ratio state array.

  • included_species (np.array of str) – List of chemical species included in the model (including bulk species).

  • bulk_species (np.array of str) – Bulk species dominating atmosphere (e.g. [‘H2’, ‘He’]).

  • param_species (np.array of str) – Chemical species with parametrised mixing ratios.

  • active_species (np.array of str) – Spectroscopically active chemical species (see supported_opac.py).

  • CIA_pairs (np.array of str) – Collisionally-induced absorption (CIA) pairs.

  • ff_pairs (np.array of str) – Free-free absorption pairs.

  • bf_species (np.array of str) – Bound-free absorption species.

  • N_sectors (int) – Number of azimuthal sectors comprising the background atmosphere.

  • N_zones (int) – Number of zenith zones comprising the background atmosphere.

  • alpha (float) – Terminator opening angle (degrees).

  • beta (float) – Day-night opening angle (degrees).

  • phi (np.array of float) – Mid-sector angles (radians).

  • theta (np.array of float) – Mid-zone angles (radians).

  • species_vert_gradient (np.array of str) – Chemical species with a vertical mixing ratio gradient.

  • He_fraction (float) – Assumed H2/He ratio (0.17 default corresponds to the solar ratio).

  • T_input (np.array of float) – Temperature profile (only if provided directly by the user).

  • X_input (2D np.array of float) – Mixing ratio profiles (only if provided directly by the user).

  • P_param_set (float) – Only used for the Madhusudhan & Seager (2009) P-T profile. Sets the pressure where the reference temperature parameter is defined (P_param_set = 1.0e-6 corresponds to that paper’s choice).

  • constant_gravity (bool) – If True, disable inverse square law gravity (only for testing).

  • chemistry_grid (dict) – For models with a pre-computed chemistry grid only, this dictionary is produced in chemistry.py.

Returns:

Temperature profile (K). n (3D np.array of float):

Number density profile (m^-3).

r (3D np.array of float):

Radial distant profile (m).

r_up (3D np.array of float):

Upper layer boundaries (m).

r_low (3D np.array of float):

Lower layer boundaries (m).

dr (3D np.array of float):

Layer thicknesses (m).

mu (3D np.array of float):

Mean molecular mass (kg).

X (4D np.array of float):

Mixing ratio profile.

X_active (4D np.array of float):

Mixing ratios of active species.

X_CIA (5D np.array of float):

Mixing ratios of CIA pairs.

X_ff (5D np.array of float):

Mixing ratios of free-free pairs.

X_bf (4D np.array of float):

Mixing ratios of bound-free species.

Bool:

True if atmosphere physical, otherwise False.

Return type:

T (3D np.array of float)