poliastro2.math.propagation

Description

This subpackage offers various propagation algorithms for orbit simulation. It implements methods to compute the evolution of orbital states over time.

Modules

  • base: Base propagation routines.

  • cowell: Implementation of the Cowell propagation method.

  • danby: Danby’s propagation algorithm.

  • farnocchia: Farnocchia propagation routine.

  • gooding: Gooding’s method.

  • markley: Markley’s propagation routine.

  • mikkola: Mikkola’s algorithm.

  • pimienta: Pimienta’s propagation routines.

  • recseries: Recursive series propagation.

  • vallado: Propagation method based on Vallado’s work.

poliastro2.math.propagation.cowell(k, r, v, tofs, rtol=1e-11, *, events=None, f=CPUDispatcher(<function func_twobody>))
poliastro2.math.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 (numpy.ndarray) – Position vector.

  • v0 (numpy.ndarray) – 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 (numpy.ndarray) – Final position vector.

  • vv (numpy.ndarray) – Final velocity vector.

Notes

This algorithm was developed by Danby in his paper The solution of Kepler Equation with DOI: https://doi.org/10.1007/BF01686811

poliastro2.math.propagation.danby_coe(k, p, ecc, inc, raan, argp, nu, tof, numiter=20, rtol=1e-08)
poliastro2.math.propagation.farnocchia(k, r0, v0, tof)

Propagates orbit using mean motion.

This algorithm depends on the geometric shape of the orbit. For the case of the strong elliptic or strong hyperbolic orbits:

\[M = M_{0} + \frac{\mu^{2}}{h^{3}}\left ( 1 -e^{2}\right )^{\frac{3}{2}}t\]

Added in version 0.9.0.

Parameters:
  • k (float) – Standar Gravitational parameter

  • r0 (numpy.ndarray) – Initial position vector wrt attractor center.

  • v0 (numpy.ndarray) – Initial velocity vector.

  • tof (float) – Time of flight (s).

Notes

This method takes initial \(\vec{r}, \vec{v}\), calculates classical orbit parameters, increases mean anomaly and performs inverse transformation to get final \(\vec{r}, \vec{v}\) The logic is based on formulae (4), (6) and (7) from http://dx.doi.org/10.1007/s10569-013-9476-9

poliastro2.math.propagation.farnocchia_coe(k, p, ecc, inc, raan, argp, nu, tof)
poliastro2.math.propagation.func_twobody(t0, u_, k)

Differential equation for the initial value two body problem.

Parameters:
  • t0 (float) – Time.

  • u (numpy.ndarray) – Six component state vector [x, y, z, vx, vy, vz] (km, km/s).

  • k (float) – Standard gravitational parameter.

poliastro2.math.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 (numpy.ndarray) – Position vector.

  • v0 (numpy.ndarray) – 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 (numpy.ndarray) – Final position vector.

  • vv (numpy.ndarray) – Final velocity vector.

Note

Original paper for the algorithm: https://doi.org/10.1007/BF01238923

poliastro2.math.propagation.gooding_coe(k, p, ecc, inc, raan, argp, nu, tof, numiter=150, rtol=1e-08)
poliastro2.math.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:
  • k (float) – Standar Gravitational parameter.

  • r0 (numpy.ndarray) – Initial position vector wrt attractor center.

  • v0 (numpy.ndarray) – Initial velocity vector.

  • tof (float) – Time of flight.

Returns:

  • rr (numpy.ndarray) – Final position vector.

  • vv (numpy.ndarray) – Final velocity vector.

Notes

The following algorithm was taken from http://dx.doi.org/10.1007/BF00691917.

poliastro2.math.propagation.markley_coe(k, p, ecc, inc, raan, argp, nu, tof)
poliastro2.math.propagation.mikkola(k, r0, v0, tof, rtol=None)

Raw algorithm for Mikkola’s Kepler solver.

