POSEIDON.core

POSEIDON CORE ROUTINE.

Copyright 2023, Ryan J. MacDonald.

Module Contents

Functions

create_star(R_s, T_eff, log_g, Met[, T_eff_error, ...])

Initialise the stellar dictionary object used by POSEIDON.

create_planet(planet_name, R_p[, mass, gravity, ...])

Initialise the planet dictionary object used by POSEIDON.

define_model(model_name, bulk_species, param_species)

Create the model dictionary defining the configuration of the user-specified

wl_grid_constant_R(wl_min, wl_max, R)

Create a wavelength array with constant spectral resolution (R = wl/dwl).

wl_grid_line_by_line(wl_min, wl_max[, line_by_line_res])

Create a wavelength array with constant spectral resolution (R = wl/dwl).

read_opacities(model, wl[, opacity_treatment, T_fine, ...])

Load the various cross sections required by a given model. When using

make_atmosphere(planet, model, P, P_ref, R_p_ref[, ...])

Generate an atmosphere from a user-specified model and parameter set. In

check_atmosphere_physical(atmosphere, opac)

Checks that a specific model atmosphere is physical.

compute_spectrum(planet, star, model, atmosphere, opac, wl)

Calculate extinction coefficients, then solve the radiative transfer

compute_spectrum_c(planet, star, model, atmosphere, ...)

Calculate extinction coefficients, then solve the radiative transfer

compute_spectrum_p(planet, star, model, atmosphere, ...)

Calculate extinction coefficients, then solve the radiative transfer

load_data(data_dir, datasets, instruments, wl_model[, ...])

Load the user provided datasets into the format expected by POSEIDON.

set_priors(planet, star, model, data[, prior_types, ...])

Initialise the priors for each free parameter for a POSEIDON retrieval.

Attributes

cp

comm

rank

block

thread

POSEIDON.core.cp
POSEIDON.core.comm
POSEIDON.core.rank
POSEIDON.core.block
POSEIDON.core.thread
POSEIDON.core.create_star(R_s, T_eff, log_g, Met, T_eff_error=100.0, log_g_error=0.1, stellar_grid='blackbody', stellar_contam=None, f_het=None, T_het=None, log_g_het=None, f_spot=None, f_fac=None, T_spot=None, T_fac=None, log_g_spot=None, log_g_fac=None, wl=[], interp_backend='pysynphot')

Initialise the stellar dictionary object used by POSEIDON.

Stellar spectra are only required to compute transmission spectra if the star has unocculted stellar heterogeneities (spots/faculae) - since the stellar intensity cancels out in the transit depth. Hence a stellar spectrum is only added to the stellar dictionary if the user requests it.

Parameters:
  • R_s (float) – Stellar radius (m).

  • T_eff (float) – Stellar effective temperature (K).

  • log_g (float) – Stellar log surface gravity (log10(cm/s^2) by convention).

  • Met (float) – Stellar metallicity [log10(Fe/H_star / Fe/H_solar)].

  • T_eff_error (float) – A priori 1-sigma error on stellar effective temperature (K).

  • T_eff_error – A priori 1-sigma error on stellar log g (log10(cm/s^2)).

  • stellar_grid (string) –

    Chosen stellar model grid (Options: blackbody / cbk04 [for pysynphot] / phoenix [for pysynphot] /

    Goettingen-HiRes [for pymsg]).

  • stellar_contam (str) – Chosen prescription for modelling unocculted stellar contamination (Options: one_spot / one_spot_free_log_g / two_spots).

  • f_het (float) – For the ‘one_spot’ model, the fraction of stellar photosphere covered by either spots or faculae.

  • T_het (float) – For the ‘one_spot’ model, the temperature of the heterogeneity (K).

  • log_g_het (float) – For the ‘one_spot’ model, the log g of the heterogeneity (log10(cm/s^2)).

  • f_spot (float) – For the ‘two_spots’ model, the fraction of stellar photosphere covered by spots.

  • f_fac (float) – For the ‘two_spots’ model, the fraction of stellar photosphere covered by faculae.

  • T_spot (float) – For the ‘two_spots’ model, the temperature of the spot (K).

  • T_fac (float) – For the ‘two_spots’ model, the temperature of the facula (K).

  • log_g_spot (float) – For the ‘two_spots’ model, the log g of the spot (log10(cm/s^2)).

  • log_g_fac (float) – For the ‘two_spots’ model, the log g of the facula (log10(cm/s^2)).

  • wl (np.array of float) – Wavelength grid on which to output the stellar spectra (μm). If not provided, a fiducial grid from 0.2 to 5.4 μm will be used.

  • interp_backend (str) – Stellar grid interpolation package for POSEIDON to use. (Options: pysynphot / pymsg).

