Note

This tutorial was generated from an IPython notebook that can be downloaded here.

Fitting a detached eclipsing binary

Note

You will need exoplanet version 0.3.1 or later to run this tutorial.

In this case study, we’ll go through the steps required to fit the light curve and radial velocity measurements for the detached eclipsing binary system HD 23642. This is a bright system that has been fit by many authors (1, 2, 3, 4, and 5 to name a few) so this is a good benchmark for our demonstration.

The light curve that we’ll use is from K2 and we’ll use the same radial velocity measurements as David+ (2016) compiled from here and here. We’ll use a somewhat simplified model for the eclipses that treats the stars as spherical and ignores the phase curve (we’ll model it using a Gaussian process instead of a more physically motivated model). But, as you’ll see, these simplifying assumptions are sufficient for this case of a detached and well behaved system. Unlike some previous studies, we will fit an eccentric orbit instead of fixing the eccentricity to zero. This probably isn’t really necessary here, but it’s useful to demonstrate how you would fit a more eccentric system. Finally, we model the phase curve and other triends in both the light curve and radial velocities using Gaussian processes. This will account for unmodeled stellar variability and residual systematics, drifts, and other effects left over from the data reduction procedure.

Data access

First, let’s define some values from the literature that will be useful below. Here we’re taking the period and eclipse time from David+ (2016) as initial guesses for these parameters in our fit. We’ll also include the same prior on the flux ratio of the two stars that was computed for the Kepler bandpass by David+ (2016).

lit_period = 2.46113408
lit_t0 = 119.522070 + 2457000 - 2454833

# Prior on the flux ratio for Kepler
lit_flux_ratio = (0.354, 0.035)

Then we’ll download the Kepler data. In this case, the pipeline aperture photometry isn’t very good (because this star is so bright!) so we’ll just download the target pixel file and co-add all the pixels.

import numpy as np
import matplotlib.pyplot as plt
import lightkurve as lk

tpf = lk.search_targetpixelfile("EPIC 211082420").download()
lc = tpf.to_lightcurve(aperture_mask="all")
lc = lc.remove_nans().normalize()

hdr = tpf.hdu[1].header
texp = hdr["FRAMETIM"] * hdr["NUM_FRM"]
texp /= 60.0 * 60.0 * 24.0

x = np.ascontiguousarray(lc.time, dtype=np.float64)
y = np.ascontiguousarray(lc.flux, dtype=np.float64)
mu = np.median(y)
y = (y / mu - 1) * 1e3

plt.plot((x - lit_t0 + 0.5 * lit_period) % lit_period - 0.5 * lit_period, y, ".k")
plt.xlim(-0.5 * lit_period, 0.5 * lit_period)
plt.xlabel("time since primary eclipse [days]")
_ = plt.ylabel("relative flux [ppt]")
../../_images/eb_7_0.png

Then we’ll enter the radial velocity data. I couldn’t find these data online anywhere so I manually transcribed the data from the referenced papers (typos are my own!).

ref1 = 2453000
ref2 = 2400000
rvs = np.array(
    [
        # https://arxiv.org/abs/astro-ph/0403444
        (39.41273 + ref1, -85.0, 134.5),
        (39.45356 + ref1, -88.0, 139.0),
        (39.50548 + ref1, -91.0, 143.0),
        (43.25049 + ref1, 105.5, -136.0),
        (46.25318 + ref1, 29.5, -24.5),
        # https://ui.adsabs.harvard.edu/abs/2007A%26A...463..579G/abstract
        (52629.6190 + ref2, 88.8, -127.0),
        (52630.6098 + ref2, -48.0, 68.0),
        (52631.6089 + ref2, -9.5, 13.1),
        (52632.6024 + ref2, 63.6, -90.9),
        (52633.6162 + ref2, -94.5, 135.0),
        (52636.6055 + ref2, 10.3, -13.9),
        (52983.6570 + ref2, 18.1, -25.1),
        (52987.6453 + ref2, -80.6, 114.5),
        (52993.6322 + ref2, 49.0, -70.7),
        (53224.9338 + ref2, 39.0, -55.7),
        (53229.9384 + ref2, 57.2, -82.0),
    ]
)
rvs[:, 0] -= 2454833
rvs = rvs[np.argsort(rvs[:, 0])]

