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, 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:
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
.
Only one of incl
and b
can be given.
If a value is given for ecc
then omega
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 be 1 * R_sun
. Otherwise, at most two of m_star
,
r_star
, and rho_star
can be defined.
Either t0
(reference transit) or t_periastron
must be given,
but not both.
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, light_delay=False)¶The planets’ positions in the barycentric frame
t – The times where the position should be evaluated.
light_delay – account for the light travel time delay? Default is False.
The components of the position vector at t
in units of
R_sun
.
get_planet_velocity
(t)¶Get the planets’ velocity vectors
t – The times where the velocity should be evaluated.
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.
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
.
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, 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.
t – The times where the position should be evaluated.
light_delay – account for the light travel time delay? Default is False.
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.
t – The times where the position should be evaluated.
light_delay – account for the light travel time delay? Default is False.
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.
t – The times where the velocity should be evaluated.
The components of the velocity vector at t
in units of
R_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.
t – The times where the position should be evaluated.
light_delay – account for the light travel time delay? Default is False.
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.
t – The times where the velocity should be evaluated.
The components of the velocity vector at t
in units of
M_sun/day
.
in_transit
(t, r=0.0, texp=None, light_delay=False)¶Get a list of timestamps that are in transit
t (vector) – A vector of timestamps to be evaluated.
r (Optional) – The radii of the planets.
texp (Optional[float]) – The exposure time.
The indices of the timestamps that are in transit.
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
.
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, light_delay=False)¶The planets’ positions in the barycentric frame
t – The times where the position should be evaluated.
light_delay – account for the light travel time delay? Default is False.
The components of the position vector at t
in units of
R_sun
.
get_planet_velocity
(t)¶Get the planets’ velocity vectors
t – The times where the velocity should be evaluated.
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.
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
.
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, 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.
t – The times where the position should be evaluated.
light_delay – account for the light travel time delay? Default is False.
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.
t – The times where the position should be evaluated.
light_delay – account for the light travel time delay? Default is False.
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.
t – The times where the velocity should be evaluated.
The components of the velocity vector at t
in units of
R_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.
t – The times where the position should be evaluated.
light_delay – account for the light travel time delay? Default is False.
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.
t – The times where the velocity should be evaluated.
The components of the velocity vector at t
in units of
M_sun/day
.
in_transit
(t, r=0.0, texp=None, light_delay=False)¶Get a list of timestamps that are in transit
t (vector) – A vector of timestamps to be evaluated.
r (Optional) – The radii of the planets.
texp (Optional[float]) – The exposure time.
The indices of the timestamps that are in transit.
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, light_delay=False)¶The planets’ positions in the barycentric frame
t – The times where the position should be evaluated.
The components of the position vector at t
in units of
R_sun
.
get_planet_velocity
(t)¶Get the planets’ velocity vectors
t – The times where the velocity should be evaluated.
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.
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
.
The reflex radial velocity evaluated at t
in units of
output_units
.
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.
t – The times where the position should be evaluated.
light_delay – account for the light travel time delay? Default is False.
The separation (arcseconds) and position angle (radians, measured east of north) of the planet relative to the star.
get_relative_position
(t, light_delay=False)¶The planets’ positions relative to the star in the X,Y,Z frame.
t – The times where the position should be evaluated.
The components of the position vector at t
in units of
R_sun
.
get_relative_velocity
(t)¶The planets’ velocity relative to the star
t – The times where the velocity should be evaluated.
The components of the velocity vector at t
in units of
R_sun/day
.
get_star_position
(t, light_delay=False)¶The star’s position in the barycentric frame
Note
Unlike the KeplerianOrbit
, this will not return
the contributions from each planet separately.
t – The times where the position should be evaluated.
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.
t – The times where the velocity should be evaluated.
The components of the velocity vector at t
in units of
M_sun/day
.
in_transit
(t, r=0.0, texp=None, light_delay=False)¶This is a no-op and all points are assumed to be in transit
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.
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
t – The times where the position should be evaluated.
The components of the position vector at t
in units of
R_sun
.
in_transit
(t, r=None, texp=None, light_delay=False)¶Get a list of timestamps that are in transit
t (vector) – A vector of timestamps to be evaluated.
r (Optional) – The radii of the planets.
texp (Optional[float]) – The exposure time.
The indices of the timestamps that are in transit.
exoplanet.
LimbDarkLightCurve
(u, model=None)¶A limb darkened light curve computed using starry
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=None, light_delay=False)¶Get the light curve for an orbit at a set of times
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
.
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.
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
The radius ratio that approximately corresponds to the depth
delta
at impact parameter b
.
ror
exoplanet.
SecondaryEclipseLightCurve
(u_primary, u_secondary, surface_brightness_ratio, model=None)¶The light curve for a secondary eclipse model computed using starry
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.
exoplanet.
estimate_semi_amplitude
(periods, x, y, yerr=None, t0s=None)¶Estimate the RV semi-amplitudes for planets in an RV series
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.
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
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.
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
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)
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.
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)
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()
.
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.distributions.
get_log_abs_det_jacobian
(in_vars, out_vars, **kwargs)¶Get the log absolute determinant of the Jacobian between parameter sets
This returns a Theano tensor representation of \(\log \left|\mathrm{det} J \right|\) where the elements of \(J\) are \(J_{nm} = \partial y_n / \partial x_m\). Here, \(y_n\) and \(x_m\) should be Theano tensors or PyMC3 variables.
in_vars – A list of input parameters \(x_m\)
out_vars – A list of output parameters \(y_n\)
exoplanet.distributions.
estimate_inverse_gamma_parameters
(lower, upper, target=0.01, initial=None, **kwargs)¶Estimate an inverse Gamma with desired tail probabilities
This method numerically solves for the parameters of an inverse Gamma distribution where the tails have a given probability. In other words \(P(x < \mathrm{lower}) = \mathrm{target}\) and similarly for the upper bound. More information can be found in part 4 of this blog post.
RuntimeError – If the solver does not converge.
A dictionary with the keys alpha
and beta
for the
parameters of the distribution.
exoplanet.distributions.
UnitUniform
(*args: Any, **kwargs: Any)¶A uniform distribution between zero and one
This can sometimes get better performance than pm.Uniform.dist(0, 1)
.
exoplanet.distributions.
UnitVector
(*args: Any, **kwargs: Any)¶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.
UnitDisk
(*args: Any, **kwargs: Any)¶Two dimensional parameters constrianed to live within the unit disk
This distribution is constrained such that the sum of squares along the zeroth axis will always be less than one. For example, in this code block:
import theano.tensor as tt
disk = UnitDisk("disk")
radius = tt.sum(disk ** 2, axis=0)
the tensor radius
will always have a value in the range [0, 1)
.
Note that the shape of this distribution must be two in the zeroth axis.
exoplanet.distributions.
Angle
(*args: Any, **kwargs: Any)¶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
(*args: Any, **kwargs: Any)¶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)
.
lower – The lower bound on the range.
upper – The upper bound on the range.
exoplanet.distributions.
QuadLimbDark
(*args: Any, **kwargs: Any)¶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.
ImpactParameter
(*args: Any, **kwargs: Any)¶The impact parameter distribution for a transiting planet
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
.
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).
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.
lower (float, optional) – Restrict the eccentricity to be larger than this value.
upper (float, optional) – Restrict the eccentricity to be smaller than this value.
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.
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.
lower (float, optional) – Restrict the eccentricity to be larger than this value.
upper (float, optional) – Restrict the eccentricity to be smaller than this value.
The eccentricity distribution.
exoplanet.
optimize
(start=None, vars=None, model=None, return_info=False, verbose=True, progress_bar=True, **kwargs)¶Maximize the log prob of a PyMC3 model using scipy
All extra arguments are passed directly to the scipy.optimize.minimize
function.
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
progress_bar – A tqdm
progress bar instance. Set to True
(default) to use tqdm.auto.tqdm()
. Set to False
to disable.
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.
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.
Depending on return_func
, either the value of var
at point
,
or this value, the Theano function, and the input arguments.
exoplanet.
get_theano_function_for_var
(var, model=None, **kwargs)¶Get a callable function to evaluate a component of a PyMC3 model
This should then be called using the arguments returned by the
get_args_for_theano_function()
function.
var – The model component to evaluate. This can be a Theano tensor, a PyMC3 variable, or a list of these
model (optional) – The PyMC3 model object. By default, the current context will be used.
exoplanet.
get_args_for_theano_function
(point=None, model=None)¶Get the arguments required to evaluate a PyMC3 model component
Use the result as the arguments for the callable returned by the
get_theano_function_for_var()
function.
point (dict, optional) – The point in parameter space where the model componenet should be evaluated
model (optional) – The PyMC3 model object. By default, the current context will be used.
exoplanet.
get_samples_from_trace
(trace, size=1)¶Generate random samples from a PyMC3 MultiTrace
trace – The MultiTrace
.
size – The number of samples to generate.
exoplanet.
get_dense_nuts_step
(start=None, adaptation_window=101, doubling=True, initial_weight=10, use_hessian=False, use_hessian_diag=False, hessian_regularization=1e-08, 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
.
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
exoplanet.units.
with_unit
(obj, unit)¶Decorate a Theano tensor with Astropy units
obj – The Theano tensor
unit (astropy.Unit) – The units for this object
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
obj – The Theano tensor
target (astropy.Unit) – The target units
A Theano tensor in the right units