models

Models and likelihood functions for use with lmfit

pycheops.models.qpower2

Fast and accurate transit light curves for the power-2 limb-darkening law

The power-2 limb-darkening law is

\[I(\mu) = 1 - c (1 - \mu^\alpha)\]

Light curves are calculated using the qpower2 approximation 2. The approximation is accurate to better than 100ppm for radius ratio k < 0.1.

N.B. qpower2 is untested/inaccurate for values of k > 0.2

2

Maxted, P.F.L. & Gill, S., 2019A&A…622A..33M

Parameters
  • z – star-planet separation on the sky cf. star radius (array)

  • k – planet-star radius ratio (scalar, k<1)

  • c – power-2 limb darkening coefficient

  • a – power-2 limb darkening exponent

Returns

light curve (observed flux)

Example

>>> from pycheops.models import qpower2
>>> from pycheops.funcs import t2z
>>> from numpy import linspace
>>> import matplotlib.pyplot as plt
>>> t = linspace(-0.025,0.025,1000)
>>> sini = 0.999
>>> rstar = 0.05
>>> ecc = 0.2
>>> om = 120
>>> tzero = 0.0
>>> P = 0.1
>>> z=t2z(t,tzero,P,sini,rstar,ecc,om)
>>> c = 0.5
>>> a = 0.7
>>> k = 0.1
>>> f = qpower2(z,k,c,a)
>>> plt.plot(t,f)
>>> plt.show()
pycheops.models.ueclipse

Eclipse light curve for a planet with uniform surface brightness by a star

Parameters
  • z – star-planet separation on the sky cf. star radius (array)

  • k – planet-star radius ratio (scalar, k<1)

Returns

light curve (observed flux from eclipsed source)

class pycheops.models.TransitModel(independent_vars=['t'], prefix='', nan_policy='raise', **kwargs)

Light curve model for the transit of a spherical star by an opaque spherical body (planet).

Limb-darkening is described by the power-2 law:

\[I(\mu) = 1 - c (1 - \mu^\alpha)\]

The transit depth, width shape are parameterised by D, W and b. These parameters are defined below in terms of the radius of the star and planet, R_s and R_p, respectively, the semi-major axis, a, and the orbital inclination, i. The eccentricy and longitude of periastron for the star’s orbit are e and omega, respectively.

The following parameters are defined for convenience:

  • k = R_p/R_s;

  • aR = a/R_s;

  • rho = 0.013418*aR**3/(P/d)**2.

N.B. the mean stellar density in solar units is rho, but only if the mass ratio q = M_planet/M_star is q << 1.

Parameters
  • t

    • independent variable (time)

  • T_0

    • time of mid-transit

  • P

    • orbital period

  • D

    • (R_p/R_s)**2 = k**2

  • W

    • (R_s/a)*sqrt((1+k)**2 - b**2)/pi

  • b

    • a*cos(i)/R_s

  • f_c

    • sqrt(ecc)*cos(omega)

  • f_s

    • sqrt(ecc)*sin(omega)

  • h_1

    • I(0.5) = 1 - c*(1-0.5**alpha)

  • h_2

    • I(0.5) - I(0) = c*0.5**alpha

The flux value outside of transit is 1. The light curve is calculated using the qpower2 algorithm, which is fast but only accurate for k < ~0.3.

If the input parameters are invalid or k>0.5 the model is returned as an array of value 1 everywhere.

class pycheops.models.EclipseModel(independent_vars=['t'], prefix='', nan_policy='raise', **kwargs)

Light curve model for the eclipse by a spherical star of a spherical body (planet) with no limb darkening.

The transit depth, width shape are parameterised by D, W and b. These parameters are defined below in terms of the radius of the star and planet, R_s and R_p, respectively, the semi-major axis, a, and the orbital inclination, i. The eccentricy and longitude of periastron for the star’s orbit are e and omega, respectively. These are the same parameters used in TransitModel. The flux level outside of eclipse is 1 and inside eclipse is (1-L). The apparent time of mid-eclipse includes the correction a_c for the light travel time across the orbit, i.e., for a circular orbit the time of mid-eclipse is (T_0 + 0.5*P) + a_c.