Returns:

Collection of stellar properties used by POSEIDON.

Return type:

star (dict)

POSEIDON.core.create_planet(planet_name, R_p, mass=None, gravity=None, log_g=None, T_eq=None, d=None, d_err=None, b_p=0.0)

Initialise the planet dictionary object used by POSEIDON.

The user only need provide one out of ‘mass’, ‘gravity’, or ‘log_g’ (depending on their unit of preference). Note that ‘gravity’ follows SI units (m/s^2), whilst ‘log_g’ follows the CGS convention (log_10 cm/s^2).

Parameters:
  • planet_name (str) – Identifier for planet object (e.g. HD209458b).

  • R_p (float) – Planetary radius (m).

  • mass (float) – Planetary mass (kg).

  • gravity (float) – Planetary gravity corresponding to observed radius (m/s^2).

  • log_g (float) – Instead of g, the user can provide log_10 (g / cm/s^2).

  • T_eq (float) – Planetary equilibrium temperature (zero albedo) (K).

  • d (float) – Distance to system (m).

  • d_err (float) – Measured error on system distance (m).

  • b_p (float) – Impact parameter of planetary orbit (m) – NOT in stellar radii!

Returns:

Collection of planetary properties used by POSEIDON.

Return type:

planet (dict)

POSEIDON.core.define_model(model_name, bulk_species, param_species, object_type='transiting', PT_profile='isotherm', X_profile='isochem', cloud_model='cloud-free', cloud_type='deck', opaque_Iceberg=False, gravity_setting='fixed', stellar_contam=None, offsets_applied=None, error_inflation=None, radius_unit='R_J', distance_unit='pc', PT_dim=1, X_dim=1, cloud_dim=1, TwoD_type=None, TwoD_param_scheme='difference', species_EM_gradient=[], species_DN_gradient=[], species_vert_gradient=[], surface=False, sharp_DN_transition=False, reference_parameter='R_p_ref', disable_atmosphere=False)

Create the model dictionary defining the configuration of the user-specified forward model or retrieval.

