Circular orbit—radial-velocity shift

The radial velocity is the radial component of the velocity of a source relative to an observer and is usually inferred spectroscopically. Light from an object being part of a binary system and orbiting the common center of mass will be subject to the Doppler effect. Given a circular orbit, its radial velocity shift will show sinusoidal variations determined by the object’s orbital elements.

class PyAstronomy.modelSuite.radVel.SinRadVel

Spectroscopic radial-velocity (RV) shift due to circular orbital motion.

Fit parameters:
  • P - float, Orbital period [d]

  • T0 - float, Central transit time [d]

  • K - float, radial velocity semi-amplitude [km/s]

  • rv0 - float, constant offset in radial velocity [km/s]

By default all parameters remain frozen.

The convention is that at phase zero, also the orbital radial velocity is zero. With increasing phase the object becomes bluer first.

Methods

MCMCautoParameters(ranges[, picky, ...])

Convenience function to generate parameters for MCMC fit.

addConditionalRestriction(*args)

Define a conditional restriction.

assignValue(specval)

Assign new values to variables.

assignValues(specval)

Assign new values to variables.

autoFitMCMC(x, y, ranges[, picky, stepsize, ...])

Convenience function to using auto-generated sampling parameters in MCMC.

availableParameters()

Provides a list of existing parameters.

delRestriction(parName)

Delete restriction

description([parenthesis])

Returns a description of the model based on the names of the individual components.

errorConfInterval(par[, dstat, statTol, ...])

Calculate confidence interval for a parameter.

evaluate(x)

Calculates and returns radial velocity shift according to current model parameters.

fit(x, y[, yerr, X0, minAlgo, mAA, ...])

Carries out a fit.

fitEMCEE([x, y, yerr, nwalker, priors, ...])

MCMC sampling using emcee package.

fitMCMC(x, y, X0, Lims, Steps[, yerr, ...])

Carry out MCMC fit/error estimation.

freeParamNames()

Get the names of the free parameters.

freeParameters()

Get names and values of free parameters.

freeze(specifiers)

Consider variables free to float.

frozenParameters()

Get names and values of frozen parameters.

getRelationsOf(specifier)

Return relations of a variable.

getRestrictions()

Get all restrictions.

hasVariable(specifier)

Determine whether the variable exists.

numberOfFreeParams()

Get number of free parameters.

parameterSummary([toScreen, prefix, sorting])

Writes a summary of the parameters in text form.

parameters()

Obtain parameter names and values.

relate(dependentVar, independentVars[, func])

Define a relation.

removeConditionalRestriction(*args)

Remove an existing conditional constraint.

renameVariable(oldName, newName)

Change name of variable.

restoreState(resource)

Restores parameter values from file or dictionary.

saveState(*args, **kwargs)

Save the state of the fitting object.

setObjectiveFunction([miniFunc])

Define the objective function.

setPenaltyFactor(penalFac)

Change the penalty factor.

setRestriction(restricts)

Define restrictions.

setRootName(root[, rename])

Define the root name of the model.

showConditionalRestrictions(**kwargs)

Show conditional restrictions.

steppar(pars, ranges[, extractFctVal, quiet])

Allows to step a parameter through a specified range.

thaw(specifiers)

Consider variables fixed.

untie(parName[, forceFree])

Remove all relations of parameter parName, i.e., the parameter is not dependend on other parameters.

updateModel()

Recalculate the model using current settings.

evaluate(x)

Calculates and returns radial velocity shift according to current model parameters.

Parameters
xarray

The time stamps at which to calculate the model RV curve.

Example code - Fit RV curve

from __future__ import print_function, division
# Import some unrelated modules
from numpy import arange, random, ones
import matplotlib.pylab as plt
# ... and now the radVel module
from PyAstronomy.modelSuite import radVel as rv

# Create Radial Velocity object
r = rv.SinRadVel()
# Set parameters
r.assignValue({"P": 1.8, "T0": 0.25, "K": 0.5, "rv0": 10.0})
# Choose some time axis and calculate model
time = arange(100)/100.0 * 3.0 - 1.5
y = r.evaluate(time)

# Create some faked data by adding noise
rvData = y + random.normal(0.0, 0.05, y.size)

# Randomize starting parameters for fit
for p, v in r.parameters().items():
    r[p] = v + (random.random() - 0.5) * v
# Show starting values
print("Starting values for fit:")
r.parameterSummary()

# Thaw all parameters
r.thaw(list(r.parameters().keys()))
# Start the fit
r.fit(time, rvData, yerr=ones(y.size)*0.05)

# Show fit results
print("Fitted values:")
r.parameterSummary()

# Let's see what happened...
plt.ylabel("Radial velocity [km/s]")
plt.xlabel("Radial velocity [d]")
plt.errorbar(time, rvData, yerr=ones(y.size)*0.05, fmt='b.')
plt.plot(time, y, 'r-')
plt.show()