Orbital Elements#

In pykep the default osculating orbital elements used are the classical set \([a, e, i, \Omega, \omega, f]\) (\(f\) is the True Anomaly) together with the Cartesian position and velocity \([\mathbf r, \mathbf v]\). Support is given also for the set \([a, e, i, \Omega, \omega, M]\) (\(M\) is the Mean Anomaly) is supported as well as the Mean Equinoctial Elements [WIO85] defined as:

\[\begin{split}\left\{ \begin{array}{l} p = a (1 - e^2) \\ f = e\cos(\omega + \Omega) \\ g = e\sin(\omega + \Omega) \\ h = \tan\left(\frac i2\right)\cos\Omega \\ k = \tan\left(\frac i2\right)\sin\Omega \\ L = \Omega + \omega + f \end{array} \right.\end{split}\]

These are avoid of singularities, except at \(i = \pi\), in which case the retrogade version of the elements is to be used.

Note

In pykep` the convention \(a<0\) for hyperbolas is enforced. The user will thus not be able to instantiate orbital elements where \(a(1-e) < 0\)

A number of functions are provided to convert to and from the various orbital parameters.


class pykep.el_type#

Members:

KEP_M : Keplerian Elements \([a,e,i,\Omega,\omega,M]\) (Mean anomaly)

KEP_F : Keplerian Elements \([a,e,i,\Omega,\omega,f]\) (True anomaly)

MEQ : Modified Equinoctial Elements \([p,f,g,h,k,L]\) (Mean Longitude)

MEQ_R : Modified Equinoctial Elements (retrograde) \([p,f,g,h,k,L]\) (Mean Longitude)

POSVEL : Position and Velocity

pykep.ic2par(posvel, mu)#

Converts Cartesian state vectors (position and velocity) to Keplerian osculating orbital elements.

Args:

posvel (list [list, list]): A list containing two 3D vectors: the position vector [x, y, z] in units L and the velocity vector [vx, vy, vz] in units L/T.

mu (float): Gravitational parameter of the central body (in units L^3/T^2).

Returns:
[float, float, float, float, float, float]:

A list of six Keplerian orbital elements:

  • a: semi-major axis (in units L, positive for ellipses, negative for hyperbolae)

  • e: eccentricity (unitless)

  • i: inclination (radians, in [0, π])

  • Ω: RAAN (radians, in [0, 2π])

  • ω: argument of periapsis (radians, in [0, 2π])

  • f: true anomaly (radians, in [0, 2π])

Examples:
>>> import pykep as pk
>>> r = [7000e3, 0, 0]
>>> v = [0, 7.5e3, 0]
>>> mu = pk.MU_EARTH
>>> elements = pk.ic2par([r, v], mu)
>>> a, e, i, Omega, omega, f = elements
>>> print(f"Semi-major axis: {a} m, Eccentricity: {e}")
pykep.par2ic(par, mu)#

Converts Keplerian osculating orbital elements to Cartesian state vectors (position and velocity).

Args:

par (list of float): A list of six Keplerian orbital elements:

  • a: semi-major axis (in units L, positive for ellipses, negative for hyperbolae)

  • e: eccentricity (unitless)

  • i: inclination (radians, in [0, π])

  • Ω: longitude of ascending node (radians, in [0, 2π])

  • ω: argument of periapsis (radians, in [0, 2π])

  • f: true anomaly (radians, in [0, 2π])

mu (float):

Gravitational parameter of the central body (in units L^3/T^2)

Returns:
[list of float, list of float]:

A list containing two 3D vectors:

  • position: Cartesian position vector [x, y, z] in L

  • velocity: Cartesian velocity vector [vx, vy, vz] in L/T

Raises:

ValueError: If the semi-major axis and eccentricity are incompatible ValueError: If the true anomaly is beyond the asymptotes for a hyperbolic trajectory

Examples:
>>> import pykep as pk
>>> a = 10000e3
>>> e = 0.1
>>> i = 0.1
>>> Omega = 0.5
>>> omega = 1.0
>>> f = 2.0
>>> mu = pk.MU_EARTH
>>> r, v = pk.par2ic([a, e, i, Omega, omega, f], mu)
>>> print("Position vector:", r)
>>> print("Velocity vector:", v)
pykep.ic2eq(posvel, mu, retrogade)#

Converts Cartesian state vectors (position and velocity) to equinoctial orbital elements.

Equinoctial elements provide a non-singular representation of orbital motion, especially useful for near-circular or near-equatorial orbits. The retrograde flag allows switching between the standard and retrograde equinoctial elements. These last are not singular for inclinations of π.

Args:

posvel (list [list, list]): A list containing two 3D vectors: the position vector [x, y, z] in units L and the velocity vector [vx, vy, vz] in units L/T.

mu (float): Gravitational parameter of the central body (in units L^3/T^2).

retrogade (bool): Whether to use the retrograde equinoctial frame.

Returns:
list of float:

A list of six equinoctial orbital elements:

  • p: semi-latus rectum (in units L)

  • f: eccentricity vector times cos(Ω+ω)

  • g: eccentricity vector times sin(Ω+ω)

  • h: tan(i/2) cos Ω

  • k: tan(i/2) sin Ω

  • L: true longitude (radians, in [0, 2π])

Examples:
>>> import pykep as pk
>>> r = [7000e3, 0.0, 0.0]
>>> v = [0.0, 7.5e3, 1.0e3]
>>> mu = pk.MU_EARTH
>>> retro = False
>>> eq = pk.ic2eq([r, v], mu, retro)
>>> print("Equinoctial elements:", eq)
pykep.eq2ic(eq_elem, mu, retrogade)#

Converts equinoctial orbital elements to Cartesian state vectors (position and velocity).

Equinoctial elements provide a non-singular representation of orbital motion, especially useful for near-circular or near-equatorial orbits. The retrograde flag allows switching between the standard and retrograde equinoctial frames. These last are not singular for inclinations of pi.

Args:
eq_elem (list [float]): A list of six equinoctial elements:
  • p: semi-latus rectum (in units L)

  • f: eccentricity vector times cos(Ω+ω)

  • g: eccentricity vector times sin(Ω+ω)

  • h: tan(i/2) cos Ω

  • k: tan(i/2) sin Ω

  • L: true longitude (radians, in [0, 2π])

mu (float): Gravitational parameter of the central body (in units L^3/T^2).

retrogade (bool): Whether to use the retrograde equinoctial frame.

Returns:
list of list:

A list containing two 3D vectors:

  • position vector [x, y, z] in units L

  • velocity vector [vx, vy, vz] in units L/T

Examples:
>>> import pykep as pk
>>> eq = [7000e3, 0.01, 0.01, 0.01, 0.01, 0.0]
>>> mu = pk.MU_EARTH
>>> retro = False
>>> r, v = pk.eq2ic(eq, mu, retro)
>>> print("Position:", r)
>>> print("Velocity:", v)
pykep.eq2par(eq_elem, retrogade)#

Converts equinoctial orbital elements to classical Keplerian elements.

This function transforms the non-singular equinoctial elements into classical orbital elements, which are more intuitive but can be singular for certain inclinations or eccentricities. The retrograde flag selects the appropriate transformation for orbits with inclination near pi.

Args:
eq_elem (list [float]): A list of six equinoctial elements:
  • p: semi-latus rectum (in units L)

  • f: eccentricity vector times cos(Ω+ω)

  • g: eccentricity vector times sin(Ω+ω)

  • h: tan(i/2) cos Ω

  • k: tan(i/2) sin Ω

  • L: true longitude (radians, in [0, 2π])

retrogade (bool): Whether to use the retrograde equinoctial frame.

Returns:
list of float:

A list of six Keplerian orbital elements:

  • a: semi-major axis (in units L, positive for ellipses, negative for hyperbolae)

  • e: eccentricity (unitless)

  • i: inclination (radians, in [0, π])

  • Ω: longitude of ascending node (radians, in [0, 2π])

  • ω: argument of periapsis (radians, in [0, 2π])

  • f: true anomaly (radians, in [0, 2π])

Examples:
>>> eq = [7000e3, 0.01, 0.02, 0.001, 0.002, 0.5]
>>> retro = False
>>> par = eq2par(eq, retro)
>>> print("Keplerian elements:", par)
pykep.par2eq(par, retrogade)#

Converts classical Keplerian orbital elements to equinoctial elements.

This function provides a non-singular representation of orbits by transforming Keplerian elements into equinoctial elements. The retrograde flag allows conversion to a frame that remains non-singular for inclinations near pi.

Args:
par (list [float]): A list of six Keplerian elements:
  • a: semi-major axis (in units L, positive for ellipses, negative for hyperbolae)

  • e: eccentricity (unitless)

  • i: inclination (radians, in [0, π])

  • Ω: longitude of ascending node (radians, in [0, 2π])

  • ω: argument of periapsis (radians, in [0, 2π])

  • f: true anomaly (radians, in [0, 2π])

retrogade (bool): Whether to use the retrograde equinoctial frame.

Returns:
list of float:

A list of six equinoctial orbital elements:

  • p: semi-latus rectum (in units L)

  • f: eccentricity vector times cos(Ω+ω)

  • g: eccentricity vector times sin(Ω+ω)

  • h: tan(i/2) cos Ω

  • k: tan(i/2) sin Ω

  • L: true longitude (radians, in [0, 2π])

Examples:
>>> par = [7000e3, 0.01, 0.1, 1.0, 0.5, 0.3]
>>> retro = False
>>> eq = par2eq(par, retro)
>>> print("Equinoctial elements:", eq)