Numerical Propagation#

The backbone of numerical propagation in pykep is based on Lagrangian coefficients for Kepler’s dynamics and Taylor numerical integration, as implemented in the Heyoka [BI21] python package, for all other cases. The state transition matrix is also available and provided, in the case of numerical integration, seamlessly via variational equations.

The main routines are listed here:

Keplerian#

pykep.propagate_lagrangian(rv=[[1, 0, 0], [0, 1, 0]], tof=pi / 2, mu=1, stm=False)#

Propagates (Keplerian) the state for an assigned time and computes the State Transition Matrix (if requested) using the Lagrangian coefficients.

Args:

rv (2D array-like): Cartesian components of the initial position vector and velocity [[x0, y0, z0], [v0, vy0, vz0]]. Defaults to [[1,0,0], [0,1,0]].

tof (float): time of flight. Defaults to \(\frac{\pi}{2}\).

mu (float): gravitational parameter. Defaults to 1.

stm (bool): requests the computations of the State Transition Matrix

Returns:

tuple (list, list): r and v, that is the final position and velocity after the propagation. (if stm is False) tuple (list [list, list], numpy.ndarray (6,6)): [r,v] and the STM. (if stm is True)

Examples:
>>> import pykep as pk
>>> import numpy as np
>>> r0 = [1,0,0]
>>> v0 = [0,1,0]
>>> tof = pi/2
>>> mu = 1
>>> [r1,v1], stm = pk.propagate_lagrangian(rv=[r0,v0], tof = tof, mu = mu, stm = True)
>>> [r1,v1] = pk.propagate_lagrangian(rv=[r0,v0], tof = tof, mu = mu, stm = False)
pykep.propagate_lagrangian_v(rv=[[1, 0, 0], [0, 1, 0]], tofs=[pi / 2], mu=1, stm=False)#

This is the vectorized version of pykep.propagate_lagrangian(). Vectorization allows to compute many different time of flights at once. Note that this is not necessarily more efficient than calling pykep.propagate_lagrangian() in a loop, since there is no parallelization nor SIMD magic implemented atm. Nevertheless we offer this interface for cenvenience as it may allow more compact code.

Args:

rv (2D array-like): Cartesian components of the initial position vector and velocity [[x0, y0, z0], [v0, vy0, vz0]]. Defaults to [[1,0,0], [0,1,0]].

tof (1D array-like): time of flight. Defaults to [\(\frac{\pi}{2}\)].

mu (float): gravitational parameter. Defaults to 1.

stm (bool): requests the computations of the State Transition Matrix

Returns:

list [tuple ( list , list ) ]: For each time of flight: [r,v], that is the final position and velocity after the propagation and the flattened stm (if requested).

Examples:
>>> import pykep as pk
>>> import numpy as np
>>> r0 = [1,0,0]
>>> v0 = [0,1,0]
>>> tofs = [pi/2, pi, 3/4pi]
>>> mu = 1
>>> res = pk.propagate_lagrangian_v(rv = [r0, v0], tofs = tofs, mu = mu, stm = True)
>>> rs = [it[0][0] for it in res]
>>> vs = [it[0][1] for it in res]
>>> stms = [it[1] for it in res]

Non-Keplerian#

class pykep.stark_problem(mu=1., veff=1., tol=1e-16)#

Class representing the Stark problem. In pykep, abusing a term well established in electrodynamics, this is the initial value problem of a fixed inertial thrust mass-varying spacecraft orbiting a main body and described by the equations:

\[\begin{split}\left\{ \begin{array}{l} \dot{\mathbf r} = \mathbf v \\ \dot{\mathbf v} = -\frac{\mu}{r^3} \mathbf r + \frac{\mathbf T}{m} \\ \dot m = - \frac{|\mathbf T|}{I_{sp} g_0} \end{array}\right.\end{split}\]

Note

Similar and connected functionality is provided by the functions stark(), stark_var(): and stark_dyn():.

Args:

mu (float): central body gravitational parameter. Defaults to 1.

veff (float): propulsion system effective velocity (Isp g0). Defaults to 1.

tol (float): tolerance of the Taylor adaptive integration. Defaults to 1e-16.

Note

Units need to be consistent upon construction as well as later on when calling the propagate methods.

Examples:
>>> import pykep as pk
>>> import numpy as np
>>> mu = pk.MU_SUN
>>> veff = 3000. * pk.G0
>>> tol = 1e-14
>>> sp = pk.stark_problem(mu, veff, tol)
>>> sp.propagate(rvm_state = [1., 0., 0., 0., 1., 0., 1], thrust = [0., 0., 1e-8], tof = 7.32)
[0.5089647068650076, 0.8607873878989034, 0.0, -0.8607873878989032, 0.5089647068650074, 0.0, 1.0]
property mu#

The central body gravity parameter.

propagate(rvm_state, thrust, tof)#

Stark problem numerical propagation. In pykep, abusing a term well established in electrodynamics, this is the initial value problem of a fixed inertial thrust mass-varying spacecraft orbiting a main body.

The propagation will be singular for vanishing masses (infinite acceleration) and raise an exception.

Args:

rvm_state (list (7,)): position, velocity and mass flattened into a 7D list.

thrust (list (3,)): thrust flattened into a 3D list.

tof (float): time of flight.

Returns:

list (7,) : position, velocity and mass after the numerical propagation, flattened into a 7D list.

Note

Units need to be consistent with the ones used upon constructing the instance.

Examples:
>>> import pykep as pk
>>> import numpy as np
>>> mu = pk.MU_SUN
>>> veff = 3000. * pk.G0
>>> tol = 1e-14
>>> sp = pk.stark_problem(mu, veff, tol)
>>> sp.propagate(rvm_state = [1., 0., 0., 0., 1., 0., 1], thrust = [0., 0., 1e-8], 7.32)
[0.5089647068650076, 0.8607873878989034, 0.0, -0.8607873878989032, 0.5089647068650074, 0.0, 1.0]
propagate_var(rvm_state, thrust, tof)#

Stark problem numerical propagation via variational equations.

In pykep, abusing a term well established in electrodynamics, this is the initial value problem of a fixed inertial thrust mass-varying spacecraft orbiting a main body.

The propagation will be singular for vanishing masses (infinite acceleration) and raise an exception.

It also computes the system State Transition Matrix:

\[\mathbf M = \frac{d\mathbf x_f}{d\mathbf x_0}\]

as well as the gradients of the final states with respect to the thrust direction.

\[\mathbf U = = \frac{d\mathbf x_f}{d\mathbf u}\]
Args:

rvm_state (list (7,)): position, velocity and mass flattened into a 7D list.

thrust (list (3,)): thrust flattened into a 3D list.

tof (float): time of flight.

Returns:

tuple (list (7,), numpy.ndarray (7,7), numpy.ndarray (7,3)): position, velocity and mass (after propagation) flattened into a 7D list, the state transition matrix and the gradient with respect to thrust.

Note

Units need to be consistent upon construction and later on when calling the propagate methods.

Examples:
>>> import pykep as pk
>>> import numpy as np
>>> mu = pk.MU_SUN
>>> veff = 3000. * pk.G0
>>> tol = 1e-14
>>> sp = pk.stark_problem(mu, veff, tol)
>>> sp.propagate(state = [1., 0., 0., 0., 1., 0., 1], thrust = [0., 0., 1e-8], 7.32)
property tol#

The Taylor integrator tolerance.

property veff#

The effective velocity (Isp g0)