Priors

class orbitize.priors.GaussianPrior(mu, sigma, no_negatives=True)[source]

Gaussian prior.

\[log(p(x|\sigma, \mu)) \propto \frac{(x - \mu)}{\sigma}\]
Parameters
  • mu (float) – mean of the distribution

  • sigma (float) – standard deviation of the distribution

  • no_negatives (bool) – if True, only positive values will be drawn from

  • prior (this) –

  • 0 (and the probability of negative values will be) –

  • (default – True).

(written) Sarah Blunt, 2018

compute_lnprob(element_array)[source]

Compute log(probability) of an array of numbers wrt a Gaussian distibution. Negative numbers return a probability of -inf.

Parameters

element_array (float or np.array of float) – array of numbers. We want the probability of drawing each of these from the appopriate Gaussian distribution

Returns

array of log(probability) values, corresponding to the probability of drawing each of the numbers in the input element_array.

Return type

numpy array of float

draw_samples(num_samples)[source]

Draw positive samples from a Gaussian distribution. Negative samples will not be returned.

Parameters

num_samples (float) – the number of samples to generate

Returns

samples drawn from the appropriate Gaussian distribution. Array has length num_samples.

Return type

numpy array of float

transform_samples(u)[source]

Transform uniform 1D samples, u, to samples drawn from a Gaussian distribution.

Parameters

u (array of floats) – list of samples with values 0 < u < 1.

Returns

1D u samples transformed to a Gaussian distribution.

Return type

numpy array of floats

class orbitize.priors.KDEPrior(gaussian_kde, total_params, bounds=[], log_scale_arr=[])[source]

Gaussian kernel density estimation (KDE) prior. This class is a wrapper for scipy.stats.gaussian_kde.

Parameters
  • gaussian_kde (scipy.stats.gaussian_kde) – scipy KDE object containing the KDE defined by the user

  • total_params (float) – number of parameters in the KDE

  • bounds (array_like of bool, optional) – bounds for the KDE out of which the prob returned is -Inf

  • bounds – if True for a parameter the parameter is fit to the KDE in log-scale

Written: Jorge LLop-Sayson, Sarah Blunt (2021)

compute_lnprob(element_array)[source]

Compute log(probability) of an array of numbers wrt a the defined KDE. Negative numbers return a probability of -inf.

Parameters

element_array (float or np.array of float) – array of numbers. We want the probability of drawing each of these from the KDE.

Returns

array of log(probability) values, corresponding to the probability of drawing each of the numbers in the input element_array.

Return type

numpy array of float

draw_samples(num_samples)[source]

Draw positive samples from the KDE. Negative samples will not be returned.

Parameters

num_samples (float) – the number of samples to generate.

Returns

samples drawn from the KDE distribution. Array has length num_samples.

Return type

numpy array of float

increment_param_num()[source]

Increment the index to evaluate the appropriate parameter.

class orbitize.priors.LinearPrior(m, b)[source]

Draw samples from the probability distribution:

\[p(x) \propto mx+b\]

where m is negative, b is positive, and the range is [0,-b/m].

Parameters
  • m (float) – slope of line. Must be negative.

  • b (float) – y intercept of line. Must be positive.

draw_samples(num_samples)[source]

Draw samples from a descending linear distribution.

Parameters

num_samples (float) – the number of samples to generate

Returns

samples ranging from [0, -b/m) as floats.

Return type

np.array

transform_samples(u)[source]

Transform uniform 1D samples, u, to samples drawn from a Linear distribution.

Parameters

u (array of floats) – list of samples with values 0 < u < 1.

Returns

1D u samples transformed to a Linear distribution.

Return type

numpy array of floats

class orbitize.priors.LogUniformPrior(minval, maxval)[source]

This is the probability distribution \(p(x) \propto 1/x\)

The __init__ method should take in a “min” and “max” value of the distribution, which correspond to the domain of the prior. (If this is not implemented, the prior has a singularity at 0 and infinite integrated probability).

Parameters
  • minval (float) – the lower bound of this distribution

  • maxval (float) – the upper bound of this distribution

compute_lnprob(element_array)[source]

Compute the prior probability of each element given that its drawn from a Log-Uniofrm prior

Parameters

element_array (float or np.array of float) – array of paramters to compute the prior probability of

Returns

array of prior probabilities

Return type

np.array

draw_samples(num_samples)[source]

Draw samples from this 1/x distribution.

Parameters

num_samples (float) – the number of samples to generate

Returns

samples ranging from [minval, maxval) as floats.

Return type

np.array

transform_samples(u)[source]

Transform uniform 1D samples, u, to samples drawn from a Log Uniform distribution.

Parameters

u (array of floats) – list of samples with values 0 < u < 1.

Returns

1D u samples transformed to a Log Uniform distribution.

Return type

numpy array of floats

class orbitize.priors.NearestNDInterpPrior(interp_fct, total_params)[source]

Nearest Neighbor interp. This class is a wrapper for scipy.interpolate.NearestNDInterpolator.

