API documentation

Orbits

class exoplanet.orbits.KeplerianOrbit(period=None, a=None, t0=None, t_periastron=None, incl=None, b=None, duration=None, ecc=None, omega=None, Omega=None, m_planet=0.0, m_star=None, r_star=None, rho_star=None, ror=None, m_planet_units=None, rho_star_units=None, model=None, contact_points_kwargs=None, **kwargs)

A system of bodies on Keplerian orbits around a common central

Given the input parameters, the values of all other parameters will be computed so a KeplerianOrbit instance will always have attributes for each argument. Note that the units of the computed attributes will all be in the standard units of this class (R_sun, M_sun, and days) except for rho_star which will be in g / cm^3.

There are only specific combinations of input parameters that can be used:

  1. First, either period or a must be given. If values are given for both parameters, then neither m_star or rho_star can be defined because the stellar density implied by each planet will be computed in rho_star.
  2. Only one of incl and b can be given.
  3. If a value is given for ecc then omega must also be given.
  4. If no stellar parameters are given, the central body is assumed to be the sun. If only rho_star is defined, the radius of the central is assumed to be 1 * R_sun. Otherwise, at most two of m_star, r_star, and rho_star can be defined.
  5. Either t0 (reference transit) or t_periastron must be given, but not both.
Parameters:
  • period – The orbital periods of the bodies in days.
  • a – The semimajor axes of the orbits in R_sun.
  • t0 – The time of a reference transit for each orbits in days.
  • t_periastron – The epoch of a reference periastron passage in days.
  • incl – The inclinations of the orbits in radians.
  • b – The impact parameters of the orbits.
  • ecc – The eccentricities of the orbits. Must be 0 <= ecc < 1.
  • omega – The arguments of periastron for the orbits in radians.
  • Omega – The position angles of the ascending nodes in radians.
  • m_planet – The masses of the planets in units of m_planet_units.
  • m_star – The mass of the star in M_sun.
  • r_star – The radius of the star in R_sun.
  • rho_star – The density of the star in units of rho_star_units.
  • m_planet_units – An astropy.units compatible unit object giving the units of the planet masses. If not given, the default is M_sun.
  • rho_star_units – An astropy.units compatible unit object giving the units of the stellar density. If not given, the default is g / cm^3.
get_planet_position(t, parallax=None)

The planets’ positions in the barycentric frame

Parameters:t – The times where the position should be evaluated.
Returns:The components of the position vector at t in units of R_sun.
get_planet_velocity(t)

Get the planets’ velocity vectors

Parameters:t – The times where the velocity should be evaluated.
Returns:The components of the velocity vector at t in units of M_sun/day.
get_radial_velocity(t, K=None, output_units=None)

Get the radial velocity of the star

Note

The convention in exoplanet is that positive z points towards the observer. However, for consistency with radial velocity literature this method returns values where positive radial velocity corresponds to a redshift as expected.

Parameters:
  • t – The times where the radial velocity should be evaluated.
  • K (Optional) – The semi-amplitudes of the orbits. If provided, the m_planet and incl parameters will be ignored and this amplitude will be used instead.
  • output_units (Optional) – An AstroPy velocity unit. If not given, the output will be evaluated in m/s. This is ignored if a value is given for K.
Returns:

The reflex radial velocity evaluated at t in units of output_units. For multiple planets, this will have one row for each planet.

get_relative_angles(t, parallax=None)

The planets’ relative position to the star in the sky plane, in separation, position angle coordinates.

Note

This treats each planet independently and does not take the other planets into account when computing the position of the star. This is fine as long as the planet masses are small.

Parameters:t – The times where the position should be evaluated.
Returns:The separation (arcseconds) and position angle (radians, measured east of north) of the planet relative to the star.
get_relative_position(t, parallax=None)

The planets’ positions relative to the star in the X,Y,Z frame.

Note

This treats each planet independently and does not take the other planets into account when computing the position of the star. This is fine as long as the planet masses are small. In other words, the reflex motion of the star caused by the other planets is neglected when computing the relative coordinates of one of the planets.

Parameters:t – The times where the position should be evaluated.
Returns:The components of the position vector at t in units of R_sun.
get_relative_velocity(t)