Parameters:
  • model_name (str) – Identifier for model in output files and plots.

  • bulk_species (list of str) – The chemical species (or two for H2+He) filling most of the atmosphere.

  • param_species (list of str) – Chemical species with parametrised mixing ratios (trace species).

  • object_type (str) – Type of planet / brown dwarf the user wishes to model (Options: transiting / directly_imaged).

  • 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 / chem_eq).

  • cloud_model (str) – Chosen cloud parametrisation (Options: cloud-free / MacMad17 / Iceberg).

  • cloud_type (str) – Cloud extinction type to consider (Options: deck / haze / deck_haze).

  • opaque_Iceberg (bool) – If using the Iceberg cloud model, True disables the kappa parameter.

  • gravity_setting (str) – Whether log_g is fixed or a free parameter. (Options: fixed / free).

  • stellar_contam (str) –

    Chosen prescription for modelling unocculted stellar contamination (Options: one_spot / one_spot_free_log_g / two_spots /

    two_spots_free_log_g).

  • offsets_applied (str) – Whether a relative offset should be applied to a dataset (Options: single_dataset).

  • error_inflation (str) – Whether to consider inflation of error bars in a retrieval (Options: Line15).

  • radius_unit (str) – Planet radius unit used to report retrieval results (Options: R_J / R_E)

  • distance_unit (str) – Distance to system unit used to report retrieval results (Options: pc)

  • PT_dim (int) – Dimensionality of the pressure-temperature field (uniform -> 1, a day-night or evening-morning gradient -> 2, both day-night and evening-morning gradients -> 3) (Options: 1 / 2 / 3).

  • X_dim (int) – Max dimensionality of the mixing ratio field (not all species need have gradients, this just specifies the highest dimensionality of chemical gradients – see the species_XX_gradient arguments) (Options: 1 / 2 / 3).

  • cloud_dim (int) – Dimensionality of the cloud model prescription (only the Iceberg cloud model supports 3D clouds) (Options: 1 / 2 / 3).

  • TwoD_type (str) – For 2D models, specifies whether the model considers day-night gradients or evening-morning gradients (Options: D-N / E-M).

  • TwoD_param_scheme (str) – For 2D models, specifies which quantities should be consider as free parameters (e.g. day & night vs. terminator & day-night diff.) (Options: absolute / difference).

  • species_EM_gradient (list of str) – List of chemical species with an evening-morning mixing ratio gradient.

  • species_DN_gradient (list of str) – List of chemical species with a day-night mixing ratio gradient.

  • species_vert_gradient (list of str) – List of chemical species with a vertical mixing ratio gradient.

  • surface (bool) – If True, model a surface via an opaque cloud deck.

  • sharp_DN_transition (bool) – For 2D / 3D models, sets day-night transition width (beta) to 0.

  • reference_parameter (str) – For retrievals, whether R_p_ref or P_ref will be a free parameter (Options: R_p_ref / P_ref).

  • disable_atmosphere (bool) – If True, returns a flat planetary transmission spectrum @ (Rp/R*)^2

Returns:

Dictionary containing the description of the desired POSEIDON model.

Return type:

model (dict)

POSEIDON.core.wl_grid_constant_R(wl_min, wl_max, R)

Create a wavelength array with constant spectral resolution (R = wl/dwl).

Parameters:
  • wl_min (float) – Minimum wavelength of grid (μm).

  • wl_max (float) – Maximum wavelength of grid (μm).

  • R (int or float) – Spectral resolution of desired wavelength grid.

Returns:

Model wavelength grid (μm).

Return type:

wl (np.array of float)

POSEIDON.core.wl_grid_line_by_line(wl_min, wl_max, line_by_line_res=0.01)

Create a wavelength array with constant spectral resolution (R = wl/dwl).

Parameters:
  • wl_min (float) – Minimum wavelength of grid (μm).

  • wl_max (float) – Maximum wavelength of grid (μm).

  • line_by_line_res (float) – Wavenumber resolution of pre-computer opacity database (0.01 cm^-1).

Returns:

Model wavelength grid (μm).

Return type:

wl (np.array of float)

POSEIDON.core.read_opacities(model, wl, opacity_treatment='opacity_sampling', T_fine=None, log_P_fine=None, opacity_database='High-T', device='cpu', wl_interp='sample', testing=False)

Load the various cross sections required by a given model. When using opacity sampling, the native high-resolution are pre-interpolated onto ‘fine’ temperature and pressure grids, then sampled onto the desired wavelength grid, and stored in memory. This removes the need to interpolate opacities during a retrieval. For line-by-line models, this function only stores Rayleigh scattering cross sections in memory (cross section interpolation is handled in other functions later).

Parameters:
  • model (dict) – A specific description of a given POSEIDON model.

  • wl (np.array of float) – Model wavelength grid (μm).

  • opacity_treatment (str) – Choice between opacity sampling or line-by-line cross sections (Options: opacity_sampling / line_by_line).

  • T_fine (np.array of float) – Fine temperature grid for opacity pre-interpolation.

  • log_P_fine (np.array of float) – Fine pressure grid for opacity pre-interpolation.

  • opacity_database (str) – Choice between high-temperature or low-temperature opacity databases (Options: High-T / Temperate).

  • wl_interp (str) – When initialising cross sections, whether to use opacity sampling or linear interpolation over wavenumber (Options: sample / linear) .

  • testing (bool) – For GitHub Actions automated tests. If True, disables reading the full opacity database (since GitHub Actions can’t handle downloading the full database - alas, 30+ GB is a little too large!).

Returns:

Collection of numpy arrays storing cross sections for the molecules, atoms, and ions contained in the model. The separate arrays store standard cross sections, CIA, free-free + bound-free opacity, and Rayleigh scattering cross sections.

Return type:

opac (dict)

POSEIDON.core.make_atmosphere(planet, model, P, P_ref, R_p_ref, PT_params=[], log_X_params=[], cloud_params=[], geometry_params=[], log_g=None, T_input=[], X_input=[], P_surf=None, P_param_set=0.01, He_fraction=0.17, N_slice_EM=2, N_slice_DN=4, constant_gravity=False, chemistry_grid=None)

Generate an atmosphere from a user-specified model and parameter set. In full generality, this function generates 3D pressure-temperature and mixing ratio fields, the radial extent of atmospheric columns, geometrical properties of the atmosphere, and cloud properties.

Parameters:
  • planet (dict) – Collection of planetary properties used POSEIDON.

  • model (dict) – A specific description of a given POSEIDON model.

  • P (np.array of float) – Model pressure grid (bar).

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

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

  • PT_params (np.array of float) – Parameters defining the pressure-temperature field.

  • log_X_params (np.array of float) – Parameters defining the log-mixing ratio field.

  • cloud_params (np.array of float) – Parameters defining atmospheric aerosols.

  • geometry_params (np.array of float) – Terminator opening angle parameters.

  • log_g (float) – Gravitational field of planet - only needed for free log_g parameter.

  • 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_surf (float) – Surface pressure of the planet.

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

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

  • N_slice_EM (even int) – Number of azimuthal slices in the evening-morning transition region.

  • N_slice_DN (even int) – Number of zenith slices in the day-night transition region.

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

Collection of atmospheric properties required to compute the resultant spectrum of the planet.

Return type:

atmosphere (dict)

POSEIDON.core.check_atmosphere_physical(atmosphere, opac)

Checks that a specific model atmosphere is physical.

Parameters:
  • atmosphere (dict) – Collection of atmospheric properties.

  • opac (dict) – Collection of cross sections and other opacity sources.

Returns:

True if atmosphere physical, otherwise False.

Return type:

Bool

POSEIDON.core.compute_spectrum(planet, star, model, atmosphere, opac, wl, spectrum_type='transmission', save_spectrum=False, disable_continuum=False, suppress_print=False, Gauss_quad=2, use_photosphere_radius=True, device='cpu', y_p=np.array([0.0]))

Calculate extinction coefficients, then solve the radiative transfer equation to compute the spectrum of the model atmosphere.

Parameters:
  • planet (dict) – Collection of planetary properties used by POSEIDON.

  • star (dict) – Collection of stellar properties used by POSEIDON.

  • model (dict) – A specific description of a given POSEIDON model.

  • atmosphere (dict) – Collection of atmospheric properties.

  • opac (dict) – Collection of cross sections and other opacity sources.

  • wl (np.array of float) – Model wavelength grid (μm).

  • spectrum_type (str) –

    The type of spectrum for POSEIDON to compute (Options: transmission / emission / direct_emission /

    transmission_time_average).

  • save_spectrum (bool) – If True, writes the spectrum to ‘./POSEIDON_output/PLANET/spectra/’.

  • disable_continuum (bool) – If True, turns off CIA and Rayleigh scattering opacities.

  • suppress_print (bool) – if True, turn off opacity print statements (in line-by-line mode).

  • Gauss_quad (int) – Gaussian quadrature order for integration over emitting surface * Only for emission spectra * (Options: 2 / 3).

  • use_photosphere_radius (bool) – If True, use R_p at tau = 2/3 for emission spectra prefactor.

  • device (str) – Experimental: use CPU or GPU (only for emission spectra) (Options: cpu / gpu)

  • y_p (np.array of float) – Coordinate of planet centre along orbit at the time the spectrum is computed (y_p = 0, the default, corresponds to mid-transit). For non-grazing transits of uniform stellar disks, the spectrum is identical at all times due to translational symmetry, so y_p = 0 is good for all times post second contact and pre third contact. Units are in m, not in stellar radii.

Returns:

The spectrum of the atmosphere (transmission or emission).

Return type:

spectrum (np.array of float)

POSEIDON.core.compute_spectrum_c(planet, star, model, atmosphere, opac, wl, spectrum_type='transmission', save_spectrum=False, disable_continuum=False, suppress_print=False, Gauss_quad=2, use_photosphere_radius=True, device='cpu', y_p=np.array([0.0]), contribution_molecule_list=[], bulk=False)

Calculate extinction coefficients, then solve the radiative transfer equation to compute the spectrum of the model atmosphere.

Parameters:
  • planet (dict) – Collection of planetary properties used by POSEIDON.

  • star (dict) – Collection of stellar properties used by POSEIDON.

  • model (dict) – A specific description of a given POSEIDON model.

  • atmosphere (dict) – Collection of atmospheric properties.

  • opac (dict) – Collection of cross sections and other opacity sources.

  • wl (np.array of float) – Model wavelength grid (μm).

  • spectrum_type (str) –

    The type of spectrum for POSEIDON to compute (Options: transmission / emission / direct_emission /

    transmission_time_average).

  • save_spectrum (bool) – If True, writes the spectrum to ‘./POSEIDON_output/PLANET/spectra/’.

  • disable_continuum (bool) – If True, turns off CIA and Rayleigh scattering opacities.

  • suppress_print (bool) – if True, turn off opacity print statements (in line-by-line mode).

  • Gauss_quad (int) – Gaussian quadrature order for integration over emitting surface * Only for emission spectra * (Options: 2 / 3).

  • use_photosphere_radius (bool) – If True, use R_p at tau = 2/3 for emission spectra prefactor.

  • device (str) – Experimental: use CPU or GPU (only for emission spectra) (Options: cpu / gpu)

  • y_p (np.array of float) – Coordinate of planet centre along orbit at the time the spectrum is computed (y_p = 0, the default, corresponds to mid-transit). For non-grazing transits of uniform stellar disks, the spectrum is identical at all times due to translational symmetry, so y_p = 0 is good for all times post second contact and pre third contact. Units are in m, not in stellar radii.

  • [] (contribution_molecule_list =) – Has a list of (active) molecules you want the spectral contribution function of

  • bulk (bool) – If true, also returns the spectral contribution of the bulk species

Returns:

The spectrum of the atmosphere (transmission or emission). If contribution_molecule_list != []

Returns spectrum (np.array of float), [‘Molecule’,contribution_spectrum] The second is a list where the first element is a string of the molecule being considered And the second element is the contribution_spectrum [np.array of float]

Return type:

spectrum (np.array of float)

POSEIDON.core.compute_spectrum_p(planet, star, model, atmosphere, opac, wl, spectrum_type='transmission', save_spectrum=False, disable_continuum=False, suppress_print=False, Gauss_quad=2, use_photosphere_radius=True, device='cpu', y_p=np.array([0.0]), contribution_molecule='', layer_to_ignore=0, total=False)

Calculate extinction coefficients, then solve the radiative transfer equation to compute the spectrum of the model atmosphere.

Parameters:
  • planet (dict) – Collection of planetary properties used by POSEIDON.

  • star (dict) – Collection of stellar properties used by POSEIDON.

  • model (dict) – A specific description of a given POSEIDON model.

  • atmosphere (dict) – Collection of atmospheric properties.

  • opac (dict) – Collection of cross sections and other opacity sources.

  • wl (np.array of float) – Model wavelength grid (μm).

  • spectrum_type (str) –

    The type of spectrum for POSEIDON to compute (Options: transmission / emission / direct_emission /

    transmission_time_average).

  • save_spectrum (bool) – If True, writes the spectrum to ‘./POSEIDON_output/PLANET/spectra/’.

  • disable_continuum (bool) – If True, turns off CIA and Rayleigh scattering opacities.

  • suppress_print (bool) – if True, turn off opacity print statements (in line-by-line mode).

  • Gauss_quad (int) – Gaussian quadrature order for integration over emitting surface * Only for emission spectra * (Options: 2 / 3).

  • use_photosphere_radius (bool) – If True, use R_p at tau = 2/3 for emission spectra prefactor.

  • device (str) – Experimental: use CPU or GPU (only for emission spectra) (Options: cpu / gpu)

  • y_p (np.array of float) – Coordinate of planet centre along orbit at the time the spectrum is computed (y_p = 0, the default, corresponds to mid-transit). For non-grazing transits of uniform stellar disks, the spectrum is identical at all times due to translational symmetry, so y_p = 0 is good for all times post second contact and pre third contact. Units are in m, not in stellar radii.

  • contribution_molecule (string) – The molecule that will be turned off in the layer to ignore

  • layer_to_ignore (int) – This is looped over in the contribution_func.py Indicates the pressure layer in which the molecule will be turned off

  • total (bool) – If true, the total pressure contribution spectrum is returned