x_rv = np.ascontiguousarray(rvs[:, 0], dtype=np.float64)
y1_rv = np.ascontiguousarray(rvs[:, 1], dtype=np.float64)
y2_rv = np.ascontiguousarray(rvs[:, 2], dtype=np.float64)

fold = (rvs[:, 0] - lit_t0 + 0.5 * lit_period) % lit_period - 0.5 * lit_period
plt.plot(fold, rvs[:, 1], ".", label="primary")
plt.plot(fold, rvs[:, 2], ".", label="secondary")
plt.legend(fontsize=10)
plt.xlim(-0.5 * lit_period, 0.5 * lit_period)
plt.ylabel("radial velocity [km / s]")
_ = plt.xlabel("time since primary eclipse [days]")
../../_images/eb_9_0.png

Probabilistic model

Then we define the probabilistic model using PyMC3 and exoplanet. This is similar to the other tutorials and case studies, but here we’re using a exoplanet.SecondaryEclipseLightCurve to generate the model light curve and we’re modeling the radial velocity trends using a Gaussian process instead of a polynomial. Otherwise, things should look pretty familiar!

After defining the model, we iteratively clip outliers in the light curve using sigma clipping and then estimate the maximum a posteriori parameters.

import pymc3 as pm
import theano.tensor as tt

import exoplanet as xo


