API¶
radial.dataset module¶
-
class
radial.dataset.
RVDataSet
(file, t_col=0, rv_col=1, rv_unc_col=2, skiprows=0, delimiter=None, t_offset=None, rv_offset=None, t_unit=None, rv_unit=None, instrument_name=None, target_name=None, other_meta=None)¶ Bases:
object
Read and store the data and information of the radial velocities dataset in an intelligible manner. This class utilizes the power and convenience of
astropy
units and tables.Parameters: - file (
str
) – Name of the file that contains the radial velocities data. - t_col (
int
, optional) – Column number of the data file that corresponds to time. Default is 0. - rv_col (
int
, optional) – Column number of the data file that corresponds to radial velocities. Default is 1. - rv_unc_col (
int
, optional) – Column number of the data file that corresponds to radial velocity uncertainties. Default is 2. - skiprows (
int
, optional) – Number of rows to skip from the data file. Default is 0. - delimiter (
str
orNone
, optional) – String that is used to separate the columns in the data file. IfNone
, uses the default value fromnumpy.loadtxt
. Default isNone
. - t_offset (
float
,astropy.units.Quantity
orNone
, optional) – Numerical offset to be summed to the time array. IfNone
, no offset is applied. Default isNone
. - rv_offset (
str
,float
,astropy.units.Quantity
orNone
,) – optional Numerical offset to be summed to the radial velocities array. IfNone
, no offset is applied.str
options are ‘subtract_median’ and ‘subtract_mean’ (self-explanatory). Default isNone
. - t_unit (
astropy.units
orNone
, optional) – The unit of the time array, inastropy.units
. IfNone
, uses days. Default isNone
. - rv_unit (
astropy.units
orNone
, optional) – The unit of the radial velocities and uncertainties arrays, inastropy.units
. IfNone
, uses km/s. Default isNone
. - instrument_name (
str
orNone
, optional) – Name of the instrument, which will be saved in the metadata. Default isNone
. - target_name (
str
orNone
, optional) – Name of the observed target, which will be saved in the metadata. Default isNone
. - other_meta (
dict
orNone
, optional) – Other metadata to be saved in the table. IfNone
, no addition is made. Default isNone
.
-
plot
()¶ Plot the data set.
- file (
radial.estimate module¶
-
class
radial.estimate.
FullOrbit
(datasets, guess, use_add_sigma=False, parametrization=u'mc10')¶ Bases:
object
A class that computes the orbital parameters of a binary system given its radial velocities (and their uncertainties) in function of time. This class is optimized for time series that contain at least one full or almost full orbital period.
Parameters: - datasets (sequence or
radial.dataset.RVDataSet
) – A list ofRVDataSet
objects or oneRVDataSet
object that contains the data to be fit. If a sequence is passed, the order that the data sets in the sequence will dictate which instrumental parameter (gamma, sigma) index correspond to each data set. - guess (
dict
) – First guess of the orbital parameters. The keywords must match to the names of the parameters to be fit. These names are:'k'
,'period'
,'t0'
,'omega'
,'ecc'
,'gamma_X'
,'sigma_X'
(and so forth), where ‘X’ is the index of the data set. - use_add_sigma (
bool
, optional) – IfTrue
, the code will use additional parameter to estimate an extra uncertainty term for each RV data set. Default isFalse
. - parametrization (
str
, optional) – The parametrization for the orbital parameter search. Currently available options are'mc10'
and'exofast'
. Default is'mc10'
.
-
compute_dynamics
(main_body_mass=1.0)¶ Compute the mass and semi-major axis of the companion defined by the orbital parameters.
Parameters: main_body_mass ( float
, optional) – The mass of the main body which the companion orbits, in units of solar masses. Default is 1.0.
-
emcee_orbit
(nwalkers=20, nsteps=1000, p_scale=2.0, nthreads=1, ballsizes=None)¶ Calculates samples of parameters that best fit the signal rv.
Parameters: - nwalkers (
int
) – Number of walkers - nsteps (
int
) – Number of burning-in steps - p_scale (
float
, optional) – The proposal scale parameter. Default is 2.0. - nthreads (
int
) – Number of threads in your machine - ballsizes (
dict
) – The one-dimensional size of the volume from which to generate a first position to start the chain.
Returns: sampler – The resulting sampler object.
Return type: emcee.EnsembleSampler
- nwalkers (
-
lmfit_orbit
(vary=None, verbose=True, update_guess=False, minimize_mode=u'Nelder')¶ Perform a fit to the radial velocities datasets using
lmfit.minimize
.Parameters: - vary (
dict
orNone
, optional) – Dictionary with keywords corresponding to each parameter, and entries that areTrue
if the parameter is to be left to vary, orFalse
if the parameter is to be fixed in the value provided by the guess. IfNone
, all parameters will vary. Default isNone
. - verbose (
bool
, optional) – IfTrue
, print output in the screen. Default isFalse
. - update_guess (
bool
, optional) – IfTrue
, updates~estimate.FullOrbit.guess
with the estimated values from the minimization. Default isFalse
. - minimize_mode (
str
, optional) – The minimization algorithm string. See the documentation oflmfit.minimize
for a list of options available. Default is'Nelder'
.
Returns: result – The resulting
MinimizerResult
object.Return type: lmfit.MinimizerResult
- vary (
-
lnlike
(theta)¶ Log-likelihood of a given set of parameters to adequately describe the observed data.
Parameters: theta ( dict
orlmfit.Parameters
) – The orbital parameters.Returns: sum_res – The log-likelihood. Return type: scalar
-
lnprob
(theta_list)¶ This function calculates the ln of the probabilities to be used in the MCMC estimation.
Parameters: theta_list (sequence) – Returns: - The probability of the signal rv being the result of a model with the
- parameters theta
-
lomb_scargle
(dset_index, freqs)¶ Compute a Lomb-Scargle periodogram for a given data set using
scipy.signal.lombscargle
.Parameters: - dset_index (
int
) – Index of the data set to have the periodogram calculated for. - freqs (
array_like
) – Angular frequencies for output periodogram.
Returns: - pgram (
array_like
) – Lomb-Scargle periodogram. - fig (
matplotlib.pyplot.figure
) – Periodogram plot.
- dset_index (
-
make_chains
(ncut, outfile=None)¶ Make a chains object that represent the posterior distributions of the orbital parameters.
Parameters: - ncut (
int
) – Number of points of the burn-in phase to be ignored. - outfile (
str
orNone
) – A string containing the path to the file where the chains will be saved. This is useful when you do not want to keep runningemcee
frequently. IfNone
, no output file is produced. Default isNone
.
Returns: emcee_chains – The chains of the
emcee
run, with the burn-in phase removed.Return type: numpy.ndarray
- ncut (
-
plot_corner
()¶ Produce a corner (a.k.a. triangle) plot of the posterior distributions of the orbital parameters estimated with
emcee
.Returns: Return type: fig
-
plot_emcee_sampler
(outfile=None, n_cols=2, fig_size=(12, 12))¶ Plot the
emcee
sampler so that the user can check for convergence.Parameters: - outfile (
str
orNone
, optional) – Name of the output image file to be saved. IfNone
, then no output file is produced, and the plot is displayed on screen. Default isNone
. - n_cols (
int
, optional) – Number of columns of the plot. Default is 2. - fig_size (tuple, optional) – Sizes of each panel of the plot, where the first element of the tuple corresponds to the x-direction size, and the second element corresponds to the y-direction size. Default is (12, 12).
- outfile (
-
plot_rvs
(legend_loc=None, symbols=None, plot_guess=False, plot_samples=False, fold=False, numpoints=1000)¶ Plot the data sets.
Parameters: - legend_loc (
int
orNone
, optional) – Location of the legend. IfNone
, use the default frommatplotlib
. Default isNone
. - symbols (sequence or
None
, optional) – List of symbols for each data set in the plot. IfNone
, use the default list frommatplotlib
markers. Default isNone
. - plot_guess (
bool
, optional) – IfTrue
, also plots the guess as a black curve, and the RVs of each data set is shifted by its respective gamma value. - plot_samples (
bool
, optional) – IfTrue
, also plots the samples obtained with theemcee
estimation. Default isFalse
. - fold (
bool
, optional) – IfTrue
, plot the radial velocities by folding them around the estimated orbital period. Default isFalse
. - numpoints (
int
, optional) – Number of points to compute the radial velocities curve. Default is1000
.
- legend_loc (
-
prepare_params
(theta, vary_param=None)¶ Prepare a
lmfit.Parameters
object to be used in thelmfit
estimation.Parameters: - theta (
dict
) – The current orbital parameters. - vary_param (
dict
) – Dictionary that says which parameters will vary in the estimation. By default, all parameters vary. A parameter can be fixed by setting its key toFalse
.
Returns: params – The
lmfit.Parameters
object for the estimation.Return type: lmfit.Parameters
- theta (
-
print_emcee_result
(main_star_mass=None, mass_sigma=None)¶
- datasets (sequence or
radial.body module¶
-
class
radial.body.
Companion
(k=None, period_orb=None, t_0=None, omega=None, ecc=None, msini=None, semi_a=None, name=None, main_star=None, mass=None, sini=None)¶ Bases:
object
The massive companion class. It can be either a binary star, an exoplanet, maybe even a black hole! It can be anything that has a mass and orbits another massive object. General relativity effects are not implemented yet.
Parameters: - k (
astropy.units.Quantity
orNone
, optional) – The radial velocity semi-amplitude. Default isNone
. - period_orb (
astropy.units.Quantity
orNone
, optional) – The orbital period. Default isNone
. - t_0 (
astropy.units.Quantity
orNone
, optional) – The time of periastron passage. Default isNone
. - omega (
astropy.units.Quantity
orNone
, optional) – Argument of periapse. Default isNone
. - ecc (
float
orNone
, optional) – Eccentricity of the orbit. Default isNone
. - msini (
astropy.units.Quantity
orNone
, optional) – Mass of the companion multiplied by sine of the inclination of the orbital plane in relation to the line of sight. Default isNone
. - semi_a (
astropy.units.Quantity
orNone
, optional) – Semi-major axis of the orbit. Default isNone
. - name (
str
orNone
, optional) – Name of the companion. Default isNone
.
- k (
-
class
radial.body.
MainStar
(mass, name=None)¶ Bases:
object
The main star of a given system.
- mass :
astropy.units.Quantity
- The mass of the star.
- name :
str
orNone
- The name of the star. Default is
None
.
- mass :
-
class
radial.body.
System
(main_star, companion, time=None, name=None, dataset=None)¶ Bases:
object
The star-companions system class.
Parameters: - main_star (
radial.object.MainStar
) – The main star of the system. - companion (list) – Python list containing the all the
radial.object.Companion
of the system. - time (
astropy.units.Quantity
orNone
) – A scalar ornumpy.ndarray
containing the times in which the radial velocities are measured. Default isNone
. - name (
str
orNone
, optional) – Name of the system. Default isNone
. - dataset (sequence,
radial.dataset.RVDataSet
orNone
, optional) – A list ofRVDataSet
objects or oneRVDataSet
object that contains the data to be fit. Default isNone
.
-
compute_rv
()¶ Compute the radial velocities of the main star, both the individual RVs (corresponding to each companion) and the total RVs.
-
mass_func
()¶ Compute the mass functions of all the companions of the system. This method will also compute the msini and semi_a of the companions and save the values in their respective parameters.
-
plot_rv
(companion_index=None, plot_title=None)¶ Parameters: companion_index ( int
orNone
) – The companion index indicates which set of radial velocities will be plotted. IfNone
, then the total radial velocities are plotted. Default isNone
.Returns: - fig
- ax
- main_star (
radial.orbit module¶
-
class
radial.orbit.
BinarySystem
(k, period, t0, omega=None, ecc=None, sqe_cosw=None, sqe_sinw=None, gamma=0)¶ Bases:
object
A class that computes the radial velocities given the orbital parameters of the binary system.
Parameters: - k (scalar) – The radial velocity semi-amplitude K in m / s.
- period (scalar) – The orbital period in days.
- t0 (scalar) – Time of pariastron passage in days.
- omega (scalar or
None
, optional) – Argument of periapse in radians. IfNone
, bothsqe_cosw
andsqe_sinw
will be required. Default isNone
. - ecc (scalar or
None
, optional) – Eccentricity of the orbit. IfNone
, bothsqe_cosw
andsqe_sinw
will be required. Default isNone
. - sqe_cosw (scalar or
None
, optional) – The square root of the eccentricity multiplied by the cosine of the argument of periapse. IfNone
, bothomega
andecc
will be required. Default isNone
. - sqe_sinw (scalar or
None
, optional) – The square root of the eccentricity multiplied by the sine of the argument of periapse. IfNone
, bothomega
andecc
will be required. Default isNone
. - gamma (scalar or
None
, optional) – Proper motion of the barycenter in m / s. Default is 0.
-
get_rvs
(ts)¶ Computes the radial velocity given the orbital parameters.
Parameters: ts (scalar or numpy.ndarray
) – Time in days.Returns: rvs – Radial velocities Return type: scalar or numpy.ndarray
-
kep_eq
(e_ano, m_ano)¶ The Kepler equation.
Parameters: - e_ano (scalar) – Eccentric anomaly in radians.
- m_ano (scalar) – Mean anomaly in radians.
Returns: kep – Value of E-e*sin(E)-M
Return type: scalar
-
rv_eq
(f)¶ The radial velocities equation.
Parameters: f (scalar or numpy.ndarray
) – True anomaly in radians.Returns: rvs – Radial velocity Return type: scalar or numpy.ndarray
radial.prior module¶
-
radial.prior.
flat
(theta, parametrization)¶ Computes a flat prior probability for a given set of parameters theta.
Parameters: theta ( dict
) – The orbital and instrumental parameters. This dictionary must have keywords identical to the parameter names (which depend on the number of data sets, the parametrization option and if the additional uncertainties are also being estimated.Returns: prob – The prior probability for a given set of orbital and instrumental parameters. Return type: float
radial.rv_model module¶
-
radial.rv_model.
exofast
(t, log_k, log_period, t0, sqe_cosw, sqe_sinw, gamma)¶ The radial velocities model from EXOFAST (Eastman et al. 2013).
Parameters: - t (scalar) – Time in days.
- log_k (scalar) – Base-10 logarithm of the radial velocity semi-amplitude in dex(m / s).
- log_period (scalar) – Base-10 logarithm of the orbital period in dex(d).
- t0 (scalar) – Time of pariastron passage in days.
- sqe_cosw (scalar) – sqrt(ecc) * cos(omega).
- sqe_sinw (scalar) – sqrt(ecc) * sin(omega).
- gamma (scalar) – Instrumental radial velocity offset in m / s.
Returns: rvs – Radial velocity in m / s.
Return type: scalar
-
radial.rv_model.
mc10
(t, log_k, log_period, t0, omega, ecc, gamma)¶ The radial velocities model from Murray & Correia 2010.
Parameters: - t (scalar) – Time in days.
- log_k (scalar) – Base-10 logarithm of the radial velocity semi-amplitude in dex(m / s).
- log_period (scalar) – Base-10 logarithm of the orbital period in dex(d).
- t0 (scalar) – Time of pariastron passage in days.
- omega (scalar) – Argument of periapse in radians.
- ecc (scalar) – Base-10 logarithm of the eccentricity of the orbit.
- gamma (scalar) – Instrumental radial velocity offset in m / s.
Returns: rvs – Radial velocity in m / s.
Return type: scalar