The orbit plotting module

Name

Type

Description

pykep.orbit_plots.plot_planet()

function

Plots a planet and its orbit

pykep.orbit_plots.plot_lambert()

function

Plots a lambert arc (also multirev)

pykep.orbit_plots.plot_kepler()

function

Plots a keplerian propagated arc

pykep.orbit_plots.plot_taylor()

function

Plots a constant thrust propagated arc

pykep.orbit_plots.plot_sf_leg()

function

Plots a Sims_Flanagan leg

Detailed Documentation

pykep.orbit_plots.plot_planet(*args)

ax = plot_planet(plnt, t0=0, tf=None, N=60, units=1.0, color=’k’, alpha=1.0, s=40, legend=(False, False), axes=None):

  • axes: 3D axis object created using fig.gca(projection=’3d’)

  • plnt: pykep.planet object we want to plot

  • t0: a pykep.epoch or float (mjd2000) indicating the first date we want to plot the planet position

  • tf: a pykep.epoch or float (mjd2000) indicating the final date we want to plot the planet position.

    if None this is computed automatically from the orbital period (prone to error for non periodic orbits)

  • units: the length unit to be used in the plot

  • color: color to use to plot the orbit (passed to matplotlib)

  • s: planet size (passed to matplotlib)

  • legend 2-D tuple of bool or string: The first element activates the planet scatter plot,

    the second to the actual orbit. If a bool value is used, then an automated legend label is generated (if True), if a string is used, the string is the legend. Its also possible but deprecated to use a single bool value. In which case that value is used for both the tuple components.

Plots the planet position and its orbit.

Example:

import pykep as pk
import matplotlib.pyplot as plt

fig = plt.figure()
ax = fig.gca(projection = '3d')
pl = pk.planet.jpl_lp('earth')
t_plot = pk.epoch(219)
ax = pk.orbit_plots.plot_planet(pl, ax = ax, color='b')

pykep.orbit_plots.plot_lambert(*args)

ax = plot_lambert(l, N=60, sol=0, units=’pykep.AU’, legend=’False’, axes=None, alpha=1.)

  • axes: 3D axis object created using fig.gca(projection=’3d’)

  • l: pykep.lambert_problem object

  • N: number of points to be plotted along one arc

  • sol: solution to the Lambert’s problem we want to plot (must be in 0..Nmax*2)

    where Nmax is the maximum number of revolutions for which there exist a solution.

  • units: the length unit to be used in the plot

  • color: matplotlib color to use to plot the line

  • legend: when True it plots also the legend with info on the Lambert’s solution chosen

Plots a particular solution to a Lambert’s problem

Example:

import pykep as pk
import matplotlib.pyplot as plt

fig = plt.figure()
ax = fig.gca(projection='3d')

t1 = pk.epoch(0)
t2 = pk.epoch(640)
dt = (t2.mjd2000 - t1.mjd2000) * pk.DAY2SEC

pl = pk.planet.jpl_lp('earth')
pk.orbit_plots.plot_planet(pl, t0=t1, axes=ax, color='k')
rE,vE = pl.eph(t1)

pl = pk.planet.jpl_lp('mars')
pk.orbit_plots.plot_planet(pl, t0=t2, axes=ax, color='r')
rM, vM = pl.eph(t2)

l = lambert_problem(rE,rM,dt,pk.MU_SUN)
pk.orbit_plots.plot_lambert(l, ax=ax, color='b')
pk.orbit_plots.plot_lambert(l, sol=1, axes=ax, color='g')
pk.orbit_plots.plot_lambert(l, sol=2, axes=ax, color='g', legend = True)

plt.show()

pykep.orbit_plots.plot_kepler(*args)

ax = plot_kepler(r0, v0, tof, mu, N=60, units=1, color=’b’, label=None, axes=None):

  • axes: 3D axis object created using fig.gca(projection=’3d’)

  • r0: initial position (cartesian coordinates)

  • v0: initial velocity (cartesian coordinates)

  • tof: propagation time

  • mu: gravitational parameter

  • N: number of points to be plotted along one arc

  • units: the length unit to be used in the plot

  • color: matplotlib color to use to plot the line

  • label adds a label to the plotted arc.

Plots the result of a keplerian propagation

Example:

import pykep as pk
pi = 3.14
pk.orbit_plots.plot_kepler(r0 = [1,0,0], v0 = [0,1,0], tof = pi/3, mu = 1)

pykep.orbit_plots.plot_taylor(*args)

ax = plot_taylor(r0, v0, m0, thrust, tof, mu, veff, N=60, units=1, color=’b’, legend=False, axes=None):

  • axes: 3D axis object created using fig.gca(projection=’3d’)

  • r0: initial position (cartesian coordinates)

  • v0: initial velocity (cartesian coordinates)

  • m0: initial mass

  • u: cartesian components for the constant thrust

  • tof: propagation time

  • mu: gravitational parameter

  • veff: the product Isp * g0

  • N: number of points to be plotted along one arc

  • units: the length unit to be used in the plot

  • color: matplotlib color to use to plot the line

  • legend: when True it plots also the legend

Plots the result of a taylor propagation of constant thrust

Example:

import pykep as pk
import matplotlib.pyplot as plt
pi = 3.14

fig = plt.figure()
ax = fig.gca(projection = '3d')
pk.orbit_plots.plot_taylor([1,0,0],[0,1,0],100,[1,1,0],40, 1, 1, N = 1000, axes = ax)
plt.show()

pykep.orbit_plots.plot_sf_leg(*args)

ax = plot_sf_leg(leg, N=5, units=1, color=’b’, legend=False, no_trajectory=False, axes=None):

  • axes: 3D axis object created using fig.gca(projection=’3d’)

  • leg: a pykep.sims_flanagan.leg

  • N: number of points to be plotted along one arc

  • units: the length unit to be used in the plot

  • color: matplotlib color to use to plot the trajectory and the grid points

  • legend when True it plots also the legend

  • plot_line: when True plots also the trajectory (between mid-points and grid points)

Plots a Sims-Flanagan leg

Example:

from mpl_toolkits.mplot3d import Axes3D
import matplotlib.pyplot as plt

fig = plt.figure()
ax = fig.gca(projection='3d')
t1 = epoch(0)
pl = planet_ss('earth')
rE,vE = pl.eph(t1)
plot_planet(pl,t0=t1, units=AU, axes=ax)

t2 = epoch(440)
pl = planet_ss('mars')
rM, vM = pl.eph(t2)
plot_planet(pl,t0=t2, units=AU, axes=ax)

sc = sims_flanagan.spacecraft(4500,0.5,2500)
x0 = sims_flanagan.sc_state(rE,vE,sc.mass)
xe = sims_flanagan.sc_state(rM,vM,sc.mass)
l = sims_flanagan.leg(t1,x0,[1,0,0]*5,t2,xe,sc,MU_SUN)

plot_sf_leg(l, units=AU, axes=ax)