def build_model(mask):

    with pm.Model() as model:

        # Systemic parameters
        mean_lc = pm.Normal("mean_lc", mu=0.0, sd=5.0)
        mean_rv = pm.Normal("mean_rv", mu=0.0, sd=50.0)
        u1 = xo.QuadLimbDark("u1")
        u2 = xo.QuadLimbDark("u2")

        # Parameters describing the primary
        M1 = pm.Lognormal("M1", mu=0.0, sigma=10.0)
        R1 = pm.Lognormal("R1", mu=0.0, sigma=10.0)

        # Secondary ratios
        k = pm.Lognormal("k", mu=0.0, sigma=10.0)  # radius ratio
        q = pm.Lognormal("q", mu=0.0, sigma=10.0)  # mass ratio
        s = pm.Lognormal("s", mu=np.log(0.5), sigma=10.0)  # surface brightness ratio

        # Prior on flux ratio
        pm.Normal(
            "flux_prior",
            mu=lit_flux_ratio[0],
            sigma=lit_flux_ratio[1],
            observed=k ** 2 * s,
        )

        # Parameters describing the orbit
        b = xo.ImpactParameter("b", ror=k, testval=1.5)
        period = pm.Lognormal("period", mu=np.log(lit_period), sigma=1.0)
        t0 = pm.Normal("t0", mu=lit_t0, sigma=1.0)

        # Parameters describing the eccentricity: ecs = [e * cos(w), e * sin(w)]
        ecs = xo.UnitDisk("ecs", testval=np.array([1e-5, 0.0]))
        ecc = pm.Deterministic("ecc", tt.sqrt(tt.sum(ecs ** 2)))
        omega = pm.Deterministic("omega", tt.arctan2(ecs[1], ecs[0]))

        # Build the orbit
        R2 = pm.Deterministic("R2", k * R1)
        M2 = pm.Deterministic("M2", q * M1)
        orbit = xo.orbits.KeplerianOrbit(
            period=period,
            t0=t0,
            ecc=ecc,
            omega=omega,
            b=b,
            r_star=R1,
            m_star=M1,
            m_planet=M2,
        )

        # Track some other orbital elements
        pm.Deterministic("incl", orbit.incl)
        pm.Deterministic("a", orbit.a)

        # Noise model for the light curve
        sigma_lc = pm.InverseGamma(
            "sigma_lc", testval=1.0, **xo.estimate_inverse_gamma_parameters(0.1, 2.0)
        )
        S_tot_lc = pm.InverseGamma(
            "S_tot_lc", testval=2.5, **xo.estimate_inverse_gamma_parameters(1.0, 5.0)
        )
        ell_lc = pm.InverseGamma(
            "ell_lc", testval=2.0, **xo.estimate_inverse_gamma_parameters(1.0, 5.0)
        )
        kernel_lc = xo.gp.terms.SHOTerm(
            S_tot=S_tot_lc, w0=2 * np.pi / ell_lc, Q=1.0 / 3
        )

        # Noise model for the radial velocities
        sigma_rv1 = pm.InverseGamma(
            "sigma_rv1", testval=1.0, **xo.estimate_inverse_gamma_parameters(0.5, 5.0)
        )
        sigma_rv2 = pm.InverseGamma(
            "sigma_rv2", testval=1.0, **xo.estimate_inverse_gamma_parameters(0.5, 5.0)
        )
        S_tot_rv = pm.InverseGamma(
            "S_tot_rv", testval=2.5, **xo.estimate_inverse_gamma_parameters(1.0, 5.0)
        )
        ell_rv = pm.InverseGamma(
            "ell_rv", testval=2.0, **xo.estimate_inverse_gamma_parameters(1.0, 5.0)
        )
        kernel_rv = xo.gp.terms.SHOTerm(
            S_tot=S_tot_rv, w0=2 * np.pi / ell_rv, Q=1.0 / 3
        )

        # Set up the light curve model
        lc = xo.SecondaryEclipseLightCurve(u1, u2, s)

        def model_lc(t):
            return (
                mean_lc
                + 1e3 * lc.get_light_curve(orbit=orbit, r=R2, t=t, texp=texp)[:, 0]
            )

        # Condition the light curve model on the data
        gp_lc = xo.gp.GP(
            kernel_lc, x[mask], tt.zeros(mask.sum()) ** 2 + sigma_lc ** 2, mean=model_lc
        )
        gp_lc.marginal("obs_lc", observed=y[mask])

        # Set up the radial velocity model
        def model_rv1(t):
            return mean_rv + 1e-3 * orbit.get_radial_velocity(t)

        def model_rv2(t):
            return mean_rv - 1e-3 * orbit.get_radial_velocity(t) / q

        # Condition the radial velocity model on the data
        gp_rv1 = xo.gp.GP(
            kernel_rv, x_rv, tt.zeros(len(x_rv)) ** 2 + sigma_rv1 ** 2, mean=model_rv1
        )
        gp_rv1.marginal("obs_rv1", observed=y1_rv)
        gp_rv2 = xo.gp.GP(
            kernel_rv, x_rv, tt.zeros(len(x_rv)) ** 2 + sigma_rv2 ** 2, mean=model_rv2
        )
        gp_rv2.marginal("obs_rv2", observed=y2_rv)

        # Optimize the logp
        map_soln = model.test_point

        # First the RV parameters
        map_soln = xo.optimize(map_soln, [mean_rv, q])
        map_soln = xo.optimize(
            map_soln, [mean_rv, sigma_rv1, sigma_rv2, S_tot_rv, ell_rv]
        )

        # Then the LC parameters
        map_soln = xo.optimize(map_soln, [mean_lc, R1, k, s, b])
        map_soln = xo.optimize(map_soln, [mean_lc, R1, k, s, b, u1, u2])
        map_soln = xo.optimize(map_soln, [mean_lc, sigma_lc, S_tot_lc, ell_lc])
        map_soln = xo.optimize(map_soln, [t0, period])

        # Then all the parameters together
        map_soln = xo.optimize(map_soln)

        model.gp_lc = gp_lc
        model.model_lc = model_lc
        model.gp_rv1 = gp_rv1
        model.model_rv1 = model_rv1
        model.gp_rv2 = gp_rv2
        model.model_rv2 = model_rv2

        model.x = x[mask]
        model.y = y[mask]

    return model, map_soln


def sigma_clip():
    mask = np.ones(len(x), dtype=bool)
    num = len(mask)

    for i in range(10):
        model, map_soln = build_model(mask)

        with model:
            mdl = xo.eval_in_model(
                model.model_lc(x[mask]) + model.gp_lc.predict(), map_soln
            )

        resid = y[mask] - mdl
        sigma = np.sqrt(np.median((resid - np.median(resid)) ** 2))
        mask[mask] = np.abs(resid - np.median(resid)) < 7 * sigma
        print("Sigma clipped {0} light curve points".format(num - mask.sum()))
        if num == mask.sum():
            break
        num = mask.sum()

    return model, map_soln


