poliastro.core.propagation
¶
Low level propagation algorithms
Submodules¶
Package Contents¶
Functions¶
|
|
|
Differential equation for the initial value two body problem. |
|
Solves Kepler's Equation by applying a Newton-Raphson method. |
|
|
|
Raw algorithm for Mikkola's Kepler solver. |
|
|
|
Solves the kepler problem by a non-iterative method. Relative error is |
|
|
|
Raw algorithm for Adonis' Pimienta and John L. Crassidis 15th order |
|
|
|
Solves the Elliptic Kepler Equation with a cubic convergence and |
|
|
|
Kepler solver for both elliptic and parabolic orbits based on Danby's |
- poliastro.core.propagation.func_twobody(t0, u_, k)¶
Differential equation for the initial value two body problem.
This function follows Cowell’s formulation.
- poliastro.core.propagation.vallado(k, r0, v0, tof, numiter)¶
Solves Kepler’s Equation by applying a Newton-Raphson method.
If the position of a body along its orbit wants to be computed for a specific time, it can be solved by terms of the Kepler’s Equation:
\[E = M + e\sin{E}\]In this case, the equation is written in terms of the Universal Anomaly:
\[\sqrt{\mu}\Delta t = \frac{r_{o}v_{o}}{\sqrt{\mu}}\chi^{2}C(\alpha \chi^{2}) + (1 - \alpha r_{o})\chi^{3}S(\alpha \chi^{2}) + r_{0}\chi\]This equation is solved for the universal anomaly by applying a Newton-Raphson numerical method. Once it is solved, the Lagrange coefficients are returned:
\[\begin{split}\begin{align} f &= 1 \frac{\chi^{2}}{r_{o}}C(\alpha \chi^{2}) \\ g &= \Delta t - \frac{1}{\sqrt{\mu}}\chi^{3}S(\alpha \chi^{2}) \\ \dot{f} &= \frac{\sqrt{\mu}}{rr_{o}}(\alpha \chi^{3}S(\alpha \chi^{2}) - \chi) \\ \dot{g} &= 1 - \frac{\chi^{2}}{r}C(\alpha \chi^{2}) \\ \end{align}\end{split}\]Lagrange coefficients can be related then with the position and velocity vectors:
\[\begin{split}\begin{align} \vec{r} &= f\vec{r_{o}} + g\vec{v_{o}} \\ \vec{v} &= \dot{f}\vec{r_{o}} + \dot{g}\vec{v_{o}} \\ \end{align}\end{split}\]- Parameters
- Returns
f (float) – First Lagrange coefficient
g (float) – Second Lagrange coefficient
fdot (float) – Derivative of the first coefficient
gdot (float) – Derivative of the second coefficient
Note
The theoretical procedure is explained in section 3.7 of Curtis in really deep detail. For analytical example, check in the same book for example 3.6.
- poliastro.core.propagation.mikkola_coe(k, p, ecc, inc, raan, argp, nu, tof)¶
- poliastro.core.propagation.mikkola(k, r0, v0, tof, rtol=None)¶
Raw algorithm for Mikkola’s Kepler solver.
- Parameters
- Returns
rr (~np.array) – Final velocity vector.
vv (~np.array) – Final velocity vector.
Note
Original paper: https://doi.org/10.1007/BF01235850
- poliastro.core.propagation.markley_coe(k, p, ecc, inc, raan, argp, nu, tof)¶
- poliastro.core.propagation.markley(k, r0, v0, tof)¶
Solves the kepler problem by a non-iterative method. Relative error is around 1e-18, only limited by machine double-precision errors.
- Parameters
- Returns
rr (~np.array) – Final position vector.
vv (~np.array) – Final velocity vector.
Note
The following algorithm was taken from http://dx.doi.org/10.1007/BF00691917.
- poliastro.core.propagation.pimienta_coe(k, p, ecc, inc, raan, argp, nu, tof)¶
- poliastro.core.propagation.pimienta(k, r0, v0, tof)¶
Raw algorithm for Adonis’ Pimienta and John L. Crassidis 15th order polynomial Kepler solver.
- Parameters
- Returns
rr (~np.array) – Final position vector.
vv (~np.array) – Final velocity vector.
Note
This algorithm was derived from the original paper: Pimienta-Peñalver, A. & Crassidis, John. (2013). Accurate Kepler equation solver without transcendental function evaluations. Advances in the Astronautical Sciences. 147. 233-247.
- poliastro.core.propagation.gooding_coe(k, p, ecc, inc, raan, argp, nu, tof, numiter=150, rtol=1e-08)¶
- poliastro.core.propagation.gooding(k, r0, v0, tof, numiter=150, rtol=1e-08)¶
Solves the Elliptic Kepler Equation with a cubic convergence and accuracy better than 10e-12 rad is normally achieved. It is not valid for eccentricities equal or higher than 1.0.
- Parameters
k (float) – Standard gravitational parameter of the attractor.
r0 (array) – Position vector.
v0 (array) – Velocity vector.
tof (float) – Time of flight.
numiter (int, optional) – Number of iterations, defaults to 150.
rtol (float, optional) – Relative error for accuracy of the method, defaults to 1e-8.
- Returns
rr (~np.array) – Final position vector.
vv (~np.array) – Final velocity vector.
Note
Original paper for the algorithm: https://doi.org/10.1007/BF01238923
- poliastro.core.propagation.danby_coe(k, p, ecc, inc, raan, argp, nu, tof, numiter=20, rtol=1e-08)¶
- poliastro.core.propagation.danby(k, r0, v0, tof, numiter=20, rtol=1e-08)¶
Kepler solver for both elliptic and parabolic orbits based on Danby’s algorithm.
- Parameters
k (float) – Standard gravitational parameter of the attractor.
r0 (array) – Position vector.
v0 (array) – Velocity vector.
tof (float) – Time of flight.
numiter (int, optional) – Number of iterations, defaults to 20.
rtol (float, optional) – Relative error for accuracy of the method, defaults to 1e-8.
- Returns
rr (~np.array) – Final position vector.
vv (~np.array) – Final velocity vector.
Note
This algorithm was developed by Danby in his paper The solution of Kepler Equation with DOI: https://doi.org/10.1007/BF01686811