Quick fits for TESS light curves

.. note:: You will need exoplanet version 0.2.6 or later to run this tutorial.

In this tutorial, we will fit the TESS light curve for a known transiting planet. While the Fitting TESS data case study goes through the full details of an end-to-end fit, this tutorial is significantly faster to run and it can give pretty excellent results depending on your goals. Some of the main differences are:

  1. We start from the light curve rather than doing the photometry ourselves. This should pretty much always be fine unless you have a very bright, faint, or crowded target.
  2. We assume a circluar orbit, but as you’ll see later, we can approximately relax this assumption later.
  3. We only fit the data near transit. In many cases this will be just fine, but if you have predictable stellar variability (like coherent rotation) then you might do better fitting more data.

We’ll fit the planet in the HD 118203 (TIC 286923464) system that was found to transit by Pepper et al. (2019) because it is on an eccentric orbit so assumption #2 above is not valid.

First, let’s download the TESS light curve using lightkurve:

[3]:
import numpy as np
import lightkurve as lk
import matplotlib.pyplot as plt

lcfs = lk.search_lightcurvefile("TIC 286923464", mission="TESS").download_all()
lc = lcfs.PDCSAP_FLUX.stitch()
lc = lc.remove_nans().remove_outliers(sigma=7)

x = np.ascontiguousarray(lc.time, dtype=np.float64)
y = np.ascontiguousarray(1e3 * (lc.flux - 1), dtype=np.float64)
yerr = np.ascontiguousarray(1e3 * lc.flux_err, dtype=np.float64)

texp = np.min(np.diff(x))

plt.plot(x, y, "k", linewidth=0.5)
plt.xlabel("time [days]")
_ = plt.ylabel("relative flux [ppt]")
../../_images/tutorials_quick-tess_6_0.png

Then, find the period, phase and depth of the transit using box least squares:

[4]:
import exoplanet as xo

pg = xo.estimators.bls_estimator(x, y, yerr, min_period=2, max_period=20)

peak = pg["peak_info"]
period_guess = peak["period"]
t0_guess = peak["transit_time"]
depth_guess = peak["depth"]

plt.plot(pg["bls"].period, pg["bls"].power, "k", linewidth=0.5)
plt.axvline(period_guess, alpha=0.3, linewidth=5)
plt.xlabel("period [days]")
plt.ylabel("bls power")
plt.yticks([])
_ = plt.xlim(pg["bls"].period.min(), pg["bls"].period.max())
../../_images/tutorials_quick-tess_8_0.png

Then, for efficiency purposes, let’s extract just the data within 0.25 days of the transits:

[5]:
transit_mask = (
    np.abs(
        (x - t0_guess + 0.5 * period_guess) % period_guess - 0.5 * period_guess
    )
    < 0.25
)
x = np.ascontiguousarray(x[transit_mask])
y = np.ascontiguousarray(y[transit_mask])
yerr = np.ascontiguousarray(yerr[transit_mask])

plt.figure(figsize=(8, 4))
x_fold = (
    x - t0_guess + 0.5 * period_guess
) % period_guess - 0.5 * period_guess
plt.scatter(x_fold, y, c=x, s=3)
plt.xlabel("time since transit [days]")
plt.ylabel("relative flux [ppt]")
plt.colorbar(label="time [days]")
_ = plt.xlim(-0.25, 0.25)
../../_images/tutorials_quick-tess_10_0.png

That looks a little janky, but it’s good enough for now.

The probabilistic model

Here’s how we set up the PyMC3 model in this case:

[6]:
import pymc3 as pm
import theano.tensor as tt

import pymc3_ext as pmx
from celerite2.theano import terms, GaussianProcess