N.B. a_c must have the same units as P.

The following parameters are defined for convenience:

  • k = R_p/R_s;

  • aR = a/R_s;

  • rho = 0.013418*aR**3/(P/d)**2.

N.B. the mean stellar density in solar units is rho, but only if the mass ratio q = M_planet/M_star is q << 1.

Parameters
  • t

    • independent variable (time)

  • T_0

    • time of mid-transit

  • P

    • orbital period

  • D

    • (R_p/R_s)**2 = k**2

  • W

    • (R_s/a)*sqrt((1+k)**2 - b**2)/pi

  • b

    • a*cos(i)/R_s

  • L

    • Depth of eclipse

  • f_c

    • sqrt(ecc).cos(omega)

  • f_s

    • sqrt(ecc).sin(omega)

  • a_c

    • correction for light travel time across the orbit

class pycheops.models.FactorModel(independent_vars=['t'], prefix='', nan_policy='raise', dx=None, dy=None, sinphi=None, cosphi=None, bg=None, **kwargs)

Flux scaling and trend factor model

f = c*(1 + dfdt*dt + d2fdt2*dt**2 + dfdbg*bg(t) +

dfdx*dx(t) + dfdy*dy(t) + d2fdx2*dx(t)**2 + d2f2y2*dy(t)**2 + d2fdxdy*x(t)*dy(t) + dfdsinphi*sin(phi(t)) + dfdcosphi*cos(phi(t)) + dfdsin2phi*sin(2.phi(t)) + dfdcos2phi*cos(2.phi(t)) + dfdsin3phi*sin(3.phi(t)) + dfdcos3phi*cos(3.phi(t)) )

The detrending coefficients dfdx, etc. are 0 and fixed by default. If any of the coefficients dfdx, d2fdxdy or d2f2x2 is not 0, a function to calculate the x-position offset as a function of time, dx(t), must be passed as a keyword argument, and similarly for the y-position offset, dy(t). For detrending against the spacecraft roll angle, phi(t), the functions to be provided as keywords arguments are sinphi(t) and cosphi(t). The linear trend dfdbg is proportional to the estimated background flux in the aperture, bg(t). The time trend decribed by dfdt and d2fdt2 is calculated using the variable dt = t - median(t).

guess(data, **kwargs)

Estimate initial model parameter values from data.

class pycheops.models.ThermalPhaseModel(independent_vars=['t'], prefix='', nan_policy='raise', **kwargs)

Thermal phase model for a tidally-locked planet

