{ "cells": [ { "cell_type": "markdown", "id": "2666a73a-1300-4d7e-a5bf-fd3b54e221b7", "metadata": {}, "source": [ "(lc-multi)=\n", "\n", "# Light curves from multiple instruments" ] }, { "cell_type": "code", "execution_count": null, "id": "be2f2c90-ef8c-4bc7-a0bb-47c6549f0bdf", "metadata": {}, "outputs": [], "source": [ "import exoplanet\n", "\n", "exoplanet.utils.docs_setup()\n", "print(f\"exoplanet.__version__ = '{exoplanet.__version__}'\")" ] }, { "cell_type": "markdown", "id": "a4a5f70c-f401-4648-858b-434704ecff6f", "metadata": {}, "source": [ "In the {ref}`rv-multi` case study, we discussed fitting the radial velocity curve for a planetary system observed using multiple instruments.\n", "You might also want to fit data from multiple instruments when fitting the light curve of a transiting planet and that's what we work through in this example.\n", "This is a somewhat more complicated example than the radial velocity case because some of the physical properties of the system can vary as as function of the instrument.\n", "Specifically, the transit depth (or the effective radius of the planet) will be a function of the filter or effective wavelength of the observations.\n", "This is the idea behind transit spectroscopy and the method used in this case study could (and should!) be extended to that use case.\n", "In this case, we'll combine the light curves from the Kepler and TESS missions for the planet host HAT-P-11.\n", "\n", "## A brief aside on dataset \"weighting\"\n", "\n", "Before getting into the details of this case study, let's spend a minute talking about a topic that comes up a lot when discussing combining observations from different instruments or techniques.\n", "To many people, it seems intuitive that one should (and perhaps must) \"weight\" how much each dataset contributes to the likelihood based on how much they \"trust\" those data.\n", "For example, you might be worried that a dataset with more datapoints will have a larger effect on the the results than you would like.\n", "While this might seem intuitive, it's wrong: **the only way to combine datasets is to multiply their likelihood functions**.\n", "Instead, it is useful to understand what you actually mean when you say that you don't \"trust\" a dataset as much as another.\n", "**What you're really saying is that you don't believe the observation model that you wrote down**.\n", "For example, you might think that the quoted error bars are underestimated or there might be correlated noise that an uncorrelated normal observation model can't capture.\n", "The benefit of thinking about it this way is that it suggests a solution to the problem: incorporate a more flexible observation model that can capture these issues.\n", "In this case study, the 4 years of (long-cadence) Kepler observations only include about two times as many data points as one month of TESS observations.\n", "But, as you can see in the figure below, these two datasets have different noise properties (both in terms of photon noise and correlated noise) so we will fit using a different flexible Gaussian process noise model for each data set that will take these different properties into account." ] }, { "cell_type": "code", "execution_count": null, "id": "efef00be", "metadata": {}, "outputs": [], "source": [ "import numpy as np\n", "import lightkurve as lk\n", "from astropy.io import fits\n", "from collections import OrderedDict\n", "import matplotlib.pyplot as plt\n", "\n", "# Period and reference transit time from the literature for initialization\n", "lit_period = 4.887803076\n", "lit_t0 = 124.8130808\n", "\n", "\n", "def transit_mask(t):\n", " return (\n", " np.abs((t - lit_t0 + 0.5 * lit_period) % lit_period - 0.5 * lit_period)\n", " < 0.2\n", " )\n", "\n", "\n", "kepler_lcfs = lk.search_lightcurve(\n", " \"HAT-P-11\", mission=\"Kepler\", cadence=\"long\"\n", ").download_all(flux_column=\"pdcsap_flux\")\n", "kepler_lc = kepler_lcfs.stitch().remove_nans()\n", "\n", "# For the purposes of this example, let's only clip out data near the transits\n", "m = transit_mask(kepler_lc.time.value)\n", "kepler_t = np.ascontiguousarray(kepler_lc.time.value[m], dtype=np.float64)\n", "kepler_y = np.ascontiguousarray(\n", " 1e3 * (kepler_lc.flux[m] - 1), dtype=np.float64\n", ")\n", "kepler_yerr = np.ascontiguousarray(\n", " 1e3 * kepler_lc.flux_err[m], dtype=np.float64\n", ")\n", "\n", "with fits.open(kepler_lcfs[0].filename) as hdu:\n", " hdr = hdu[1].header\n", "kepler_texp = hdr[\"FRAMETIM\"] * hdr[\"NUM_FRM\"]\n", "kepler_texp /= 60.0 * 60.0 * 24.0\n", "\n", "tess_lcfs = lk.search_lightcurve(\n", " \"HAT-P-11\", mission=\"TESS\", author=\"SPOC\", sector=(14, 15)\n", ").download_all(flux_column=\"pdcsap_flux\")\n", "tess_lc = tess_lcfs.stitch().remove_nans()\n", "\n", "# For the purposes of this example, let's only clip out data near the transits\n", "t = tess_lc.time.value + 2457000 - 2454833\n", "m = transit_mask(t)\n", "tess_t = np.ascontiguousarray(t[m], dtype=np.float64)\n", "tess_y = np.ascontiguousarray(1e3 * (tess_lc.flux[m] - 1), dtype=np.float64)\n", "tess_yerr = np.ascontiguousarray(1e3 * tess_lc.flux_err[m], dtype=np.float64)\n", "\n", "with fits.open(tess_lcfs[0].filename) as hdu:\n", " hdr = hdu[1].header\n", "tess_texp = hdr[\"FRAMETIM\"] * hdr[\"NUM_FRM\"]\n", "tess_texp /= 60.0 * 60.0 * 24.0\n", "\n", "datasets = OrderedDict(\n", " [\n", " (\"Kepler\", [kepler_t, kepler_y, kepler_yerr, kepler_texp]),\n", " (\"TESS\", [tess_t, tess_y, tess_yerr, tess_texp]),\n", " ]\n", ")\n", "\n", "fig, axes = plt.subplots(1, len(datasets), sharey=True, figsize=(10, 5))\n", "\n", "for i, (name, (t, y, _, _)) in enumerate(datasets.items()):\n", " ax = axes[i]\n", " ax.plot(t, y, \"k\", lw=0.75, label=name)\n", " ax.set_xlabel(\"time [KBJD]\")\n", " ax.set_title(name, fontsize=14)\n", "\n", " x_mid = 0.5 * (t.min() + t.max())\n", " ax.set_xlim(x_mid - 10, x_mid + 10)\n", "axes[0].set_ylim(-10, 10)\n", "fig.subplots_adjust(wspace=0.05)\n", "_ = axes[0].set_ylabel(\"relative flux [ppt]\")" ] }, { "cell_type": "markdown", "id": "0f97992e", "metadata": {}, "source": [ "## The probabilistic model\n", "\n", "This model is mostly the same as the one used in [Quick fits for TESS light curves](./quick-tess.ipynb), but we're allowing for different noise variances (both the white noise component and the GP amplitude), effective planet radii, and limb-darkening coefficients for each dataset.\n", "For the purposes of demonstration, we're sharing the length scale of the GP between the two datasets, but this could just have well been a different parameter for each dataset without changing the results.\n", "The final change that we're using is to use the approximate transit depth `approx_depth` (the depth of the transit at minimum assuming the limb-darkening profile is constant under the disk of the planet) as a parameter instead of the radius ratio.\n", "This does not have a large effect on the performance or the results, but it can sometimes be a useful parameterization when dealing with high signal-to-noise transits because it reduces the covariance between the radius parameter and the limb darkening coefficients.\n", "As usual, we run a few iterations of sigma clipping and then find the maximum a posteriori parameters to check to make sure that everything is working:" ] }, { "cell_type": "code", "execution_count": null, "id": "1d6498d1", "metadata": {}, "outputs": [], "source": [ "import pymc3 as pm\n", "import pymc3_ext as pmx\n", "import exoplanet as xo\n", "import aesara_theano_fallback.tensor as tt\n", "from functools import partial\n", "from celerite2.theano import terms, GaussianProcess\n", "\n", "# Find a reference transit time near the middle of the observations to avoid\n", "# strong covariances between period and t0\n", "x_min = min(np.min(x) for x, _, _, _ in datasets.values())\n", "x_max = max(np.max(x) for x, _, _, _ in datasets.values())\n", "x_mid = 0.5 * (x_min + x_max)\n", "t0_ref = lit_t0 + lit_period * np.round((x_mid - lit_t0) / lit_period)\n", "\n", "# Do several rounds of sigma clipping\n", "for i in range(10):\n", " with pm.Model() as model:\n", " # Shared orbital parameters\n", " log_period = pm.Normal(\"log_period\", mu=np.log(lit_period), sigma=1.0)\n", " period = pm.Deterministic(\"period\", tt.exp(log_period))\n", " t0 = pm.Normal(\"t0\", mu=t0_ref, sigma=1.0)\n", " log_dur = pm.Normal(\"log_dur\", mu=np.log(0.1), sigma=10.0)\n", " dur = pm.Deterministic(\"dur\", tt.exp(log_dur))\n", " b = pmx.UnitUniform(\"b\")\n", " ld_arg = 1 - tt.sqrt(1 - b**2)\n", " orbit = xo.orbits.KeplerianOrbit(\n", " period=period, duration=dur, t0=t0, b=b\n", " )\n", "\n", " # We'll also say that the timescale of the GP will be shared\n", " rho_gp = pm.InverseGamma(\n", " \"rho_gp\",\n", " testval=2.0,\n", " **pmx.estimate_inverse_gamma_parameters(1.0, 5.0),\n", " )\n", "\n", " # Loop over the instruments\n", " parameters = dict()\n", " lc_models = dict()\n", " gp_preds = dict()\n", " gp_preds_with_mean = dict()\n", " for n, (name, (x, y, yerr, texp)) in enumerate(datasets.items()):\n", " # We define the per-instrument parameters in a submodel so that we\n", " # don't have to prefix the names manually\n", " with pm.Model(name=name, model=model):\n", " # The flux zero point\n", " mean = pm.Normal(\"mean\", mu=0.0, sigma=10.0)\n", "\n", " # The limb darkening\n", " u = xo.QuadLimbDark(\"u\")\n", " star = xo.LimbDarkLightCurve(u)\n", "\n", " # The radius ratio\n", " log_approx_depth = pm.Normal(\n", " \"log_approx_depth\", mu=np.log(4e-3), sigma=10\n", " )\n", " ld = 1 - u[0] * ld_arg - u[1] * ld_arg**2\n", " ror = pm.Deterministic(\n", " \"ror\", tt.exp(0.5 * log_approx_depth) / tt.sqrt(ld)\n", " )\n", "\n", " # Noise parameters\n", " med_yerr = np.median(yerr)\n", " std = np.std(y)\n", " sigma = pm.InverseGamma(\n", " \"sigma\",\n", " testval=med_yerr,\n", " **pmx.estimate_inverse_gamma_parameters(\n", " med_yerr, 0.5 * std\n", " ),\n", " )\n", " sigma_gp = pm.InverseGamma(\n", " \"sigma_gp\",\n", " testval=0.5 * std,\n", " **pmx.estimate_inverse_gamma_parameters(\n", " med_yerr, 0.5 * std\n", " ),\n", " )\n", "\n", " # Keep track of the parameters for optimization\n", " parameters[name] = [mean, u, log_approx_depth]\n", " parameters[f\"{name}_noise\"] = [sigma, sigma_gp]\n", "\n", " # The light curve model\n", " def lc_model(mean, star, ror, texp, t):\n", " return mean + 1e3 * tt.sum(\n", " star.get_light_curve(orbit=orbit, r=ror, t=t, texp=texp),\n", " axis=-1,\n", " )\n", "\n", " lc_model = partial(lc_model, mean, star, ror, texp)\n", " lc_models[name] = lc_model\n", "\n", " # The Gaussian Process noise model\n", " kernel = terms.SHOTerm(sigma=sigma_gp, rho=rho_gp, Q=1.0 / 3)\n", " gp = GaussianProcess(\n", " kernel, t=x, diag=yerr**2 + sigma**2, mean=lc_model\n", " )\n", " gp.marginal(f\"{name}_obs\", observed=y)\n", " gp_preds[name] = gp.predict(y, include_mean=False)\n", " gp_preds_with_mean[name] = gp_preds[name] + gp.mean_value\n", "\n", " # Optimize the model\n", " map_soln = model.test_point\n", " for name in datasets:\n", " map_soln = pmx.optimize(map_soln, parameters[name])\n", " for name in datasets:\n", " map_soln = pmx.optimize(map_soln, parameters[f\"{name}_noise\"])\n", " map_soln = pmx.optimize(map_soln, parameters[name] + [log_dur, b])\n", " map_soln = pmx.optimize(map_soln)\n", "\n", " # Do some sigma clipping\n", " num = dict((name, len(datasets[name][0])) for name in datasets)\n", " clipped = dict()\n", " masks = dict()\n", " for name in datasets:\n", " mdl = pmx.eval_in_model(gp_preds_with_mean[name], map_soln)\n", " resid = datasets[name][1] - mdl\n", " sigma = np.sqrt(np.median((resid - np.median(resid)) ** 2))\n", " masks[name] = np.abs(resid - np.median(resid)) < 7 * sigma\n", " clipped[name] = num[name] - masks[name].sum()\n", " print(f\"Sigma clipped {clipped[name]} {name} light curve points\")\n", "\n", " if all(c < 10 for c in clipped.values()):\n", " break\n", "\n", " else:\n", " for name in datasets:\n", " datasets[name][0] = datasets[name][0][masks[name]]\n", " datasets[name][1] = datasets[name][1][masks[name]]\n", " datasets[name][2] = datasets[name][2][masks[name]]" ] }, { "cell_type": "markdown", "id": "faf98001", "metadata": {}, "source": [ "Here are the two phased light curves (with the Gaussian process model removed).\n", "We can see the effect of exposure time integration and the difference in photometric precision, but everything should be looking good!" ] }, { "cell_type": "code", "execution_count": null, "id": "1d32bb53", "metadata": {}, "outputs": [], "source": [ "dt = np.linspace(-0.2, 0.2, 500)\n", "\n", "with model:\n", " trends = pmx.eval_in_model([gp_preds[k] for k in datasets], map_soln)\n", " phase_curves = pmx.eval_in_model(\n", " [lc_models[k](t0 + dt) for k in datasets], map_soln\n", " )\n", "\n", "fig, axes = plt.subplots(2, sharex=True, sharey=True, figsize=(8, 6))\n", "\n", "for n, name in enumerate(datasets):\n", " ax = axes[n]\n", "\n", " x, y = datasets[name][:2]\n", "\n", " period = map_soln[\"period\"]\n", " folded = (x - map_soln[\"t0\"] + 0.5 * period) % period - 0.5 * period\n", " m = np.abs(folded) < 0.2\n", " ax.plot(\n", " folded[m],\n", " (y - trends[n] - map_soln[f\"{name}_mean\"])[m],\n", " \".k\",\n", " alpha=0.3,\n", " mec=\"none\",\n", " )\n", " ax.plot(\n", " dt, phase_curves[n] - map_soln[f\"{name}_mean\"], f\"C{n}\", label=name\n", " )\n", " ax.annotate(\n", " name,\n", " xy=(1, 0),\n", " xycoords=\"axes fraction\",\n", " va=\"bottom\",\n", " ha=\"right\",\n", " xytext=(-3, 3),\n", " textcoords=\"offset points\",\n", " fontsize=14,\n", " )\n", "\n", "axes[-1].set_xlim(-0.15, 0.15)\n", "axes[-1].set_xlabel(\"time since transit [days]\")\n", "for ax in axes:\n", " ax.set_ylabel(\"relative flux [ppt]\")" ] }, { "cell_type": "markdown", "id": "e1691f03", "metadata": {}, "source": [ "Then we run the MCMC:" ] }, { "cell_type": "code", "execution_count": null, "id": "b05c7cea", "metadata": {}, "outputs": [], "source": [ "import platform\n", "\n", "with model:\n", " trace = pmx.sample(\n", " tune=1500,\n", " draws=1500,\n", " start=map_soln,\n", " # Parallel sampling runs poorly or crashes on macos\n", " cores=1 if platform.system() == \"Darwin\" else 2,\n", " chains=2,\n", " target_accept=0.95,\n", " return_inferencedata=True,\n", " random_seed=[11, 9],\n", " )" ] }, { "cell_type": "markdown", "id": "ca22b57e", "metadata": {}, "source": [ "And check the convergence diagnostics:" ] }, { "cell_type": "code", "execution_count": null, "id": "d1cbc411", "metadata": {}, "outputs": [], "source": [ "import arviz as az\n", "\n", "az.summary(trace)" ] }, { "cell_type": "markdown", "id": "f4c2fd5f", "metadata": {}, "source": [ "Since we fit for a radius ratio in each band, we can see if the transit depth is different in Kepler compared to TESS.\n", "The plot below demonstrates that there is no statistically significant difference between the radii measured in these two bands:" ] }, { "cell_type": "code", "execution_count": null, "id": "ba50274e", "metadata": {}, "outputs": [], "source": [ "plt.hist(\n", " trace.posterior[\"Kepler_ror\"].values.flatten(),\n", " 20,\n", " density=True,\n", " histtype=\"step\",\n", " label=\"Kepler\",\n", ")\n", "plt.hist(\n", " trace.posterior[\"TESS_ror\"].values.flatten(),\n", " 20,\n", " density=True,\n", " histtype=\"step\",\n", " label=\"TESS\",\n", ")\n", "plt.yticks([])\n", "plt.xlabel(\"effective radius ratio\")\n", "_ = plt.legend(fontsize=12)" ] }, { "cell_type": "markdown", "id": "54bd6aab", "metadata": {}, "source": [ "We can also compare the inferred limb-darkening coefficients:" ] }, { "cell_type": "code", "execution_count": null, "id": "bb21e376", "metadata": {}, "outputs": [], "source": [ "import corner\n", "\n", "fig = corner.corner(\n", " trace,\n", " var_names=[\"TESS_u\"],\n", " bins=40,\n", " color=\"C1\",\n", " range=((0.5, 0.9), (-0.5, 0.1)),\n", ")\n", "corner.corner(\n", " trace,\n", " var_names=[\"Kepler_u\"],\n", " bins=40,\n", " color=\"C0\",\n", " fig=fig,\n", " labels=[\"$u_1$\", \"$u_2$\"],\n", " range=((0.5, 0.9), (-0.5, 0.1)),\n", ")\n", "fig.axes[0].axvline(-1.0, color=\"C0\", label=\"Kepler\")\n", "fig.axes[0].axvline(-1.0, color=\"C1\", label=\"TESS\")\n", "_ = fig.axes[0].legend(\n", " fontsize=12, loc=\"center left\", bbox_to_anchor=(1.1, 0.5)\n", ")" ] }, { "cell_type": "markdown", "id": "c34de6a5", "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": "1a6bc209", "metadata": {}, "outputs": [], "source": [ "with model:\n", " txt, bib = xo.citations.get_citations_for_model()\n", "print(txt)" ] }, { "cell_type": "code", "execution_count": null, "id": "f1018631", "metadata": {}, "outputs": [], "source": [ "print(bib.split(\"\\n\\n\")[0] + \"\\n\\n...\")" ] }, { "cell_type": "code", "execution_count": null, "id": "a567800c-90b7-4a31-b98c-e253a2873d7f", "metadata": {}, "outputs": [], "source": [] } ], "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 }