model, map_soln = sigma_clip()
optimizing logp for variables: [q, mean_rv]
13it [00:03,  3.74it/s, logp=-1.307282e+04]
message: Optimization terminated successfully.
logp: -23380.5842277239 -> -13072.816791759677
optimizing logp for variables: [ell_rv, S_tot_rv, sigma_rv2, sigma_rv1, mean_rv]
40it [00:00, 121.20it/s, logp=-8.482576e+03]
message: Optimization terminated successfully.
logp: -13072.816791759677 -> -8482.576126912147
optimizing logp for variables: [b, k, s, R1, mean_lc]
124it [00:00, 173.44it/s, logp=-4.696334e+03]
message: Desired error not necessarily achieved due to precision loss.
logp: -8482.576126912147 -> -4696.333901223958
optimizing logp for variables: [u2, u1, b, k, s, R1, mean_lc]
37it [00:00, 109.56it/s, logp=-4.694047e+03]
message: Desired error not necessarily achieved due to precision loss.
logp: -4696.333901223958 -> -4694.047160949973
optimizing logp for variables: [ell_lc, S_tot_lc, sigma_lc, mean_lc]
20it [00:00, 76.01it/s, logp=-3.640456e+03]
message: Optimization terminated successfully.
logp: -4694.047160949972 -> -3640.455503590622
optimizing logp for variables: [period, t0]
83it [00:00, 159.32it/s, logp=-3.636741e+03]
message: Desired error not necessarily achieved due to precision loss.
logp: -3640.455503590622 -> -3636.7412339176194
optimizing logp for variables: [ell_rv, S_tot_rv, sigma_rv2, sigma_rv1, ell_lc, S_tot_lc, sigma_lc, ecs, t0, period, b, s, q, k, R1, M1, u2, u1, mean_rv, mean_lc]
178it [00:01, 164.38it/s, logp=-3.539431e+03]
message: Desired error not necessarily achieved due to precision loss.
logp: -3636.7412339176194 -> -3539.430818773797
Sigma clipped 34 light curve points
optimizing logp for variables: [q, mean_rv]
13it [00:00, 58.20it/s, logp=-1.271094e+04]
message: Optimization terminated successfully.
logp: -22965.237266602366 -> -12710.941974601861
optimizing logp for variables: [ell_rv, S_tot_rv, sigma_rv2, sigma_rv1, mean_rv]
40it [00:00, 127.97it/s, logp=-8.141907e+03]
message: Optimization terminated successfully.
logp: -12710.941974601868 -> -8141.907237100719
optimizing logp for variables: [b, k, s, R1, mean_lc]
33it [00:00, 46.02it/s, logp=-4.368968e+03]
message: Optimization terminated successfully.
logp: -8141.907237100712 -> -4368.968064606436
optimizing logp for variables: [u2, u1, b, k, s, R1, mean_lc]
41it [00:00, 118.12it/s, logp=-4.366608e+03]
message: Optimization terminated successfully.
logp: -4368.968064606436 -> -4366.608037969302
optimizing logp for variables: [ell_lc, S_tot_lc, sigma_lc, mean_lc]
105it [00:00, 173.15it/s, logp=-1.789111e+03]
message: Desired error not necessarily achieved due to precision loss.
logp: -4366.608037969302 -> -1789.1106169136733
optimizing logp for variables: [period, t0]
65it [00:00, 148.98it/s, logp=-1.783905e+03]
message: Desired error not necessarily achieved due to precision loss.
logp: -1789.1106169136738 -> -1783.9053609029884
optimizing logp for variables: [ell_rv, S_tot_rv, sigma_rv2, sigma_rv1, ell_lc, S_tot_lc, sigma_lc, ecs, t0, period, b, s, q, k, R1, M1, u2, u1, mean_rv, mean_lc]
197it [00:01, 169.56it/s, logp=-1.624412e+03]
message: Desired error not necessarily achieved due to precision loss.
logp: -1783.9053609029884 -> -1624.411649047605
Sigma clipped 1 light curve points
optimizing logp for variables: [q, mean_rv]
13it [00:00, 54.19it/s, logp=-1.270829e+04]
message: Optimization terminated successfully.
logp: -22960.077690259175 -> -12708.286451975617
optimizing logp for variables: [ell_rv, S_tot_rv, sigma_rv2, sigma_rv1, mean_rv]
40it [00:00, 120.35it/s, logp=-8.140053e+03]
message: Optimization terminated successfully.
logp: -12708.286451975624 -> -8140.052575477539
optimizing logp for variables: [b, k, s, R1, mean_lc]
33it [00:00, 107.14it/s, logp=-4.367498e+03]
message: Optimization terminated successfully.
logp: -8140.052575477532 -> -4367.497544998031
optimizing logp for variables: [u2, u1, b, k, s, R1, mean_lc]
111it [00:00, 182.67it/s, logp=-4.365109e+03]
message: Desired error not necessarily achieved due to precision loss.
logp: -4367.497544998031 -> -4365.1089412302845
optimizing logp for variables: [ell_lc, S_tot_lc, sigma_lc, mean_lc]
122it [00:00, 175.80it/s, logp=-1.779796e+03]
message: Desired error not necessarily achieved due to precision loss.
logp: -4365.1089412302845 -> -1779.7958190026554
optimizing logp for variables: [period, t0]
72it [00:00, 154.59it/s, logp=-1.774645e+03]
message: Desired error not necessarily achieved due to precision loss.
logp: -1779.7958190026554 -> -1774.6449775292624
optimizing logp for variables: [ell_rv, S_tot_rv, sigma_rv2, sigma_rv1, ell_lc, S_tot_lc, sigma_lc, ecs, t0, period, b, s, q, k, R1, M1, u2, u1, mean_rv, mean_lc]
233it [00:01, 169.16it/s, logp=-1.614652e+03]
message: Desired error not necessarily achieved due to precision loss.
logp: -1774.6449775292624 -> -1614.6521192609016
Sigma clipped 0 light curve points