\[a_{th}[1-\cos(\phi))/2 + b_{th}*(1+\sin(\phi)/2 + c_{th},\]

where \(\phi = 2\pi(t-T_0)/P\)

Parameters
  • t

    • independent variable (time)

  • T_0

    • time of inferior conjunction (mid-transit)

  • P

    • orbital period

  • a_th

    • coefficient of cosine-like term

  • b_th

    • coefficient of sine-like term

  • c_th

    • constant term (minimum flux)

The following parameters are defined for convenience.

  • A = sqrt(a_th**2 + b_th**2), peak-to-trough amplitude of the phase curve

  • F = c_th + (a_th + b_th + A)/2, flux at the maximum of the phase curve

  • ph_max = arctan2(b_th,-a_th)/(2*pi) = phase at maximum flux

class pycheops.models.ReflectionModel(independent_vars=['t'], prefix='', nan_policy='raise', **kwargs)

Reflected stellar light from a planet with a Lambertian phase function.

The fraction of the stellar flux reflected from the planet of radius \(R_p\) at a distance \(r\) from the star and viewed at phase angle \(\beta\) is

\[A_g(R_p/r)^2 \times [\sin(\beta) + (\pi-\beta)*\cos(\beta) ]/\pi\]

The eccentricity and longitude of periastron for the planet’s orbit are ecc and omega, respectively.

Parameters
  • t

    • independent variable (time)

  • T_0

    • time of inferior conjunction (mid-transit)

  • P

    • orbital period

  • A_g

    • geometric albedo

  • r_p

    • R_p/a, where a is the semi-major axis.

  • f_c

    • sqrt(ecc).cos(omega)

  • f_s

    • sqrt(ecc).sin(omega)

  • sini

    • sin(inclination)

class pycheops.models.RVModel(independent_vars=['t'], prefix='', nan_policy='raise', **kwargs)

Radial velocity in a Keplerian orbit

Parameters
  • t

    • independent variable (time)

  • T_0

    • time of inferior conjunction for the companion (mid-transit)

  • P

    • orbital period

  • V_0

    • radial velocity of the centre-of-mass

  • K

    • semi-amplitude of spectroscopic orbit

  • f_c

    • sqrt(ecc).cos(omega)

  • f_s

    • sqrt(ecc).sin(omega)

  • sini

    • sine of the orbital inclination

class pycheops.models.RVCompanion(independent_vars=['t'], prefix='', nan_policy='raise', **kwargs)

Radial velocity in a Keplerian orbit for the companion

Parameters
  • t

    • independent variable (time)

  • T_0

    • time of inferior conjunction for the companion (mid-transit)

  • P

    • orbital period

  • V_0

    • radial velocity of the centre-of-mass

  • K

    • semi-amplitude of spectroscopic orbit

  • f_c

    • sqrt(ecc).cos(omega)

  • f_s

    • sqrt(ecc).sin(omega)

  • sini

    • sine of the orbital inclination

class pycheops.models.EBLMModel(independent_vars=['t'], prefix='', nan_policy='raise', **kwargs)

Light curve model for the mutual eclipses by spherical stars in an eclipsing binary with one low-mass companion, e.g., F/G-star + M-dwarf.

The transit depth, width shape are parameterised by D, W and b. These parameters are defined below in terms of the radii of the stars, R_1 and R_2, the semi-major axis, a, and the orbital inclination, i. This model assumes R_1 >> R_2, i.e., k=R_2/R_1 <~0.2. The eccentricy and longitude of periastron for the star’s orbit are e and omega, respectively. These are the same parameters used in TransitModel. The flux level outside of eclipse is 1 and inside eclipse is (1-L). The apparent time of mid-eclipse includes the correction a_c for the light travel time across the orbit, i.e., for a circular orbit the time of mid-eclipse is (T_0 + 0.5*P) + a_c.

N.B. a_c must have the same units as P. The power-2 law is used to model the limb-darkening of star 1. Limb-darkening on star 2 is ignored.

The following parameters are defined for convenience:

  • k = R_2/R_1;

  • aR = a/R_1;

  • J = L/D (surface brightness ratio).

Parameters
  • t

    • independent variable (time)

  • T_0

    • time of mid-transit

  • P

    • orbital period

  • D

    • (R_2/R_1)**2 = k**2

  • W

    • (R_1/a)*sqrt((1+k)**2 - b**2)/pi

  • b

    • a*cos(i)/R_1

  • L

    • Depth of eclipse

  • f_c

    • sqrt(ecc).cos(omega)

  • f_s

    • sqrt(ecc).sin(omega)

  • h_1

    • I(0.5) = 1 - c*(1-0.5**alpha)

  • h_2

    • I(0.5) - I(0) = c*0.5**alpha

  • a_c

    • correction for light travel time across the orbit

class pycheops.models.PlanetModel(independent_vars=['t'], prefix='', nan_policy='raise', **kwargs)

Light curve model for a transiting exoplanet including transits, eclipses, and a thermal phase curve for the planet with an offset.

The flux level from the star is 1 and is assumed to be constant.

The thermal phase curve from the planet is approximated by a cosine function with amplitude A=F_max-F_min plus the minimum flux, F_min, i.e., the maximum flux is F_max = F_min+A, and this occurs at phase (ph_off+0.5) relative to the time of mid-transit, i.e.,

\[f_{\rm th} = F_{\rm min} + A[1-\cos(\phi-\phi_{\rm off})]/2\]

where \(\phi = 2\pi(t-T_0)/P\) and \(\phi_{\rm off} = 2\pi\,{\rm ph\_off}\).

The transit depth, width shape are parameterised by D, W and b. These parameters are defined below in terms of the radius of the star, R_1 and R_2, the semi-major axis, a, and the orbital inclination, i. This model assumes R_1 >> R_2, i.e., k=R_2/R_1 <~0.2. The eccentricy and longitude of periastron for the star’s orbit are e and omega, respectively. These are the same parameters used in TransitModel. The eclipse of the planet assumes a uniform flux distribution.

The apparent time of mid-eclipse includes the correction a_c for the light travel time across the orbit, i.e., for a circular orbit the time of mid-eclipse is (T_0 + 0.5*P) + a_c.

N.B. a_c must have the same units as P.

Stellar limb-darkening is described by the power-2 law:

\[I(\mu) = 1 - c (1 - \mu^\alpha)\]

The following parameters are defined for convenience:

  • k = R_2/R_1;

  • aR = a/R_1;

  • A = F_max - F_min = amplitude of thermal phase effect.

  • rho = 0.013418*aR**3/(P/d)**2.

N.B. the mean stellar density in solar units is rho, but only if the mass ratio q = M_planet/M_star is q << 1.

Parameters
  • t

    • independent variable (time)

  • T_0

    • time of mid-transit

  • P

    • orbital period

  • D

    • (R_2/R_1)**2 = k**2

  • W

    • (R_1/a)*sqrt((1+k)**2 - b**2)/pi

  • b

    • a*cos(i)/R_1

  • F_min

    • minimum flux in the thermal phase model

  • F_max

    • maximum flux in the thermal phase model

  • ph_off

    • offset phase in the thermal phase model

  • f_c

    • sqrt(ecc).cos(omega)

  • f_s

    • sqrt(ecc).sin(omega)

  • h_1

    • I(0.5) = 1 - c*(1-0.5**alpha)

  • h_2

    • I(0.5) - I(0) = c*0.5**alpha

  • a_c

    • correction for light travel time across the orbit

pycheops.models.scaled_transit_fit

Optimum scaled transit depth for data with scaled errors

Find the value of the scaling factor s that provides the best fit of the model m = 1 + s*(model-1) to the normalised input fluxes. It is assumed that the true standard errors on the flux measurements are a factor f times the nominal standard error(s) provided in sigma. Also returns standard error estimates for s and f, sigma_s and sigma_f, respectively.

param flux

Array of normalised flux measurements

param sigma

Standard error estimate(s) for flux - array or scalar

param model

Transit model to be scaled

returns

s, b, sigma_s, sigma_b

pycheops.models.minerr_transit_fit(flux, sigma, model)

Optimum scaled transit depth for data with lower bounds on errors

Find the value of the scaling factor s that provides the best fit of the model m = 1 + s*(model-1) to the normalised input fluxes. It is assumed that the nominal standard error(s) provided in sigma are lower bounds to the true standard errors on the flux measurements. 1 The probability distribution for the true standard errors is assumed to be

\[P(\sigma_{\rm true} | \sigma) = \sigma/\sigma_{\rm true}^2 \]
param flux

Array of normalised flux measurements

param sigma

Lower bound(s) on standard error for flux - array or scalar

param model

Transit model to be scaled

returns

s, sigma_s

1

Sivia, D.S. & Skilling, J., Data Analysis - A Bayesian Tutorial, 2nd ed., section 8.3.1