Fly-by routines

Fly-by routines#

The gravity assist technique, widely used in the Pioneer and Voyager missions, was developed by Michael Andrew Minovitch when he was a UCLA graduate student working as intern at NASA’s Jet Propulsion Laboratory. Many authors later refined the technique which, avoiding to solve directly the full three body problem, is based on patching the planetocentric hyperbola with incoming and outgoing (elliptic) trajectories. In pykep we offer the basic routines that allow to use this technqiue in the context of trajectory optimization and patched conics propagations.


Minovitch fly-by#

pykep.fb_con(v_rel_in, v_rel_out, mu, safe_radius)#

Alternative signature: fb_con(v_rel_in, v_rel_out, planet)

Computes the constraint violation during a fly-by modelled as an instantaneous rotation of the incoming and outgoing relative velocities (Mivovitch). The two must be identical in magnitude (equality constraint) and the angle \(\alpha\) between them must be less than the \(\alpha_{max}\): the maximum value allowed for that particular planet (inequality constraint), as computed from its gravitational parameter and safe radius using the formula:

\[\alpha_{max} = - 2 \arcsin\left(\frac{1}{e_{min}}\right)\]

where:

\[e_{min} = 1 + V_{\infty}^2 \frac{R_{safe}}{\mu}\]

Note

This function is often used in the multiple gravity assist low-thrust (MGA-LT) encoding of an interplanetary trajectory where multiple low-thrust legs are patched at the fly-by planets by forcing satisfaction for these constraints.

Args:

v_rel_in (list (3,)): Cartesian components of the incoming relative velocity.

v_rel_out (list (3,)): Cartesian components of the outgoing relative velocity.

mu (float): planet gravitational parameter

safe_radius (float): planet safe radius

planet (planet): planet (in which case mu and safe_radius will be extracted from this object). Note: this signature is slower and to be avoided.

Returns:

tuple [float, float]: The equality constraint violation (defined as the difference between the squared velocities) and the inequality constraint violation (negative if satisfied).

Examples:
>>> import pykep as pk
>>> eq, ineq = pk.fb_con(v_rel_in = [10.,1.,-4.], v_rel_out = [10.,1.,-4.], mu = 1., safe_radius = 1.)
pykep.fb_dv(v_rel_in, v_rel_out, mu, safe_radius)#

Alternative signature: fb_dv(v_rel_in, v_rel_out, planet)

Computes the \(\Delta V\) necessary to perform a fly-by modelled as an instantaneous rotation of the incoming and outgoing relative velocities (Mivovitch). If planetary gravity is not enough to patch the incoming and outcoming conditions (i.e. the fb_con() returns some constraint violation) a \(\Delta V\) is assumed at the end of the planetocentric hyperbola.

Note

This function is often used in the multiple gravity assist (MGA) encoding of an interplanetary trajectory where multiple Lambert arcs are patched at the fly-by planets by applying a hopefully vanishing \(\Delta V\).

Args:

v_rel_in (list (3,)): Cartesian components of the incoming relative velocity.

v_rel_out (list (3,)): Cartesian components of the outgoing relative velocity.

mu (float): planet gravitational parameter

safe_radius (float): planet safe radius

planet (planet): planet (in which case mu and safe_radius will be extracted from this object). Note: this signature is slower and to be avoided.

Returns:

float: The magnitude of the required \(\Delta V\)

Examples:
>>> import pykep as pk
>>> DV = pk.fb_dv(v_rel_in = [10.,1.,-4.], v_rel_out = [10.,1.,-4.], mu = 1., safe_radius = 1.)
pykep.fb_vout(v_in, v_pla, rp, beta, mu)#

Propagates incoming conditions through a planetary encounter (fly-by) assuming an instantaneous rotation of magnitude \(\delta\) of the incoming and outgoing relative velocities (Mivovitch). The planetocentric hyperbola (hence the angle \(\delta\)) is fully determined by the incoming condition \(v_{\infty} = \mathbf v_{in} - \mathbf v_{pla}\) as well as by its pericenter \(r_p\) and an angle \(\beta\) defining the orientation of the orbital plane. Eventually the outgoing conditions are computed as:

\[\mathbf v_{out} = \mathbf v_{in} + \mathbf v_{\infty}^{out}\]
\[\mathbf v_{\infty}^{out} = |\mathbf v_{\infty}^{in}| \left( \cos\delta\hat{\mathbf b}_1 +\sin\delta\cos\beta\hat{\mathbf b}_2 +\sin\delta\sin\beta\hat{\mathbf b}_3 \right)\]

where the \([\hat{\mathbf b}_1, \hat{\mathbf b}_2, \hat{\mathbf b}_3]\) frame is defined by the incoming relative velocity (normalized), \(\hat{\mathbf b}_1 \times \mathbf v_{pla}\) (normalized) and completing the right handed frame.

Note

This function is often used in the multiple gravity assist with DSM (MGA-DSM) encoding of an interplanetary trajectory where multiple impulsive manouvres are allowed betwwen planetary fly-bys.

Args:

v_in (list (3,)): Cartesian components of the incoming velocity (note: this is NOT the relative velocity, rather the absolute and in the same frame as v_pla).

v_pla (list (3,)): Cartesian components of the planet velocity.

rp (float): planetocentric hyperbola pericenter radius.

beta (float): planetocentric hyperbola plane angle.

mu (float): planet gravitational parameter.

Returns:

list (3,): The outgoing velocity in the frame of v_pla (inertial).

Examples:
>>> import pykep as pk
>>> DV = pk.fb_vout(v_in = [1.,1.,1], v_pla = [10.,1.,-4.], rp = 1., mu = 1., beta=1.2)