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, 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:
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
.incl
and b
can be given.ecc
then omega
must also be given.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.t0
(reference transit) or t_periastron
must be given,
but not both.Parameters: |
|
---|
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 vector
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: |
|
---|---|
Returns: | The reflex radial velocity evaluated at |
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.
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 . |
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.
Parameters: |
|
---|
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 vector
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: |
|
---|---|
Returns: | The reflex radial velocity evaluated at |
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.
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 . |
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: |
|
---|
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 . |
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: |
|
---|
exoplanet.gp.
GP
(kernel, x, diag, 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: |
|
---|
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.
exoplanet.gp.terms.
RealTerm
(**kwargs)¶The simplest celerite term
This term has the form
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: |
|
---|
exoplanet.gp.terms.
ComplexTerm
(**kwargs)¶A general celerite term
This term has the form
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: |
|
---|
exoplanet.gp.terms.
SHOTerm
(*args, **kwargs)¶A term representing a stochastically-driven, damped harmonic oscillator
The PSD of this term is
with the parameters S0
, Q
, and w0
.
Parameters: |
|
---|
exoplanet.gp.terms.
Matern32Term
(**kwargs)¶A term that approximates a Matern-3/2 function
This term is defined as
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
Parameters: |
|
---|
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: |
|
---|
exoplanet.
estimate_semi_amplitude
(periods, x, y, yerr=None, t0s=None)¶Estimate the RV semi-amplitudes for planets in an RV series
Parameters: |
|
---|---|
Returns: | An estimate of the semi-amplitude of each planet in units of |
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: |
|
---|---|
Returns: | An estimate of the minimum mass of each planet as an AstroPy Quantity
with units of |
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: |
|
---|---|
Returns: | A dictionary with the computed |
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: |
|
---|---|
Returns: | A dictionary with the computed autocorrelation function and the
estimated period. For compatibility with the
|
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.
exoplanet.distributions.
UnitUniform
(*args, **kwargs)¶A uniform distribution between zero and one
This can sometimes get better performance than pm.Uniform.dist(0, 1)
.
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.
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: |
|
---|
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.
exoplanet.distributions.
RadiusImpact
(*args, **kwargs)¶The Espinoza (2018) distribution over radius and impact parameter
This is an implementation of Espinoza (2018) The first axis of the shape of the parameter should be exactly 2. The radius ratio will be in the zeroth entry in the first dimension and the impact parameter will be in the first.
Parameters: |
|
---|
exoplanet.distributions.
get_joint_radius_impact
(name='', N_planets=None, min_radius=0, max_radius=1, testval_r=None, testval_b=None, **kwargs)¶Get the joint distribution over radius and impact parameter
This uses the Espinoza (2018) parameterization of the distribution (see
distributions.RadiusImpact
for more details).
Parameters: |
|
---|---|
Returns: | Two |
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: |
|
---|
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: |
|
---|---|
Returns: | Depending on |
exoplanet.
get_samples_from_trace
(trace, size=1)¶Generate random samples from a PyMC3 MultiTrace
Parameters: |
|
---|
exoplanet.
PyMC3Sampler
(start=75, finish=50, window=25, dense=True, **kwargs)¶A sampling wrapper for PyMC3 with support for dense mass matrices
This schedule is based on the method used by as described in Section 34.2 of the Stan Manual.
All extra keyword arguments are passed to the pymc3.sample
function.
Parameters: |
|
---|
extend_tune
(steps, start=None, step_kwargs=None, trace=None, step=None, **kwargs)¶Extend the tuning phase by a given number of steps
After running the sampling, the mass matrix is re-estimated based on this run.
Parameters: | steps (int) – The number of steps to run. |
---|
get_step_for_trace
(trace=None, model=None, regular_window=0, regular_variance=0.001, **kwargs)¶Get a PyMC3 NUTS step tuned for a given burn-in trace
Parameters: |
|
---|
sample
(trace=None, step=None, start=None, step_kwargs=None, **kwargs)¶Run the production sampling using the tuned mass matrix
This is a light wrapper around pymc3.sample
and any arguments used
there (for example draws
) can be used as input to this method too.
tune
(tune=1000, start=None, step_kwargs=None, **kwargs)¶Run the full tuning run for the mass matrix
This will run start
steps of warmup followed by chains with
exponentially increasing chains to tune the mass matrix.
Parameters: | tune (int) – The total number of steps to run. |
---|
warmup
(start=None, step_kwargs=None, **kwargs)¶Run an initial warmup phase to find the typical set