poliastro2.math.iod.base_iod

Functions

izzo(k, r1, r2, tof, M, prograde, lowpath, ...)

Aplies izzo algorithm to solve Lambert's problem.

vallado(k, r0, r, tof, M, prograde, lowpath, ...)

Solves the Lambert's problem.

poliastro2.math.iod.base_iod._compute_T_min(ll, M, numiter, rtol)

Compute minimum T.

poliastro2.math.iod.base_iod._compute_psi(x, y, ll)

Computes psi.

“The auxiliary angle psi is computed using Eq.(17) by the appropriate inverse function”

poliastro2.math.iod.base_iod._compute_y(x, ll)

Computes y.

poliastro2.math.iod.base_iod._find_xy(ll, T, M, numiter, lowpath, rtol)

Computes all x, y for given number of revolutions.

poliastro2.math.iod.base_iod._halley(p0, T0, ll, tol, maxiter)

Find a minimum of time of flight equation using the Halley method.

Notes

This function is private because it assumes a calling convention specific to this module and is not really reusable.

poliastro2.math.iod.base_iod._householder(p0, T0, ll, M, tol, maxiter)

Find a zero of time of flight equation using the Householder method.

Notes

This function is private because it assumes a calling convention specific to this module and is not really reusable.

poliastro2.math.iod.base_iod._initial_guess(T, ll, M, lowpath)

Initial guess.

poliastro2.math.iod.base_iod._reconstruct(x, y, r1, r2, ll, gamma, rho, sigma)

Reconstruct solution velocity vectors.

poliastro2.math.iod.base_iod._tof_equation(x, T0, ll, M)

Time of flight equation.

poliastro2.math.iod.base_iod._tof_equation_p(x, y, T, ll)
poliastro2.math.iod.base_iod._tof_equation_p2(x, y, T, dT, ll)
poliastro2.math.iod.base_iod._tof_equation_p3(x, y, _, dT, ddT, ll)
poliastro2.math.iod.base_iod._tof_equation_y(x, y, T0, ll, M)

Time of flight equation with externally computated y.

poliastro2.math.iod.base_iod.izzo(k, r1, r2, tof, M, prograde, lowpath, numiter, rtol)

Aplies izzo algorithm to solve Lambert’s problem.

Parameters:
  • k (float) – Gravitational Constant

  • r1 (numpy.ndarray) – Initial position vector

  • r2 (numpy.ndarray) – Final position vector

  • tof (float) – Time of flight between both positions

  • M (int) – Number of revolutions

  • prograde (boolean) – Controls the desired inclination of the transfer orbit.

  • lowpath (boolean) – If True or False, gets the transfer orbit whose vacant focus is below or above the chord line, respectively.

  • numiter (int) – Number of iterations

  • rtol (float) – Error tolerance

Returns:

  • v1 (numpy.ndarray) – Initial velocity vector

  • v2 (numpy.ndarray) – Final velocity vector

poliastro2.math.iod.base_iod.vallado(k, r0, r, tof, M, prograde, lowpath, numiter, rtol)

Solves the Lambert’s problem.

The algorithm returns the initial velocity vector and the final one, these are computed by the following expresions:

\[\begin{split}\vec{v_{o}} &= \frac{1}{g}(\vec{r} - f\vec{r_{0}}) \\ \vec{v} &= \frac{1}{g}(\dot{g}\vec{r} - \vec{r_{0}})\end{split}\]

Therefore, the lagrange coefficients need to be computed. For the case of Lamber’s problem, they can be expressed by terms of the initial and final vector:

\[\begin{split}\begin{align} f = 1 -\frac{y}{r_{o}} \\ g = A\sqrt{\frac{y}{\mu}} \\ \dot{g} = 1 - \frac{y}{r} \\ \end{align}\end{split}\]

Where y(z) is a function that depends on the poliastro.core.stumpff coefficients:

\[\begin{split}y = r_{o} + r + A\frac{zS(z)-1}{\sqrt{C(z)}} \\ A = \sin{(\Delta \nu)}\sqrt{\frac{rr_{o}}{1 - \cos{(\Delta \nu)}}}\end{split}\]

The value of z to evaluate the stump functions is solved by applying a Numerical method to the following equation:

\[z_{i+1} = z_{i} - \frac{F(z_{i})}{F{}'(z_{i})}\]

Function F(z) to the expression:

\[F(z) = \left [\frac{y(z)}{C(z)} \right ]^{\frac{3}{2}}S(z) + A\sqrt{y(z)} - \sqrt{\mu}\Delta t\]
Parameters:
  • k (float) – Gravitational Parameter

  • r0 (numpy.ndarray) – Initial position vector

  • r (numpy.ndarray) – Final position vector

  • tof (float) – Time of flight

  • M (int) – Number of revolutions

  • prograde (boolean) – Controls the desired inclination of the transfer orbit.

  • lowpath (boolean) – If True or False, gets the transfer orbit whose vacant focus is below or above the chord line, respectively.

  • numiter (int) – Number of iterations to

  • rtol (int) – Number of revolutions

Returns:

  • v0 (numpy.ndarray) – Initial velocity vector

  • v (numpy.ndarray) – Final velocity vector

Examples

>>> from poliastro.core.iod import vallado
>>> from astropy import units as u
>>> import numpy as np
>>> from poliastro.bodies import Earth
>>> k = Earth.k.to(u.km ** 3 / u.s ** 2)
>>> r1 = np.array([5000, 10000, 2100]) * u.km # Initial position vector
>>> r2 = np.array([-14600, 2500, 7000]) * u.km # Final position vector
>>> tof = 3600 * u.s # Time of flight
>>> v1, v2 = vallado(k.value, r1.value, r2.value, tof.value, M=0, prograde=True, lowpath=True, numiter=35, rtol=1e-8)
>>> v1 = v1 * u.km / u.s
>>> v2 = v2 * u.km / u.s
>>> print(v1, v2)
[-5.99249503  1.92536671  3.24563805] km / s [-3.31245851 -4.19661901 -0.38528906] km / s

Notes

This procedure can be found in section 5.3 of Curtis, with all the theoretical description of the problem. Analytical example can be found in the same book under name Example 5.2.