The planets’ velocity relative to the star

Note

This treats each planet independently and does not take the other planets into account when computing the position of the star. This is fine as long as the planet masses are small.

Parameters:t – The times where the velocity should be evaluated.
Returns:The components of the velocity vector at t in units of R_sun/day.
get_star_position(t, parallax=None)

The star’s position in the barycentric frame

Note

If there are multiple planets in the system, this will return one column per planet with each planet’s contribution to the motion. The star’s full position can be computed by summing over the last axis.

Parameters:t – The times where the position should be evaluated.
Returns:The components of the position vector at t in units of R_sun.
get_star_velocity(t)

Get the star’s velocity vector

Note

For a system with multiple planets, this will return one column per planet with the contributions from each planet. The total velocity can be found by summing along the last axis.

Parameters:t – The times where the velocity should be evaluated.
Returns:The components of the velocity vector at t in units of M_sun/day.
in_transit(t, r=0.0, texp=None)

Get a list of timestamps that are in transit

Parameters:
  • t (vector) – A vector of timestamps to be evaluated.
  • r (Optional) – The radii of the planets.
  • texp (Optional[float]) – The exposure time.
Returns:

The indices of the timestamps that are in transit.

class exoplanet.orbits.TTVOrbit(*args, **kwargs)

A generalization of a Keplerian orbit with transit timing variations

Only one of the arguments ttvs or transit_times can be given and the other will be computed from the one that was provided.

In practice the way this works is that the time axis is shifted to account for the TTVs before solving for a standard Keplerian orbit. To identify which times corrorspond to which transits, this will find the closest labelled transit to the timestamp and adjust the timestamp accordingly. This means that things will go very badly if the time between neighboring transits is larger than 2*period.

Parameters:
  • ttvs – A list (with on entry for each planet) of “O-C” vectors for each transit of each planet in units of days. “O-C” means the difference between the observed transit time and the transit time expected for a regular periodic orbit.
  • transit_times – A list (with on entry for each planet) of transit times for each transit of each planet in units of days. These times will be used to compute the implied (least squares) ttv_period and t0. It is possible to supply a separate period parameter that will set the shape of the transits, but care should be taken to make sure that period and ttv_period don’t diverge because things will break if the time between neighboring transits is larger than 2*period.
  • transit_inds – A list of integer value tensors giving the transit number for each transit in transit_times'' or ``ttvs. This is useful when not all transits are observed. This should always be zero indexed.
  • delta_log_period – If using the transit_times argument, this parameter specifies the difference (in natural log) between the leqast squares period and the effective period of the transit.
get_planet_position(t, parallax=None)

The planets’ positions in the barycentric frame

Parameters:t – The times where the position should be evaluated.
Returns:The components of the position vector at t in units of R_sun.
get_planet_velocity(t)

Get the planets’ velocity vectors

Parameters:t – The times where the velocity should be evaluated.
Returns:The components of the velocity vector at t in units of M_sun/day.
get_radial_velocity(t, K=None, output_units=None)

Get the radial velocity of the star

Note

The convention in exoplanet is that positive z points towards the observer. However, for consistency with radial velocity literature this method returns values where positive radial velocity corresponds to a redshift as expected.

Parameters:
  • t – The times where the radial velocity should be evaluated.
  • K (Optional) – The semi-amplitudes of the orbits. If provided, the m_planet and incl parameters will be ignored and this amplitude will be used instead.
  • output_units (Optional) – An AstroPy velocity unit. If not given, the output will be evaluated in m/s. This is ignored if a value is given for K.
Returns:

The reflex radial velocity evaluated at t in units of output_units. For multiple planets, this will have one row for each planet.

get_relative_angles(t, parallax=None)

The planets’ relative position to the star in the sky plane, in separation, position angle coordinates.

Note

This treats each planet independently and does not take the other planets into account when computing the position of the star. This is fine as long as the planet masses are small.

Parameters:t – The times where the position should be evaluated.
Returns:The separation (arcseconds) and position angle (radians, measured east of north) of the planet relative to the star.
get_relative_position(t, parallax=None)