At these best fit parameters, let’s make some plots of the model predictions compared to the observations to make sure that things look reasonable. First the phase-folded radial velocities:

period = map_soln["period"]
t0 = map_soln["t0"]
mean = map_soln["mean_rv"]

x_fold = (x_rv - t0 + 0.5 * period) % period - 0.5 * period
plt.plot(fold, y1_rv - mean, ".", label="primary")
plt.plot(fold, y2_rv - mean, ".", label="secondary")

x_phase = np.linspace(-0.5 * period, 0.5 * period, 500)
with model:
    y1_mod, y2_mod = xo.eval_in_model(
        [model.model_rv1(x_phase + t0), model.model_rv2(x_phase + t0)], map_soln
    )
plt.plot(x_phase, y1_mod - mean, "C0")
plt.plot(x_phase, y2_mod - mean, "C1")

plt.legend(fontsize=10)
plt.xlim(-0.5 * period, 0.5 * period)
plt.ylabel("radial velocity [km / s]")
plt.xlabel("time since primary eclipse [days]")
_ = plt.title("HD 23642; map model", fontsize=14)
../../_images/eb_13_0.png

And then the light curve. In the top panel, we show the Gaussian process model for the phase curve. It’s clear that there’s a lot of information there that we could take advantage of, but that’s a topic for another day. In the bottom panel, we’re plotting the phase folded light curve and we can see the ridiculous signal to noise that we’re getting on the eclipses.

with model:
    gp_pred = xo.eval_in_model(model.gp_lc.predict(), map_soln) + map_soln["mean_lc"]
    lc = xo.eval_in_model(model.model_lc(model.x), map_soln) - map_soln["mean_lc"]

fig, (ax1, ax2) = plt.subplots(2, sharex=True, figsize=(12, 7))

ax1.plot(model.x, model.y, "k.", alpha=0.2)
ax1.plot(model.x, gp_pred, color="C1", lw=1)

ax2.plot(model.x, model.y - gp_pred, "k.", alpha=0.2)
ax2.plot(model.x, lc, color="C2", lw=1)
ax2.set_xlim(model.x.min(), model.x.max())

ax1.set_ylabel("raw flux [ppt]")
ax2.set_ylabel("de-trended flux [ppt]")
ax2.set_xlabel("time [KBJD]")
ax1.set_title("HD 23642; map model", fontsize=14)

fig.subplots_adjust(hspace=0.05)

fig, ax1 = plt.subplots(1, figsize=(12, 3.5))

x_fold = (model.x - map_soln["t0"]) % map_soln["period"] / map_soln["period"]
inds = np.argsort(x_fold)

ax1.plot(x_fold[inds], model.y[inds] - gp_pred[inds], "k.", alpha=0.2)
ax1.plot(x_fold[inds] - 1, model.y[inds] - gp_pred[inds], "k.", alpha=0.2)
ax2.plot(x_fold[inds], model.y[inds] - gp_pred[inds], "k.", alpha=0.2, label="data!")
ax2.plot(x_fold[inds] - 1, model.y[inds] - gp_pred, "k.", alpha=0.2)

