API for Kepler Orbit

Class API defined in the module.

The solver for Kepler’s Equation

class PyAstronomy.pyasl.MarkleyKESolver

Solve Kepler’s Equation.

This class implements a solver for Kepler’s Equation: .. math:

M = E - sin(E),

where M is the “mean anomaly” and E is the “eccentric anomaly”. The implementation follows the prescription given by Markley (Markley 1995, CeMDA, 63, 101).

Methods

getE(M, e)

Solve Kepler's Equation for the "eccentric anomaly", E.

precisionTest([trials])

Carry out a test of the achieved precision.

getE(M, e)

Solve Kepler’s Equation for the “eccentric anomaly”, E.

Parameters
Mfloat

Mean anomaly.

efloat

Eccentricity

Returns
Eccentric anomaly: float

The solution of Kepler’s Equation

precisionTest(trials=10000)

Carry out a test of the achieved precision.

Generate random numbers for the mean anomaly and the eccentricity and calculate the eccentric anomaly. Use the number in Kepler’s Equation and compare the resulting mean anomaly with the input.

Parameters
trialsint, optional

The number of trial evaluations.

Returns
Deviationfloat

The largest determined deviation between input and output mean anomaly.

The KeplerEllipse class - the orbit

class PyAstronomy.pyasl.KeplerEllipse(a, per, e=0, tau=0, Omega=0, w=0, i=0, ks=<class 'PyAstronomy.pyasl.asl.keplerOrbit.MarkleyKESolver'>)

Calculate a Kepler orbit.

The definitions and most of the formulae used in this class derive from the book “Orbital Motion” by A.E. Roy.

Orientation of the ellipse in the coordinate system

One of the foci of the ellipse is in the origin of the coordinate system. For zero inclination the ellipse is located in the x-y plane. If the eccentricity is increased, the periastron will lie in +x direction. If the inclination is increased, the ellipse will be rotating around the x-axis, so that +y is rotated toward +z. An increase in Omega corresponds to a rotation around the z-axis so that +x is rotated toward +y. Changing w, i.e., the argument of the periastron, will not change the plane of the orbit, but rather represent a rotation of the orbit in the plane. In particular, the periapsis is shifted in the direction of motion.

Orbital angular momentum

For all parameters but semi-major axis and orbital period set to zero, the (orbital) angular momentum points into the +z direction. For an inclination of 90 deg (the remaining parameters remaining zero), it points in the -y direction.

Orientation of the ellipse in the sky

To project the ellipse onto the sky, the coordinate system should be oriented so that the +x direction points North and the +y direction points East (direction of increasing right ascension). The +z axis must be chosen so that the coordinate system becomes right handed. If the line of sight (LOS) points in the +z direction, i.e., the observer is located on the negative z axis, the parameters assume their conventional meaning.

The ascending and descending nodes

For systems outside the Solar System, the ascending node is the point where the body “crosses” the plane of the sky away from the observer. Likewise, the descending node is the point where the plane is crossed with the body approaching the observer. For the coordinate system described above and a value of zero for the longitude of the ascending node, the latter is in the North and rotates toward East (i.e., +y) when the longitude of the ascending node is increased.

The argument and longitude of periapsis

The argument of periapsis is the angle between the ascending node and the periapsis of the body measured in the direction of motion. For exoplanets with circular orbits, for which no well-defined periapsis exists, the argument of periapsis is often chosen so that time of periapsis and central transit time coincide. For the planet, this is the case if the argument of periapsis is -90 deg. However, in the exoplanet literature, the argument of periapsis often refers to the stellar orbit (see, e.g., Pollacco et al. 2008, MNRAS 385, 1576-1584, Sect. 3.2.1). In this case, the corresponding value is +90 deg. The so-called longitude of the periapsis is given by the sum of the longitude of the ascending node and the argument of periapsis. Note that for example in the literature on eclipsing binaries the terms are not always clearly distinguished and what is referred to the longitude of periapsis may actually be the argument of periapsis in the above nomenclature.

Parameters
afloat

Semi-major axis. The units of this parameter define the length scale of the orbit. Common choices are AU, solar, or stellar radii.

perfloat

Orbital period

efloat, optional

Orbital eccentricity (0-1).

taufloat, optional

Time of periapsis passage.

Omegafloat, optional

Longitude of the ascending node [deg].

wfloat, optional

Argument of periapsis [deg]. Note that the longitude if periapsis is given by Omega+w.

ifloat, optional

Orbit inclination [deg].

ksClass object, optional

The solver for Kepler’s Equation. Default is the MarkleyKESolver. Each solver must have a getE method, which takes either a float or array of float of mean anomalies and the eccentricity, and returns the associated eccentric anomalies.

Attributes
ifloat

Orbit inclination [deg].

wfloat

Argument of periapsis [deg]

Omegafloat

Longitude of the ascending node [deg]

efloat

Eccentricity

afloat

Semi-major axis

