poliastro.core.elements
¶
This module contains a set of functions that can be used to convert between different elements that define the orbit of a body.
Module Contents¶
Functions¶
|
Eccentricity vector. |
|
Compute circular velocity for a given body given thegravitational parameter and the semimajor axis. |
|
Returns r and v vectors in perifocal frame. |
|
Create a rotation matrix for coe transformation |
|
Converts from classical orbital to state vectors. |
|
|
|
Converts from classical orbital elements to modified equinoctial orbital elements. |
|
Converts from vectors to classical orbital elements. |
|
Converts from modified equinoctial orbital elements to classical |
|
Calculates position and velocity vector from modified equinoctial elements. |
- poliastro.core.elements.eccentricity_vector(k, r, v)¶
Eccentricity vector.
\[\vec{e} = \frac{1}{\mu} \left( \left( v^{2} - \frac{\mu}{r}\right ) \vec{r} - (\vec{r} \cdot \vec{v})\vec{v} \right)\]- Parameters
k (float) – Standard gravitational parameter (km^3 / s^2).
r (numpy.ndarray) – Position vector (km)
v (numpy.ndarray) – Velocity vector (km / s)
- poliastro.core.elements.circular_velocity(k, a)¶
Compute circular velocity for a given body given thegravitational parameter and the semimajor axis.
\[v = \sqrt{\frac{\mu}{a}}\]
- poliastro.core.elements.rv_pqw(k, p, ecc, nu)¶
Returns r and v vectors in perifocal frame.
- Parameters
- Returns
r (numpy.ndarray) – Position. Dimension 3 vector
v (numpy.ndarray) – Velocity. Dimension 3 vector
Notes
These formulas can be checked at Curtis 3rd. Edition, page 110. Also the example proposed is 2.11 of Curtis 3rd Edition book.
\[ \begin{align}\begin{aligned}\begin{split}\vec{r} = \frac{h^2}{\mu}\frac{1}{1 + e\cos(\theta)}\begin{bmatrix} \cos(\theta)\\ \sin(\theta)\\ 0 \end{bmatrix} \\\\\\\end{split}\\\begin{split}\vec{v} = \frac{h^2}{\mu}\begin{bmatrix} -\sin(\theta)\\ e+\cos(\theta)\\ 0 \end{bmatrix}\end{split}\end{aligned}\end{align} \]Examples
>>> from poliastro.constants import GM_earth >>> k = GM_earth.value # Earth gravitational parameter >>> ecc = 0.3 # Eccentricity >>> h = 60000e6 # Angular momentum of the orbit (m**2 / s) >>> nu = np.deg2rad(120) # True Anomaly (rad) >>> p = h**2 / k # Parameter of the orbit >>> r, v = rv_pqw(k, p, ecc, nu) >>> # Printing the results r = [-5312706.25105345 9201877.15251336 0] [m] v = [-5753.30180931 -1328.66813933 0] [m]/[s]
- poliastro.core.elements.coe_rotation_matrix(inc, raan, argp)¶
Create a rotation matrix for coe transformation
- poliastro.core.elements.coe2rv(k, p, ecc, inc, raan, argp, nu)¶
Converts from classical orbital to state vectors.
Classical orbital elements are converted into position and velocity vectors by rv_pqw algorithm. A rotation matrix is applied to position and velocity vectors to get them expressed in terms of an IJK basis.
- Parameters
k (float) – Standard gravitational parameter (km^3 / s^2).
p (float) – Semi-latus rectum or parameter (km).
ecc (float) – Eccentricity.
inc (float) – Inclination (rad).
raan (float) – Longitude of ascending node, omega (rad).
argp (float) – Argument of perigee (rad).
nu (float) – True anomaly (rad).
- Returns
r_ijk (numpy.ndarray) – Position vector in basis ijk.
v_ijk (numpy.ndarray) – Velocity vector in basis ijk.
Notes
\[\begin{split}\begin{align} \vec{r}_{IJK} &= [ROT3(-\Omega)][ROT1(-i)][ROT3(-\omega)]\vec{r}_{PQW} = \left [ \frac{IJK}{PQW} \right ]\vec{r}_{PQW}\\ \vec{v}_{IJK} &= [ROT3(-\Omega)][ROT1(-i)][ROT3(-\omega)]\vec{v}_{PQW} = \left [ \frac{IJK}{PQW} \right ]\vec{v}_{PQW}\\ \end{align}\end{split}\]Previous rotations (3-1-3) can be expressed in terms of a single rotation matrix:
\[\left [ \frac{IJK}{PQW} \right ]\]\[\begin{split}\begin{bmatrix} \cos(\Omega)\cos(\omega) - \sin(\Omega)\sin(\omega)\cos(i) & -\cos(\Omega)\sin(\omega) - \sin(\Omega)\cos(\omega)\cos(i) & \sin(\Omega)\sin(i)\\ \sin(\Omega)\cos(\omega) + \cos(\Omega)\sin(\omega)\cos(i) & -\sin(\Omega)\sin(\omega) + \cos(\Omega)\cos(\omega)\cos(i) & -\cos(\Omega)\sin(i)\\ \sin(\omega)\sin(i) & \cos(\omega)\sin(i) & \cos(i) \end{bmatrix}\end{split}\]
- poliastro.core.elements.coe2rv_many(k, p, ecc, inc, raan, argp, nu)¶
- poliastro.core.elements.coe2mee(p, ecc, inc, raan, argp, nu)¶
Converts from classical orbital elements to modified equinoctial orbital elements.
The definition of the modified equinoctial orbital elements is taken from [Walker, 1985].
The modified equinoctial orbital elements are a set of orbital elements that are useful for trajectory analysis and optimization. They are valid for circular, elliptic, and hyperbolic orbits. These direct modified equinoctial equations exhibit no singularity for zero eccentricity and orbital inclinations equal to 0 and 90 degrees. However, the h and k components are singular for an orbital inclination of 180 degrees since the retrograde factor is not yet supported.
- Parameters
- Returns
p (float) – Semi-latus rectum or parameter
f (float) – Equinoctial parameter f
g (float) – Equinoctial parameter g
h (float) – Equinoctial parameter h
k (float) – Equinoctial parameter k
L (float) – Longitude
Notes
The conversion equations are taken directly from the original paper:
\[\begin{split}\begin{align} p &= a(1-e^2) \\ f &= e\cos(\omega + \Omega) \\ g &= e\sin(\omega + \Omega) \\ h &= \tan(\frac{i}{2})\cos(\Omega) \\ k &= \tan(\frac{i}{2})\sin(\Omega) \\ L &= \Omega + \omega + \theta \\ \end{align}\end{split}\]
- poliastro.core.elements.rv2coe(k, r, v, tol=1e-08)¶
Converts from vectors to classical orbital elements.
- Parameters
k (float) – Standard gravitational parameter (km^3 / s^2)
r (numpy.ndarray) – Position vector (km)
v (numpy.ndarray) – Velocity vector (km / s)
tol (float, optional) – Tolerance for eccentricity and inclination checks, default to 1e-8
- Returns
p (float) – Semi-latus rectum of parameter (km)
ecc (float) – Eccentricity
inc (float) – Inclination (rad)
raan (float) – Right ascension of the ascending nod (rad)
argp (float) – Argument of Perigee (rad)
nu (float) – True Anomaly (rad)
Notes
This example is a real exercise from Orbital Mechanics for Engineering students by Howard D.Curtis. This exercise is 4.3 of 3rd. Edition, page 200.
First the angular momentum is computed:
\[\vec{h} = \vec{r} \times \vec{v}\]With it the eccentricity can be solved:
\[\begin{split}\begin{align} \vec{e} &= \frac{1}{\mu}\left [ \left ( v^{2} - \frac{\mu}{r}\right ) \vec{r} - (\vec{r} \cdot \vec{v})\vec{v} \right ] \\ e &= \sqrt{\vec{e}\cdot\vec{e}} \\ \end{align}\end{split}\]The node vector line is solved:
\[\begin{split}\begin{align} \vec{N} &= \vec{k} \times \vec{h} \\ N &= \sqrt{\vec{N}\cdot\vec{N}} \end{align}\end{split}\]The rigth ascension node is computed:
\[\begin{split}\Omega = \left\{ \begin{array}{lcc} cos^{-1}{\left ( \frac{N_{x}}{N} \right )} & if & N_{y} \geq 0 \\ \\ 360^{o} -cos^{-1}{\left ( \frac{N_{x}}{N} \right )} & if & N_{y} < 0 \\ \end{array} \right.\end{split}\]The argument of perigee:
\[\begin{split}\omega = \left\{ \begin{array}{lcc} cos^{-1}{\left ( \frac{\vec{N}\vec{e}}{Ne} \right )} & if & e_{z} \geq 0 \\ \\ 360^{o} -cos^{-1}{\left ( \frac{\vec{N}\vec{e}}{Ne} \right )} & if & e_{z} < 0 \\ \end{array} \right.\end{split}\]And finally the true anomaly:
\[\begin{split}\nu = \left\{ \begin{array}{lcc} cos^{-1}{\left ( \frac{\vec{e}\vec{r}}{er} \right )} & if & v_{r} \geq 0 \\ \\ 360^{o} -cos^{-1}{\left ( \frac{\vec{e}\vec{r}}{er} \right )} & if & v_{r} < 0 \\ \end{array} \right.\end{split}\]Examples
>>> from poliastro.bodies import Earth >>> from astropy import units as u >>> k = Earth.k.to_value(u.km ** 3 / u.s ** 2) >>> r = np.array([-6045., -3490., 2500.]) >>> v = np.array([-3.457, 6.618, 2.533]) >>> p, ecc, inc, raan, argp, nu = rv2coe(k, r, v) >>> print("p:", p, "[km]") p: 8530.47436396927 [km] >>> print("ecc:", ecc) ecc: 0.17121118195416898 >>> print("inc:", np.rad2deg(inc), "[deg]") inc: 153.2492285182475 [deg] >>> print("raan:", np.rad2deg(raan), "[deg]") raan: 255.27928533439618 [deg] >>> print("argp:", np.rad2deg(argp), "[deg]") argp: 20.068139973005362 [deg] >>> print("nu:", np.rad2deg(nu), "[deg]") nu: 28.445804984192122 [deg]
- poliastro.core.elements.mee2coe(p, f, g, h, k, L)¶
Converts from modified equinoctial orbital elements to classical orbital elements.
The definition of the modified equinoctial orbital elements is taken from [Walker, 1985].
\[\begin{split}\begin{align} p &= a(1 - e^{2})\\ e &= \sqrt{f^{2} + g^{2}}\\ i &= 2\arctan{(\sqrt{h^{2} + k^{2}})}\\ raan &= atan2(k, h) \pmod{2\pi}\\ argp &= (atan2(g, f) - raan) \pmod{2\pi}\\ nu &= (L - atan2(g, f)) \pmod{2\pi}\\ \end{align}\end{split}\]- Parameters
- Returns
p (float) – Semi-latus rectum
ecc (float) – Eccentricity of the orbit
inc (float) – Inclination of the orbit
raan (float) – RAAN of orbit
argp (float) – Argument of the periapsis
nu (float) – True anomaly
Notes
The conversion is always safe because arctan2 works also for 0, 0 arguments.
- poliastro.core.elements.mee2rv(p, f, g, h, k, L)¶
Calculates position and velocity vector from modified equinoctial elements.
- Parameters
- Returns
r (numpy.ndarray) – Position vector.
v (numpy.ndarray) – Velocity vector.
Note
The definition of r and v is taken from https://spsweb.fltops.jpl.nasa.gov/portaldataops/mpg/MPG_Docs/Source%20Docs/EquinoctalElements-modified.pdf, Equation 3a and 3b.