yval = model.y[inds] - gp_pred
bins = np.linspace(0, 1, 75)
num, _ = np.histogram(x_fold[inds], bins, weights=yval)
denom, _ = np.histogram(x_fold[inds], bins)
ax2.plot(0.5 * (bins[:-1] + bins[1:]) - 1, num / denom, ".w")

args = dict(lw=1)

ax1.plot(x_fold[inds], lc[inds], "C2", **args)
ax1.plot(x_fold[inds] - 1, lc[inds], "C2", **args)

ax1.set_xlim(-1, 1)
ax1.set_ylabel("de-trended flux [ppt]")
ax1.set_xlabel("phase")
_ = ax1.set_title("HD 23642; map model", fontsize=14)
../../_images/eb_15_0.png ../../_images/eb_15_1.png

Sampling

Finally we can run the MCMC:

np.random.seed(23642)
with model:
    trace = xo.sample(
        tune=3500,
        draws=3000,
        start=map_soln,
        chains=4,
        initial_accept=0.8,
        target_accept=0.95,
    )
Multiprocess sampling (4 chains in 4 jobs)
NUTS: [ell_rv, S_tot_rv, sigma_rv2, sigma_rv1, ell_lc, S_tot_lc, sigma_lc, ecs, t0, period, b, s, q, k, R1, M1, u2, u1, mean_rv, mean_lc]
Sampling 4 chains, 0 divergences: 100%|██████████| 26000/26000 [22:49<00:00, 18.98draws/s]
The acceptance probability does not match the target. It is 0.901791582425816, but should be close to 0.95. Try to increase the number of tuning steps.

As usual, we can check the convergence diagnostics for some of the key parameters.

pm.summary(trace, var_names=["M1", "M2", "R1", "R2", "ecs", "incl", "s"])
mean sd hpd_3% hpd_97% mcse_mean mcse_sd ess_mean ess_sd ess_bulk ess_tail r_hat
M1 2.246 0.029 2.190 2.301 0.000 0.000 11897.0 11897.0 12178.0 8150.0 1.0
M2 1.570 0.022 1.528 1.614 0.000 0.000 11175.0 11175.0 11714.0 7879.0 1.0
R1 1.759 0.033 1.698 1.821 0.001 0.000 4107.0 4099.0 4131.0 4982.0 1.0
R2 1.513 0.052 1.418 1.611 0.001 0.001 4016.0 4016.0 4070.0 4180.0 1.0
ecs[0] 0.000 0.000 0.000 0.000 0.000 0.000 15446.0 14246.0 15454.0 8898.0 1.0
ecs[1] -0.002 0.005 -0.011 0.006 0.000 0.000 8834.0 6580.0 8883.0 7961.0 1.0
incl 1.365 0.002 1.361 1.369 0.000 0.000 4299.0 4296.0 4491.0 4116.0 1.0
s 0.478 0.021 0.438 0.518 0.000 0.000 7815.0 7815.0 7832.0 8020.0 1.0

Results

It can be useful to take a look at some diagnostic corner plots to see how the sampling went. First, let’s look at some observables:

import corner

samples = pm.trace_to_dataframe(trace, varnames=["k", "q", "ecs"])
_ = corner.corner(
    samples,
    labels=["$k = R_2 / R_1$", "$q = M_2 / M_1$", "$e\,\cos\omega$", "$e\,\sin\omega$"],
)
../../_images/eb_21_0.png

And then we can look at the physical properties of the stars in the system. In this figure, we’re comparing to the results from David+ (2016) (shown as blue crosshairs). The orange contours in this figure show the results transformed to a uniform prior on eccentricity as discussed below. These contours are provided to demonstrate (qualitatively) that these inferences are not sensitive to the choice of prior.

samples = pm.trace_to_dataframe(trace, varnames=["R1", "R2", "M1", "M2"])
weights = 1.0 / trace["ecc"]
weights *= len(weights) / np.sum(weights)
fig = corner.corner(samples, weights=weights, plot_datapoints=False, color="C1")
_ = corner.corner(samples, truths=[1.727, 1.503, 2.203, 1.5488], fig=fig)
../../_images/eb_23_0.png

A note about eccentricities

