Anomalies Conversions#
In pykep we adopt the following naming for the various anomalies: \(M\) is the Mean Anomaly, \(E\) is the Eccentric Anomaly, \(L\) is the True Longitude, \(\lambda\) is the Mean Longitude, \(H\) is the Hyperbolic Anomaly, \(N\) is the Mean Hyperbolic Anomaly, \(\zeta\) is the Gudermannian, in functions or variable names these symbols can be spelled out (like zeta for \(\zeta\)) or made lowercase like m for \(M\).
Below we list various conversion functions that allow to convert from one anomaly to another, and their vectorized versions.
Normal#
- pykep.m2e(M, ecc)#
 Converts from Mean to Eccentric anomaly. Requires ecc < 1.
- pykep.e2m(E, ecc)#
 Converts from Eccentric to Mean anomaly. Requires ecc < 1.
- pykep.m2f(M, ecc)#
 Converts from Mean to True anomaly. Requires ecc < 1.
- pykep.f2m(f, ecc)#
 Converts from True to Mean anomaly. Requires ecc < 1.
- pykep.e2f(E, ecc)#
 Converts from eccentric to true anomaly. Requires ecc < 1.
- pykep.f2e(f, ecc)#
 Converts from True to Eccentric anomaly. Requires ecc < 1.
- pykep.n2h(N, ecc)#
 Converts from Hyperbolic Mean to Hyperbolic anomaly. Requires ecc > 1.
- pykep.h2n(H, ecc)#
 Converts from Hyperbolic to Hyperbolic Mean anomaly. Requires ecc > 1.
- pykep.n2f(N, ecc)#
 Converts from Hyperbolic Mean to True anomaly. Requires ecc > 1.
- pykep.f2n(f, ecc)#
 Converts from True to Hyperbolic Mean anomaly. Requires ecc > 1.
- pykep.h2f(H, ecc)#
 Converts from Hyperbolic to True anomaly. Requires ecc > 1.
- pykep.f2h(f, ecc)#
 Converts from True to Hyperbolic anomaly. Requires ecc > 1.
- pykep.zeta2f(zeta, ecc)#
 Converts from Gudermannian to True anomaly. Requires ecc > 1.
See Battin: “An Introduction to the Mathematics and Methods of Astrodynamics” for a definition of zeta and the treatment of the resulting equations.
Vectorized#
- pykep.m2e_v(Ms, eccs)#
 Converts from Mean to Eccentric anomaly (vectorized version). Requires ecc < 1.