with pm.Model() as model:

    # Stellar parameters
    mean = pm.Normal("mean", mu=0.0, sigma=10.0)
    u = xo.QuadLimbDark("u")
    star_params = [mean, u]

    # Gaussian process noise model
    sigma = pm.InverseGamma("sigma", alpha=3.0, beta=2 * np.median(yerr))
    sigma_gp = pm.Lognormal("sigma_gp", mu=0.0, sigma=10.0)
    rho_gp = pm.Lognormal("rho_gp", mu=np.log(10.0), sigma=10.0)
    kernel = terms.SHOTerm(sigma=sigma_gp, rho=rho_gp, Q=1.0 / 3)
    noise_params = [sigma, sigma_gp, rho_gp]

    # Planet parameters
    ror = pm.Lognormal("ror", mu=0.5 * np.log(depth_guess * 1e-3), sigma=10.0)

    # Orbital parameters
    period = pm.Lognormal("period", mu=np.log(period_guess), sigma=1.0)
    t0 = pm.Normal("t0", mu=t0_guess, sigma=1.0)
    dur = pm.Lognormal("dur", mu=np.log(0.1), sigma=10.0)
    b = xo.distributions.ImpactParameter("b", ror=ror)

    # Set up the orbit
    orbit = xo.orbits.KeplerianOrbit(period=period, duration=dur, t0=t0, b=b)

    # We're going to track the implied density for reasons that will become clear later
    pm.Deterministic("rho_circ", orbit.rho_star)

    # Set up the mean transit model
    star = xo.LimbDarkLightCurve(u)

    def lc_model(t):
        return mean + 1e3 * tt.sum(
            star.get_light_curve(orbit=orbit, r=ror, t=t), axis=-1
        )

    # Finally the GP observation model
    gp = GaussianProcess(
        kernel, t=x, diag=yerr ** 2 + sigma ** 2, mean=lc_model
    )
    gp.marginal("obs", observed=y)

    # Double check that everything looks good - we shouldn't see any NaNs!
    print(model.check_test_point())

    # Optimize the model
    map_soln = model.test_point
    map_soln = pmx.optimize(map_soln, [sigma])
    map_soln = pmx.optimize(map_soln, [ror, b, dur])
    map_soln = pmx.optimize(map_soln, noise_params)
    map_soln = pmx.optimize(map_soln, star_params)
    map_soln = pmx.optimize(map_soln)
mean                   -3.22
u_quadlimbdark__       -2.77
sigma_log__            -0.53
sigma_gp_log__         -3.22
rho_gp_log__           -3.22
ror_log__              -3.22
period_log__           -0.92
t0                     -0.92
dur_log__              -3.22
b_impact__             -1.39
obs                -26673.92
Name: Log-probability of test_point, dtype: float64
optimizing logp for variables: [sigma]
100.00% [15/15 00:00<00:00 logp = -6.366e+03]

message: Optimization terminated successfully.
logp: -26696.557571189896 -> -6365.810162508956
optimizing logp for variables: [dur, b, ror]
100.00% [32/32 00:00<00:00 logp = -4.857e+03]

message: Optimization terminated successfully.
logp: -6365.810162508956 -> -4857.110705140416
optimizing logp for variables: [rho_gp, sigma_gp, sigma]
100.00% [36/36 00:00<00:00 logp = -1.841e+03]

message: Optimization terminated successfully.
logp: -4857.110705140416 -> -1840.821686020982
optimizing logp for variables: [u, mean]
100.00% [71/71 00:00<00:00 logp = -1.835e+03]

message: Desired error not necessarily achieved due to precision loss.
logp: -1840.821686020982 -> -1834.5100241857647
optimizing logp for variables: [b, dur, t0, period, ror, rho_gp, sigma_gp, sigma, u, mean]
100.00% [112/112 00:00<00:00 logp = -1.457e+03]

message: Desired error not necessarily achieved due to precision loss.
logp: -1834.5100241857647 -> -1457.074984558115

Now we can plot our initial model:

[7]:
with model:
    gp_pred, lc_pred = xo.eval_in_model(
        [gp.predict(y, include_mean=False), lc_model(x)], map_soln
    )