The planets’ positions relative to the star in the X,Y,Z frame.

Note

This treats each planet independently and does not take the other planets into account when computing the position of the star. This is fine as long as the planet masses are small. In other words, the reflex motion of the star caused by the other planets is neglected when computing the relative coordinates of one of the planets.

Parameters:t – The times where the position should be evaluated.
Returns:The components of the position vector at t in units of R_sun.
get_relative_velocity(t)

The planets’ velocity relative to the star

Note

This treats each planet independently and does not take the other planets into account when computing the position of the star. This is fine as long as the planet masses are small.

Parameters:t – The times where the velocity should be evaluated.
Returns:The components of the velocity vector at t in units of R_sun/day.
get_star_position(t, parallax=None)

The star’s position in the barycentric frame

Note

If there are multiple planets in the system, this will return one column per planet with each planet’s contribution to the motion. The star’s full position can be computed by summing over the last axis.

Parameters:t – The times where the position should be evaluated.
Returns:The components of the position vector at t in units of R_sun.
get_star_velocity(t)

Get the star’s velocity vector

Note

For a system with multiple planets, this will return one column per planet with the contributions from each planet. The total velocity can be found by summing along the last axis.

Parameters:t – The times where the velocity should be evaluated.
Returns:The components of the velocity vector at t in units of M_sun/day.
in_transit(t, r=0.0, texp=None)

Get a list of timestamps that are in transit

Parameters:
  • t (vector) – A vector of timestamps to be evaluated.
  • r (Optional) – The radii of the planets.
  • texp (Optional[float]) – The exposure time.
Returns:

The indices of the timestamps that are in transit.

class exoplanet.orbits.ReboundOrbit(*args, **kwargs)

An N-body system powered by the rebound integrator

This takes all the same arguments as the KeplerianOrbit, but these arguments define the orbital elements at some reference time (given by the rebound_t parameter). The positions and velocities of the bodies are then computed by numerically integrating the gravitational N-body system.

rebound-specific parameters can be provided as keyword arguments prefixed by rebound_. These will then be applied to the rebound.Simulation object as properties. Therefore, if you want to change the integrator, you could use: rebound_integrator = "whfast", for example. All of these parameters are passed directly through to rebound except rebound_t (the reference time) which is converted from days to years over two pi (the default time units in rebound).

Note

exoplanet and rebound use different base units, but all of the unit conversions are handled behind the scenes in this object so that means that you should mostly use the exoplanet units when interacting with this class and you should be very cautious about setting the rebound_G argument. One example of a case where you’ll need to use the rebound units is when you want to set the integrator step size using the rebound_dt parameter.

get_planet_position(t)

The planets’ positions in the barycentric frame

Parameters:t – The times where the position should be evaluated.
Returns:The components of the position vector at t in units of R_sun.
get_planet_velocity(t)

Get the planets’ velocity vectors

Parameters:t – The times where the velocity should be evaluated.
Returns:The components of the velocity vector at t in units of M_sun/day.
get_radial_velocity(t, output_units=None)

Get the radial velocity of the star

Note

The convention in exoplanet is that positive z points towards the observer. However, for consistency with radial velocity literature this method returns values where positive radial velocity corresponds to a redshift as expected.

Note

Unlike the KeplerianOrbit implementation, the semi-amplitude K cannot be used with the ReboundOrbit. Also, the contributions of each planet are not returned separately; this will always return a single time series.

Parameters:
  • t – The times where the radial velocity should be evaluated.
  • output_units (Optional) – An AstroPy velocity unit. If not given, the output will be evaluated in m/s. This is ignored if a value is given for K.
Returns:

The reflex radial velocity evaluated at t in units of output_units.

get_relative_angles(t, parallax=None)

The planets’ relative position to the star in the sky plane, in separation, position angle coordinates.

Note

This treats each planet independently and does not take the other planets into account when computing the position of the star. This is fine as long as the planet masses are small.

Parameters:t – The times where the position should be evaluated.
Returns:The separation (arcseconds) and position angle (radians, measured east of north) of the planet relative to the star.
get_relative_position(t)

The planets’ positions relative to the star in the X,Y,Z frame.