Returns:

The spectrum of the atmosphere (transmission or emission). If contribution_molecule_list != [] or total = True

Returns a spectrum where the molecule indicated is turned off for the pressure layer to ignore If total = true, then a whole layer is turned off

Return type:

spectrum (np.array of float)

POSEIDON.core.load_data(data_dir, datasets, instruments, wl_model, offset_datasets=None, wl_unit='micron', bin_width='half', spectrum_unit='(Rp/Rs)^2', skiprows=None)

Load the user provided datasets into the format expected by POSEIDON. Also generate the functions required for POSEIDON to later calculate the binned data for each instrument (e.g. the PSFs for each instrument) corresponding to model spectra.

Parameters:
  • data_dir (str) – Path to the directory containing the user’s data files.

  • datasets (list of str) – List containing file names of the user’s data files.

  • instruments (list of str) – List containing the instrument names corresponding to each data file (e.g. WFC3_G141, JWST_NIRSpec_PRISM, JWST_NIRISS_SOSS_Ord2).

  • wl_model (np.array of float) – Model wavelength grid (μm).

  • offset_datasets (list of str) – If applying a relative offset to one or more datasets, this list gives the file names of the datasets that will have free offsets applied (note: currently only supports one offset dataset).

  • wl_unit (str) – Unit of wavelength column (first column in file) (Options: micron (or equivalent) / nm / A / m)

  • bin_width (str) – Whether bin width (second column) is half or full width (Options: half / full).

  • spectrum_unit (str) – Unit of spectrum (third column) and spectrum errors (fourth column) (Options: (Rp/Rs)^2 / Rp/Rs / Fp/Fs / Fp (or equivalent units)).

  • skiprows (int) – The number of rows to skip (e.g. use 1 if file has a header line).

Returns:

Collection of data properties required for POSEIDON’s instrument simulator (i.e. to create simulated binned data during retrievals).

Return type:

data (dict)

POSEIDON.core.set_priors(planet, star, model, data, prior_types={}, prior_ranges={})

Initialise the priors for each free parameter for a POSEIDON retrieval.

If the user does not provide a prior type or prior range for one or more of the parameters, this function will prescribe a default prior with a wide range. Thus the user can choose the degree to which they would like to ‘micromanage’ the assignment of priors.

Disclaimer: while using default priors can be good for exploratory retrievals, for a publication we strongly suggest you explicitly specify your priors - you’ll need to give your priors in a Table somewhere anyway, so it’s generally a good idea to know what they are ;)

Parameters:
  • planet (dict) – Collection of planetary properties used by POSEIDON.

  • star (dict) – Collection of stellar properties used by POSEIDON.

  • model (dict) – A specific description of a given POSEIDON model.

  • data (dict) – Collection of data properties in POSEIDON format.

  • prior_types (dict) – User-provided dictionary containing the prior type for each free parameter in the retrieval moel (Options: uniform, gaussian, sine, CLR).

  • prior_ranges (dict) –

    User-provided dictionary containing numbers defining the prior range for each free parameter in the retrieval model (Options: for ‘uniform’ [low, high], for ‘gaussian’ [mean, std],

    for ‘sine’ [high] - only for 2D/3D angle parameters, for ‘CLR’ [low] - only for mixing ratios).

Returns:

Collection of the prior types and ranges used by POSEIDON’s retrieval module.

Return type:

priors (dict)