Trajectory Optimization#

In pykep both direct and indirect optimization transcriptions as well as evolutionary encodings are provided to perform spacecraft trajectory optimization. Most direct techniques provided are based on the “Sims” transcription [SF97], while indirect techniques are based on some version of Pontryagin Maximum Principle. The evolutionary encoding are mostly based on the work performed at ESA’ Advanced Concepts Team ([Izz10], [ISM+13]).

Most of the classes in this Trajectory Optimization module are provided as optimization problems (NLP) compatible to the pagmo software suite.

Direct#

When solving optimal control problems (OCPs), the unknown is a function (i.e. the control), hence one is dealing with an infinite dimensional search space. A common approach to overcome this difficulty is to use direct methods. These are based on the idea of discretizing the problem introducing some time grid and thus transforming the OCP into a NLP (Non Linear Programming) problem to then solve it using numerical solvers available for the task.


class pykep.trajopt.direct_point2point(rvs=[array([1.49597871e+11, 1.49597871e+10, -1.49597871e+10]), array([5956.93836686, 29784.69183431, -5956.93836686])], rvf=[array([-1.79517445e+11, -1.49597871e+10, 1.49597871e+10]), array([5956.93836686, -30469.7397465, 13105.2644071])], ms=1000, mu=1.3271244004127942e+20, max_thrust=0.12, isp=3000, tof_bounds=[80, 400], mf_bounds=[200.0, 1000.0], nseg=10, cut=0.6, mass_scaling=1000, r_scaling=149597870700.0, v_scaling=29784.69183430911, with_gradient=True)#

Represents the optimal low-thrust transfer between two fixed points using a direct method.

This problem works internally using the sims_flanagan and manipulates its transfer time T, final mass mf and the controls as to link two fixed points in space with a low-thrust trajectory.

It can be used to better profile and understand performances of optimizers on this type of direct approach, but has a limited use in the design of interplanetary trajectories as per the fixed point limitation.

The decision vector is:

z = [mf, throttles, tof]

where throttles is a vector of throttles structures as [u0x, u0y,u0z, …]. By throttles we intend non dimensiona thrust levels in [0,1].

Initializes the direct_point2point instance with given parameters.

Args:

rvs (list): Initial position and velocity vectors. Defaults to two vectors scaled by AU and Earth’s velocity.

rvf (list): Final position and velocity vectors. Defaults to two vectors scaled by AU and Earth’s velocity.

ms (float): Initial spacecraft mass in kg. Defaults to 1000 kg.

mu (float): Gravitational parameter, default is for the Sun (MU_SUN).

max_thrust (float): Maximum thrust in Newtons. Defaults to 0.12 N.

isp (float): Specific impulse in seconds. Defaults to 3000 s.

tof_bounds (list): Bounds for time of flight in days. Defaults to [0, 400] days.

mf_bounds (list): Bounds for final mass in kg. Defaults to [200.0, 1000.0] kg.

nseg (int): Number of segments for the trajectory. Defaults to 10.

cut (float): Cut parameter for the sims_flanagan. Defaults to 0.6.

mass_scaling (float): Scaling factor for mass (used to scale constraints). Defaults to 1000.

r_scaling (float): Scaling factor for distance, (used to scale constraints). Defaults AU (AU).

v_scaling (float): Scaling factor for velocity (used to scale constraints). Defaults the Earth’s velocity (EARTH_VELOCITY).

with_gradient (bool): Indicates if gradient information should be used. Defaults True.

plot(x, ax=None, units=149597870700.0, show_midpoints=False, show_gridpoints=False, show_throttles=False, length=0.1, arrow_length_ratio=0.05, **kwargs)#

Plots the trajectory leg 3D axes.

Args:

x (list): The decision vector containing final mass, thrust direction, and time of flight.

ax (mpl_toolkits.mplot3d.axes3d.Axes3D, optional): The 3D axis to plot on. Defaults to None.

units (float, optional): The unit scale for the plot. Defaults to _pk.AU.

show_midpoints (bool, optional): Whether to show midpoints on the trajectory. Defaults to False.

show_gridpoints (bool, optional): Whether to show grid points on the trajectory. Defaults to False.

show_throttles (bool, optional): Whether to show throttle vectors. Defaults to False.

length (float, optional): Length of the throttle vectors. Defaults to 0.1.

arrow_length_ratio (float, optional): Arrow length ratio for the throttle vectors. Defaults to 0.05.

**kwargs: Additional keyword arguments for the plot.

Returns:

mpl_toolkits.mplot3d.axes3d.Axes3D: The 3D axis with the plotted trajectory.