Parameters:t – The times where the position should be evaluated.
Returns:The components of the position vector at t in units of R_sun.
get_relative_velocity(t)

The planets’ velocity relative to the star

Parameters:t – The times where the velocity should be evaluated.
Returns:The components of the velocity vector at t in units of R_sun/day.
get_star_position(t)

The star’s position in the barycentric frame

Note

Unlike the KeplerianOrbit, this will not return the contributions from each planet separately.

Parameters:t – The times where the position should be evaluated.
Returns:The components of the position vector at t in units of R_sun.
get_star_velocity(t)

Get the star’s velocity vector

Note

Unlike the KeplerianOrbit, this will not return the contributions from each planet separately.

Parameters:t – The times where the velocity should be evaluated.
Returns:The components of the velocity vector at t in units of M_sun/day.
in_transit(t, r=0.0, texp=None)

This is a no-op and all points are assumed to be in transit

class exoplanet.orbits.SimpleTransitOrbit(period=None, t0=0.0, b=0.0, duration=None, r_star=1.0)

An orbit representing a set of planets transiting a common central

This orbit is parameterized by the observables of a transiting system, period, phase, duration, and impact parameter.

Parameters:
  • period – The orbital period of the planets in days.
  • t0 – The midpoint time of a reference transit for each planet in days.
  • b – The impact parameters of the orbits.
  • duration – The durations of the transits in days.
  • r_star – The radius of the star in R_sun.
get_relative_position(t)

The planets’ positions relative to the star

Parameters:t – The times where the position should be evaluated.
Returns:The components of the position vector at t in units of R_sun.
in_transit(t, r=None, texp=None)

Get a list of timestamps that are in transit

Parameters:
  • t (vector) – A vector of timestamps to be evaluated.
  • r (Optional) – The radii of the planets.
  • texp (Optional[float]) – The exposure time.
Returns:

The indices of the timestamps that are in transit.

Light curve models

class exoplanet.LimbDarkLightCurve(u, model=None)

A limb darkened light curve computed using starry

Parameters:u (vector) – A vector of limb darkening coefficients.
get_light_curve(orbit=None, r=None, t=None, texp=None, oversample=7, order=0, use_in_transit=True)

Get the light curve for an orbit at a set of times

Parameters:
  • orbit – An object with a get_relative_position method that takes a tensor of times and returns a list of Cartesian coordinates of a set of bodies relative to the central source. This method should return three tensors (one for each coordinate dimension) and each tensor should have the shape append(t.shape, r.shape) or append(t.shape, oversample, r.shape) when texp is given. The first two coordinate dimensions are treated as being in the plane of the sky and the third coordinate is the line of sight with positive values pointing away from the observer. For an example, take a look at orbits.KeplerianOrbit.
  • r (tensor) – The radius of the transiting body in the same units as r_star. This should have a shape that is consistent with the coordinates returned by orbit. In general, this means that it should probably be a scalar or a vector with one entry for each body in orbit.
  • t (tensor) – The times where the light curve should be evaluated.
  • texp (Optional[tensor]) – The exposure time of each observation. This can be a scalar or a tensor with the same shape as t. If texp is provided, t is assumed to indicate the timestamp at the middle of an exposure of length texp.
  • oversample (Optional[int]) – The number of function evaluations to use when numerically integrating the exposure time.
  • order (Optional[int]) – The order of the numerical integration scheme. This must be one of the following: 0 for a centered Riemann sum (equivalent to the “resampling” procedure suggested by Kipping 2010), 1 for the trapezoid rule, or 2 for Simpson’s rule.
  • use_in_transit (Optional[bool]) – If True, the model will only be evaluated for the data points expected to be in transit as computed using the in_transit method on orbit.

Scalable Gaussian processes

class exoplanet.gp.GP(kernel, x, diag=None, mean=<exoplanet.gp.means.Zero object>, J=-1, model=None)

The interface for computing Gaussian Process models with celerite

This class implements the method described in Foreman-Mackey et al. (2017) and Foreman-Mackey (2018) for scalable evaluation of Gaussian Process (GP) models in 1D.

Note

