{ "cells": [ { "cell_type": "markdown", "id": "92509a83", "metadata": {}, "source": [ "(astrometric)=\n", "\n", "# Astrometric fitting" ] }, { "cell_type": "code", "execution_count": null, "id": "1343498d-bf19-420d-bc0e-631c76ea9175", "metadata": {}, "outputs": [], "source": [ "import exoplanet\n", "\n", "exoplanet.utils.docs_setup()\n", "print(f\"exoplanet.__version__ = '{exoplanet.__version__}'\")" ] }, { "cell_type": "markdown", "id": "8a599053", "metadata": {}, "source": [ "In this case study we'll walk through the simplest astrometric example with `exoplanet` and then explain how to build up a more complicated example with parallax measurements. For our dataset, we'll use astrometric and radial velocity observations of a binary star system.\n", "\n", "Astrometric observations usually consist of measurements of the separation and position angle of the secondary star (or directly imaged exoplanet), relative to the primary star as a function of time. The simplest astrometric orbit (in terms of number of parameters), describes the orbit using a semi-major axis `a_ang` measured in *arcseconds*, since the distance to the system is assumed to be unknown. We'll work through this example first, then introduce the extra constraints provided by parallax information.\n", "\n", "## Data\n", "First, let's load and examine the data. We'll use the astrometric measurements of HR 466 (HD 10009) as compiled by [Pourbaix 1998](https://ui.adsabs.harvard.edu/#abs/1998A&AS..131..377P/abstract). The speckle observations are originally from [Hartkopf et al. 1996](https://ui.adsabs.harvard.edu/#abs/1996AJ....111..370H/abstract)." ] }, { "cell_type": "code", "execution_count": null, "id": "c08172cd", "metadata": {}, "outputs": [], "source": [ "from astropy.io import ascii\n", "from astropy.time import Time\n", "\n", "# grab the formatted data and do some munging\n", "dirname = \"https://gist.github.com/iancze/262aba2429cb9aee3fd5b5e1a4582d4d/raw/c5fa5bc39fec90d2cc2e736eed479099e3e598e3/\"\n", "\n", "astro_data_full = ascii.read(\n", " dirname + \"astro.txt\", format=\"csv\", fill_values=[(\".\", \"0\")]\n", ")\n", "\n", "# convert UT date to JD\n", "astro_dates = Time(astro_data_full[\"date\"].data, format=\"decimalyear\")\n", "\n", "# Following the Pourbaix et al. 1998 analysis, we'll limit ourselves to the highest quality data\n", "# since the raw collection of data outside of these ranges has some ambiguities in swapping\n", "# the primary and secondary star\n", "ind = (\n", " (astro_dates.value > 1975.0)\n", " & (astro_dates.value < 1999.73)\n", " & (~astro_data_full[\"rho\"].mask)\n", " & (~astro_data_full[\"PA\"].mask)\n", ")\n", "\n", "astro_data = astro_data_full[ind]\n", "\n", "astro_yrs = astro_data[\"date\"]\n", "astro_dates.format = \"jd\"\n", "astro_jds = astro_dates[ind].value" ] }, { "cell_type": "markdown", "id": "36310478", "metadata": {}, "source": [ "Many of these measurements in this heterogeneous dataset do not have reported error measurements. For these, we assume a modest uncertainty of $1^\\circ$ in position angle and $0.01^{\\prime\\prime}$ in separation for the sake of specifying something, but we'll include a jitter term for both of these measurements as well. The scatter in points around the final solution will be a decent guide of what the measurement uncertainties actually were." ] }, { "cell_type": "code", "execution_count": null, "id": "4438fcad", "metadata": {}, "outputs": [], "source": [ "import numpy as np\n", "\n", "astro_data[\"rho_err\"][astro_data[\"rho_err\"].mask == True] = 0.01\n", "astro_data[\"PA_err\"][astro_data[\"PA_err\"].mask == True] = 1.0\n", "\n", "# Convert all masked frames to be raw np arrays, since theano has issues with astropy masked columns\n", "rho_data = np.ascontiguousarray(astro_data[\"rho\"], dtype=float) # arcsec\n", "rho_err = np.ascontiguousarray(astro_data[\"rho_err\"], dtype=float)\n", "\n", "# The position angle measurements come in degrees in the range [0, 360].\n", "# We'll convert this to radians in the range [-pi, pi]\n", "deg = np.pi / 180.0\n", "theta_data = np.ascontiguousarray(astro_data[\"PA\"] * deg, dtype=float)\n", "theta_data[theta_data > np.pi] -= 2 * np.pi\n", "\n", "theta_err = np.ascontiguousarray(astro_data[\"PA_err\"] * deg) # radians" ] }, { "cell_type": "markdown", "id": "b085f46c", "metadata": {}, "source": [ "## Astrometric conventions\n", "\n", "The conventions describing the orientation of the orbits are described in detail in the *exoplanet* paper; we summarize them briefly here. Generally, we follow the conventions from Pourbaix et al. 1998, which are a consistent set conforming to the right-hand-rule and the conventions of the visual binary field, where the ascending node is that where the secondary is *receeding* from the observer (without radial velocity information, there is a $\\pi$ degeneracy in which node is ascending, and so common practice in the literature is to report a value in the range $[0,\\pi]$). The orbital inclination ranges from $[0, \\pi$]. $i = 0$ describes a face-on orbit rotating counter-clockwise on the sky plane, while $i=\\pi$ describes a face-on orbit rotating clockwise on the sky. $i = \\pi/2$ is an edge-on orbit.\n", "\n", "The observer frame $X$, $Y$, $Z$ is oriented on the sky such that $+Z$ points towards the observer, $X$ is the north axis, and $Y$ is the east axis. *All* angles are measured in radians, and the position angle is returned in the range $[-\\pi, \\pi]$, which is the degrees east of north (be sure to check your data is in this format too!) The radial velocity is still defined such that a positive radial velocity corresponds to motion away from the observer.\n", "\n", "In an astrometric-only orbit, it is common practice in the field to report $\\omega = \\omega_\\mathrm{secondary}$, whereas with an RV orbit it is generally common practice to report $\\omega = \\omega_\\mathrm{primary}$. The result is that unless the authors specify what they're using, in a joint astrometric-RV orbit there is an ambiguity to which $\\omega$ the authors mean, since $\\omega_\\mathrm{primary} = \\omega_\\mathrm{secondary} + \\pi$. To standardize this across the *exoplanet* package, in all orbits (including astrometric-only) $\\omega = \\omega_\\mathrm{primary}$." ] }, { "cell_type": "code", "execution_count": null, "id": "cdd2fe55", "metadata": {}, "outputs": [], "source": [ "import matplotlib.pyplot as plt\n", "\n", "# Make a plot of the astrometric data on the sky\n", "# The convention is that North is up and East is left\n", "fig, ax = plt.subplots(nrows=1, figsize=(4, 4))\n", "\n", "xs = rho_data * np.cos(theta_data) # X is north\n", "ys = rho_data * np.sin(theta_data) # Y is east\n", "ax.plot(ys, xs, \".k\")\n", "ax.set_ylabel(r\"$\\Delta \\delta$ ['']\")\n", "ax.set_xlabel(r\"$\\Delta \\alpha \\cos \\delta$ ['']\")\n", "ax.invert_xaxis()\n", "ax.plot(0, 0, \"k*\")\n", "ax.set_aspect(\"equal\", \"datalim\")" ] }, { "cell_type": "markdown", "id": "4ba61c0e", "metadata": {}, "source": [ "The plot on the sky is helpful to look at, but the \"raw\" measurements are the values of $\\rho$ (separation) and $\\theta$ (also called P.A., position angle) that we listed in our data table, and that the measurement uncertainties live on these values as nice Gaussians. So, to visualize this space more clearly, we can plot $\\rho$ vs. time and P.A. vs. time." ] }, { "cell_type": "code", "execution_count": null, "id": "6b5eb6da", "metadata": {}, "outputs": [], "source": [ "fig, ax = plt.subplots(nrows=2, sharex=True)\n", "ax[0].errorbar(astro_yrs, rho_data, yerr=rho_err, fmt=\".k\", lw=1, ms=5)\n", "ax[0].set_ylabel(r'$\\rho\\,$ [\"]')\n", "\n", "ax[1].errorbar(astro_yrs, theta_data, yerr=theta_err, fmt=\".k\", lw=1, ms=5)\n", "ax[1].set_ylabel(r\"P.A. [radians]\")\n", "_ = ax[1].set_xlabel(\"time [years]\")" ] }, { "cell_type": "markdown", "id": "3e783255", "metadata": {}, "source": [ "## Fitting the astrometric orbit with *exoplanet*\n", "\n", "To get started, let's import the relative packages from *exoplanet*, plot up a preliminary orbit from the literature, and then sample." ] }, { "cell_type": "raw", "id": "ad195029", "metadata": {}, "source": [ ".. note:: Orbits in *exoplanet* generally specify the semi-major axis in units of solar radii `R_sun`. For transits and RV orbits, you usually have enough external information (e.g., estimate of stellar mass from spectral type) to put a physical scale onto the orbit. For the most basic of astrometric orbits without parallax information, however, this information can be lacking and thus it makes sense to fit for the semi-major axis in units of `arcseconds`. But, `exoplanet` is modeling a real orbit (where semi-major axis is in units of `R_sun`), so we do need to at least provide a fake parallax to convert from arcseconds to `R_sun.`" ] }, { "cell_type": "code", "execution_count": null, "id": "292d128d", "metadata": {}, "outputs": [], "source": [ "import pymc3 as pm\n", "import pymc3_ext as pmx\n", "\n", "from aesara_theano_fallback import aesara as theano\n", "import aesara_theano_fallback.tensor as tt\n", "\n", "import exoplanet as xo\n", "\n", "from astropy import constants\n", "\n", "# conversion constant from au to R_sun\n", "au_to_R_sun = (constants.au / constants.R_sun).value\n", "\n", "# Just to get started, let's take a look at the orbit using the parameter estimates from Pourbaix et al. 1998\n", "\n", "# Orbital elements from Pourbaix et al. 1998\n", "# For the relative astrometric fit, we only need the following parameters\n", "a_ang = 0.324 # arcsec\n", "parallax = 1 # arcsec (meaningless choice for now)\n", "a = a_ang * au_to_R_sun / parallax\n", "e = 0.798\n", "i = 96.0 * deg # [rad]\n", "omega = 251.6 * deg - np.pi # Pourbaix reports omega_2, but we want omega_1\n", "Omega = 159.6 * deg\n", "P = 28.8 * 365.25 # days\n", "\n", "T0 = Time(1989.92, format=\"decimalyear\")\n", "T0.format = \"jd\"\n", "T0 = T0.value # [Julian Date]\n", "\n", "# instantiate the orbit\n", "orbit = xo.orbits.KeplerianOrbit(\n", " a=a, t_periastron=T0, period=P, incl=i, ecc=e, omega=omega, Omega=Omega\n", ")\n", "\n", "# The position functions take an optional argument parallax to convert from\n", "# physical units back to arcseconds\n", "t = np.linspace(T0 - P, T0 + P, num=200) # days\n", "rho, theta = theano.function([], orbit.get_relative_angles(t, parallax))()\n", "\n", "# Plot the orbit\n", "fig, ax = plt.subplots(nrows=1, figsize=(4, 4))\n", "\n", "xs = rho * np.cos(theta) # X is north\n", "ys = rho * np.sin(theta) # Y is east\n", "ax.plot(ys, xs, color=\"C0\", lw=1)\n", "\n", "# plot the data\n", "xs = rho_data * np.cos(theta_data) # X is north\n", "ys = rho_data * np.sin(theta_data) # Y is east\n", "ax.plot(ys, xs, \".k\")\n", "\n", "ax.set_ylabel(r\"$\\Delta \\delta$ ['']\")\n", "ax.set_xlabel(r\"$\\Delta \\alpha \\cos \\delta$ ['']\")\n", "ax.invert_xaxis()\n", "ax.plot(0, 0, \"k*\")\n", "ax.set_aspect(\"equal\", \"datalim\")\n", "ax.set_title(\"initial orbit\")\n", "\n", "fig, ax = plt.subplots(nrows=2, sharex=True, figsize=(6, 6))\n", "ax[0].errorbar(astro_jds, rho_data, yerr=rho_err, fmt=\".k\", lw=1, ms=5)\n", "ax[0].plot(t, rho, color=\"C0\", lw=1)\n", "ax[0].set_ylabel(r'$\\rho\\,$ [\"]')\n", "ax[0].set_title(\"initial orbit\")\n", "\n", "ax[1].errorbar(astro_jds, theta_data, yerr=theta_err, fmt=\".k\", lw=1, ms=5)\n", "ax[1].plot(t, theta, color=\"C0\", lw=1)\n", "ax[1].set_ylabel(r\"P.A. [radians]\")\n", "_ = ax[1].set_xlabel(\"time [JD]\")" ] }, { "cell_type": "markdown", "id": "ecd7ec33", "metadata": {}, "source": [ "Now that we have an initial orbit, we can set the model up using PyMC3 to do inference." ] }, { "cell_type": "code", "execution_count": null, "id": "5c4fdff3", "metadata": {}, "outputs": [], "source": [ "yr = 365.25\n", "\n", "# for predicted orbits\n", "t_fine = np.linspace(astro_jds.min() - 500, astro_jds.max() + 500, num=1000)\n", "\n", "\n", "def get_model(parallax=None):\n", " with pm.Model() as model:\n", " if parallax is None:\n", " # Without an actual parallax measurement, we can model the orbit in units of arcseconds\n", " # by providing a fake_parallax and conversion constant\n", " plx = 1 # arcsec\n", " else:\n", " # Below we will run a version of this model where a measurement of parallax is provided\n", " # The measurement is in milliarcsec\n", " m_plx = pm.Bound(pm.Normal, lower=0, upper=100)(\n", " \"m_plx\", mu=parallax[0], sd=parallax[1], testval=parallax[0]\n", " )\n", " plx = pm.Deterministic(\"plx\", 1e-3 * m_plx)\n", "\n", " a_ang = pm.Uniform(\"a_ang\", 0.1, 1.0, testval=0.324)\n", " a = pm.Deterministic(\"a\", a_ang / plx)\n", "\n", " # We expect the period to be somewhere in the range of 25 years,\n", " # so we'll set a broad prior on logP\n", " logP = pm.Normal(\n", " \"logP\", mu=np.log(25 * yr), sd=10.0, testval=np.log(28.8 * yr)\n", " )\n", " P = pm.Deterministic(\"P\", tt.exp(logP))\n", "\n", " # For astrometric-only fits, it's generally better to fit in\n", " # p = (Omega + omega)/2 and m = (Omega - omega)/2 instead of omega and Omega\n", " # directly\n", " omega0 = 251.6 * deg - np.pi\n", " Omega0 = 159.6 * deg\n", " p = pmx.Angle(\"p\", testval=0.5 * (Omega0 + omega0))\n", " m = pmx.Angle(\"m\", testval=0.5 * (Omega0 - omega0))\n", " omega = pm.Deterministic(\"omega\", p - m)\n", " Omega = pm.Deterministic(\"Omega\", p + m)\n", "\n", " # For these orbits, it can also be better to fit for a phase angle\n", " # (relative to a reference time) instead of the time of periastron\n", " # passage directly\n", " phase = pmx.Angle(\"phase\", testval=0.0)\n", " tperi = pm.Deterministic(\"tperi\", T0 + P * phase / (2 * np.pi))\n", "\n", " # Geometric uniform prior on cos(incl)\n", " cos_incl = pm.Uniform(\n", " \"cos_incl\", lower=-1, upper=1, testval=np.cos(96.0 * deg)\n", " )\n", " incl = pm.Deterministic(\"incl\", tt.arccos(cos_incl))\n", " ecc = pm.Uniform(\"ecc\", lower=0.0, upper=1.0, testval=0.798)\n", "\n", " # Set up the orbit\n", " orbit = xo.orbits.KeplerianOrbit(\n", " a=a * au_to_R_sun,\n", " t_periastron=tperi,\n", " period=P,\n", " incl=incl,\n", " ecc=ecc,\n", " omega=omega,\n", " Omega=Omega,\n", " )\n", " if parallax is not None:\n", " pm.Deterministic(\"M_tot\", orbit.m_total)\n", "\n", " # Compute the model in rho and theta\n", " rho_model, theta_model = orbit.get_relative_angles(astro_jds, plx)\n", " pm.Deterministic(\"rho_model\", rho_model)\n", " pm.Deterministic(\"theta_model\", theta_model)\n", "\n", " # Add jitter terms to both separation and position angle\n", " log_rho_s = pm.Normal(\n", " \"log_rho_s\", mu=np.log(np.median(rho_err)), sd=2.0\n", " )\n", " log_theta_s = pm.Normal(\n", " \"log_theta_s\", mu=np.log(np.median(theta_err)), sd=2.0\n", " )\n", " rho_tot_err = tt.sqrt(rho_err**2 + tt.exp(2 * log_rho_s))\n", " theta_tot_err = tt.sqrt(theta_err**2 + tt.exp(2 * log_theta_s))\n", "\n", " # define the likelihood function, e.g., a Gaussian on both rho and theta\n", " pm.Normal(\"rho_obs\", mu=rho_model, sd=rho_tot_err, observed=rho_data)\n", "\n", " # We want to be cognizant of the fact that theta wraps so the following is equivalent to\n", " # pm.Normal(\"obs_theta\", mu=theta_model, observed=theta_data, sd=theta_tot_err)\n", " # but takes into account the wrapping. Thanks to Rob de Rosa for the tip.\n", " theta_diff = tt.arctan2(\n", " tt.sin(theta_model - theta_data), tt.cos(theta_model - theta_data)\n", " )\n", " pm.Normal(\"theta_obs\", mu=theta_diff, sd=theta_tot_err, observed=0.0)\n", "\n", " # Set up predicted orbits for later plotting\n", " rho_dense, theta_dense = orbit.get_relative_angles(t_fine, plx)\n", " rho_save = pm.Deterministic(\"rho_save\", rho_dense)\n", " theta_save = pm.Deterministic(\"theta_save\", theta_dense)\n", "\n", " # Optimize to find the initial parameters\n", " map_soln = model.test_point\n", " map_soln = pmx.optimize(map_soln, vars=[log_rho_s, log_theta_s])\n", " map_soln = pmx.optimize(map_soln, vars=[phase])\n", " map_soln = pmx.optimize(map_soln, vars=[p, m, ecc])\n", " map_soln = pmx.optimize(map_soln, vars=[logP, a_ang, phase])\n", " map_soln = pmx.optimize(map_soln)\n", "\n", " return model, map_soln\n", "\n", "\n", "model, map_soln = get_model()" ] }, { "cell_type": "markdown", "id": "dfddbce6", "metadata": {}, "source": [ "Now that we have a maximum a posteriori estimate of the parameters, let's take a look at the results to make sure that they seem reasonable." ] }, { "cell_type": "code", "execution_count": null, "id": "434ee0d1", "metadata": {}, "outputs": [], "source": [ "ekw = dict(fmt=\".k\", lw=0.5)\n", "\n", "fig, ax = plt.subplots(nrows=4, sharex=True, figsize=(6, 8))\n", "ax[0].set_ylabel(r'$\\rho\\,$ [\"]')\n", "ax[1].set_ylabel(r\"$\\rho$ residuals\")\n", "ax[2].set_ylabel(r\"P.A. [radians]\")\n", "ax[3].set_ylabel(r\"P.A. residuals\")\n", "\n", "tot_rho_err = np.sqrt(rho_err**2 + np.exp(2 * map_soln[\"log_rho_s\"]))\n", "tot_theta_err = np.sqrt(theta_err**2 + np.exp(2 * map_soln[\"log_theta_s\"]))\n", "\n", "ax[0].errorbar(astro_jds, rho_data, yerr=tot_rho_err, **ekw)\n", "ax[0].plot(t_fine, map_soln[\"rho_save\"], \"C1\")\n", "\n", "ax[1].axhline(0.0, color=\"0.5\")\n", "ax[1].errorbar(\n", " astro_jds, rho_data - map_soln[\"rho_model\"], yerr=tot_rho_err, **ekw\n", ")\n", "\n", "\n", "ax[2].plot(t_fine, map_soln[\"theta_save\"], \"C1\")\n", "ax[2].errorbar(astro_jds, theta_data, yerr=tot_theta_err, **ekw)\n", "\n", "ax[3].axhline(0.0, color=\"0.5\")\n", "ax[3].errorbar(\n", " astro_jds, theta_data - map_soln[\"theta_model\"], yerr=tot_theta_err, **ekw\n", ")\n", "\n", "ax[3].set_xlim(t_fine[0], t_fine[-1])\n", "_ = ax[0].set_title(\"map orbit\")" ] }, { "cell_type": "markdown", "id": "9d74dc5b", "metadata": {}, "source": [ "Now let's sample the posterior." ] }, { "cell_type": "code", "execution_count": null, "id": "b4a6779c", "metadata": {}, "outputs": [], "source": [ "np.random.seed(1234)\n", "with model:\n", " trace = pmx.sample(\n", " tune=1000,\n", " draws=1000,\n", " start=map_soln,\n", " cores=2,\n", " chains=2,\n", " target_accept=0.9,\n", " return_inferencedata=True,\n", " )" ] }, { "cell_type": "markdown", "id": "a5e747e3", "metadata": {}, "source": [ "First we can check the convergence for some of the key parameters." ] }, { "cell_type": "code", "execution_count": null, "id": "2f665945", "metadata": {}, "outputs": [], "source": [ "import arviz as az\n", "\n", "az.summary(\n", " trace,\n", " var_names=[\"P\", \"tperi\", \"a_ang\", \"omega\", \"Omega\", \"incl\", \"ecc\"],\n", ")" ] }, { "cell_type": "markdown", "id": "c3712f75", "metadata": {}, "source": [ "That looks pretty good.\n", "Now here's a corner plot showing the covariances between parameters." ] }, { "cell_type": "code", "execution_count": null, "id": "46980c98", "metadata": {}, "outputs": [], "source": [ "import corner\n", "\n", "_ = corner.corner(\n", " trace, var_names=[\"P\", \"tperi\", \"a_ang\", \"omega\", \"Omega\", \"incl\", \"ecc\"]\n", ")" ] }, { "cell_type": "markdown", "id": "dee00e74", "metadata": {}, "source": [ "Finally, we can plot the posterior constraints on $\\rho$ and $\\theta$.\n", "This figure is much like the one for the MAP solution above, but this time the orange is a contour (not a line) showing the 68% credible region for the model." ] }, { "cell_type": "code", "execution_count": null, "id": "049df130", "metadata": {}, "outputs": [], "source": [ "ekw = dict(fmt=\".k\", lw=0.5)\n", "\n", "fig, ax = plt.subplots(nrows=2, sharex=True, figsize=(6, 6))\n", "ax[0].set_ylabel(r'$\\rho\\,$ [\"]')\n", "ax[1].set_ylabel(r\"P.A. [radians]\")\n", "\n", "tot_rho_err = np.sqrt(\n", " rho_err**2\n", " + np.exp(2 * np.median(trace.posterior[\"log_rho_s\"].values, axis=(0, 1)))\n", ")\n", "tot_theta_err = np.sqrt(\n", " theta_err**2\n", " + np.exp(2 * np.median(trace.posterior[\"log_theta_s\"].values, axis=(0, 1)))\n", ")\n", "\n", "ax[0].errorbar(astro_jds, rho_data, yerr=tot_rho_err, **ekw)\n", "q = np.percentile(trace.posterior[\"rho_save\"].values, [16, 84], axis=(0, 1))\n", "ax[0].fill_between(t_fine, q[0], q[1], color=\"C1\", alpha=0.8, lw=0)\n", "\n", "ax[1].errorbar(astro_jds, theta_data, yerr=tot_theta_err, **ekw)\n", "q = np.percentile(trace.posterior[\"theta_save\"].values, [16, 84], axis=(0, 1))\n", "ax[1].fill_between(t_fine, q[0], q[1], color=\"C1\", alpha=0.8, lw=0)\n", "\n", "ax[-1].set_xlim(t_fine[0], t_fine[-1])\n", "_ = ax[0].set_title(\"posterior inferences\")" ] }, { "cell_type": "markdown", "id": "e2decdf6", "metadata": {}, "source": [ "As we can see from the narrow range of orbits (the orange swath appears like a thin line), the orbit is actually highly constrained by the astrometry.\n", "We also see two outlier epochs in the vicinity of 2445000 - 2447000, since adjacent epochs seem to be right on the orbit.\n", "It's likely the uncertainties were not estimated correctly for these, and the simplistic jitter model we implemented isn't sophisticated enough to apply more weight to only these discrepant points." ] }, { "cell_type": "markdown", "id": "82791242", "metadata": {}, "source": [ "## Including parallax\n", "\n", "While this is encouraging that we fit an astrometric orbit, a simple astrometric fit to just $\\rho$ and $\\theta$ isn't actually that physically satisfying, since many of the orbital parameters simply have to do with the orientation relative to us ($i$, $\\omega$, and $\\Omega$). The only truly intrinsic parameters are $P$ and $e$. To learn more about some of the physical parameters, such as the total mass of the system, we'd like to incorporate distance information to put a physical scale to the problem.\n", "\n", "The *Gaia* DR2 parallax is $\\varpi = 24.05 \\pm 0.45$ mas.\n", "\n", "We can use exactly the same model as above with only an added parallax constraint:" ] }, { "cell_type": "code", "execution_count": null, "id": "30088e53", "metadata": {}, "outputs": [], "source": [ "plx_model, plx_map_soln = get_model(parallax=[24.05, 0.45])" ] }, { "cell_type": "code", "execution_count": null, "id": "5a634a4a", "metadata": {}, "outputs": [], "source": [ "np.random.seed(5432)\n", "with plx_model:\n", " plx_trace = pmx.sample(\n", " tune=1000,\n", " draws=1000,\n", " start=plx_map_soln,\n", " cores=2,\n", " chains=2,\n", " target_accept=0.9,\n", " return_inferencedata=True,\n", " )" ] }, { "cell_type": "markdown", "id": "1c1107f4", "metadata": {}, "source": [ "Check the convergence diagnostics." ] }, { "cell_type": "code", "execution_count": null, "id": "500ca9c6", "metadata": {}, "outputs": [], "source": [ "az.summary(\n", " plx_trace,\n", " var_names=[\n", " \"P\",\n", " \"tperi\",\n", " \"a_ang\",\n", " \"omega\",\n", " \"Omega\",\n", " \"incl\",\n", " \"ecc\",\n", " \"M_tot\",\n", " \"plx\",\n", " ],\n", ")" ] }, { "cell_type": "markdown", "id": "40ae60b6", "metadata": {}, "source": [ "And make the corner plot for the physical parameters." ] }, { "cell_type": "code", "execution_count": null, "id": "c7a478ad", "metadata": {}, "outputs": [], "source": [ "_ = corner.corner(plx_trace, var_names=[\"P\", \"tperi\", \"a\", \"ecc\", \"M_tot\"])" ] }, { "cell_type": "markdown", "id": "c0445f0f", "metadata": {}, "source": [ "## Citations\n", "\n", "As described in the [citation tutorial](https://docs.exoplanet.codes/en/stable/tutorials/citation/), we can use [citations.get_citations_for_model](https://docs.exoplanet.codes/en/stable/user/api/#exoplanet.citations.get_citations_for_model) to construct an acknowledgement and BibTeX listing that includes the relevant citations for this model." ] }, { "cell_type": "code", "execution_count": null, "id": "dd90ab44", "metadata": {}, "outputs": [], "source": [ "with model:\n", " txt, bib = xo.citations.get_citations_for_model()\n", "print(txt)" ] }, { "cell_type": "code", "execution_count": null, "id": "5c903782", "metadata": {}, "outputs": [], "source": [ "print(bib.split(\"\\n\\n\")[0] + \"\\n\\n...\")" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.9.4" } }, "nbformat": 4, "nbformat_minor": 5 }