pretty(x)#

Prints a detailed representation of the Point to point problem.

Args:

x (list): The decision vector containing final mass, thrust direction, and time of flight.


class pykep.trajopt.direct_pl2pl(pls, plf, ms=1500, mu=_pk.MU_SUN, max_thrust=0.12, isp=3000, t0_bounds=[6700.0, 6800.0], tof_bounds=[200.0, 300.0], mf_bounds=[1300.0, 1500.0], vinfs=3.0, vinff=0.0, nseg=10, cut=0.6, mass_scaling=1500, r_scaling=pk.AU, v_scaling=pk.EARTH_VELOCITY, with_gradient=True)#

Represents the optimal low-thrust transfer between two planet using a direct method.

This problem works internally using the sims_flanagan and manipulates its initial and final states, as well as its transfer time T, final mass mf and the controls as to link the two planets with a low-thrust trajectory.

The particular transcription used is suitable only for few revolutions, after which convergence will start to be problematic.

The decision vector is:

z = [t0, mf, Vsx, Vsy, Vsz, Vfx, Vfy, Vfz, throttles, tof] - all in S.I. units except t0 and tof in days

where throttles is a vector of throttles structured as [u0x, u0y,u0z, …]. By throttles we intend non dimensiona thrust levels in [0,1].

Args:

pls (planet): Initial planet. Defaults to jpl_lp Earth.

plf (planet): Final planet. Defaults to jpl_lp Mars.

ms (float): Initial spacecraft mass in kg. Defaults to 1000 kg.

mu (float): Gravitational parameter, default is for the Sun (MU_SUN).

max_thrust (float): Maximum thrust in Newtons. Defaults to 0.12 N.

isp (float): Specific impulse in seconds. Defaults to 3000 s.

t0_bounds (list): Bounds for departure epoch in MJD2000. Defaults to [6700.0, 6800.0].

tof_bounds (list): Bounds for time of flight in days. Defaults to [200, 300] days.

mf_bounds (list): Bounds for final mass in kg. Defaults to [1300.0, 1500.0] kg.

vinfs (float): Allowed magnitude for the departure’s relative velocity in km/s. Defaults to 3.

vinff (float): Allowed magnitude for the arrival’s relative velocity in km/s. Defaults to 0.

nseg (int): Number of segments for the trajectory. Defaults to 10.

cut (float): Cut parameter for the sims_flanagan. Defaults to 0.6.

mass_scaling (float): Scaling factor for mass (used to scale constraints). Defaults to 1500.

r_scaling (float): Scaling factor for distance, (used to scale constraints). Defaults AU (AU).

v_scaling (float): Scaling factor for velocity (used to scale constraints). Defaults the Earth’s velocity (EARTH_VELOCITY).

with_gradient (bool): Indicates if gradient information should be used. Defaults True.

plot(x, ax=None, units=149597870700.0, show_midpoints=False, show_gridpoints=False, show_throttles=False, length=0.1, arrow_length_ratio=0.05, **kwargs)#

Plots the trajectory leg 3D axes.

Args:

x (list): The decision vector containing final mass, thrust direction, and time of flight.

ax (mpl_toolkits.mplot3d.axes3d.Axes3D, optional): The 3D axis to plot on. Defaults to None.

units (float, optional): The unit scale for the plot. Defaults to _pk.AU.

show_midpoints (bool, optional): Whether to show midpoints on the trajectory. Defaults to False.

show_gridpoints (bool, optional): Whether to show grid points on the trajectory. Defaults to False.

show_throttles (bool, optional): Whether to show throttle vectors. Defaults to False.

length (float, optional): Length of the throttle vectors. Defaults to 0.1.

arrow_length_ratio (float, optional): Arrow length ratio for the throttle vectors. Defaults to 0.05.

**kwargs: Additional keyword arguments for the plot.

Returns:

mpl_toolkits.mplot3d.axes3d.Axes3D: The 3D axis with the plotted trajectory.

pretty(x)#

Prints a human readable representation of the transfer.

Args:

x (list): The decision vector containing final mass, thrust direction, and time of flight.


Evolutionary Encodings#

Some type of interplanetary trajectories can be evolved: shocking?. This AI approach to trajectory design resulted in two silver ([TAF+05], [PGJ+18]) and one gold ([ISM+13]) Humies award for human-competitive results that were produced by any form of genetic and evolutionary computation. In pykep we offer some classes that help instantiating these type of transfers amenable to evolutionary techniques.