Parameters:
  • k (float) – Standard gravitational parameter of the attractor.

  • r0 (numpy.ndarray) – Position vector.

  • v0 (numpy.ndarray) – Velocity vector.

  • tof (float) – Time of flight.

  • rtol (float) – This method does not require tolerance since it is non-iterative.

Returns:

  • rr (numpy.ndarray) – Final velocity vector.

  • vv (numpy.ndarray) – Final velocity vector.

Note

Original paper: https://doi.org/10.1007/BF01235850

poliastro2.math.propagation.mikkola_coe(k, p, ecc, inc, raan, argp, nu, tof)
poliastro2.math.propagation.pimienta(k, r0, v0, tof)

Raw algorithm for Adonis’ Pimienta and John L. Crassidis 15th order polynomial Kepler solver.

Parameters:
  • k (float) – Standar Gravitational parameter.

  • r0 (numpy.ndarray) – Initial position vector wrt attractor center.

  • v0 (numpy.ndarray) – Initial velocity vector.

  • tof (float) – Time of flight.

Returns:

  • rr (numpy.ndarray) – Final position vector.

  • vv (numpy.ndarray) – Final velocity vector.

Notes

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.

poliastro2.math.propagation.pimienta_coe(k, p, ecc, inc, raan, argp, nu, tof)
poliastro2.math.propagation.recseries(k, r0, v0, tof, method='rtol', order=8, numiter=100, rtol=1e-08)

Kepler solver for elliptical orbits with recursive series approximation method. The order of the series is a user defined parameter.

Parameters:
  • k (float) – Standard gravitational parameter of the attractor.

  • r0 (numpy.ndarray) – Position vector.

  • v0 (numpy.ndarray) – Velocity vector.

  • tof (float) – Time of flight.

  • method (str) – Type of termination method (‘rtol’,’order’)

  • order (int, optional) – Order of recursion, defaults to 8.

  • numiter (int, optional) – Number of iterations, defaults to 100.

  • rtol (float, optional) – Relative error for accuracy of the method, defaults to 1e-8.

Returns:

  • rr (numpy.ndarray) – Final position vector.

  • vv (numpy.ndarray) – Final velocity vector.

Notes

This algorithm uses series discussed in the paper Recursive solution to Kepler’s problem for elliptical orbits - application in robust Newton-Raphson and co-planar closest approach estimation with DOI: http://dx.doi.org/10.13140/RG.2.2.18578.58563/1

poliastro2.math.propagation.recseries_coe(k, p, ecc, inc, raan, argp, nu, tof, method='rtol', order=8, numiter=100, rtol=1e-08)
poliastro2.math.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:
  • k (float) – Standard gravitational parameter.

  • r0 (numpy.ndarray) – Initial position vector.

  • v0 (numpy.ndarray) – Initial velocity vector.

  • tof (float) – Time of flight.

  • numiter (int) – Number of iterations.

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

Notes

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.

Modules

base

cowell(k, r, v, tofs[, rtol, events, f])

danby(k, r0, v0, tof[, numiter, rtol])

Kepler solver for both elliptic and parabolic orbits based on Danby's algorithm.

farnocchia(k, r0, v0, tof)

Propagates orbit using mean motion.

gooding(k, r0, v0, tof[, numiter, rtol])

Solves the Elliptic Kepler Equation with a cubic convergence and accuracy better than 10e-12 rad is normally achieved.

markley(k, r0, v0, tof)

Solves the kepler problem by a non-iterative method.

mikkola(k, r0, v0, tof[, rtol])

Raw algorithm for Mikkola's Kepler solver.

pimienta(k, r0, v0, tof)

Raw algorithm for Adonis' Pimienta and John L.

recseries(k, r0, v0, tof[, method, order, ...])

Kepler solver for elliptical orbits with recursive series approximation method.

vallado(k, r0, v0, tof, numiter)

Solves Kepler's Equation by applying a Newton-Raphson method.