This module contains classes that help instantiating interplanetary trajectory optimization problems. While the resulting problems are compatible with the optimization package pygmo, which can thus be used to assemble a solution strategy, they can be also interfaced with any other optimization package / method.
Name |
Type |
Description |
---|---|---|
class |
A Multiple Gravity Assist Trajectory with no deep space manouvres |
|
class |
A multiple Gravity Assist Trajectory with one deep space manouvre at each leg |
|
class |
A single leg transfer with N impulses |
|
class |
A cubesat mission to near Earth asteroids. Solar Electric Propulsion and Earth gravity are modelled. |
|
class |
A multiple randezvous low-thrust optimization problem (e.g. for asteroids in the main belt) |
|
class |
A low-thrust transfer between planets using a direct transcription. |
|
class |
A low-thrust transfer between Cartesian states using an indirect transcription. |
|
class |
A low-thrust transfer between orbits using an indirect transcription. |
|
class |
Contains functors delivering the performance of launchers in terms of mass delivered |
Some instances of the classes above are provided in an Interplanetary Trajectory Optimization Gym. A problem test set we use to test and tune optimization approaches able to tackle a wide variety of problems with little / no tuning.
Note
Most problem have a pretty and a plot method that can be used to visualize decision vectors (chromosomes).
Name |
Type |
Description |
---|---|---|
instance |
Cassini inspired MGA problem. Time of flights are encoded directly in the chromosome. |
|
instance |
Cassini inspired MGA problem. Time of flights use alpha encoding. |
|
instance |
Cassini inspired MGA problem. Time of flights use eta encoding. |
|
instance |
Earth-Mars transfer with 3 impulses |
|
instance |
Earth-Mars transfer with 5 impulses |
|
instance |
Earth-Mars transfer with 7 impulses |
|
instance |
Earth-Venus-Earth transfer. Time of flights are encoded directly in the chromosome. |
|
instance |
Earth-Venus-Earth transfer. Time of flights use alpha encoding. |
|
instance |
Earth-Venus-Earth transfer. Time of flights use eta encoding. |
|
instance |
Cassini inspired MGA-1DSM problem. |
|
instance |
Rosetta inspired MGA-1DSM problem. |
|
class |
TandEM inspired MGA-1DSM problem. |
|
instance |
JUICE inspired MGA-1DSM problem. |
|
instance |
JUICE inspired MGA-1DSM multiobjective problem. |
|
instance |
Messenger inspired MGA-1DSM problem. |
pykep.trajopt.
mga
(*args)¶This class transcribes a Multiple Gravity Assist (MGA) trajectory with no deep space manouvres into an optimisation problem. It may be used as a User Defined Problem (UDP) for the pygmo (http://esa.github.io/pygmo/) optimisation suite.
Izzo, Dario. “Global optimization and space pruning for spacecraft trajectory design.” Spacecraft Trajectory Optimization 1 (2010): 178-200.
The decision vector (chromosome) is:
direct encoding: [t0, T1, T2 ... ] in [mjd2000, days, days ... ]
alpha encoding: [t0, T, a1, a2 ...] in [mjd2000, days, nd, nd ... ]
eta encoding: [t0, n1, n2, n3 ...] in [mjd2000, nd, nd ...]
Note
The time of flights of a MGA trajectory (and in general) can be encoded in different ways. When they are directly present in the decision vector, we have the direct encoding. This is the most ‘evolvable’ encoding but also the one that requires the most problem knowledge (e.g. to define the bounds on each leg) and is not very flexible in dealing with constraints on the total time of flight. The alpha and eta encodings, instead, allow to only specify bounds on the time of flight of the entire trajectory, and not on the single legs: a property that is attractive for multi-objective optimization, for example.
In the alpha encoding each leg time-of-flight is decoded as follows, T_i = T log(alpha_i) / sum_n(log(alpha_n)). In the eta encoding each leg time-of-flight is decoded as follows, T_i = (tof_max - sum_0^(i-1)(T_j)) * eta_i
The chromosome dimension for the direct and eta encoding is the same, while the alpha encoding requires one more gene.
Note
The resulting problem is box-bounded (unconstrained).
__init__
(*args)¶mga(seq=[pk.planet.jpl_lp(‘earth’), pk.planet.jpl_lp(‘venus’), pk.planet.jpl_lp(‘earth’)], t0=[0, 1000], tof=[100, 500], vinf=2.5, multi_objective=False, alpha_encoding=False, orbit_insertion=False, e_target=None, rp_target=None)
seq (list of pk.planet
): sequence of body encounters including the starting object
t0 (list of pk.epoch
): the launch window
tof (list
or float
): bounds on the time of flight. If tof_encoding is ‘direct’, this contains a list
of 2D lists defining the upper and lower bounds on each leg. If tof_encoding is ‘alpha’,
this contains a 2D list with the lower and upper bounds on the time-of-flight. If tof_encoding
is ‘eta’ tof is a float defining the upper bound on the time-of-flight
vinf (float
): the vinf provided at launch for free
multi_objective (bool
): when True constructs a multiobjective problem (dv, T). In this case, ‘alpha’ or eta encodings are reccomended
tof_encoding (str
): one of ‘direct’, ‘alpha’ or ‘eta’. Selects the encoding for the time of flights
orbit_insertion (bool
): when True the arrival dv is computed as that required to acquire a target orbit defined by e_target and rp_target
e_target (float
): if orbit_insertion is True this defines the target orbit eccentricity around the final planet
rp_target (float
): if orbit_insertion is True this defines the target orbit pericenter around the final planet (in m)
ValueError: if planets do not share the same central body (checked on the mu_central_body attribute)
ValueError: if t0 does not contain objects able to construct a epoch (e.g. pk. epoch or floats)
ValueError: if tof is badly defined
ValueError: it the target orbit is not defined and orbit_insertion is True
alpha2direct
(*args)¶alpha2direct(x)
x (array-like
): a chromosome encoding an MGA trajectory in the alpha encoding
numpy.array
: a chromosome encoding the MGA trajectory using the direct encoding
direct2alpha
(*args)¶direct2alpha(x)
x (array-like
): a chromosome encoding an MGA trajectory in the direct encoding
numpy.array
: a chromosome encoding the MGA trajectory using the alpha encoding
eta2direct
(*args)¶eta2direct(x)
x (array-like
): a chromosome encoding an MGA trajectory in the eta encoding
numpy.array
: a chromosome encoding the MGA trajectory using the direct encoding
ValueError: when the tof_encoding is not ‘eta’
direct2eta
(*args)¶direct2eta(x)
x (array-like
): a chromosome encoding an MGA trajectory in the direct encoding
numpy.array
: a chromosome encoding the MGA trajectory using the eta encoding
ValueError: when the tof_encoding is not ‘eta’
pretty
(*args)¶pretty(x)
x (list
, tuple
, numpy.ndarray
): Decision chromosome, e.g. (pygmo.population.champion_x
).
Prints human readable information on the trajectory represented by the decision vector x
plot
(*args)¶plot(self, x, axes=None, units=pk.AU, N=60)
Plots the spacecraft trajectory.
x (tuple
, list
, numpy.ndarray
): Decision chromosome.
axes (matplotlib.axes._subplots.Axes3DSubplot
): 3D axes to use for the plot
units (float
, int
): Length unit by which to normalise data.
N (float
): Number of points to plot per leg
pykep.trajopt.
mga_1dsm
(*args)¶This class transcribes a Multiple Gravity Assist trajectory with one deep space manouvre per leg (MGA-1DSM) into an optimisation problem. It may be used as a User Defined Problem (UDP) for the pygmo (http://esa.github.io/pygmo/) optimisation suite.
Izzo, Dario. “Global optimization and space pruning for spacecraft trajectory design.” Spacecraft Trajectory Optimization 1 (2010): 178-200.
The decision vector (chromosome) is:
direct encoding: [t0] + [u, v, Vinf, eta1, T1] + [beta, rp/rV, eta2, T2] + ...
alpha encoding: [t0] + [u, v, Vinf, eta1, a1] + [beta, rp/rV, eta2, a2] + ... + [T]
eta encoding: [t0] + [u, v, Vinf, eta1, n1] + [beta, rp/rV, eta2, n2] + ...
where t0 is a mjd2000, Vinf is in km/s, T in days, beta in radians and the rest non dimensional.
Note
The time of flights of a MGA-1DSM trajectory (and in general) can be encoded in different ways. When they are directly present in the decision vector, we have a direct encoding. This is the most ‘evolvable’ encoding but also the one that requires the most problem knowledge (e.g. to define the bounds on each leg) and is not very flexible in dealing with constraints on the total time of flight. The alpha and eta encodings, instead, allow to only specify bounds on the time of flight of the entire trajectory, and not on the single legs: a property that is attractive for multi-objective optimization, for example.
In the alpha encoding each leg time-of-flight is decoded as follows, T_i = T log(alpha_i) / sum_n(log(alpha_n)). In the eta encoding each leg time-of-flight is decoded as follows, T_i = (tof_max - sum_0^(i-1)(T_j)) * eta_i
The chromosome dimension for the direct and eta encoding is the same, while the alpha encoding requires one more gene.
Note
The resulting problem is box-bounded (unconstrained).
__init__
(*args)¶pykep.trajopt.mga_1dsm(seq = [jpl_lp(‘earth’), jpl_lp(‘venus’), jpl_lp(‘earth’)], t0 = [epoch(0),epoch(1000)], tof = [1.0,5.0], vinf = [0.5, 2.5], multi_objective = False, add_vinf_dep = False, add_vinf_arr=True)
seq (list
of pykep.planet
): the encounter sequence (including the starting launch)
t0 (list
of pykep.epoch
or floats
): the launch window (in mjd2000 if floats)
list
or float
): bounds on the time of flight (days). If tof_encoding is ‘direct’, this contains a listof 2D lists defining the upper and lower bounds on each leg. If tof_encoding is ‘alpha’, this contains a list of two floats containing the lower and upper bounds on the time-of-flight. If tof_encoding is ‘eta’ tof is a float defining the upper bound on the time-of-flight
vinf (list
): the minimum and maximum allowed initial hyperbolic velocity (at launch), in km/sec
add_vinf_dep (bool
): when True the computed Dv includes the initial hyperbolic velocity (at launch)
add_vinf_arr (bool
): when True the computed Dv includes the final hyperbolic velocity (at the last planet)
tof_encoding (str
): one of ‘direct’, ‘alpha’ or ‘eta’. Selects the encoding for the time of flights
multi_objective (bool
): when True constructs a multiobjective problem (dv, T)
orbit_insertion (bool
): when True the arrival dv is computed as that required to acquire a target orbit defined by e_target and rp_target
e_target (float
): if orbit_insertion is True this defines the target orbit eccentricity around the final planet
rp_target (float
): if orbit_insertion is True this defines the target orbit pericenter around the final planet (in m)
pretty
(*args)¶prob.plot(x)
x: encoded trajectory
Prints human readable information on the trajectory represented by the decision vector x
Example:
print(prob.pretty(x))
plot
(*args)¶ax = prob.plot(x, ax=None)
x: encoded trajectory
ax: matplotlib axis where to plot. If None figure and axis will be created
[out] ax: matplotlib axis where to plot
Plots the trajectory represented by a decision vector x on the 3d axis ax
Example:
ax = prob.plot(x)
pykep.trajopt.
pl2pl_N_impulses
(*args)¶This class is a pygmo (http://esa.github.io/pygmo/) problem representing a single leg transfer between two planets allowing up to a maximum number of impulsive Deep Space Manouvres.
The decision vector is:
[t0,T] + [alpha,u,v,V_inf] * (N-2) + [alpha] + ([tf])
… in the units: [mjd2000, days] + [nd, nd, m/sec, nd] + [nd] + [mjd2000]
Each time-of-flight can be decoded as follows, T_n = T log(alpha_n) / sum_i(log(alpha_i))
Note
The resulting problem is box-bounded (unconstrained). The resulting trajectory is time-bounded.
__init__
(*args)¶prob = pykep.trajopt.pl2pl_N_impulses(start=jpl_lp(‘earth’), target=jpl_lp(‘venus’), N_max=3, tof=[20., 400.], vinf=[0., 4.], phase_free=True, multi_objective=False, t0=None)
start (pykep.planet
): the starting planet
target (pykep.planet
): the target planet
N_max (int
): maximum number of impulses
tof (list
): the box bounds [lower,upper] for the time of flight (days)
vinf (list
): the box bounds [lower,upper] for each DV magnitude (km/sec)
phase_free (bool
): when True, no randezvous condition are enforced and start and arrival anomalies will be free
multi_objective (bool
): when True, a multi-objective problem is constructed with DV and time of flight as objectives
t0 (list
): the box bounds on the launch window containing two pykep.epoch. This is not needed if phase_free is True.
pretty
(*args)¶plot
(*args)¶ax = prob.plot_trajectory(x, axes=None)
x: encoded trajectory
axes: matplotlib axis where to plot. If None figure and axis will be created
[out] ax: matplotlib axis where to plot
Plots the trajectory represented by a decision vector x on the 3d axis ax
Example:
ax = prob.plot(x)
pykep.trajopt.
lt_margo
(*args)¶This class can be used as a User Defined Problem (UDP) in the pygmo2 software and if successfully solved, represents a low-thrust interplanetary trajectory from the Earth (or from the Sun-Earth L1 or L2 Lagrangian point) to a target NEO. The trajectory is modeled using the Sims-Flanagan model, extended to include the Earth’s gravity (assumed constant along each segment). The propulsion model can be both nuclear (NEP) or solar (SEP).
This problem was developed and used during the European Space Agency interplanetary cubesat M-ARGO study.
The decision vector (chromosome) is:
[t0, tof, mf] + [throttles1] + [throttles2] + ...
__init__
(*args)¶prob = lt_margo(target = pk.planet.mpcorb(mpcorbline), n_seg = 30, grid_type = “uniform”, t0 = [epoch(8000), epoch(9000)], tof = [200, 365.25*3], m0 = 20.0, Tmax = 0.0017, Isp = 3000.0, earth_gravity = False, sep = False, start = “earth”)
target (pykep.planet
): target planet
n_seg (int
): number of segments to use in the problem transcription (time grid)
grid_type (string
): “uniform” for uniform segments, “nonuniform” toi use a denser grid in the first part of the trajectory
t0 (list
): list of two pykep.epoch defining the bounds on the launch epoch
tof (list
): list of two floats defining the bounds on the time of flight (days)
m0 (float
): initial mass of the spacecraft (kg)
Tmax (float
): maximum thrust at 1 AU (N)
Isp (float
): engine specific impulse (s)
earth_gravity (bool
): activates the Earth gravity effect in the dynamics
sep (bool
): Activates a Solar Electric Propulsion model for the thrust - distance dependency.
start(string
): Starting point (“earth”, “l1”, or “l2”).
Note
L1 and L2 are approximated as the points on the line connecting the Sun and the Earth at a distance of, respectively, 0.99 and 1.01 AU from the Sun.
Note
If the Earth’s gravity is enabled, the starting point cannot be the Earth
pretty
(*args)¶prob.pretty(x)
x (list
, tuple
, numpy.ndarray
): Decision chromosome, e.g. (pygmo.population.champion_x
).
Prints human readable information on the trajectory represented by the decision vector x
plot_traj
(*args)¶ax = prob.plot_traj(self, x, units=AU, plot_segments=True, plot_thrusts=False, axes=None)
x (list
, tuple
, numpy.ndarray
): Decision chromosome, e.g. (pygmo.population.champion_x
).
units (float
): the length unit to be used in the plot
plot_segments (bool
): when True plots also the segments boundaries
plot_thrusts (bool
): when True plots also the thrust vectors
matplotlib.axes: axes where to plot
Visualizes the the trajectory in a 3D plot
plot_dists_thrust
(*args)¶axes = prob.plot_dists_thrust(self, x, axes=None)
x (list
, tuple
, numpy.ndarray
): Decision chromosome, e.g. (pygmo.population.champion_x
).
matplotlib.axes: axes where to plot
Plots the distance of the spacecraft from the Earth/Sun and the thrust profile
double_segments
(*args)¶new_prob, new_x = prob.double_segments(self,x)
x (list
, tuple
, numpy.ndarray
): Decision chromosome, e.g. (pygmo.population.champion_x
).
the new udp having twice the segments
list: the new chromosome to be used as initial guess
Returns the decision vector encoding a low trust trajectory having double the number of segments with respect to x and a ‘similar’ throttle history. In case high fidelity is True, and x is a feasible trajectory, the returned decision vector also encodes a feasible trajectory that can be further optimized Returns the problem and the decision vector encoding a low-thrust trajectory having double the number of segments with respect to x and the same throttle history. If x is a feasible trajectory, the new chromosome is “slightly unfeasible”, due to the new refined Earth’s gravity and/or SEP thrust that are now computed in the 2 halves of each segment.
pykep.trajopt.
direct_pl2pl
(*args)¶Represents a direct transcription transfer between solar system planets.
This problem works by manipulating the starting epoch t0, the transfer time T the final mass mf and the controls The dicision vector is:
z = [t0, T, mf, Vxi, Vyi, Vzi, Vxf, Vyf, Vzf, controls]
__init__
(*args)¶Initialises a direct transcription orbit to orbit problem.
p0 (str
): Departure planet name. (will be used to construct a planet.jpl_lp object)
pf (str
): Arrival planet name. (will be used to construct a planet.jpl_lp object)
mass (float
, int
): Spacecraft wet mass [kg].
thrust (float
, int
): Spacecraft maximum thrust [N].
isp (float
, int
): Spacecraft specific impulse [s].
nseg (int
): Number of colocation nodes.
t0 (list
): Launch epochs bounds [mjd2000].
tof (list
): Transfer time bounds [days].
vinf_dep (float
): allowed launch DV [km/s]
vinf_arr (float
): allowed arrival DV [km/s]
hf (bool
): High-fidelity. Activates a continuous representation for the thrust.
pretty
(*args)¶pretty(x)
x (list
, tuple
, numpy.ndarray
): Decision chromosome, e.g. (pygmo.population.champion_x
).
Prints human readable information on the trajectory represented by the decision vector x
plot_traj
(*args)¶This function plots the 3 dimensional spacecraft trajectory, given a solution chromosome.
z (list
, tuple
, numpy.ndarray
): Decision chromosome, e.g. (pygmo.population.champion_x
).
units (float
, int
): units by which to scale the trajectory dimensions.
N (int
): Number of points to be plotted along one arc.
plot_control
(*args)¶Plots the control profile of the trajectory, as a function of time.
z (list
, tuple
, numpy.ndarray
): Decision chromosome, e.g. (pygmo.population.champion_x
).
mark (string
): matplotlib marker.
time (bool
): If True
, x-axis is time in mjd2000
. If False
, x-axis is node index.
get_traj
(*args)¶Retrieves the trajectory information. The returned np.array contains:
traj = [[t0, x0, y0, z0, vx0, vy0, vz0, m0, u0, ux0, uy0, uz0]
...
[tf, xf, yf, zf, vxf, vyf, vzf, mf, uf, uxf, uyf, uzf]]
z (list
, tuple
, numpy.ndarray
): Decision chromosome, e.g. (pygmo.population.champion_x
).
np.array full information on states and controls along the trajectory nodes
pykep.trajopt.
mr_lt_nep
(*args)¶This class represents, as a global optimization problem (linearly constrained, high diemensional), a Multiple Randezvous trajectory of a low-thrust spacecraft equipped with a nuclear electric propulsion engine.
Izzo, D. et al., GTOC7 - Solution Returned by the ACT/ISAS team
The decision vector (chromosome) is: [t1, tof1, rest1, m_f1] + [throttles1] + [t2, tof2, rest2, m_f2] + [throttles2] + …. …. + [total_tof]
where the units are [mjd2000, days, days, kg] + [n/a] + …. + [days]
Note
The resulting problem is non linearly constrained. The resulting trajectory is not time-bounded.
__init__
(*args)¶leg_tof=[1, 365.25 * 3], rest=[30., 365.25], mass=[800, 2000], Tmax=0.3, Isp=3000., traj_tof=365.25 * 6, objective=’mass’, c_tol=1e-05)
seq: list of pykep.planet defining the encounter sequence for the trajectoty (including the initial planet)
n_seg: list of integers containing the number of segments to be used for each leg (len(n_seg) = len(seq)-1)
t0: list of two pykep epochs defining the launch window
leg_tof: list of two floats defining the minimum and maximum time of each leg (days)
rest: list of two floats defining the minimum and maximum time the spacecraft can rest at one planet (days)
mass: list of two floats defining the minimum final spacecraft mass and the starting spacecraft mass (kg)
Tmax: maximum thrust (N)
Isp: engine specific impulse (sec)
traj_tof maximum total mission duration (days)
c_tol: tolerance on the constraints
plot
(*args)¶ax = prob.plot(x, ax=None)
x: encoded trajectory
ax: matplotlib axis where to plot. If None figure and axis will be created
[out] ax: matplotlib axis where to plot
Plots the trajectory represented by a decision vector x on the 3d axis ax
Example:
ax = prob.plot(x)
pykep.trajopt.
indirect_pt2pt
(_indirect_base)¶Represents an indirect trajectory optimisation problem between two Cartesian states with heliocentric dynamics. The class can be used as UDP in pagmo.
The decision chromosome is
z = [T, l0]
__init__
(*args)¶Constructs an instance of the pykep.trajopt.indirect_pt2pt
problem.
x0 (list
, tuple
, numpy.ndarray
): Departure state [m, m, m, m/s, m/s, m/s, kg].
xf (list
, tuple
, numpy.ndarray
): Arrival state [m, m, m, m/s, m/s, m/s, kg].
tof (list
): Transfer time bounds [days].
thrust (float
, int
): Spacecraft maximum thrust [N].
isp (float
, int
): Spacecraft specific impulse [s].
mu (float
): Gravitational parametre of primary body [m^3/s^2].
freetime (bool
): Activates final time transversality condition. Allows final time to vary.
alpha (float
, int
): Homotopy parameter (0 -quadratic control, 1 - mass optimal)
bound (bool
): Activates bounded control, in which the control throttle is bounded between 0 and 1, otherwise the control throttle is allowed to unbounded.
atol (float
, int
): Absolute integration solution tolerance.
rtol (float
, int
): Relative integration solution tolerance.
pykep.trajopt.
indirect_or2or
(_indirect_base)¶Represents an indirect trajectory optimisation problem between two orbits.
Decision chromosome is
z = [T, M0, Mf, l0]
__init__
(*args)¶Initialises pykep.trajopt.indirect_or2or
problem.
elem0 (list
, tuple
, numpy.ndarray
): Departure Keplerian elements (mutable eccentric anomaly).
elemf (list
, tuple
, numpy.ndarray
): Arrival Keplerian elements (mutable eccentric anomaly).
mass (float
, int
): Spacecraft wet mass [kg].
thrust (float
, int
): Spacecraft maximum thrust [N].
isp (float
, int
): Spacecraft specific impulse [s].
atol (float
, int
): Absolute integration solution tolerance.
rtol (float
, int
): Relative integration solution tolerance.
tof (list
): Transfer time bounds [days].
freetime (bool
): Activates final time transversality condition. Allows final time to vary.
alpha (float
, int
): Homotopy parametre, governing the degree to which the theoretical control law is intended to reduce propellant expenditure or energy.
bound (bool
): Activates bounded control, in which the control throttle is bounded between 0 and 1, otherwise the control throttle is allowed to unbounded.
mu (float
): Gravitational parametre of primary body [m^3/s^2].
pykep.trajopt.
indirect_pt2or
(_indirect_base)¶Represents an indirect trajectory optimisation problem between a Cartesian state and an orbit.
Decision chromosome is
z = [T, Mf, l0]
__init__
(*args)¶Initialises pykep.trajopt.indirect_pt2or
problem.
x0 (list
, tuple
, numpy.ndarray
): Departure state [m, m, m, m/s, m/s, m/s, kg].
elemf (list
, tuple
, numpy.ndarray
): Arrival Keplerian elements SI units. (mean anomaly will be changed).
mass (float
, int
): Spacecraft wet mass [kg].
thrust (float
, int
): Spacecraft maximum thrust [N].
isp (float
, int
): Spacecraft specific impulse [s].
atol (float
, int
): Absolute integration solution tolerance.
rtol (float
, int
): Relative integration solution tolerance.
tof (list
): Transfer time bounds [days].
freetime (bool
): Activates final time transversality condition. Allows final time to vary.
alpha (float
, int
): Homotopy parametre, governing the degree to which the theoretical control law is intended to reduce propellant expenditure or energy.
bound (bool
): Activates bounded control, in which the control throttle is bounded between 0 and 1, otherwise the control throttle is allowed to unbounded.
mu (float
): Gravitational parametre of primary body [m^3/s^2].
pykep.trajopt.
_launchers
¶This class contains a few functions that return the mass launchers can deliver to a certain declination / vinf.
Note
In pykep the object pykep.trajopt.launchers is already an instance of this class and is to be used as it has all the data preallocated upon import.
>>> import pykep as pk
>>> mass = pk.trajopt.launchers.soyuzf(4.5, 33.21)
>>> mass2 = pk.trajopt.launchers.atlas501(4.5, 33.21)
atlas501
(*args)¶atlas501(vinfs, decls)
Computes the mass that the Atlas 501 launcher can deliver to a certain vinf and declination. If the inputs are arrays, then a mesh is considered and the mass is returned on points of the mesh
vinfs (float
or array-like): the hyperbolic escape velocity in km/s
decls (float
or array-like): the declination in degrees
Numpy array containg the mass delivered to escape with said declinations and magnitudes.
soyuzf
(*args)¶soyuzf(vinfs, decls)
Computes the mass that the Soyutz-Fregat launcher can deliver to a certain vinf and declination. If the inputs are arrays, then a mesh is considered and the mass is returned on points of the mesh
vinfs (float
or array-like): the hyperbolic escape velocity in km/s
decls (float
or array-like): the declination in degrees
Numpy array containg the mass delivered to escape with said declinations and magnitudes.
ariane5
(*args)¶ariane5(vinfs, decls)
Computes the mass that the Ariane5 launcher can deliver to a certain vinf and declination, assuming a launch from Kourou. Data provided to ESOC by Arianespace when Ariane launch for ExoMars was an option. If the inputs are arrays, then a mesh is considered and the mass is returned on points of the mesh.
vinfs (float
or array-like): the hyperbolic escape velocity in km/s
decls (float
or array-like): the declination in degrees
Numpy array containg the mass delivered to escape with said declinations and magnitudes.
gym.
cassini1
= Cassini MGA direct tof encoding (Trajectory Optimisation Gym P1)¶This is an MGA problem inspired to the Cassini spacecraft interplanetary transfer to Saturn. The objective of this mission is to reach Saturn and to be captured by its gravity into an orbit having pericenter radius r_p=108950 km, and eccentricity e=0.98. The planetary fly-by sequence considered is E-VVEJ-S (as the one used by Cassini spacecraft). As objective function we use the total deltaV accumulated during the mission, including the launch deltaV and the various deltaV one needs to give at the planets and upon arrival to perform the final orbit injection.
Note
A more complex representation of this interplanetary trajectory transfer can be found
in pykep.trajopt.gym.cassini2
gym.
cassini1_a
= Cassini MGA alpha tof encoding (Trajectory Optimisation Gym P2)¶This is the same as pykep.trajopt.gym.cassini1
,
but the time of flights are encoded using the alpha encoding technique.
gym.
cassini1_n
= Cassini1 MGA eta tof encoding (Trajectory Optimisation Gym P3)¶This is the same as pykep.trajopt.gym.cassini1
,
but the time of flights are encoded using the eta encoding technique.
gym.
em3imp
= Earth-Mars 3 impulses (Trajectory Optimisation Gym P4)¶This is a simple Earth Mars transfer that has the possibility to use multiple impulses. While the actual solution may not need more than a few, the problem is used to study the alpha encoding for manouvre timings and its evolvability.
Note
The problem is rather simple and is used to study the multiple impulse transription (using alpha encoding)
comparing it to the folowing problems pykep.trajopt.gym.em5imp
and pykep.trajopt.gym.em7imp
having increasing number of impulses.
gym.
em5imp
= Earth-Mars 5 impulses (Trajectory Optimisation Gym P5)¶This is the same as pykep.trajopt.gym.em3imp
,
but the number of impulses is increased.
gym.
em7imp
= Earth-Mars 7 impulses (Trajectory Optimisation Gym P6)¶This is the same as pykep.trajopt.gym.em3imp
,
but the number of impulses is increased.
gym.
eve_mga1dsm
= Earth-Venus-Earth mga-1dsm (Trajectory Optimisation Gym P7-9)¶This is an Earth - Venus - Earth transfer where 1 deep space manouvre is allowed at each leg.
Note
Together with the problems pykep.trajopt.gym.eve_mga1dsm_a
and pykep.trajopt.gym.eve_mga1dsm_n
this simple problem is used to study the alpha, eta and direct encoding for the time of flights of an mga-1DSM
type of problem.
gym.
eve_mga1dsm_a
= Earth-Venus-Earth mga-1dsm (Trajectory Optimisation Gym P7-9)¶This is the same as pykep.trajopt.gym.eve_mga1dsm
,
but the time of flights are encoded using the alpha encoding technique.
gym.
eve_mga1dsm_n
= Earth-Venus-Earth mga-1dsm (Trajectory Optimisation Gym P7-9)¶This is the same as pykep.trajopt.gym.eve_mga1dsm
,
but the time of flights are encoded using the eta encoding technique.
gym.
cassini2
= Cassini 2 (Trajectory Optimisation Gym P11)¶This is an MGA-1DSM problem inspired to the Cassini spacecraft interplanetary transfer to Saturn. The objective of this mission is to reach Saturn and to be captured by its gravity into an orbit having pericenter radius r_p=108950 km, and eccentricity e=0.98. The planetary fly-by sequence considered is E-VVEJ-S (as the one used by Cassini spacecraft). As objective function we use the total deltaV accumulated during the mission by the spacecraft, including the one necessary for the final orbit injection.
Note
A similar problem was also part of the ESA’s GTOP database with the same name, but different implementation details and mission definition. They should not be comapred.
gym.
rosetta
= Rosetta (Trajectory Optimisation Gym P10)¶This represents a rendezvous mission to the comet 67P/Churyumov-Gerasimenko modelled as an MGA-1DSM transfer. The fly-by sequence selected (i.e. E-EMEE-C) is similar to the one planned for the spacecraft Rosetta. The objective function considered is the total mission delta V. No launcher model is employed and a final randezvous with the comet is included in the delta V computations.
Note
A similar problem was also part of the ESA’s GTOP database with the same name, but different implementation details and mission definition. They should not be comapred.
pykep.trajopt.gym.
tandem
¶alias of pykep.trajopt.gym._tandem._tandem_udp
This class represents a rendezvous mission to Saturn modelled as an MGA-1DSM transfer. Mission parameters are inspired to the TandEM mission. A launcher model (i.e. Atlas 501) is also used, so that the final mass delivered at Saturn is the main objective of this optimization problem.
The problem draws inspiration from the work performed in April 2008 by the European Space Agency working group on mission analysis on the mission named TandEM. TandEM is an interplanetary mission aimed at reaching Titan and Enceladus (two moons of Saturn).
Note
The Titan and Enceladus Mission (TandEM), an ambitious scientific mission to study the Saturnian system with particular emphasis on the moons Titan and Enceladus, was selected in October 2007 as a candidate mission within the ESA Cosmic Vision plan. In February 2009, TandEM exited the Cosmic Vision programme when ESA and NASA chose EJSM-Laplace as the L-class outer Solar System mission candidate.
Note
A significantly similar version of this problem was part of the no longer maintained GTOP database, https://www.esa.int/gsp/ACT/projects/gtop/gtop.html. The exact definition is, though, different and results cannot thus not be compared to those posted in GTOP.
_tandem_udp.
__init__
(prob_id=1, constrained=True)¶The TandEM problem of the trajectory gym consists in 48 different instances varying in fly-by sequence and the presence of a time constraint.
prob_id (int
): The problem id defines the fly-by sequence.
constrained (bool
): Activates the constraint on the time of flight
(fitness will thus return two numbers, the objectove function and the inequality constraint violation)
gym.
juice
= Juice (Trajectory Optimization Gym P13-14)¶This class represents a rendezvous mission to Jupiter modelled as an MGA-1DSM transfer. The selected fly-by sequence, E-EVEME-J, and other parameters are inspired to the ESA JUICE mission. A launcher model (i.e. Ariane 5) is also used, so that the final mass delivered at Saturn is the main objective of this optimization problem.
Note
JUICE - JUpiter ICy moons Explorer - is the first large-class mission in ESA’s Cosmic Vision 2015-2025 programme. Planned for launch in 2022 and arrival at Jupiter in 2029, it will spend at least three years making detailed observations of the giant gaseous planet Jupiter and three of its largest moons, Ganymede, Callisto and Europa.
gym.
juice_mo
= Juice (Trajectory Optimization Gym P13-14)¶This is the multiobjective version of pykep.trajopt.gym.juice
.
Time of flights are encoded using the alpha encoding.
gym.
messenger
= <pykep.trajopt.gym._messenger._messenger_udp object>¶This class represents a rendezvous mission to Mercury modelled as an MGA-1DSM transfer. The selected fly-by sequence, E-VVMeMeMe-Me, and other parameters are inspired to the Messenger mission. We have only omitted the first Earth fly-by that was used to correct for launcher performances, since we here do not make use of a luncher model. As far as chemical propelled interplanetary trajectories go, this particular one is particularly complex and difficult to design. The time of flights among successive Mercury fly-bys allow for multiple rvolutions and resonances, making optimization techniques struggle to find the correct combination. The amount of specialistic knowledge that needs to be used to obtain a successfull design is significant. Finding a global optimization approach able to find a good trajectory in complete autonomy without making use of additional problem knowledge is possible, but limiting the number of fitness call is difficult.
Note
A similar problem was also part of the ESA’s GTOP database with the same name, but different implementation details and mission definition. They should not be comapred.