Retrievals: Basic Retrieval Tutorial#
Written by Evert Nasedkin. Please cite pRT’s retrieval package (Nasedkin et al. 2024) in addition to pRT (Mollière et al. 2019) if you make use of the retrieval package for your work.
What is a retrieval?
Let’s start with the basics! An atmospheric retrieval is the process through which we fit a model of a spectrum to data. In the exoplanet community this usually means Bayesian model fitting, for which you can find an introduction here. The model usually incorporates elements of chemistry and radiative transfer in order to compute a spectrum - see the getting started page for an example.
If you’re here, we’re going to assume that you already have an exoplanet emission or transmission spectrum, and you want to learn something about the planet’s properties. We’ll pRT spectral modeling capabilities combined with a sampler which will allow us to statistically measure the properties of the planet of interest. As inputs for this, you’ll need a 1D spectrum with error bars (or a covariance matrix), and a model function that generates a synthetic spectrum. We have a few example datasets included in petitRADTRANS/retrieval/examples/, and several forward models implemented in petitRADTRANS/retrieval/models.py, which you can use directly, or as templates to build your own model. These files are included in the pRT package folder. Alternatively you can access them through the git links (just click on the folder / file names above).
In this example we will make use of an HST transmission spectrum, located locally in petitRADTRANS/retrieval/examples/transmission/observations/HST/hst_example_clear_spec.txt, or, again, on git. This dataset is included when you install pRT, and we’ll automatically locate the data and output files throughout the notebook. The outputs of
the retrieval example will be a best-fit spectrum, an atmospheric pressure-temperature profile, and the posterior distributions of all of the parameters of your model atmosphere. Of course, you do not have to use pRT’s retrieval package to run retrievals with pRT-generated spectra, but we believe that what is described below may be sufficient for most people, and will thus save you from most of the tedious implementation work.
Getting started#
You should already have an installation of pRT on your machine; please see the installation manual if you have not installed it yet. As mentioned on the installation page, pyMultiNest - and therefore Multinest, which implements a special form of nested sampling, are also required. Multinest allows for easy parallelization of the retrievals when run on large clusters, but can also be installed locally on
your laptop to run small retrievals, such as our example here (you can also parallelize your runs then). Ultranest is required if you want to use that as your sampler rather than pyMultiNest. See the Ultranest documentation for why you might want to choose this method. Using nested sampling rather than, for example,
MCMC is often faster, handles multimodal cases better, and directly provides the Bayesian evidence, which allows for easier model comparison. Again, see here if you need an introduction to Bayesian model fitting and its terminology. We’ve also found that in pyMultinest is often faster than Ultranest, and generally recommend its use.
In this tutorial, we will outline the process of setting up a RetrievalConfig object, which is the class used to set up a pRT retrieval. The basic process will always be to set up the configuration, and then pass it to the Retrieval class to run the retrieval using pyMultiNest. Several standard plotting outputs will also be produced by the retrieval class. Most of the classes and functions used in this tutorial have more advanced features than what will be explained here, so it’s highly recommended to take a look at the code and API documentation. There should be enough flexibility built in to cover most typical retrieval studies, but if you have feature requests please get in touch, or open an issue on gitlab.
[ ]:
# Let's start by importing everything we need.
from petitRADTRANS.retrieval import Retrieval, RetrievalConfig
import tensorflow_probability.substrates.jax as tfp
from petitRADTRANS.retrieval.models import isothermal_transmission
from petitRADTRANS import physical_constants as cst
from jax import config
import os
# To not have numpy start parallelizing on its own when using MPI, OpenMP, mpi4py
os.environ["OMP_NUM_THREADS"] = "1"
os.environ["DYLD_LIBRARY_PATH"] = "/Users/nasedkin/software/MultiNest/lib/:$DYLD_LIBRARY_PATH"
os.environ["JAX_CHECK_TRACER_LEAKS"] = "False"
os.environ["XLA_FLAGS"] = "--xla_force_host_platform_device_count=2"
config.update("jax_enable_x64", True)
# Import the class used to set up the retrieval.
# Import an atmospheric model function
# Physical constants
# Import prior distributions
tfpd = tfp.distributions
Let’s start with a simple model. We’ll use an isothermal temperature profile, and use “free chemistry” - allowing the chemical abundances to vary freely within fixed boundaries, but requiring them to be vertically constant throughout the atmosphere. The basic process for every retrieval will be the same. We need to set up a RetrievalConfig object, adding in the Data objects that we want to fit together with the Parameter objects and their priors, which we can use to customize our retrieval.
Retrieval configuration object and parameters#
Let’s start out by setting up a simple run definition. We’ll add the data AFTER we define the model function below. We start out by initialising the RetrievalConfig object (which we call retrieval_config here), and then adding each of the parameters of interest to it.
Each parameter must be given a name, which matches the name used in the model function defined below. Parameters can be set to fixed or free. Fixed parameters must be given a value, while free parameters are given a function that transforms the unit hypercube generated by multinest into physical prior space. Various prior functions are stored in retrieval.utils, see the API documentation. Full details of the parameters can be found
in the API documentation.
[ ]:
# Lets start out by setting up a simple run definition
# We'll add the data AFTER we define the model function below
# Full details of the parameters can be found in retrieval_config.py
# If our retrieval has already run before, you can set the run_mode to 'evaluate' to only make some plots.
retrieval_config = RetrievalConfig(
retrieval_name="hst_example_clear_spec",
run_mode="retrieval", # this can be 'evaluate' to get the results of a retrieval
sampler_type="jaxnestedsampler", # The sampler to use. Here we use the JAXNS sampler.
adaptive_mesh_refinement=False, # We won't be using adaptive mesh refinement for the pressure grid
# This would turn on scattering when calculating emission spectra (always activated in transmission).
scattering_in_emission=False,
)
# Let's start with the parameters that will remain fixed during the retrieval
retrieval_config.add_parameter(
name="stellar_radius", # name
is_free_parameter=False, # is_free_parameter, so stellar_radius is not retrieved here.
value=0.651 * cst.r_sun,
plot_in_corner=False, # we won't be plotting this parameter in the corner plot
)
# Log of the surface gravity in cgs units.
# The transform_prior_cube_coordinate function transforms a uniform random sample x, drawn between 0 and 1.
# to the physical units we're using in pRT. We use Python's lambda functions here.
retrieval_config.add_parameter(
name="log_g",
is_free_parameter=True, # is_free_parameter: we are retrieving log(g) here
distribution=tfpd.Uniform(low=2.0, high=5.5), # this means that log g can vary between 2 and 5.5
plot_in_corner=True, # we will plot this parameter in the corner plot
corner_ranges=(2.0, 5.5), # this sets the range of the axis in the corner plot
corner_label=r"$\log(g)$", # this sets the label of the axis in the corner plot
)
# Note that logg is usually not a free parameter in retrievals of
# transmission spectra, at least it is much better constrained as being
# assumed here, as the planetary mass and radius are usually known.
# Planet radius in cm
retrieval_config.add_parameter(
name="planet_radius",
is_free_parameter=True,
distribution=tfpd.Uniform(low=0.2 * cst.r_jup, high=0.4 * cst.r_jup), # radius varies between 0.2 and 0.4 R_jup
plot_in_corner=True, # we will plot this parameter in the corner plot
corner_ranges=(0.2, 0.4), # this sets the range of the axis in the corner plot
corner_label=r"R$_{pl}$", # this sets the label of the axis in the corner plot
corner_transform=lambda x: x / cst.r_jup_mean, # this transforms the units in the corner plot to R_jup
)
# Atmospheric temperature in Kelvin
retrieval_config.add_parameter(
name="temperature",
is_free_parameter=True,
distribution=tfpd.Uniform(low=300, high=2000),
plot_in_corner=True, # we will plot this parameter in the corner plot
corner_ranges=(300, 2000), # this sets the range of the axis in the corner plot
corner_label=r"Temp", # this sets the label of the axis in the corner plot
)
# Let's include a grey cloud as well to see what happens, even though we assumed a clear atmosphere.
retrieval_config.add_parameter(
"log_Pcloud",
is_free_parameter=True,
distribution=tfpd.Uniform(low=-6, high=2), # the atmosphere can have an opaque cloud between 1e-6 and 100 bar
plot_in_corner=True, # we will plot this parameter in the corner plot
corner_ranges=(-6, 2), # this sets the range of the axis in the corner plot
corner_label=r"$\log P_{cl}$", # this sets the label of the axis in the corner plot
)
Opacities#
petitRADTRANS can be setup to use line-by-line opacities at \(\lambda/\Delta \lambda = 10^6\), or a lower resolution correlated-k mode (which we use here). The names of each species must follow the naming convention explained here, which also contains a list of available species. ‘Getting Started’ also has an example for how to load opacities in pRT.
In addition to line opacities, pRT can also include, for example, Rayleigh scattering, collision-induced absorption (CIA) opacities, or opacities for various condensate cloud species. Again, see here for a list. If a cloud species is used, it will have a set of parameters added in order to retrieve the cloud abundance. We give such an example below. The parameters added depend on the model chosen.
If you aren’t sure which species are available, the RetrievalConfig class has a few functions that will list the currently available line, CIA and cloud opacities, see the API documentation.
Opacity resolution#
Do not append resolution identifiers, for example .R100, to the species name, if you want to use a correlated-k species at lower resolution of \(\lambda/\Delta \lambda = 100\), instead of the nominal \(\lambda/\Delta \lambda = 1000\). Instead, you will specify a model resolution when adding data (see below). If you set the model resolution for your data, and do not have existing c-k tables at that resolution,
Exo_k will be used to bin the higher resolution tables down to a lower resolution. Please cite Exo_k and the corresponding paper (Leconte et al. 2021) if you make use of this functionality in pRT.
You can use retrieval_config.list_available_line_species() to see which line species are available on your machine, but note that pRT will automatically download those listed here if you request some that you do not have yet.
Let’s add some opacities below now:
[3]:
retrieval_config.set_rayleigh_species(["H2", "He"])
retrieval_config.set_continuum_opacities(["H2--H2", "H2--He"])
# Here we setup the line species for a free retrieval,
# setting the prior bounds of the log10 abundance with the abund_lim parameter
# So the retrieved value is the log mass fraction.
retrieval_config.set_line_species(
["1H2-16O__POKAZATEL", "12C-1H4__HITEMP", "12C-16O2__UCL-4000", "C-O-NatAbund__HITEMP"],
use_equilibrium_chemistry=False,
free_mass_fraction_limits=(-6.0, 0.0),
plot_in_corner=True,
)
Model Functions#
In addition to setting up a retrieval with a standard model, you may also want to write your own model to be used in the retrieval setup, as shown in the basic retrieval tutorial.. For this example we will use the built-in isothermal_transmission model.
Data#
With the model set up, we can read in the data. The data must be a 1D spectrum with error bars or a covariance matrix. The input format for ASCII data is
# Wavelength [micron], Flux [W/m2/micron or (Rp/Rstar)^2], Flux Error [W/m2/micron or (Rp/Rstar)^2]
Optionally, there can also be an additional column to specify the wavelength bins, which would then have the format:
# Wavelength [micron], Bins [micron], Flux [W/m2/micron or (Rp/Rstar)^2], Flux Error [W/m2/micron or (Rp/Rstar)^2]
The file can be comma or space separated, and by default ‘#’ is the comment character. The wavelength column must be sorted and is in ascending order. If you write your own model function, the units of your data can vary from this template, but the units of your data and the output of your model function must match. Throughout our documentation we will use W/m^2/micron for emission spectra and (Rp/Rstar)^2 for consistency between our model functions and the example data.
The fits file requirements are more strict. An extension labelled ‘SPECTRUM’ is required, with three fields in a binary table: ‘WAVELENGTH’, ‘FLUX’ and ‘COVARIANCE’. An example of this is the HR8799e_spectra.fits file included in petitRADTRANS/retrieval/examples/emission/, see pRT’s gitlab page. If your file is not of that format, but you want to
add a covariance matrix, you could also fill in these attributes of the data object that gets added to the retrieval configuration after calling add_data():
covariance: the covariance matrixinv_cov: the inverse of the covariance matrixlog_covariance_determinant: the natural logarithm of the derterminant of the covariance matrix
The Data class is arguably the most important part of setting up the retrieval. Not only do you input your data here, but you also choose your model function and resolution. This means that you can design a retrieval around different datatypes and retrieve simultaneously on both - for example, if you want the day and nightside of a planet, or want to combine the eastward and westward limbs of a transmission spectrum with different models.
By adjusting the resolution you can also speed up the retrieval. By default, pRT uses the correlated-k method for computing opacities, which works up to \(\lambda/\Delta \lambda = 1000\). Here you need to pick 'c-k' as line_opacity_mode, see below. Often low resolution data does not require this resolution, and the computation can be sped up using the model_resolution parameter, which will use the Exo-k package to bin down the correlated-k tables.
If you use the 'lbl' (line-by-line) mode, pRT can uses opacities up to \(\lambda/\Delta \lambda = 10^6\), though this is quite slow and memory intensive. By setting the model_resolution parameter, pRT will sample the high resolution line list at the specified resolution (e.g., model_resolution=10**5. will down-sample by a factor of 10). In this case you must make sure that the down-sampled calculations match the binned-down results of the \(\lambda/\Delta \lambda = 10^6\)
mode at the required data resolution (i.e., sampling can introduce noise and you have to test for that).
The Data class can also handle photometric data, which is described in the emission tutorial.
[ ]:
# Finally, we associate the model with the data, and we can run our retrieval!
# Here we import petitRADTRANS to find the path of the example files on your machine.
# In general this is not required, you just put the files in the folder that you are running
# Your script in, for example.
path_to_data = "./" # Default location for the example
retrieval_config.add_data(
"HST", # simulated HST data
f"{path_to_data}retrievals/transmission/observations/HST/hst_example_clear_spec.txt", # where is the data
model_generating_function=isothermal_transmission, # the atmospheric model to use
line_opacity_mode="c-k",
data_resolution=60, # the spectral resolution of the data
model_resolution=120, # the resolution of the c-k tables for the lines
)
# This model is a noise-free, clear atmosphere transmission spectrum for a sub-neptune type planet
# The model function used to calculate it was slightly different, and used different opacity tables
# than what we'll use in the retrieval, but the results don't significantly change.
# You could even see it as a simulation of reality since we never have the perfect model or opacities
# that "Nature" is using.
#
# In general, it is also necessary to use the data_resolution and model_resolution arguments.
# The data_resolution argument will convolve the generated model spectrum by the instrumental
# resolution before to calculating the chi squared value. As mentioned at the bottom of
# "Rebinning opacities", the model resolution should be >= 2 times larger than the data resolution.
#
# The model_resolution parameter will cause pRT to use exo-k to generate low-resolution correlated-k opacity tables
# in order to speed up the radiative transfer calculations (if this has not been done yet).
Plotting#
Let’s also set up some plotting details so that we can generate nice outputs after the retrieval has finished.
Each parameter can be added to the corner plot, its label changed, and the values transformed to more digestible units (e.g., the planet radius in jupiter radii, rather than cm). We can also set the bounds and scaling on the best-fit spectrum plot and the limits for the P-T profile plot. With this complete, our retrieval is ready to go.
Most parameters include a default setting, so the plots will be created even if you don’t set any plotting parameters, but the outputs might not be very informative. In general, the possible arguments to plot_kwargs follow the naming conventions of matplotlib arguments and functions, with some additions. Full details of the plotting can be found in the Retrieval class, see the API documentation.
[ ]:
# Define axis properties of spectral plot if run_mode == 'evaluate'
def configure_retrieval_plotter(retrieval):
retrieval.plotter.spec_xlabel = "Wavelength [micron]"
retrieval.plotter.spec_ylabel = "Transit Depth [ppm]"
retrieval.plotter.y_axis_scaling = 1e6
retrieval.plotter.xscale = "linear"
retrieval.plotter.yscale = "linear"
retrieval.plotter.resolution = None # maximum resolution, will rebin the data
# Define from which observation object to take P-T in evaluation mode.
retrieval.plotter.reference_data_name = "HST"
retrieval.plotter.temp_limits = [150, 3000]
retrieval.plotter.press_limits = [1e2, 1e-5]
# If in doubt, define all of the plot_kwargs used here.
Running the retrieval#
At this point, we are almost ready to run the retrieval! We only need to do is pass the RunDefinition object “retrieval_config” we just created to an object of the Retrieval class first.
[6]:
output_directory = "./retrievals/runs"
retrieval = Retrieval(
configuration=retrieval_config,
output_directory=output_directory,
use_mpi=False,
evaluate_sample_spectra=False,
use_prt_plot_style=True,
)
Setting up Radtrans object for data 'HST'...
Loading Radtrans opacities...
Loading line opacities of species '1H2-16O__POKAZATEL.R120' from file '/Users/nasedkin/python-packages/petitRADTRANS/input_data/opacities/lines/correlated_k/H2O/1H2-16O/1H2-16O__POKAZATEL.R120_0.3-50mu.ktable.petitRADTRANS.h5'... Done.
Loading line opacities of species '12C-1H4__HITEMP.R120' from file '/Users/nasedkin/python-packages/petitRADTRANS/input_data/opacities/lines/correlated_k/CH4/12C-1H4/12C-1H4__HITEMP.R120_0.1-250mu.ktable.petitRADTRANS.h5'...
---------------------------------------------------------------------------
KeyboardInterrupt Traceback (most recent call last)
Cell In[6], line 3
1 output_directory = f"./retrievals/runs"
----> 3 retrieval = Retrieval(
4 configuration=retrieval_config,
5 output_directory=output_directory,
6 use_mpi=False,
7 evaluate_sample_spectra=False,
8 use_prt_plot_style=True,
9 )
File ~/python-packages/petitRADTRANS/petitRADTRANS/retrieval/retrieval.py:178, in Retrieval.__init__(self, configuration, output_directory, use_mpi, evaluate_sample_spectra, corner_plot_names, use_prt_plot_style, test_plotting, uncertainties_mode, print_log_likelihood_for_debugging, generate_mock_data)
174 os.makedirs(os.path.join(self.output_directory,
175 'evaluate_' + self.configuration.retrieval_name), exist_ok=True)
177 # Setup pRT Objects for each data structure.
--> 178 self.setup_data()
179 self._n_free_params_total_retrieval = 0
180 free_param_names_list = []
File ~/python-packages/petitRADTRANS/petitRADTRANS/retrieval/retrieval.py:3157, in Retrieval.setup_data(self, scaling, width)
3154 pressures = self.configuration.pressures
3156 # Set up the pRT objects for the given dataset
-> 3157 radtrans = Radtrans(
3158 pressures=pressures,
3159 line_species=copy.copy(tuple(species)),
3160 rayleigh_species=copy.copy(self.configuration.rayleigh_species),
3161 gas_continuum_contributors=copy.copy(self.configuration.continuum_opacities),
3162 cloud_species=copy.copy(self.configuration.cloud_species),
3163 line_opacity_mode=data.line_opacity_mode,
3164 wavelength_boundaries=data.wavelength_boundaries,
3165 scattering_in_emission=self.configuration.scattering_in_emission,
3166 line_by_line_opacity_sampling=line_by_line_sampling
3167 )
3168 data.radtrans_object = radtrans
3170 # Setup for non-uniform spectral resolution
File ~/python-packages/petitRADTRANS/petitRADTRANS/radtrans.py:258, in Radtrans.__init__(self, pressures, wavelength_boundaries, line_species, gas_continuum_contributors, rayleigh_species, cloud_species, line_opacity_mode, line_by_line_opacity_sampling, scattering_in_emission, emission_angle_grid, anisotropic_cloud_scattering, path_input_data)
254 raise ValueError(f"emission_angle_grid must be a dictionary or an iterable of shape (2, n_angles), but "
255 f"is of type {type(emission_angle_grid)}")
257 # Load all opacities
--> 258 self.load_all_opacities()
File ~/python-packages/petitRADTRANS/petitRADTRANS/radtrans.py:4225, in Radtrans.load_all_opacities(self)
4222 print("Loading Radtrans opacities...")
4224 # Load line opacities
-> 4225 self.load_line_opacities(self._path_input_data)
4227 # Read continuum opacities
4228 # CIA
4229 if len(self._gas_continuum_contributors) > 0:
File ~/python-packages/petitRADTRANS/petitRADTRANS/radtrans.py:4632, in Radtrans.load_line_opacities(self, path_input_data)
4629 print(f" Loading line opacities of species '{species}' from file '{hdf5_file}'...", end='')
4631 if self._line_opacity_mode == 'c-k':
-> 4632 self._lines_loaded_opacities['opacity_grid'][species] = self.load_hdf5_ktables(
4633 file_path_hdf5=hdf5_file,
4634 frequencies=self._frequencies,
4635 g_size=self._lines_loaded_opacities['g_gauss'].size,
4636 temperature_pressure_grid_size=self._lines_loaded_opacities['temperature_pressure_grid'][
4637 species].shape[0]
4638 )
4639 elif self._line_opacity_mode == 'lbl':
4640 self._lines_loaded_opacities['opacity_grid'][species] = self.load_hdf5_line_opacity_table(
4641 file_path_hdf5=hdf5_file,
4642 frequencies=self._frequencies,
4643 line_by_line_opacity_sampling=self._line_by_line_opacity_sampling
4644 )
File ~/python-packages/petitRADTRANS/petitRADTRANS/radtrans.py:4425, in Radtrans.load_hdf5_ktables(cls, file_path_hdf5, frequencies, g_size, temperature_pressure_grid_size, return_radtrans_opacities)
4423 """Load k-coefficient tables in HDF5 format, based on the ExoMol setup."""
4424 with h5py.File(file_path_hdf5, 'r') as f:
-> 4425 k_table = np.array(f['kcoeff'])
4427 if not return_radtrans_opacities:
4428 return (
4429 f['DOI'][()][0].decode('utf-8'),
4430 f['bin_centers'][:],
(...)
4438 f['weights'][:]
4439 )
File h5py/_objects.pyx:54, in h5py._objects.with_phil.wrapper()
File h5py/_objects.pyx:55, in h5py._objects.with_phil.wrapper()
File /opt/anaconda3/envs/jaxprt/lib/python3.12/site-packages/h5py/_hl/dataset.py:1091, in Dataset.__array__(self, dtype, copy)
1088 if self.size == 0:
1089 return arr
-> 1091 self.read_direct(arr)
1092 return arr
File /opt/anaconda3/envs/jaxprt/lib/python3.12/site-packages/h5py/_hl/dataset.py:1047, in Dataset.read_direct(self, dest, source_sel, dest_sel)
1044 dest_sel = sel.select(dest.shape, dest_sel)
1046 for mspace in dest_sel.broadcast(source_sel.array_shape):
-> 1047 self.id.read(mspace, fspace, dest, dxpl=self._dxpl)
KeyboardInterrupt:
[ ]:
configure_retrieval_plotter(retrieval)
Now, we can call the retrieval.run() method and run the retrieval!
There are a few additional parameters that can be adjusted, but the defaults should work well for almost all use cases. In general it may not be wise to run a retrieval on your laptop, as it can be quite computationally expensive. However, the included HST example should be able to be run fairly easily! (minutes!). We only use 40 live points in the example so that it runs in a reasonable length of time, but we highly recommend using more live points for your own retrievals! If you want to run
the retrieval on a cluster, or just on multiple cores on your local machine, set use_MPI=True. You should copy the contents of this notebook into a single .py file, which then you can run using mpiexec -n $NUMBER_OF_CORES python my_retrieval_file.py, substituting in the number of cores you want to use to run the retrieval. Most of the various parameters used to control pyMultiNest or Ultranest can be set in the retrieval.run() function, see the API
documentation.
[ ]:
# %%time
retrieval.run(
num_live_points=40,
max_samples=1e6,
verbose=True,
)
Starting retrieval hst_example_clear_spec with jaxnestedsampler
Testing data 'HST':
wavelengths:
OK (no NaN, infinite, or negative value detected)
spectrum:
OK (no NaN, infinite, or negative value detected)
uncertainties:
OK (no NaN, infinite, or negative value detected)
Testing model function for data 'HST'...
No errors detected in the model functions!
/Users/nasedkin/python-packages/jaxns/src/jaxns/framework/jaxify.py:29: UserWarning: You're using a non-JAX log-likelihood function. This may be slower than a JAX log-likelihood function. Also, you are responsible for ensuring that the function is deterministic. Also, you cannot use learnable parameters in the likelihood call.
warnings.warn(
/opt/anaconda3/envs/jaxprt/lib/python3.12/site-packages/jax/_src/numpy/lax_numpy.py:5820: FutureWarning: None encountered in jnp.array(); this is currently treated as NaN. In the future this will result in an error.
return array(a, dtype=dtype, copy=bool(copy), order=order, device=device)
/Users/nasedkin/python-packages/petitRADTRANS/petitRADTRANS/radtrans.py:1367: UserWarning: temperature array ('temperatures') contains negative or null values ([0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0.
0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0.
0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0.
0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0. 0.
0. 0. 0. 0.])
warnings.warn(f"temperature array ('temperatures') contains negative or null values ({temperatures})")
INFO:jaxns:Number of Markov-chains set to: 40
Starting JAXNS retrieval: hst_example_clear_spec
Creating initial state with 40 live points.
Once the retrieval is complete, we can use several functions to generate plots of the best fit spectrum, the pressure-temperature profile and the corner plots. A standard output file will be produced when the retrieval is run, describing the prior boundaries, the data used in the retrieval, and if the retrieval is finished it will include the best fit and median model parameters. Check out this file, together with the other outputs in ./retrievals/runs/evaluate_hst_example_clear_spec.
Let’s make some plots to see how our retrieval went! We need to start by reading in the results. If we want to read results from multiple retrievals, use the ret_names argument. Currently, only the corner plots make use of this feature.
[ ]:
# These are dictionaries in case we want to look at multiple retrievals.
# The keys for the dictionaries are the retrieval_names
sample_dict, parameter_dict = retrieval.get_samples(output_directory=output_directory)
# Pick the current retrieval to look at.
samples_use = sample_dict[retrieval.configuration.retrieval_name]
parameters_read = parameter_dict[retrieval.configuration.retrieval_name]
Best Fit Spectrum#
Here we will plot the best fit spectrum, together with the data used in the retrieval and the residuals between the two. Using the mode parameter, it is also possible to plot the spectrum using the 'median' fit instead of the maximum likelihood sample. The chi squared value presented in the figure is the reduced chi square, that is the chi square divided by the degrees of freedom (number of data points - number of parameters).
[ ]:
# Plotting the best fit spectrum
# This will generate a few warnings, but it's fine.
fig, ax, ax_r = retrieval.plot_spectra(samples_use, parameters_read, refresh=True, mode="bestfit")
Not in evaluate mode. Changing run mode to evaluate.
Plotting Best-fit spectrum
Best fit likelihood = 172.32
Loading Radtrans opacities...
Loading line opacities of species 'H2O__POKAZATEL' from file '/home/dblain/petitRADTRANS/input_data/opacities/lines/correlated_k/H2O/1H2-16O/1H2-16O__POKAZATEL.R1000_0.3-50mu.ktable.petitRADTRANS.h5'... Done.
Loading line opacities of species 'CH4__HITEMP' from file '/home/dblain/petitRADTRANS/input_data/opacities/lines/correlated_k/CH4/12C-1H4/12C-1H4__HITEMP.R1000_0.1-250mu.ktable.petitRADTRANS.h5'... Done.
Loading line opacities of species 'CO2' from file '/home/dblain/petitRADTRANS/input_data/opacities/lines/correlated_k/CO2/12C-16O2/12C-16O2__UCL-4000.R1000_0.3-50mu.ktable.petitRADTRANS.h5'... Done.
Loading line opacities of species 'CO-NatAbund__HITEMP' from file '/home/dblain/petitRADTRANS/input_data/opacities/lines/correlated_k/CO/C-O-NatAbund/C-O-NatAbund__HITEMP.R1000_0.1-250mu.ktable.petitRADTRANS.h5'... Done.
Successfully loaded all line opacities
Loading CIA opacities for H2-H2 from file '/home/dblain/petitRADTRANS/input_data/opacities/continuum/collision_induced_absorptions/H2--H2/H2--H2-NatAbund/H2--H2-NatAbund__BoRi.R831_0.6-250mu.ciatable.petitRADTRANS.h5'... Done.
Loading CIA opacities for H2-He from file '/home/dblain/petitRADTRANS/input_data/opacities/continuum/collision_induced_absorptions/H2--He/H2--He-NatAbund/H2--He-NatAbund__BoRi.DeltaWavenumber2_0.5-500mu.ciatable.petitRADTRANS.h5'... Done.
Successfully loaded all CIA opacities
Successfully loaded all opacities
Best fit 𝛘^2 = 1.31
Best fit 𝛘^2/DoF = 0.11
/mnt/c/Users/doria/OneDrive/Documents/programs/Python/petitRADTRANS/petitRADTRANS/retrieval/retrieval.py:4374: UserWarning: This figure includes Axes that are not compatible with tight_layout, so results might be incorrect.
plt.tight_layout()
Pressure-Temperature profile#
Here we show the retrieved pressure-temperature profile throughout the atmosphere. The contours show the 1, 2 and 3 sigma confidence intervals in temperature around the median retrieved profile. If contribution=True, the opacity shading and dashed line show the region of the atmosphere that contributes to the observed spectrum, highlighting which part of the pressure profile is being probed.
[ ]:
# Plotting the PT profile
fig, ax = retrieval.plot_pt(sample_dict, parameters_read, contribution=False)
Plotting PT profiles
Corner plot#
Here we show the corner plot for the retrieval, which shows the marginal posterior probability distributions for each parameter (though it is also possibe to plot only a subset of the parameters). The corner plot produces 1, 2 and 3 sigma contours for the 2D plots, and highlights the 1 sigma confidence interval in the histograms for each parameter. Optional, there are several kwarg parameters that can be used to customise the plotting aesthetics.
[ ]:
# Corner plot
# The corner plot produces 1,2 and 3 sigma contours for the 2D plots
retrieval.plot_corner(sample_dict, parameter_dict, parameters_read, title_kwargs={"fontsize": 10}, smooth=1)
Making corner plot
Note that all plots, but also additional files (including the best-fit spectrum at high-resolution, so \(\lambda/\Delta\lambda=1000\), regardless of what your model resolution for the retrieval was), will also be put into the evaluate_hst_example_clear_spec folder.
Plot everything at once#
It is possible to plot every figures at once using the plot_all function. However, control on the figures is limited.
[ ]:
retrieval.plot_all(contribution=True)
Best fit parameters
Best fit likelihood = 172.32
log_g 2.74209718736825
planet_radius 1635752360.22101
temperature 852.6696022849248
log_Pcloud 1.8713243922076401
H2O__POKAZATEL -2.4379260627288124
CH4__HITEMP -2.672435452971601
CO2 -5.615758586571792
CO-NatAbund__HITEMP -3.1393823614914034
Plotting Best-fit spectrum
Best fit likelihood = 172.32
Loading Radtrans opacities...
Loading line opacities of species 'H2O__POKAZATEL' from file '/home/dblain/petitRADTRANS/input_data/opacities/lines/correlated_k/H2O/1H2-16O/1H2-16O__POKAZATEL.R1000_0.3-50mu.ktable.petitRADTRANS.h5'... Done.
Loading line opacities of species 'CH4__HITEMP' from file '/home/dblain/petitRADTRANS/input_data/opacities/lines/correlated_k/CH4/12C-1H4/12C-1H4__HITEMP.R1000_0.1-250mu.ktable.petitRADTRANS.h5'... Done.
Loading line opacities of species 'CO2' from file '/home/dblain/petitRADTRANS/input_data/opacities/lines/correlated_k/CO2/12C-16O2/12C-16O2__UCL-4000.R1000_0.3-50mu.ktable.petitRADTRANS.h5'... Done.
Loading line opacities of species 'CO-NatAbund__HITEMP' from file '/home/dblain/petitRADTRANS/input_data/opacities/lines/correlated_k/CO/C-O-NatAbund/C-O-NatAbund__HITEMP.R1000_0.1-250mu.ktable.petitRADTRANS.h5'... Done.
Successfully loaded all line opacities
Loading CIA opacities for H2-H2 from file '/home/dblain/petitRADTRANS/input_data/opacities/continuum/collision_induced_absorptions/H2--H2/H2--H2-NatAbund/H2--H2-NatAbund__BoRi.R831_0.6-250mu.ciatable.petitRADTRANS.h5'... Done.
Loading CIA opacities for H2-He from file '/home/dblain/petitRADTRANS/input_data/opacities/continuum/collision_induced_absorptions/H2--He/H2--He-NatAbund/H2--He-NatAbund__BoRi.DeltaWavenumber2_0.5-500mu.ciatable.petitRADTRANS.h5'... Done.
Successfully loaded all CIA opacities
Successfully loaded all opacities
Best fit 𝛘^2 = 1.31
Best fit 𝛘^2/DoF = 0.11
/mnt/c/Users/doria/OneDrive/Documents/programs/Python/petitRADTRANS/petitRADTRANS/retrieval/retrieval.py:4374: UserWarning: This figure includes Axes that are not compatible with tight_layout, so results might be incorrect.
plt.tight_layout()
Plotting PT profiles
Best fit likelihood = 172.32
Loading Radtrans opacities...
Loading line opacities of species 'H2O__POKAZATEL' from file '/home/dblain/petitRADTRANS/input_data/opacities/lines/correlated_k/H2O/1H2-16O/1H2-16O__POKAZATEL.R1000_0.3-50mu.ktable.petitRADTRANS.h5'... Done.
Loading line opacities of species 'CH4__HITEMP' from file '/home/dblain/petitRADTRANS/input_data/opacities/lines/correlated_k/CH4/12C-1H4/12C-1H4__HITEMP.R1000_0.1-250mu.ktable.petitRADTRANS.h5'... Done.
Loading line opacities of species 'CO2' from file '/home/dblain/petitRADTRANS/input_data/opacities/lines/correlated_k/CO2/12C-16O2/12C-16O2__UCL-4000.R1000_0.3-50mu.ktable.petitRADTRANS.h5'... Done.
Loading line opacities of species 'CO-NatAbund__HITEMP' from file '/home/dblain/petitRADTRANS/input_data/opacities/lines/correlated_k/CO/C-O-NatAbund/C-O-NatAbund__HITEMP.R1000_0.1-250mu.ktable.petitRADTRANS.h5'... Done.
Successfully loaded all line opacities
Loading CIA opacities for H2-H2 from file '/home/dblain/petitRADTRANS/input_data/opacities/continuum/collision_induced_absorptions/H2--H2/H2--H2-NatAbund/H2--H2-NatAbund__BoRi.R831_0.6-250mu.ciatable.petitRADTRANS.h5'... Done.
Loading CIA opacities for H2-He from file '/home/dblain/petitRADTRANS/input_data/opacities/continuum/collision_induced_absorptions/H2--He/H2--He-NatAbund/H2--He-NatAbund__BoRi.DeltaWavenumber2_0.5-500mu.ciatable.petitRADTRANS.h5'... Done.
Successfully loaded all CIA opacities
Successfully loaded all opacities
Making corner plot
WARNING:root:Too few points to create valid contours
WARNING:root:Too few points to create valid contours
WARNING:root:Too few points to create valid contours
WARNING:root:Too few points to create valid contours
WARNING:root:Too few points to create valid contours
WARNING:root:Too few points to create valid contours
WARNING:root:Too few points to create valid contours
WARNING:root:Too few points to create valid contours
WARNING:root:Too few points to create valid contours
Plotting Best-fit contribution function
Best fit likelihood = 172.32
Loading best fit spectrum and contribution from file
Plotting Abundances profiles
Best fit likelihood = 172.32
Loading best fit spectrum and contribution from file
Finished generating all plots!
Contact
If you need any additional help, don’t hesitate to contact Evert Nasedkin.