Parameters
  • interp_fct (scipy.interpolate.NearestNDInterpolator) – scipy Interpolator object containing the NDInterpolator defined by the user

  • total_params (float) – number of parameters

Written: Jorge LLop-Sayson (2021)

compute_lnprob(element_array)[source]

Compute log(probability) of an array of numbers wrt a the defined ND interpolator. Negative numbers return a probability of -inf.

Parameters

element_array (float or np.array of float) – array of numbers. We want the probability of drawing each of these from the ND interpolator.

Returns

array of log(probability) values, corresponding to the probability of drawing each of the numbers in the input element_array.

Return type

numpy array of float

draw_samples(num_samples)[source]

Draw positive samples from the ND interpolator. Negative samples will not be returned.

Parameters

num_samples (float) – the number of samples to generate.

Returns

samples drawn from the ND interpolator distribution. Array has length num_samples.

Return type

numpy array of float

increment_param_num()[source]

Increment the index to evaluate the appropriate parameter.

class orbitize.priors.ObsPrior(epochs, ra_err, dec_err, mtot, period_lims=(0, inf), tp_lims=(-inf, inf), tau_ref_epoch=58849)[source]

Implements the observation-based priors described in O’Neil+ 2018 (https://ui.adsabs.harvard.edu/abs/2019AJ….158….4O/abstract)

Parameters
  • epochs (np.array of float) – array of epochs at which observations are taken [mjd]

  • ra_err (np.array of float) – RA errors of observations [mas]

  • dec_err (np.array of float) – decl errors of observations [mas]

  • mtot (float) – total mass of system [Msol]

  • period_lims (2-tuple of float) – optional lower and upper prior limits for the orbital period [yr]

  • tp_lims (2-tuple of float) – optional lower and upper prior limits for the time of periastron passage [mjd]

  • tau_ref_epoch (float) – epoch [mjd] tau is defined relative to.

Note

This implementation is designed to be mathematically identical to the implementation in O’Neil+ 2018. There are several limitations of our implementation, in particular:

  1. ObsPrior only works with MCMC (not OFTI)

  2. ObsPrior only works with relative astrometry (i.e. you can’t use RVs or other data types)

  3. ObsPrior only works when the input astrometry is given in RA/decl. format (i.e. not sep/PA)

  4. ObsPrior assumes total mass (mtot) and parallax (plx) are fixed.

  5. ObsPrior only works for systems with one secondary object (no multi-planet systems)

  6. You must use ObsPrior with the orbitize.basis.ObsPriors orbital basis.

None of these are inherent limitations of the observation-based technique, so let us know if you have a science case that would benefit from implementing one or more of these things!

draw_samples(num_samples)[source]

Draws num_samples samples from uniform distributions in log(per), ecc, and tp. This is used for initializing the MCMC walkers.

Warning

The behavior of orbitize.priors.ObsPrior.draw_samples() is different from the draw_samples() methods of other Prior objects, which draws random samples from the prior itself.

class orbitize.priors.Prior[source]

Abstract base class for prior objects. All prior objects should inherit from this class.

Written: Sarah Blunt, 2018

class orbitize.priors.SinPrior[source]

This is the probability distribution \(p(x) \propto sin(x)\)

The domain of this prior is [0,pi].

compute_lnprob(element_array)[source]

Compute the prior probability of each element given that its drawn from a sine prior

Parameters

element_array (float or np.array of float) – array of paramters to compute the prior probability of

Returns

array of prior probabilities

Return type

np.array

draw_samples(num_samples)[source]

Draw samples from a Sine distribution.

Parameters

num_samples (float) – the number of samples to generate

Returns

samples ranging from [0, pi) as floats.

Return type

np.array

transform_samples(u)[source]

Transform uniform 1D samples, u, to samples drawn from a Sine distribution.

Parameters

u (array of floats) – list of samples with values 0 < u < 1.

Returns

1D u samples transformed to a Sine distribution.

Return type

numpy array of floats

class orbitize.priors.UniformPrior(minval, maxval)[source]

This is the probability distribution p(x) propto constant.

Parameters
  • minval (float) – the lower bound of the uniform prior

  • maxval (float) – the upper bound of the uniform prior

compute_lnprob(element_array)[source]

Compute the prior probability of each element given that its drawn from this uniform prior

Parameters

element_array (float or np.array of float) – array of paramters to compute the prior probability of

Returns

array of prior probabilities

Return type

np.array

draw_samples(num_samples)[source]

Draw samples from this uniform distribution.

Parameters

num_samples (float) – the number of samples to generate

Returns

samples ranging from [0, pi) as floats.

Return type

np.array

transform_samples(u)[source]

Transform uniform 1D samples, u, to samples drawn from a uniform distribution.

Parameters

u (array of floats) – list of samples with values 0 < u < 1.

Returns

1D u samples transformed to a uniform distribution.

Return type

numpy array of floats

orbitize.priors.all_lnpriors(params, priors)[source]

Calculates log(prior probability) of a set of parameters and a list of priors

Parameters
  • params (np.array) – size of N parameters

  • priors (list) – list of N prior objects corresponding to each parameter

Returns

prior probability of this set of parameters

Return type

float