API Documentation

DirectOrbits.KeplerianElementsType

KeplerianElements( a, # semi-major axis [AU] e, # eccentricity i, # inclination [rad] ω, # argument of periapsis [rad] Ω, # longitude of ascending node [rad] τ, # epoch of periastron passage at MJD=0 M, # mass of primary [M⊙] plx, # parallax [mas]; defines the distance to the primary )

Represents the Keplerian elements of a secondary body orbiting a primary. Values can be specified by keyword argument or named tuple for convenience.

See also KeplerianElementsDeg for a convenience constructor accepting units of degrees instead of radians for i, ω, and Ω.

source
DirectOrbits.KeplerianElementsDegFunction

KeplerianElementsDeg(a, e, i, ω, Ω, τ, M, plx)

A convenience function for constructing KeplerianElements where i, ω, and Ω are provided in units of degrees instead of radians.

source
DirectOrbits.orbitsolveFunction
orbitsolve(elements, t)

Given a set of orbital elements with a time t in days, get the position and velocity of the secondary body (e.g. planet around a star).

This will output an AbstractOrbitSolution struct with the following properties:

  • x: δ right ascension [mas]
  • y: δ declination [mas]
  • : right ascension proper motion anomaly [mas/year]
  • : declination proper motion anomaly [mas/year]
  • : radial velocity of the secondary [m/s]
  • : right ascension acceleration [mas/year^2]
  • : declination acceleration [mas/year^2]

You can access the properties by name .x. There are helper functions to calculate each of these properties individually, but if you need more than one it is most efficient to calculate them in one go.

radvel can optionally accept the mass of the primary to calculate the impact of the secondary body on radial velocity of the primary, instead of the radial velocity of the secondary body itself.

Note: these calculations use the small angle approximation, so are only accurate when the star is much further way from the observer than the secondary is from the primary.

See also: orbitsolve_ν, projectedseparation, raoff, decoff, radvel, propmotionanom.

source
DirectOrbits.orbitsolve_νFunction
orbitsolve_ν(elem, ν)

Solve a keplerian orbit from a given true anomaly [rad]. See orbitsolve for the same function accepting a given time.

source
orbitsolve_ν(elem, ν)

Solve a keplerian orbit from a given true anomaly [rad]. See orbitsolve for the same function accepting a given time.

source
Missing docstring.

Missing docstring for OrbitSolution. Check Documenter's build log for details.

DirectOrbits.periastronFunction

periastron(elements, tref=58849)

Compute the MJD of periastron passage most recently after the reference epoch tref. N.B. mjd of 58849 = 2020-01-01

source
DirectOrbits.raoffFunction
raoff(elem, t)

Get the offset [mas] from the primary body in Right Ascension at the time t [days].

raoff(o)

Get the offset [mas] from the primary body in Right Ascension from an instance of AbstractOrbitSolution.

source
DirectOrbits.decoffFunction
decoff(elem, t)

Get the offset [mas] from the primary body in Declination at the time t [days].

decoff(elem, t)

Get the offset [mas] from the primary body in Declination from an instance of AbstractOrbitSolution.

source
DirectOrbits.posangleFunction
posangle(elem, t)

Calculate the position angle [rad] of the secondary about its primary from our perspective at the time t [days].

posangle(o)

Calculate the position angle [rad] of the secondary about its primary from our perspective from an instance of AbstractOrbitSolution.

source
DirectOrbits.projectedseparationFunction
projectedseparation(elem, t)

Calculate the projected separation [mas] of the secondary from its primary at the time t [days].

projectedseparation(o)

Calculate the projected separation [mas] of the secondary from its primary from an instance of AbstractOrbitSolution.

source
DirectOrbits.propmotionanomFunction
propmotionanom(elem, t)

Get the instantaneous proper motion anomaly [mas/year] of the secondary at the time t [days].

propmotionanom(o)

Get the instantaneous proper motion anomaly [mas/year] of the secondary from an instance of AbstractOrbitSolution.

source
propmotionanom(elem, t, M_planet)

Get the instantaneous proper motion anomaly [mas/year] of the primary in at the time t [days]. The units of M_planet and elem.M must match.

source
DirectOrbits.radvelFunction
radvel(elem, t)

Get the radial velocity [m/s] of the secondary along the line of sight at the time t [days].

radvel(o)

Get the radial velocity [m/s] of the secondary along the line of sight from an instance of AbstractOrbitSolution.

source
DirectOrbits.accelerationFunction
acceleration(elem, t)

Get the instantaneous acceleration [mas/year^2] of the secondary at the time t [days].

acceleration(o)

Get the instantaneous acceleration [mas/year^2] of the secondary from an instance of AbstractOrbitSolution.

source
acceleration(elem, t, M_planet)

Get the instantaneous acceleration [mas/year^2] of the primary in at the time t [days]. The units of M_planet and elem.M must match.

acceleration(o)

Get the instantaneous acceleration [mas/year^2] of the primary from an instance of AbstractOrbitSolution. The units of M_planet and elem.M must match.

source