The input coordinates x must be sorted in ascending order, but this is not checked in the code. If the values are not sorted, the behavior of the algorithm is undefined.

Parameters:
  • kernel – A exoplanet.gp.terms.Term object the specifies the GP kernel.
  • x – The input coordinates. This should be a 1D array and the elements must be sorted. Otherwise the results are undefined.
  • diag (Optional) – The extra diagonal to add to the covariance matrix. This should have the same length as x and correspond to the excess variance for each data point. Note: this is different from the usage in the celerite package where the standard deviation (instead of variance) is provided.
  • mean (Optional) – The mean function for the GP. This can be a constant scalar value or a callable that will be called with a single, one dimensional tensor argument specifying the input coordinates where the mean should be evaluated.
  • J (Optional) –

    The width of the system. This is the J parameter from Foreman-Mackey (2018) (not the original paper) so a real term contributes J += 1 and a complex term contributes J += 2. If you know this value in advance, you can provide it. Otherwise, the code will try to work it out.

marginal(name, **kwargs)

The marginal likelihood

Parameters:
  • name – The name of the node
  • observed – The observed data
Returns:

A pymc3.DensityDist with the likelihood

predict(t=None, return_var=False, return_cov=False, predict_mean=False, kernel=None, _fast_mean=True)

Compute the conditional distribution

Parameters:
  • t (optional) – The independent coordinates where the prediction should be evaluated. If not provided, this will be evaluated at the observations.
  • return_var (bool, optional) – Return the variance of the conditional distribution.
  • return_cov (bool, optional) – Return the full covariance matrix of the conditional distribution.
  • predict_mean (bool, optional) – Include the mean function in the prediction.
  • kernel (optional) – If provided, compute the conditional distribution using a different kernel. This is generally used to separate the contributions from different model components.
Raises:

RuntimeError – if GP.condition() or marginal() are not called first.

class exoplanet.gp.terms.Term(**kwargs)

The abstract base “term” that is the superclass of all other terms

Subclasses should overload the terms.Term.get_real_coefficients() and terms.Term.get_complex_coefficients() methods.

class exoplanet.gp.terms.RealTerm(**kwargs)

The simplest celerite term

This term has the form

\[k(\tau) = a_j\,e^{-c_j\,\tau}\]

with the parameters a and c.

Strictly speaking, for a sum of terms, the parameter a could be allowed to go negative but since it is somewhat subtle to ensure positive definiteness, we recommend keeping both parameters strictly positive. Advanced users can build a custom term that has negative coefficients but care should be taken to ensure positivity.

Parameters:
  • a or log_a (tensor) – The amplitude of the term.
  • c or log_c (tensor) – The exponent of the term.
class exoplanet.gp.terms.ComplexTerm(**kwargs)

A general celerite term

This term has the form

\[k(\tau) = \frac{1}{2}\,\left[(a_j + b_j)\,e^{-(c_j+d_j)\,\tau} + (a_j - b_j)\,e^{-(c_j-d_j)\,\tau}\right]\]

with the parameters a, b, c, and d.

This term will only correspond to a positive definite kernel (on its own) if \(a_j\,c_j \ge b_j\,d_j\).

Parameters:
  • a or log_a (tensor) – The real part of amplitude.
  • b or log_b (tensor) – The imaginary part of amplitude.
  • c or log_c (tensor) – The real part of the exponent.
  • d or log_d (tensor) – The imaginary part of exponent.
class exoplanet.gp.terms.SHOTerm(*args, **kwargs)

A term representing a stochastically-driven, damped harmonic oscillator

The PSD of this term is

\[S(\omega) = \sqrt{\frac{2}{\pi}} \frac{S_0\,\omega_0^4} {(\omega^2-{\omega_0}^2)^2 + {\omega_0}^2\,\omega^2/Q^2}\]

with the parameters S0, Q, and w0.

Parameters:
  • S0 or log_S0 (tensor) – The parameter \(S_0\).
  • Q or log_Q (tensor) – The parameter \(Q\).
  • w0 or log_w0 (tensor) – The parameter \(\omega_0\).
  • Sw4 or log_Sw4 (tensor) – It can sometimes be more efficient to parameterize the amplitude of a SHO kernel using \(S_0\,{\omega_0}^4\) instead of \(S_0\) directly since \(S_0\) and \(\omega_0\) are strongly correlated. If provided, S0 will be computed from Sw4 and w0.