perfloat

The orbital period.

taufloat

Time of periapsis passage

ksClass object

Solver for Kepler’s equation

_nfloat

Circular frequency.

Methods

eccentricAnomaly(t)

Calculates the eccentric anomaly.

meanAnomaly(t)

Calculate the mean anomaly.

orbAngMomentum([t])

The specific orbital angular momentum of a body in the current orbit.

projPlaStDist(t)

Calculate the sky-projected planet-star separation.

radius(t[, E])

Calculate the orbit radius.

trueAnomaly(t[, E])

Calculate the true anomaly.

xyzApastron()

The position of the apastron.

xyzCenter()

Center of the ellipse

xyzFoci()

Calculate the foci of the ellipse

xyzNodes()

Calculate the nodes of the orbit.

xyzNodes_LOSZ([los, getTimes])

Calculate the nodes of the orbit for LOS in +/-z direction.

xyzPeriastron()

The position of the periastron.

xyzPos(t[, getTA])

Calculate orbit position.

xyzVel(t)

Calculate orbit velocity.

xzCrossingTime()

Calculate times of crossing the xz-plane.

yzCrossingTime([ordering])

Calculate times of crossing the yz-plane.

eccentricAnomaly(t)

Calculates the eccentric anomaly.

Parameters
tarray of float

The times at which to calculate the eccentric anomaly, E.

Returns
EArray of float

The eccentric anomaly.

meanAnomaly(t)

Calculate the mean anomaly.

Parameters
tfloat or array

The time axis.

Returns
Mean anomalyfloat or array

The mean anomaly (whether float or array depends on input).

orbAngMomentum(t=0.0)

The specific orbital angular momentum of a body in the current orbit.

Parameters
tfloat, optional

The time used to calculate the angular momentum. As is it constant for a Kepler orbit, this parameter is of no relevance.

property per

The orbital period.

projPlaStDist(t)

Calculate the sky-projected planet-star separation.

Parameters
tfloat or array

The time axis.

Returns
Positionarray

The sky-projected planet-star separation at the given time. If the input was an array, the output will be an array, holding the separations at the given times.

radius(t, E=None)

Calculate the orbit radius.

Parameters
tfloat or array

The time axis.

Efloat or array, optional

If known, the eccentric anomaly corresponding to the time points. If not given, the numbers will be calculated.

Returns
Radiusfloat or array

The orbit radius at the given points in time. Type depends on input type.

trueAnomaly(t, E=None)

Calculate the true anomaly.

Parameters
tfloat or array

The time axis.

Efloat or array, optional

If known, the eccentric anomaly corresponding to the time points (note that t will be ignored then). If not given, the numbers will be calculated.

Returns
Radiusfloat or array

The orbit radius at the given points in time. Type depends on input type.

xyzApastron()

The position of the apastron.

The apastron is the point of greatest distance.

Returns
Apastronarray of float

The x, y, and z coordinates of the apastron

xyzCenter()

Center of the ellipse

Returns
Centerarray of float

x, y, and z coordinates of the center of the Ellipse.

xyzFoci()

Calculate the foci of the ellipse

Returns
FociTuple of array of float

A tuple containing two arrays, which hold the x, y, and z coordinates of the foci.

xyzNodes()

Calculate the nodes of the orbit.

The nodes of the orbit are the points at which the orbit cuts the observing plane. In this case, these are the points at which the z-coordinate vanishes, i.e., the x-y plane is regarded the plane of observation.

Returns
NodesTuple of two coordinate arrays

Returns the xyz coordinates of both nodes.

xyzNodes_LOSZ(los='+z', getTimes=False)

Calculate the nodes of the orbit for LOS in +/-z direction.

The nodes of the orbit are the points at which the orbit cuts the plane of the sky. In this case, these are the points at which the z-coordinate vanishes, i.e., the x-y plane is regarded the plane of the sky.

Parameters
losstring, {+z,-z}, optional

Line of sight points either in +z direction (observer at -z) or vice versa. Changing the direction interchanges the ascending and descending node.

getTimesboolean, optional

If True, also the times will be returned at which the nodes are reached. Default is False.

Returns
NodesTuple of two coordinate arrays

Returns the xyz coordinates of both nodes. The first is the ascending node and the second is the descending node. If getTimes is True, tuples of the form (xyz-pos, time) will be returned.

xyzPeriastron()

The position of the periastron.

Returns
Periastronarray of float

The x, y, and z coordinates of the periastron.

xyzPos(t, getTA=False)

Calculate orbit position.

Parameters
tfloat or array

The time axis.

getTAboolean, optional

If True, also returns the “true anomaly” as a function of time (default is False).

Returns
Positionarray

The x, y, and z coordinates of the body at the given time. If the input was an array, the output will be an array of shape (input-length, 3), holding the positions at the given times.

True anomalyfloat or array

Is returned only if getTA is set to True. The true anomaly at the specified times.

xyzVel(t)

