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, sin_omega=None, cos_omega=None, Omega=None, m_planet=0.0, m_star=None, r_star=None, rho_star=None, ror=None, model=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
, anddays
) except forrho_star
which will be ing / cm^3
.There are only specific combinations of input parameters that can be used:
First, either
period
ora
must be given. If values are given for both parameters, then neitherm_star
orrho_star
can be defined because the stellar density implied by each planet will be computed inrho_star
.Only one of
incl
andb
can be given.If a value is given for
ecc
thenomega
must also be given.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 be1 * R_sun
. Otherwise, at most two ofm_star
,r_star
, andrho_star
can be defined.Either
t0
(reference transit) ort_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 isM_sun
.rho_star_units – An
astropy.units
compatible unit object giving the units of the stellar density. If not given, the default isg / cm^3
.
- get_planet_position(t, parallax=None, light_delay=False)#
The planets’ positions in the barycentric frame
- Parameters:
t – The times where the position should be evaluated.
light_delay – account for the light travel time delay? Default is False.
- Returns:
The components of the position vector at
t
in units ofR_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 ofM_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
andincl
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 forK
.
- Returns:
The reflex radial velocity evaluated at
t
in units ofoutput_units
. For multiple planets, this will have one row for each planet.
- get_relative_angles(t, parallax=None, light_delay=False)#
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.
light_delay – account for the light travel time delay? Default is False.
- 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, light_delay=False)#
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.
light_delay – account for the light travel time delay? Default is False.
- Returns:
The components of the position vector at
t
in units ofR_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 ofR_sun/day
.
- get_star_position(t, parallax=None, light_delay=False)#
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.
light_delay – account for the light travel time delay? Default is False.
- Returns:
The components of the position vector at
t
in units ofR_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 ofM_sun/day
.
- in_transit(t, r=0.0, texp=None, light_delay=False)#
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
ortransit_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
andt0
. It is possible to supply a separateperiod
parameter that will set the shape of the transits, but care should be taken to make sure thatperiod
andttv_period
don’t diverge because things will break if the time between neighboring transits is larger than2*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, light_delay=False)#
The planets’ positions in the barycentric frame
- Parameters:
t – The times where the position should be evaluated.
light_delay – account for the light travel time delay? Default is False.
- Returns:
The components of the position vector at
t
in units ofR_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 ofM_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
andincl
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 forK
.
- Returns:
The reflex radial velocity evaluated at
t
in units ofoutput_units
. For multiple planets, this will have one row for each planet.
- get_relative_angles(t, parallax=None, light_delay=False)#
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.
light_delay – account for the light travel time delay? Default is False.
- 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, light_delay=False)#
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.
light_delay – account for the light travel time delay? Default is False.
- Returns:
The components of the position vector at
t
in units ofR_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 ofR_sun/day
.
- get_star_position(t, parallax=None, light_delay=False)#
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.
light_delay – account for the light travel time delay? Default is False.
- Returns:
The components of the position vector at
t
in units ofR_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 ofM_sun/day
.
- in_transit(t, r=0.0, texp=None, light_delay=False)#
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.SimpleTransitOrbit(period, duration, t0=0.0, b=0.0, r_star=1.0, ror=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, light_delay=False)#
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 ofR_sun
.
- in_transit(t, r=None, texp=None, light_delay=False)#
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(u1, u2=None, model=None)#
A quadratically limb darkened light curve
Note
Earlier versions of exoplanet supported higher (and lower) order limb darkening, but this support was removed in v0.5.0. For higher order limb darkening profiles, use starry directly.
- Parameters:
u1 (scalar) – The first limb darkening coefficient
u2 (scalar) – The second limb darkening coefficient
- get_light_curve(orbit=None, r=None, t=None, texp=None, oversample=7, order=0, use_in_transit=None, light_delay=False)#
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 shapeappend(t.shape, r.shape)
orappend(t.shape, oversample, r.shape)
whentexp
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 atorbits.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 byorbit
. In general, this means that it should probably be a scalar or a vector with one entry for each body inorbit
. Note that this is a different quantity than the planet-to-star radius ratio; do not confuse the two!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
. Iftexp
is provided,t
is assumed to indicate the timestamp at the middle of an exposure of lengthtexp
.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, or2
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 thein_transit
method onorbit
.
- get_ror_from_approx_transit_depth(delta, b, jac=False)#
Get the radius ratio corresponding to a particular transit depth
This result will be approximate and it requires
|b| < 1
because it relies on the small planet approximation.- Parameters:
delta (tensor) – The approximate transit depth in relative units
b (tensor) – The impact parameter
jac (bool) – If true, the Jacobian
d ror / d delta
is also returned
- Returns:
The radius ratio that approximately corresponds to the depth
delta
at impact parameterb
.- Return type:
ror
- class exoplanet.SecondaryEclipseLightCurve(u_primary, u_secondary, surface_brightness_ratio, model=None)#
The light curve for a secondary eclipse model computed using starry
- Parameters:
u_primary (vector) – A vector of limb darkening coefficients for the primary body.
u_secondary (vector) – A vector of limb darkening coefficients for the secondary body.
surface_brightness_ratio (scalar) – The surface brightness ratio of the secondary with respect to the primary.
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 tomax_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 keypeaks
.
- 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 keypeaks
.
Distributions#
Physical distributions#
- class exoplanet.distributions.QuadLimbDark(name, *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(name, *args, **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 between0
and1+ror
.
Eccentricity distributions#
- exoplanet.distributions.eccentricity.kipping13(name, fixed=False, long=None, lower=None, upper=None, model=None, **kwargs)#
The beta eccentricity distribution fit by Kipping (2013)
The beta distribution parameters fit by Kipping (2013b).
- 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. IfFalse
, use the parameters for the short period fit. If not given, the parameters fit using the full dataset are used.lower (float, optional) – Restrict the eccentricity to be larger than this value.
upper (float, optional) – Restrict the eccentricity to be smaller than this value.
- Returns:
The eccentricity distribution.
- exoplanet.distributions.eccentricity.vaneylen19(name, fixed=False, multi=False, lower=None, upper=None, 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. IfFalse
(default), use the distribution for systems with only one detected transiting planet.lower (float, optional) – Restrict the eccentricity to be larger than this value.
upper (float, optional) – Restrict the eccentricity to be smaller than this value.
- Returns:
The eccentricity distribution.
Miscellaneous#
- 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
- 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
- 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.