class exoplanet.gp.terms.Matern32Term(**kwargs)

A term that approximates a Matern-3/2 function

This term is defined as

\[k(\tau) = \sigma^2\,\left[ \left(1+1/\epsilon\right)\,e^{-(1-\epsilon)\sqrt{3}\,\tau/\rho} \left(1-1/\epsilon\right)\,e^{-(1+\epsilon)\sqrt{3}\,\tau/\rho} \right]\]

with the parameters sigma and rho. The parameter eps controls the quality of the approximation since, in the limit \(\epsilon \to 0\) this becomes the Matern-3/2 function

\[\lim_{\epsilon \to 0} k(\tau) = \sigma^2\,\left(1+ \frac{\sqrt{3}\,\tau}{\rho}\right)\, \exp\left(-\frac{\sqrt{3}\,\tau}{\rho}\right)\]
Parameters:
  • sigma or log_sigma (tensor) – The parameter \(\sigma\).
  • rho or log_rho (tensor) – The parameter \(\rho\).
  • eps (Optional[float]) – The value of the parameter \(\epsilon\). (default: 0.01)
class exoplanet.gp.terms.RotationTerm(**kwargs)

A mixture of two SHO terms that can be used to model stellar rotation

This term has two modes in Fourier space: one at period and one at 0.5 * period. This can be a good descriptive model for a wide range of stochastic variability in stellar time series from rotation to pulsations.

Parameters:
  • amp or log_amp (tensor) – The amplitude of the variability.
  • period or log_period (tensor) – The primary period of variability.
  • Q0 or log_Q0 (tensor) – The quality factor (or really the quality factor minus one half) for the secondary oscillation.
  • deltaQ or log_deltaQ (tensor) – The difference between the quality factors of the first and the second modes. This parameterization (if deltaQ > 0) ensures that the primary mode alway has higher quality.
  • mix – The fractional amplitude of the secondary mode compared to the primary. This should probably always be 0 < mix < 1.

Estimators

exoplanet.estimate_semi_amplitude(periods, x, y, yerr=None, t0s=None)

Estimate the RV semi-amplitudes for planets in an RV series

Parameters:
  • periods – The periods of the planets. Assumed to be in days if not an AstroPy Quantity.
  • x – The observation times. Assumed to be in days if not an AstroPy Quantity.
  • y – The radial velocities. Assumed to be in m/s if not an AstroPy Quantity.
  • yerr (Optional) – The uncertainty on y.
  • t0s (Optional) – The time of a reference transit for each planet, if known.
Returns:

An estimate of the semi-amplitude of each planet in units of m/s.

exoplanet.estimate_minimum_mass(periods, x, y, yerr=None, t0s=None, m_star=1)

Estimate the minimum mass(es) for planets in an RV series

Parameters:
  • periods – The periods of the planets. Assumed to be in days if not an AstroPy Quantity.
  • x – The observation times. Assumed to be in days if not an AstroPy Quantity.
  • y – The radial velocities. Assumed to be in m/s if not an AstroPy Quantity.
  • yerr (Optional) – The uncertainty on y.
  • t0s (Optional) – The time of a reference transit for each planet, if known.
  • m_star (Optional) – The mass of the star. Assumed to be in M_sun if not an AstroPy Quantity.
Returns:

An estimate of the minimum mass of each planet as an AstroPy Quantity with units of M_jupiter.

exoplanet.lomb_scargle_estimator(x, y, yerr=None, min_period=None, max_period=None, filter_period=None, max_peaks=2, **kwargs)

Estimate period of a time series using the periodogram

Parameters:
  • x (ndarray[N]) – The times of the observations
  • y (ndarray[N]) – The observations at times x
  • yerr (Optional[ndarray[N]]) – The uncertainties on y
  • min_period (Optional[float]) – The minimum period to consider
  • max_period (Optional[float]) – The maximum period to consider
  • filter_period (Optional[float]) – If given, use a high-pass filter to down-weight period longer than this
  • max_peaks (Optional[int]) – The maximum number of peaks to return (default: 2)