Calculate orbit velocity.

Parameters
tfloat or array

The time axis.

Returns
Velocityarray

The x, y, and z components of the body’s velocity at the given time. If the input was an array, the output will be an array of shape (input-length, 3), holding the velocity components at the given times. The unit is that of the semi-major axis divided by that of the period.

xzCrossingTime()

Calculate times of crossing the xz-plane.

This method calculates the times at which the xz-plane is crossed by the orbit. This is equivalent to finding the times where y=0.

Returns
Time 1float

First crossing time defined as having POSITIVE x position.

Time 2float

Second crossing time defined as having NEGATIVE x position.

yzCrossingTime(ordering='y')

Calculate times of crossing the yz-plane.

This method calculates the times at which the yz-plane is crossed by the orbit. This is equivalent to finding the times where x=0.

Parameters
orderingstring, {y,z}

Determines whether the sign of y or z axis position at crossing time is used for ordering the output.

Returns
Time 1float

First crossing time defined as having POSITIVE y/z position.

Time 2float

Second crossing time defined as having NEGATIVE y/z position.

The BinaryOrbit class

class PyAstronomy.pyasl.BinaryOrbit(m2m1, mtot, per, e=0, tau=0, Omega=0, w=0, i=0, msun=None, ks=<class 'PyAstronomy.pyasl.asl.keplerOrbit.MarkleyKESolver'>)

Keplerian orbits of binary components

The binary system consists of the primary mass (m1) and the secondary mass (m2), which orbit a common center of mass. The latter is in the origin of the coordinate system, stationary, and identical with a focus of either of the orbit ellipses. The system is parameterized by its total mass (m1+m2), the mass ratio (m2/m1), and the orbital elements of the primary orbit. The orbit of the secondary has a different semi-major axis than that of the primary and its argument of periapsis is shifted by 180 degrees. The orientation of the orbit ellipses follows the convention of KeplerEllipse

Parameters
m2m1float

Mass ratio between secondary and primary component (e.g., low for planetary systems)

mtotfloat

Total mass of the system [solar masses]

perfloat

Orbital period [d]

efloat, optional

Orbital eccentricity (0-1).

taufloat, optional

Time of periapsis passage [d].

Omegafloat, optional

Longitude of the ascending node [deg].

wfloat, optional

Argument of periapsis (primary orbit) [deg]. Note that the longitude if periapsis is given by Omega+w.

ifloat, optional

Orbit inclination [deg].

ksClass object, optional

The solver for Kepler’s Equation. Default is the MarkleyKESolver. Each solver must have a getE method, which takes either a float or array of float of mean anomalies and the eccentricity, and returns the associated eccentric anomalies.

msunfloat, optional

Solar mass [kg] to be adopted. By default the value from PyA’s constants will be used.

Attributes
a

Total semi-major axis (a1+a2) [m]

a1

Semi-major axis of primary orbit (a1) [m]

a2

Semi-major axis of secondary orbit (a2) [m]

m2m1

Mass ratio m2/m1 (>=0)

msun

Adopted solar mass [kg]

mtot

Total mass of system [MSun]

per

Orbital period [days]

tau

Time of periapsis [days]

Methods

getKeplerEllipse_primary()

Get Kepler ellipse model for the primary component

getKeplerEllipse_secondary()

Get Kepler ellipse model for the secondary component

xyzPos(t)

Calculate orbit position.

xyzVel(t)

Calculate orbit velocity.

property a

Total semi-major axis (a1+a2) [m]

property a1

Semi-major axis of primary orbit (a1) [m]

property a2

Semi-major axis of secondary orbit (a2) [m]

getKeplerEllipse_primary()

Get Kepler ellipse model for the primary component

Returns
keKeplerEllipse

KeplerEllipse object pertaining to primary object.

getKeplerEllipse_secondary()

Get Kepler ellipse model for the secondary component

Returns
keKeplerEllipse

KeplerEllipse object pertaining to secondary object.

property m2m1

Mass ratio m2/m1 (>=0)

property msun

Adopted solar mass [kg]

property mtot

Total mass of system [MSun]

property per

Orbital period [days]

property tau

Time of periapsis [days]

xyzPos(t)

Calculate orbit position.

Parameters
tfloat or array

The time axis [s]

Returns
Position primaryarray

The x, y, and z coordinates of the body at the given time [m]. If the input was an array, the output will be an array of shape (input-length, 3), holding the positions at the given times.

Position secondaryarray

The x, y, and z coordinates of the secondary [m].

xyzVel(t)

Calculate orbit velocity.

Parameters
tfloat or array

The time axis [s].

Returns
Velocity primaryarray

The x, y, and z components of the body’s velocity at the given time [m/s]. If the input was an array, the output will be an array of shape (input-length, 3), holding the velocity components at the given times. The unit is that of the semi-major axis divided by that of the period.

Velocity secondaryarray

The x, y, and z velocity components of the secondary [m/s].