Angles module

digraph {
   "poliastro.core.angles" -> "D_to_nu", "nu_to_D", "nu_to_E", "nu_to_F", "E_to_nu", "F_to_nu",
                              "M_to_E", "M_to_F", "M_to_D", "E_to_M", "F_to_M", "D_to_M", "M_to_nu",
                              "nu_to_M", "fp_angle";
}

The poliastro.core.angles module contains functions related to conversion between different angles used to define different orbital elements.

poliastro.core.angles.D_to_nu

True anomaly from parabolic eccentric anomaly.

\[\nu = 2 \cdot \arctan{(D)}\]
Parameters:D (Quantity) – Eccentric anomaly.
Returns:nu – True anomaly.
Return type:Quantity

Note

Taken from Farnocchia, Davide, Davide Bracali Cioci, and Andrea Milani. “Robust resolution of Kepler’s equation in all eccentricity regimes.” Celes

poliastro.core.angles.nu_to_D

Parabolic eccentric anomaly from true anomaly.

\[D = \tan{\frac{\nu}{2}}\]
Parameters:nu (Quantity) – True anomaly.
Returns:D – Hyperbolic eccentric anomaly.
Return type:Quantity

Note

Taken from Farnocchia, Davide, Davide Bracali Cioci, and Andrea Milani. “Robust resolution of Kepler’s equation in all eccentricity regimes.” Celestial Mechanics and Dynamical Astronomy 116, no. 1 (2013): 21-34.

poliastro.core.angles.nu_to_E

Eccentric anomaly from true anomaly.

New in version 0.4.0.

\[E = 2\arctan{\sqrt{\frac{1-e}{1+e}}\tan{\frac{\nu}{2}}}\]
Parameters:
Returns:

E – Eccentric anomaly.

Return type:

Quantity

poliastro.core.angles.nu_to_F

Hyperbolic eccentric anomaly from true anomaly.

\[F = ln{\left ( \frac{\sin{(\nu)}\sqrt{e^{2}-1} + \cos{\nu} + e}{1+e\cos{(\nu)}} \right )}\]
Parameters:
Returns:

F – Hyperbolic eccentric anomaly.

Return type:

Quantity

Note

Taken from Curtis, H. (2013). Orbital mechanics for engineering students. 167

poliastro.core.angles.E_to_nu

True anomaly from eccentric anomaly.

New in version 0.4.0.

\[\nu = 2\arctan{\left ( \sqrt{\frac{1+e}{1-e}}\tan{\frac{E}{2}} \right )}\]
Parameters:
Returns:

nu – True anomaly.

Return type:

Quantity

poliastro.core.angles.F_to_nu

True anomaly from hyperbolic eccentric anomaly.

Parameters:
  • F (Quantity) – Hyperbolic eccentric anomaly.
  • ecc (Quantity) – Eccentricity (>1).
Returns:

nu – True anomaly.

Return type:

Quantity

poliastro.core.angles.M_to_E

Eccentric anomaly from mean anomaly.

New in version 0.4.0.

Parameters:
Returns:

E – Eccentric anomaly.

Return type:

Quantity

poliastro.core.angles.M_to_F

Hyperbolic eccentric anomaly from mean anomaly.

Parameters:
Returns:

F – Hyperbolic eccentric anomaly.

Return type:

Quantity

poliastro.core.angles.M_to_D

Parabolic eccentric anomaly from mean anomaly.

Parameters:
Returns:

D – Parabolic eccentric anomaly.

Return type:

Quantity

poliastro.core.angles.E_to_M

Mean anomaly from eccentric anomaly.

New in version 0.4.0.

Parameters:
Returns:

M – Mean anomaly.

Return type:

Quantity

poliastro.core.angles.F_to_M

Mean anomaly from eccentric anomaly.

Parameters:
  • F (Quantity) – Hyperbolic eccentric anomaly.
  • ecc (Quantity) – Eccentricity (>1).
Returns:

M – Mean anomaly.

Return type:

Quantity

poliastro.core.angles.D_to_M

Mean anomaly from eccentric anomaly.

Parameters:
  • D (Quantity) – Parabolic eccentric anomaly.
  • ecc (Quantity) – Eccentricity.
Returns:

M – Mean anomaly.

Return type:

Quantity

poliastro.core.angles.M_to_nu

True anomaly from mean anomaly.

New in version 0.4.0.

Parameters:
  • M (Quantity) – Mean anomaly.
  • ecc (Quantity) – Eccentricity.
  • delta (float (optional)) – threshold of near-parabolic regime definition (from Davide Farnocchia et al)
Returns:

nu – True anomaly.

Return type:

Quantity

Examples

>>> from numpy import radians, degrees
>>> degrees(M_to_nu(radians(30.0), 0.06))
33.67328493021166
poliastro.core.angles.nu_to_M

Mean anomaly from true anomaly.

New in version 0.4.0.

Parameters:
  • nu (Quantity) – True anomaly.
  • ecc (Quantity) – Eccentricity.
  • delta (float (optional)) – threshold of near-parabolic regime definition (from Davide Farnocchia et al)
Returns:

M – Mean anomaly.

Return type:

Quantity

poliastro.core.angles.fp_angle

Returns the flight path angle.

\[\gamma = \arctan{\frac{e\sin{\theta}}{1 + e\cos{\theta}}}\]
Parameters:

Note

Algorithm taken from Vallado 2007, pp. 113.