Returns:

A dictionary with the computed periodogram and the parameters for up to max_peaks peaks in the periodogram.

exoplanet.autocorr_estimator(x, y, yerr=None, min_period=None, max_period=None, oversample=2.0, smooth=2.0, max_peaks=10)

Estimate the period of a time series using the autocorrelation function

Note

The signal is interpolated onto a uniform grid in time so that the autocorrelation function can be computed.

Parameters:
  • x (ndarray[N]) – The times of the observations
  • y (ndarray[N]) – The observations at times x
  • yerr (Optional[ndarray[N]]) – The uncertainties on y
  • min_period (Optional[float]) – The minimum period to consider
  • max_period (Optional[float]) – The maximum period to consider
  • oversample (Optional[float]) – When interpolating, oversample the times by this factor (default: 2.0)
  • smooth (Optional[float]) – Smooth the autocorrelation function by this factor times the minimum period (default: 2.0)
  • max_peaks (Optional[int]) – The maximum number of peaks to identify in the autocorrelation function (default: 10)
Returns:

A dictionary with the computed autocorrelation function and the estimated period. For compatibility with the lomb_scargle_estimator(), the period is returned as a list with the key peaks.

exoplanet.bls_estimator(x, y, yerr=None, duration=0.2, min_period=None, max_period=None, objective=None, method=None, oversample=10, **kwargs)

Estimate the period of a time series using box least squares

All extra keyword arguments are passed directly to astropy.timeseries.BoxLeastSquares.autopower().

Parameters:
  • x (ndarray[N]) – The times of the observations
  • y (ndarray[N]) – The observations at times x
  • yerr (Optional[ndarray[N]]) – The uncertainties on y
  • min_period (Optional[float]) – The minimum period to consider
  • max_period (Optional[float]) – The maximum period to consider
Returns:

A dictionary with the computed autocorrelation function and the estimated period. For compatibility with the lomb_scargle_estimator(), the period is returned as a list with the key peaks.

Distributions

Base distributions

class exoplanet.distributions.UnitUniform(*args, **kwargs)

A uniform distribution between zero and one

This can sometimes get better performance than pm.Uniform.dist(0, 1).

class exoplanet.distributions.UnitVector(*args, **kwargs)

A vector where the sum of squares is fixed to unity

For a multidimensional shape, the normalization is performed along the last dimension.

class exoplanet.distributions.Angle(*args, **kwargs)

An angle constrained to be in the range -pi to pi

The actual sampling is performed in the two dimensional vector space (sin(theta), cos(theta)) so that the sampler doesn’t see a discontinuity at pi.

class exoplanet.distributions.Periodic(lower=0, upper=1, **kwargs)

An periodic parameter in a given range

Like the Angle distribution, the actual sampling is performed in a two dimensional vector space (sin(theta), cos(theta)) and then transformed into the range [lower, upper).

Parameters:
  • lower – The lower bound on the range.
  • upper – The upper bound on the range.

Physical distributions

class exoplanet.distributions.QuadLimbDark(*args, **kwargs)

An uninformative prior for quadratic limb darkening parameters

This is an implementation of the Kipping (2013) reparameterization of the two-parameter limb darkening model to allow for efficient and uninformative sampling.

class exoplanet.distributions.ImpactParameter(ror=None, **kwargs)

The impact parameter distribution for a transiting planet

Parameters:ror – A scalar, tensor, or PyMC3 distribution representing the radius ratio between the planet and star. Conditioned on a value of ror, this will be uniformly distributed between 0 and 1+ror.

Eccentricity distributions

exoplanet.distributions.eccentricity.kipping13(name, fixed=False, long=None, model=None, **kwargs)

The beta eccentricity distribution fit by Kipping (2013)

The beta distribution parameters fit by Kipping (2013).

Parameters:
  • name (str) – The name of the eccentricity variable.
  • fixed (bool, optional) – If True, use the posterior median hyperparameters. Otherwise, marginalize over the parameters.
  • long (bool, optional) – If True, use the parameters for the long period fit. If False, use the parameters for the short period fit. If not given, the parameters fit using the full dataset are used.