plt.figure(figsize=(8, 4))
x_fold = (x - map_soln["t0"] + 0.5 * map_soln["period"]) % map_soln[
    "period"
] - 0.5 * map_soln["period"]
inds = np.argsort(x_fold)
plt.scatter(x_fold, y - gp_pred - map_soln["mean"], c=x, s=3)
plt.plot(x_fold[inds], lc_pred[inds] - map_soln["mean"], "k")
plt.xlabel("time since transit [days]")
plt.ylabel("relative flux [ppt]")
plt.colorbar(label="time [days]")
_ = plt.xlim(-0.25, 0.25)
../../_images/tutorials_quick-tess_14_0.png

That looks better!

Now on to sampling:

[8]:
np.random.seed(286923464)
with model:
    trace = pmx.sample(
        tune=2000, draws=2000, start=map_soln, chains=2, cores=2
    )
Multiprocess sampling (2 chains in 2 jobs)
NUTS: [b, dur, t0, period, ror, rho_gp, sigma_gp, sigma, u, mean]
100.00% [8000/8000 03:31<00:00 Sampling 2 chains, 0 divergences]
Sampling 2 chains for 2_000 tune and 2_000 draw iterations (4_000 + 4_000 draws total) took 211 seconds.

Then we can take a look at the summary statistics:

[9]:
with model:
    summary = pm.summary(trace)
summary
[9]:
mean sd hdi_3% hdi_97% mcse_mean mcse_sd ess_mean ess_sd ess_bulk ess_tail r_hat
mean 0.276 0.128 0.041 0.524 0.002 0.002 4147.0 3288.0 4282.0 2670.0 1.0
t0 1712.662 0.000 1712.662 1712.662 0.000 0.000 5245.0 5245.0 5239.0 3016.0 1.0
u[0] 0.288 0.083 0.130 0.435 0.002 0.001 2104.0 2104.0 2356.0 1337.0 1.0
u[1] 0.192 0.113 -0.008 0.413 0.002 0.002 2245.0 1407.0 2345.0 1420.0 1.0
sigma 0.185 0.007 0.171 0.198 0.000 0.000 4487.0 4487.0 4564.0 2528.0 1.0
sigma_gp 0.527 0.091 0.372 0.688 0.002 0.001 3472.0 3147.0 4169.0 2606.0 1.0
rho_gp 0.652 0.135 0.418 0.899 0.002 0.002 3795.0 3597.0 3884.0 3074.0 1.0
ror 0.055 0.000 0.054 0.056 0.000 0.000 3984.0 3974.0 4054.0 2215.0 1.0
period 6.135 0.000 6.135 6.135 0.000 0.000 5538.0 5538.0 5542.0 2980.0 1.0
dur 0.222 0.001 0.221 0.223 0.000 0.000 3687.0 3684.0 3815.0 2267.0 1.0
b 0.166 0.093 0.002 0.313 0.002 0.001 1952.0 1952.0 1800.0 1569.0 1.0
rho_circ 0.325 0.016 0.295 0.346 0.000 0.000 2567.0 2504.0 1931.0 2300.0 1.0

And plot the posterior covariances compared to the values from Pepper et al. (2019):

[10]:
import corner
import astropy.units as u

samples = pm.trace_to_dataframe(trace, varnames=["period", "ror", "b"])
_ = corner.corner(samples, truths=[6.134980, 0.05538, 0.125])
../../_images/tutorials_quick-tess_20_0.png

Bonus: eccentricity

As discussed above, we fit this model assuming a circular orbit which speeds things up for a few reasons. First, setting eccentricity to zero means that the orbital dynamics are much simpler and more computationally efficient, since we don’t need to solve Kepler’s equation numerically. But this isn’t actually the main effect! Instead the bigger issues come from the fact that the degeneracies between eccentricity, arrgument of periasteron, impact parameter, and planet radius are hard for the sampler to handle, causing the sampler’s performance to plummet. In this case, by fitting with a circular orbit where duration is one of the parameters, everything is well behaved and the sampler runs faster.

