Various approximations for orbital transfers#
Computing the exact optimal transfer between orbits is often more expensive than allowed in preliminary phases of the mission design. For this reason, pykep offers a number of approximations that can be used as surrogates of the actual complete computation. These can be used to quickly pre-screen a large number of possible transfers, to then focus on the most promising ones. Many of these approximations have been used and developed during the various edition of the [GTOC competition](https://sophia.estec.esa.int/gtoc_portal/) and can be found in our papers, e.g. [HIL16], [IMartensB+25].
Basic transfers#
Sometimes basic ideal transfers models can be used to get some bounds or information on the possible difficulty of a certain orbital geometry.
- pykep.hohmann(r1, r2, mu)#
Computes the delta v required for a Hohmann transfer between two circular orbits.
- Args:
r1 (
float
): radius of the first orbitr2 (
float
): radius of the second orbitmu (
float
): gravitational parameter of the central body- Returns:
[
float
, [float
,float
]]: [total delta v, [first delta v, second delta v]]- Examples:
>>> import pykep as pk >>> r1 = 7000000 >>> r2 = 9000000 >>> mu = pk.MU_EARTH >>> dv_total, [dv1, dv2] = pk.hohmann(r1, r2, mu) >>> print("Total delta v:", dv_total, "m/s") >>> print("First delta v:", dv1, "m/s") >>> print("Second delta v:", dv2, "m/s")
- pykep.bielliptic(r1, r2, rb, mu)#
Computes the delta v required for a bielliptic transfer between two circular orbits.
- Args:
r1 (
float
): radius of the first orbitr2 (
float
): radius of the second orbitrb (
float
): radius of the intermediate orbitmu (
float
): gravitational parameter of the central body- Returns:
- Examples:
>>> import pykep as pk >>> r1 = 7000000 >>> r2 = 9000000 >>> rb = 11000000 >>> mu = pk.MU_EARTH >>> dv_total, [dv1, dv2, dv3] = pk.bielliptic(r1, r2, rb, mu) >>> print("Total delta v:", dv_total, "m/s") >>> print("First delta v:", dv1, "m/s") >>> print("Second delta v:", dv2, "m/s") >>> print("Third delta v:", dv3, "m/s")
Maximum initial mass approximation (MIMA)#
If we are computing the low-thrust transfer between two arbitrary orbits and we know the starting and final time of arrival, these approximations can be used to compute the maximum initial mass (MIM) that will result in a feasible transfer. In other words if our spacecraft is fat (heavier than the MIM), the transfer will not be feasible. In general the best way to approximate the MIM is via neural method and machine learning [ABI24], but that requires learning from a vast database of precomputed trajectories and sometimes we do not have that luxury.
- pykep.mim_from_hop(pl_s, pl_t, when_s, when_t, Tmax, veff, MASS=1000.0, max_trials=10)#
The Maximum Initial Mass from a hop transfer. A hop transfer is fixed time randevouz between two planets. This function is slower then its approximations
mima_from_hop
andmima2_from_hop
, if succesfull, returns the exact solution to the problem. Behind the scene Potryagin Maximum Principle is applied to the problem and the TPBVP is solved using ipopt / heyoka.Args:
pl_s (
planet
): the source planetpl_t (
planet
): the target planetwhen_s (
epoch
): the start epochwhen_f (
epoch
): the final epochTmax (
float
, optional): Maximum spacecraft thrust.veff (
float
, optional): Isp*G0.MASS (
float
, optional): Units to use to scale the mass intrnally. It also influences the bounds for the MIM in the optimization problem.max_trials (
int
, optional): Maximum number of trials to find the MIM. If the algorithm does not converge, it will return None.- Returns:
float
: The Maximum Initial Mass
Examples:
>>> import pykep as pk >>> ... # assuming to have two planets pl_s and pl_t >>> when_s = pk.epoch(0.0) >>> when_f = pk.epoch(100.0) >>> Tmax = 0.6 >>> veff = 3000*pk.G0 >>> mima = pk.mim_from_hop(pl_s, pl_t, when_s, when_f, Tmax, veff) >>> print("Maximum initial mass:", mima, "kg")
- pykep.mima(dv1, dv2, tof, Tmax, veff)#
The Maximum Initial Mass Approximation.
Izzo, D., Hennes, D., Simões, L. F., & Märtens, M. (2016). Designing complex interplanetary trajectories for the global trajectory optimization competitions. Space Engineering: Modeling and Optimization with Case Studies, 151-176. https://www.esa.int/gsp/ACT/doc/MAD/pub/ACT-RPR-MAD-2016-SSCI-Fast_approximators.pdf
Having computed a two-impulse transfer, this approximation allows to compute the maximum initial mass that a spacecraft can have as to be able to perform that transfer in low-thrust.
- Args:
dv1 (
list
): First vectorial delta v (m/s, or Any velocity units)dv2 (
list
): Second vectorial delta v (m/s or Any velocity units)tof (
float
): Time of flight (Any time units)Tmax (
float
, optional): Maximum spacecraft thrust.veff (
float
, optional): Isp*G0.- Returns:
float
,float
: mima and magnitude of the acceleration required (units induced by the inputs)- Examples:
>>> import numpy as np >>> import pykep as pk >>> dv1 = np.array([320,-345,43]) #m/s >>> dv2 = np.array([-510,175,87]) #m/s >>> tof = 150*24*60*60 #seconds >>> mima, a_required = pk.mima(dv1, dv2, tof, Tmax = 0.6, veff=3000*pk.G0) >>> print("Maximum initial mass:", mima, "kg") >>> print("Required acceleration:", a_required*1000, "mm/s)
- pykep.mima_from_hop(pl_s, pl_t, when_s, when_t, Tmax, veff, mu)#
The Maximum Initial Mass Approximation from a hop transfer. A hop transfer is fixed time randevouz between two planets.
Izzo, D., Hennes, D., Simões, L. F., & Märtens, M. (2016). Designing complex interplanetary trajectories for the global trajectory optimization competitions. Space Engineering: Modeling and Optimization with Case Studies, 151-176. https://www.esa.int/gsp/ACT/doc/MAD/pub/ACT-RPR-MAD-2016-SSCI-Fast_approximators.pdf
- Args:
pl_s (
planet
): the source planetpl_t (
planet
): the target planetwhen_s (
epoch
): the departure epochwhen_t (
epoch
): the arrival epochTmax (
float
, optional): Maximum spacecraft thrust.veff (
float
, optional): Isp*G0.mu (
float
, optional): gravitational parameter of the central body.- Returns:
float
,float
: mima and magnitude of the acceleration required- Examples:
>>> import pykep as pk >>> ... # assuming to have two planets pl_s and pl_t >>> when_s = pk.epoch(0.0) >>> when_t = pk.epoch(100.0) >>> Tmax = 0.6 >>> veff = 3000*pk.G0 >>> mu = pk.MU_SUN >>> mima, a_required = pk.mima_from_hop(pl_s, pl_t, when_s, when_t, Tmax, veff, mu) >>> print("Maximum initial mass:", mima, "kg") >>> print("Required acceleration:", a_required*1000, "mm/s)
- pykep.mima2(posvel1, dv1, dv2, tof, Tmax, veff, mu)#
More accurate approximation of mima.
Having computed a two-impulse transfer, this approximation allows to compute the maximum initial mass that a spacecraft can have as to be able to perform that transfer in low-thrust.
Izzo, D., … & Yam, C. H. (2025). Asteroid mining: ACT&Friends’ results for the GTOC12 problem. Astrodynamics, 9(1), 19-40. (https://arxiv.org/pdf/2410.20839)
- Args:
posvel1 (
list
[list
,list
]): initial position and velocty ALONG THE LAMBERT TRANSFER.dv1 (
numpy.ndarray
): First delta v (m/s, or Any velocity units)dv2 (
numpy.ndarray
): Second delta v (m/s or Any velocity units)tof (
float
): Time of flight (Any time units)Tmax (
float
): Maximum spacecraft thrust.veff (
float
): Isp*G0.mu (
float
): gravitational parameter of the central body.- Returns:
float
,float
: mima and magnitude of the acceleration required (units induced by the inputs)
- pykep.mima2_from_hop(pl_s, pl_t, when_s, when_t, Tmax, veff, mu)#
More accurate approximation of mima from a hop transfer. A hop transfer is fixed time transfer between two planets.
Izzo, D., … & Yam, C. H. (2025). Asteroid mining: ACT&Friends’ results for the GTOC12 problem. Astrodynamics, 9(1), 19-40. (https://arxiv.org/pdf/2410.20839)
- Args:
pl_s (
planet
): the source planetpl_t (
planet
): the target planetwhen_s (
epoch
): the departure epochwhen_t (
epoch
): the arrival epochTmax (
float
, optional): Maximum spacecraft thrust.veff (
float
, optional): Isp*G0.mu (
float
, optional): gravitational parameter of the central body.- Returns:
float
,float
: mima and magnitude of the acceleration required- Examples:
>>> import pykep as pk >>> ... # assuming to have two planets pl_s and pl_t >>> when_s = pk.epoch(0.0) >>> when_t = pk.epoch(100.0) >>> Tmax = 0.6 >>> veff = 3000*pk.G0 >>> mu = pk.MU_SUN >>> mima2, a_required = pk.mima2_from_hop(pl_s, pl_t, when_s, when_t, Tmax, veff, mu) >>> print("Maximum initial mass:", mima, "kg") >>> print("Required acceleration:", a_required*1000, "mm/s)
Phasing indicators#
- class pykep.utils.knn(planet_list, when, metric='orbital', ref_r=149597870700.0, ref_v=29784.69183430911, tof=180.0)#
The class finds the k-nearest neighbours to a given planet at some epoch from a list of planets. The idea being that under some definition of “closeness” close-by planets are good candidates for orbital transfers (i.e. resulting in a low \(\Delta V\)). The use of this class is thus in preliminary mission phases of multiple randevous trajectories where target selection is to be performed efficiently.
The problem of finding who is “close-by” under a given metric can be efficiently solved using appropriate data structures. Here a kdtree is employed bringing complexity down to O(log N). The k-d-tree can then be queried efficiently for all asteroid within a given distance (‘ball’ query) or for all k closest asteroids (‘knn’ query).
The notion of distance used (metric) can be:
‘euclidean’: the simple Euclidean distance over the asteroid’s (r,v). The position and velocity vectors are scaled w.r.t. some reference values.
‘orbital’, the distance is computed with respect to \(\frac{r}{T} + v\), \(\frac{r}{T}\), corresponding to the \(\Delta V\) computed over a linear model of the orbital transfer. The distance returned will thus be in m/s.
The class is initialized with a list of planets, an epoch, and the metric to be used.
Initializes the knn class.
- Args:
planet_list (
list
ofplanet
): list of pykep planets (typically thousands).when (
epoch
): epoch.metric (
str
, optional): one of [‘euclidean’, ‘orbital’]. Defaults to ‘orbital’.ref_r (
float
, optional): reference radius (used as a scaling factor for r if the metric is ‘euclidean’). Defaults to AU.ref_v (
float
, optional): reference velocity (used as a scaling factor for v if the metric is ‘euclidean’). Defaults to EARTH_VELOCITY.tof (
float
, optional): average transfer time in days (used in the definition of the ‘orbital’ metric). Defaults to 180.0.
Example:
from pykep import * pl_list = ...... # a list of planets knn = phasing.knn(pl_list, epoch(t0), metric='orbital', tof=180) neighb, ids, dists = knn.find_neighbours(pl_list[ast_0], query_type='knn', k=10000) neighb, ids, _ = knn.find_neighbours(pl_list[ast_0], query_type='ball', r=5000)
- find_neighbours(query_planet, query_type='knn', *args, **kwargs)#
Finds the neighbours of a given planet at a given epoch. The user may query for the k-nearest neighbours or for all neighbours within a given distance.
- Args:
query_planet (
planet
orint
): the planet we want to find neighbours of. Can be an integer, in which case it refers to the idx in self.asteroid_list.query_type (
str
, optional): one of ‘knn’ or ‘ball’. Defaults to ‘knn’.*args: according to the query type (read below).
**kwargs: according to the query type (read below).
- Returns:
tuple
: (neighb, neighb_ids, dists), where dist is only computed if ‘knn’ type query is made.
The following kinds of spatial queries are currently implemented:
- query_type = ‘knn’:
The kwarg ‘k’ determines how many k-nearest neighbours are returned. For arguments, see: https://docs.scipy.org/doc/scipy/reference/generated/scipy.spatial.cKDTree.query.html
- query_type = ‘ball’:
The kwarg ‘r’ determines the distance within which all asteroids are returned. For arguments, see: https://docs.scipy.org/doc/scipy/reference/generated/scipy.spatial.cKDTree.query_ball_point.html