Returns:

The eccentricity distribution.

exoplanet.distributions.eccentricity.vaneylen19(name, fixed=False, multi=False, model=None, **kwargs)

The eccentricity distribution for small planets

The mixture distribution fit by Van Eylen et al. (2019) to a population of well-characterized small transiting planets observed by Kepler.

Parameters:
  • name (str) – The name of the eccentricity variable.
  • fixed (bool, optional) – If True, use the posterior median hyperparameters. Otherwise, marginalize over the parameters.
  • multi (bool, optional) – If True, use the distribution for systems with multiple transiting planets. If False (default), use the distribution for systems with only one detected transiting planet.
Returns:

The eccentricity distribution.

Utilities

exoplanet.optimize(start=None, vars=None, model=None, return_info=False, verbose=True, **kwargs)

Maximize the log prob of a PyMC3 model using scipy

All extra arguments are passed directly to the scipy.optimize.minimize function.

Parameters:
  • start – The PyMC3 coordinate dictionary of the starting position
  • vars – The variables to optimize
  • model – The PyMC3 model
  • return_info – Return both the coordinate dictionary and the result of scipy.optimize.minimize
  • verbose – Print the success flag and log probability to the screen
exoplanet.eval_in_model(var, point=None, return_func=False, model=None, **kwargs)

Evaluate a Theano tensor or PyMC3 variable in a PyMC3 model

This method builds a Theano function for evaluating a node in the graph given the required parameters. This will also cache the compiled Theano function in the current pymc3.Model to reduce the overhead of calling this function many times.

Parameters:
  • var – The variable or tensor to evaluate.
  • point (Optional) – A dict of input parameter values. This can be model.test_point (default), the result of pymc3.find_MAP, a point in a pymc3.MultiTrace or any other representation of the input parameters.
  • return_func (Optional[bool]) – If False (default), return the evaluated variable. If True, return the result, the Theano function and the list of arguments for that function.
Returns:

Depending on return_func, either the value of var at point, or this value, the Theano function, and the input arguments.

exoplanet.get_samples_from_trace(trace, size=1)

Generate random samples from a PyMC3 MultiTrace

Parameters:
  • trace – The MultiTrace.
  • size – The number of samples to generate.
exoplanet.get_dense_nuts_step(start=None, adaptation_window=101, doubling=True, model=None, **kwargs)

Get a NUTS step function with a dense mass matrix

The entries in the mass matrix will be tuned based on the sample covariances during tuning. All extra arguments are passed directly to pymc3.NUTS.

Parameters:
  • start (dict, optional) – A starting point in parameter space. If not provided, the model’s test_point is used.
  • adaptation_window (int, optional) – The (initial) size of the window used for sample covariance estimation.
  • doubling (bool, optional) – If True (default) the adaptation window is doubled each time the matrix is updated.
exoplanet.orbits.ttv.compute_expected_transit_times(min_time, max_time, period, t0)

Compute the expected transit times within a dataset

Parameters:
  • min_time (float) – The start time of the dataset
  • max_time (float) – The end time of the dataset
  • period (array) – The periods of the planets
  • t0 (array) – The reference transit times for the planets
Returns:

A list of arrays of expected transit times for each planet

Units

exoplanet.units.with_unit(obj, unit)

Decorate a Theano tensor with Astropy units

Parameters:
  • obj – The Theano tensor
  • unit (astropy.Unit) – The units for this object
Raises:

TypeError – If the tensor already has units

exoplanet.units.has_unit(obj)

Does an object have units as defined by exoplanet?

exoplanet.units.to_unit(obj, target)

Convert a Theano tensor with units to a target set of units

Parameters:
  • obj – The Theano tensor
  • target (astropy.Unit) – The target units
Returns:

A Theano tensor in the right units

Citations

exoplanet.citations.get_citations_for_model(model=None, width=79)

Get the citations for the components used an exoplanet PyMC3

Returns: The acknowledgement text for exoplanet and its dependencies and a string containing the BibTeX entries for the citations in the acknowledgement.