But, in this case, the planet is actually on an eccentric orbit, so that assumption isn’t justified. It has been recognized by various researchers over the years (I first learned about this from Bekki Dawson) that, to first order, the eccentricity mainly just changes the transit duration. The key realization is that this can be thought of as a change in the impled density of the star. Therefore, if you fit the transit using stellar density (or duration, in this case) as one of the parameters (note: you must have a different stellar density parameter for each planet if there are more than one), you can use an independent measurement of the stellar density to infer the eccentricity of the orbit after the fact. All the details are described in Dawson & Johnson (2012), but here’s how you can do this here using the stellar density listed in the TESS input catalog:

[11]:
from astroquery.mast import Catalogs

star = Catalogs.query_object("TIC 286923464", catalog="TIC", radius=0.001)
tic_rho_star = float(star["rho"]), float(star["e_rho"])
print("rho_star = {0} ± {1}".format(*tic_rho_star))

# Extract the implied density from the fit
rho_circ = np.repeat(trace["rho_circ"], 100)

# Sample eccentricity and omega from their priors (the math might
# be a little more subtle for more informative priors, but I leave
# that as an exercise for the reader...)
ecc = np.random.uniform(0, 1, len(rho_circ))
omega = np.random.uniform(-np.pi, np.pi, len(rho_circ))

# Compute the "g" parameter from Dawson & Johnson and what true
# density that implies
g = (1 + ecc * np.sin(omega)) / np.sqrt(1 - ecc ** 2)
rho = rho_circ / g ** 3

# Re-weight these samples to get weighted posterior samples
log_weights = -0.5 * ((rho - tic_rho_star[0]) / tic_rho_star[1]) ** 2
weights = np.exp(log_weights - np.max(log_weights))

# Estimate the expected posterior quantiles
q = corner.quantile(ecc, [0.16, 0.5, 0.84], weights=weights)
print(
    "eccentricity = {0:.2f} +{1[1]:.2f} -{1[0]:.2f}".format(q[1], np.diff(q))
)

_ = corner.corner(
    np.vstack((ecc, omega)).T,
    weights=weights,
    truths=[0.316, None],
    plot_datapoints=False,
    labels=["eccentricity", "omega"],
)
rho_star = 0.121689 ± 0.0281776
eccentricity = 0.46 +0.25 -0.13
../../_images/tutorials_quick-tess_22_1.png

As you can see, this eccentricity estimate is consistent (albeit with large uncertainties) with the value that Pepper et al. (2019) measure using radial velocities and it is definitely clear that this planet is not on a circular orbit.

Citations

As described in the citation tutorial, we can use citations.get_citations_for_model to construct an acknowledgement and BibTeX listing that includes the relevant citations for this model.

[12]:
with model:
    txt, bib = xo.citations.get_citations_for_model()
print(txt)
This research made use of \textsf{exoplanet} \citep{exoplanet} and its
dependencies \citep{celerite2:foremanmackey17, celerite2:foremanmackey18,
exoplanet:agol20, exoplanet:astropy13, exoplanet:astropy18,
exoplanet:exoplanet, exoplanet:kipping13, exoplanet:luger18, exoplanet:pymc3,
exoplanet:theano}.
[13]:
print("\n".join(bib.splitlines()[:10]) + "\n...")

@misc{exoplanet:exoplanet,
  author = {Daniel Foreman-Mackey and Rodrigo Luger and Ian Czekala and
            Eric Agol and Adrian Price-Whelan and Timothy D. Brandt and
            Tom Barclay and Luke Bouma},
   title = {exoplanet-dev/exoplanet v0.4.0},
   month = oct,
    year = 2020,
     doi = {10.5281/zenodo.1998447},
     url = {https://doi.org/10.5281/zenodo.1998447}
...