If you looked closely at the model defined above, you might have noticed that we chose a slightly odd eccentricity prior: \(p(e) \propto e\). This is implied by sampling with \(e\,\cos\omega\) and \(e\,\sin\omega\) as the parameters, as has been discussed many times in the literature. There are many options for correcting for this prior and instead assuming a uniform prior on eccentricity (for example, sampling with \(\sqrt{e}\,\cos\omega\) and \(\sqrt{e}\,\sin\omega\) as the parameters), but you’ll find much worse sampling performance for this problem if you try any of these options (trust us, we tried!) because the geometry of the posterior surface becomes much less suitable for the sampling algorithm in PyMC3. Instead, we can re-weight the samples after running the MCMC to see how the results change under the new prior. Most of the parameter inferences are unaffected by this change (because the data are very constraining!), but the inferred eccentricity (and especially \(e\,\sin\omega\)) will depend on this choice. The following plots show how these parameter inferences are affected. Note, especially, how the shape of the \(e\,\sin\omega\) density changes.

plt.hist(
    trace["ecc"] * np.sin(trace["omega"]),
    50,
    density=True,
    histtype="step",
    label="$p(e) = e / 2$",
)
plt.hist(
    trace["ecc"] * np.sin(trace["omega"]),
    50,
    density=True,
    histtype="step",
    weights=1.0 / trace["ecc"],
    label="$p(e) = 1$",
)
plt.xlabel("$e\,\sin(\omega)$")
plt.ylabel("$p(e\,\sin\omega\,|\,\mathrm{data})$")
plt.yticks([])
plt.legend(fontsize=12)

plt.figure()
plt.hist(trace["ecc"], 50, density=True, histtype="step", label="$p(e) = e / 2$")
plt.hist(
    trace["ecc"],
    50,
    density=True,
    histtype="step",
    weights=1.0 / trace["ecc"],
    label="$p(e) = 1$",
)
plt.xlabel("$e$")
plt.ylabel("$p(e\,|\,\mathrm{data})$")
plt.yticks([])
plt.xlim(0, 0.015)
_ = plt.legend(fontsize=12)
../../_images/eb_25_0.png ../../_images/eb_25_1.png

We can then use the corner.quantile function to compute summary statistics of the weighted samples as follows. For example, here how to compute the 90% posterior upper limit for the eccentricity:

weights = 1.0 / trace["ecc"]
print(
    "for p(e) = e/2: p(e < x) = 0.9 -> x = {0:.5f}".format(
        corner.quantile(trace["ecc"], [0.9])[0]
    )
)
print(
    "for p(e) = 1:   p(e < x) = 0.9 -> x = {0:.5f}".format(
        corner.quantile(trace["ecc"], [0.9], weights=weights)[0]
    )
)
for p(e) = e/2: p(e < x) = 0.9 -> x = 0.00851
for p(e) = 1:   p(e < x) = 0.9 -> x = 0.00389

Or, the posterior mean and variance for the radius of the primary:

samples = trace["R1"]

print(
    "for p(e) = e/2: R1 = {0:.3f} ± {1:.3f}".format(np.mean(samples), np.std(samples))
)

mean = np.sum(weights * samples) / np.sum(weights)
sigma = np.sqrt(np.sum(weights * (samples - mean) ** 2) / np.sum(weights))
print("for p(e) = 1:   R1 = {0:.3f} ± {1:.3f}".format(mean, sigma))
for p(e) = e/2: R1 = 1.759 ± 0.033
for p(e) = 1:   R1 = 1.766 ± 0.029

As you can see (and as one would hope) this choice of prior does not significantly change our inference of the primary radius.

Citations

As described in the Citing exoplanet & its dependencies tutorial, we can use exoplanet.citations.get_citations_for_model() to construct an acknowledgement and BibTeX listing that includes the relevant citations for this model.

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{exoplanet:agol19, exoplanet:astropy13, exoplanet:astropy18,
exoplanet:exoplanet, exoplanet:foremanmackey17, exoplanet:foremanmackey18,
exoplanet:kipping13, exoplanet:luger18, exoplanet:pymc3, exoplanet:theano}.
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 Tom Barclay},
   title = {exoplanet-dev/exoplanet v0.3.2},
   month = may,
    year = 2020,
     doi = {10.5281/zenodo.1998447},
     url = {https://doi.org/10.5281/zenodo.1998447}
}
...