- Args:
 Ms (
numpy.ndarrayorfloat): the Mean anomaly (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the Eccentric anomaly in [-pi, pi] (rad.)- Examples:
 >>> import pykep as pk >>> import numpy as np >>> Ms = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 0.375 >>> Es = pk.m2e_v(Ms, ecc) >>> np.shape(Es) (100,)
- pykep.e2m_v(Es, eccs)#
 Converts from Eccentric to Mean anomaly (vectorized version). Requires ecc < 1.
- Args:
 Es (
numpy.ndarrayorfloat): the Eccentric anomaly (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the Mean anomaly (rad.)- Examples:
 >>> import pykep as pk >>> import numpy as np >>> Es = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 0.86345 >>> Ms = pk.e2m_v(Es, ecc) >>> np.shape(Ms) (100,)
- pykep.m2f_v(Ms, eccs)#
 Converts from Mean to True anomaly (vectorized version). Requires ecc < 1.
- Args:
 Ms (
numpy.ndarrayorfloat): the Mean anomaly (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the True anomaly in [-pi, pi] (rad.)- Examples:
 >>> import pykep as pk >>> import numpy as np >>> Ms = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 0.4 >>> fs = pk.m2f_v(Ms, ecc) >>> np.shape(fs) (100,)
- pykep.f2m_v(fs, eccs)#
 Converts from True to Mean anomaly (vectorized version). Requires ecc < 1.
- Args:
 fs (
numpy.ndarrayorfloat): the True anomaly (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the Mean anomaly in [-pi, pi] (rad.)- Examples:
 >>> import pykep as pk >>> import numpy as np >>> fs = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 0.4 >>> Ms = pk.f2m_v(fs, ecc) >>> np.shape(Ms) (100,)
- pykep.e2f_v(Es, eccs)#
 Converts from eccentric to true anomaly (vectorized version). Requires ecc < 1.
- Args:
 Es (
numpy.ndarrayorfloat): the Eccentric anomaly (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the True anomaly in [-pi, pi] (rad.)- Examples:
 >>> import pykep as pk >>> import numpy as np >>> Es = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 0.0256 >>> fs = pk.e2f_v(Es, ecc) >>> np.shape(fs) (100,)
- pykep.f2e_v(fs, eccs)#
 Converts from True to Eccentric anomaly (vectorized version). Requires ecc < 1.
- Args:
 fs (
numpy.ndarrayorfloat): the True anomaly (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the Eccentric anomaly in [-pi, pi] (rad.)- Examples:
 >>> import pykep as pk >>> import numpy as np >>> fs = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 0.23 >>> Es = pk.f2e_v(fs, ecc) >>> np.shape(Es) (100,)
- pykep.n2h_v(Ns, eccs)#
 Converts from Hyperbolic Mean to Hyperbolic anomaly (vectorized version). Requires ecc > 1.
- Args:
 Ns (
numpy.ndarrayorfloat): the Hyperbolic Mean anomaly (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the Hyperbolic anomaly (rad.)- Examples:
 >>> import pykep as pk >>> import numpy as np >>> Ns = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 4.5 >>> Hs = pk.n2h_v(Ns, ecc) >>> np.shape(Hs) (100,)
- pykep.h2n_v(Hs, eccs)#
 Converts from Hyperbolic to Hyperbolic Mean anomaly (vectorized version). Requires ecc > 1.
- Args:
 Hs (
numpy.ndarrayorfloat): the Hyperbolic anomaly (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the Hyperbolic Mean anomaly (rad.)- Examples:
 >>> import pykep as pk >>> import numpy as np >>> Hs = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 4.5 >>> Ns = pk.h2n_v(Hs, ecc) >>> np.shape(Ns) (100,)
- pykep.n2f_v(Ns, eccs)#
 Converts from Hyperbolic Mean to True anomaly (vectorized version). Requires ecc > 1.
- Args:
 Ns (
numpy.ndarrayorfloat): the Hyperbolic Mean anomaly (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the True anomaly- Examples:
 >>> import pykep as pk >>> import numpy as np >>> Ns = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 13.45 >>> fs = pk.n2f_v(Ns, ecc) >>> np.shape(fs) (100,)
- pykep.f2n_v(fs, eccs)#
 Converts from True to Hyperbolic Mean anomaly (vectorized version). Requires ecc > 1.
- Args:
 fs (
numpy.ndarrayorfloat): the True anomaly (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the Hyperbolic Mean anomaly- Examples:
 >>> import pykep as pk >>> import numpy as np >>> fs = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 5.7 >>> Ns = pk.n2f_v(fs, ecc) >>> np.shape(Ns) (100,)
- pykep.h2f_v(Hs, eccs)#
 Converts from Hyperbolic to True anomaly (vectorized version). Requires ecc > 1.
- Args:
 Hs (
numpy.ndarrayorfloat): the Hyperbolic anomaly (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the True anomaly in [-pi, pi] (rad.)- Examples:
 >>> import pykep as pk >>> import numpy as np >>> Hs = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 4.5 >>> fs = pk.h2f_v(Hs, ecc) >>> np.shape(fs) (100,)
- pykep.f2h_v(fs, eccs)#
 Converts from True to Hyperbolic anomaly (vectorized version). Requires ecc > 1.
- Args:
 fs (
numpy.ndarrayorfloat): the True anomaly (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the Hyperbolic anomaly- Examples:
 >>> import pykep as pk >>> import numpy as np >>> fs = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 5.7 >>> Hs = pk.n2h_v(fs, ecc) >>> np.shape(Hs) (100,)
- pykep.zeta2f_v(zetas, eccs)#
 Converts from Gudermannian to True anomaly (vectorized version). Requires ecc > 1.
See Battin: “An Introduction to the Mathematics and Methods of Astrodynamics” for a definition of zeta and the treatment of the resulting equations.
- Args:
 zetas (
numpy.ndarrayorfloat): the Gudermannian (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the True anomaly- Examples:
 >>> import pykep as pk >>> import numpy as np >>> zetas = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 2.2 >>> fs = pk.zeta2f_v(zetas, ecc) >>> np.shape(fs) (100,)
- pykep.f2zeta_v(fs, eccs)#
 Converts from True anomaly to Gudermannian (vectorized version). Requires ecc > 1.
- Args:
 fs (
numpy.ndarrayorfloat): the True anomaly (rad.)eccs (
numpy.ndarrayorfloat): the eccentricity- Returns:
 numpy.ndarrayorfloat: the Gudermannian- Examples:
 >>> import pykep as pk >>> import numpy as np >>> fs = np.linspace(-np.pi/2, np.pi/2, 100) >>> ecc = 10.2 >>> zetas = pk.f2zeta_v(fs, ecc) >>> np.shape(zetas) (100,)