class pykep.trajopt.mga(seq, t0, tof, vinf, multi_objective=False, tof_encoding="direct"", orbit_insertion=False, e_target=None, rp_target=None)#

The Multiple Gravity Assist (MGA) encoding of an interplanetary trajectory.

This class may be used as a User Defined Problem (UDP) for the pygmo (http://esa.github.io/pygmo/) optimisation suite.

The decision vector (chromosome) is:

direct encoding: z = [t0, T1, T2 ... ] in [mjd2000, days, days ... ]
alpha encoding: z = [t0, T, a1, a2 ...] in [mjd2000, days, nd, nd ... ]
eta encoding: z = [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 talk about 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 = (T_{max} - \sum_{j=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).

Args:

seq (list [planet]): sequence of planetary encounters including the departure body.

t0 (list [float or epoch]): lower and upper bounds for the launch epoch. When floats are used MJD2000 is assumed.

tof (list or float): defines the 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 total time-of-flight. If tof_encoding is ‘eta’ tof is a float defining an upper bound for 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 recommended

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)

max_revs (int): maximal number of revolutions for lambert transfer

plot(x, ax=None, units=149597870700.0, N=60, c_orbit='dimgray', c='indianred', leg_ids=[], figsize=(5, 5), **kwargs)#

Plots the trajectory leg 3D axes.

Args:

x (list): The decision vector in the correct tof encoding.

ax (mpl_toolkits.mplot3d.axes3d.Axes3D, optional): The 3D axis to plot on. Defaults to None.

units (float, optional): The unit scale for the plot. Defaults to pk.AU.

N (int, optional): The number of points to use when plotting the trajectory. Defaults to 60.

figsize (tuple): The figure size (only used if a*ax* is None and axis ave to be created.), Defaults to (5, 5).

leg_ids (list): selects the legs to plot. Optional, defaults to all legs.

**kwargs: Additional keyword arguments to pass to the trajectory plot (all Lambert arcs)

Returns:

mpl_toolkits.mplot3d.axes3d.Axes3D: The 3D axis where the trajectory was plotted.

pretty(x)#

Prints a human readable representation of the transfer.

Args:

x (list): The decision vector in the correct tof encoding.

to_planet(x: List[float])#

Returns a planet representing the spacecraft trajectory encoded into x.

Args:

x (list): The decision vector in the correct tof encoding.

Example:

sc = mga_udp.to_planet(population.champion_x) r, v = sc.eph(7000.)

Utilities#

In order to facilitate the use of the classes in this module, some utilities are provided.


class pykep.trajopt._launchers#

This class contains an API to access the tabulated data on the mass launchers can deliver to a certain declination / vinf.

Note

In pykep the object pykep.trajopt.launchers is an instance of this class created upon import and is used to access launchers data.

Examples:
>>> import pykep as pk
>>> mass = pk.trajopt.launchers.soyuzf(4.5, 33.21)
>>> mass2 = pk.trajopt.launchers.atlas501(4.5, 33.21)
ariane5(vinfs, decls)#

Estimates 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 constructed and the mass is returned on all points of the mesh.

Args:

vinfs (float or numpy.ndarray (N)): the hyperbolic escape velocities in km/s.

decls (float or numpy.ndarray (N)):: the declinations in degrees (in the ECF frame).

Returns:

numpy.ndarray (N, N): masses delivered to escape with said declinations and magnitudes.

atlas501(vinfs, decls)#

Estimates the mass that the Atlas 501 launcher can deliver to a certain vinf and declination. If the inputs are arrays, then a mesh is constructed and the mass is returned on all points of the mesh.

Args:

vinfs (float or numpy.ndarray (N)): the hyperbolic escape velocities in km/s.

decls (float or numpy.ndarray (N)):: the declinations in degrees (in the ECF frame).

Returns:

numpy.ndarray (N, N): masses delivered to escape with said declinations and magnitudes.

atlas551(vinfs)#

Estimates the mass that the Atlas 551 launcher can deliver to a certain vinf

Args:

vinfs (float or numpy.ndarray (N)): the hyperbolic escape velocities in km/s.

Returns:

numpy.ndarray (N): masses delivered to escape with said vinf magnitudes.

soyuzf(vinfs, decls)#

Estimates the mass that the Soyutz-Fregat launcher can deliver to a certain vinf and declination. If the inputs are arrays, then a mesh is constructed and the mass is returned on all points of the mesh.

Args:

vinfs (float or numpy.ndarray (N)): the hyperbolic escape velocities in km/s.

decls (float or numpy.ndarray (N)):: the declinations in degrees (in the ECF frame).

Returns:

numpy.ndarray (N, N): masses delivered to escape with said declinations and magnitudes.