The core module, imported in the main pykep namespace, contains the classes and helper functions to perform basic space-flight mechanics computations. You can browse the list of classes and functions below. Find the detailed documentation later on this page.
Name |
Type |
Description |
---|---|---|
class |
represents an epoch (i.e. a fixed point in time) |
|
function |
helper function to construct an epoch from a string containing a date in the format YYYY-MM-DD HH:MM:SS |
|
function |
helper function to construct an epoch from a string containing a date in the ISO format YYYYMMDDTHHMMSS |
|
class |
solves the multirevolution lambert problem |
|
function |
propagates pure keplerian motion using Lagrange coefficients and universal variables |
|
function |
propagates keplerian motion disturbed by a constant inertial thrust using Taylor integration method |
|
function |
returns violation of velocity and angular constraint during a fly-by |
|
function |
returns the violation of the velocity and angular constraint during a fly-by in terms of one single DV |
|
function |
propoagates forward a fly-by hyperbola returning the new inetrial velocity of a spacecraft after the planetary encounter |
|
function |
computes the (parabolic) time-of-flight from the Barker equation |
|
function |
Transforms r and v into the osculating orbital elements |
|
function |
Transforms osculating orbital elements into r and v |
|
function |
Transforms osculating orbital elements into modified equinoctial |
|
function |
Transforms modified equinoctial into osculating orbital elements |
|
function |
Transforms modified equinoctial elements into r and v |
|
function |
Transforms r and v into modified equinoctial elements |
pykep.
lambert_problem
(*args)¶Represents a multiple revolution Lambert’s problem
__init__
(*args)¶lambert_problem(r1 = [1,0,0], r2 = [0,1,0], tof = pi/2, mu = 1., cw = False, max_revs = 0)
r1: starting position (x1,y1,z1)
r2: final position (x2,y2,z2)
tof: time of flight
mu: gravitational parameter (default is 1)
cw: True for retrograde motion (clockwise), False if counter-clock wise
max_revs: Maximum number of multirevs to be computed
Note
Units need to be consistent. The multirev Lambert’s problem will be solved upon construction and its solution stored in data members.
Example:
l = lambert_problem(r1 = [1,0,0], r2 = [0,1,0], tof = 5 * pi / 2.)
get_v1
()¶Returns a sequence of vectors containing the velocities at r1 of all computed solutions to the Lambert’s Problem
Solutions are stored in order 0 rev, 1rev, 1rev, 2rev, 2rev, …
Example (extracts v1 for the 0 revs solution):
v10 = l.get_v1()[0]
get_v2
()¶Returns a sequence of vectors containing the velocities at r2 of all computed solutions to the Lambert’s Problem
Solutions are stored in order 0 rev, 1rev, 1rev, 2rev, 2rev, …
Example (extracts v2 for the 0 revs solution):
v20 = l.get_v2()[0]
get_x
()¶Returns a sequence containing the x values of all computed solutions to the Lambert’s Problem
Solutions are stored in order 0 rev, 1rev, 1rev, 2rev, 2rev, …
Example (extracts x for the 0 revs solution):
x0 = l.get_x()[0]
get_Nmax
()¶Returns the maximum number of revolutions allowed. The total number of solution to the Lambert’s problem will thus be n_sol = Nmax*2 + 1
Example:
Nmax = l.get_Nmax()
n_sol = Nmax*2+1
get_iters
()¶Returns a sequence containing the number of iterations employed to compute each solution to the Lambert’s Problem
Solutions are stored in order 0 rev, 1rev, 1rev, 2rev, 2rev, …
Example (extracts the number of iterations employed for the 0 revs solution):
p0 = l.get_iters()[0]
pykep.
epoch
(*args)¶A precise point in time. Julian dates are supported.
__init__
(*args)¶pykep.epoch(julian_date=0, julian_date_type=”mjd2000”)
julian_date: a julian date, defaults to 0
julian_date_type: julian date type, one of “jd”, “mjd” and “mjd2000”, defaults to “mjd2000”
Examples:
e1 = epoch(0)
e1 = epoch(0,"mjd2000")
e2 = epoch(54333, "mjd")
jd
¶Returns the Julian Date
Example:
jd = e.jd
mjd
¶Returns the Modified Julian Date
Example:
jd = e.mjd
mjd2000
¶Returns the Modified Julian Date 2000
Example:
jd = e.mjd2000
pykep.
epoch_from_string
(s)¶pykep.epoch_from_string(s)
s: a string containing the date in the format ‘YYYY-MM-DD HH:MM:SS’
Returns a pykep.epoch
Note:
Excess digits in fractional seconds will be dropped.
Ex: '1:02:03.123456999' => '1:02:03.123456'. This behavior depends on the precision
defined in astro_constant.hpp and used to compile the keplerian toolbox library.
The function is based on the corresponding boost date_time library function
Example:
e = pykep.epoch_from_string('2002-01-20 23:59:54.003')
pykep.
epoch_from_iso_string
(s)¶pykep.epoch_from_iso_string(s)
s: string containing a date in the ISO format ‘YYYYMMDDTHHMMSS’
Returns a pykep.epoch
object constructed from a from a non delimited iso form string containing a date.
The function is based on the corresponding boost date_time library function
Example:
e = pykep.epoch_from_iso_string('20020120T235954')
pykep.
propagate_lagrangian
(*args)¶propagate_lagrangian(r0 = [1,0,0], v0 = [0,1,0], tof = pi/2, mu = 1)
r: start position, x,y,z
v: start velocity, vx,vy,vz
tof: propagation time
mu: central body gravity constant
Returns a tuple (rf, vf) containing the final position and velocity after the propagation.
Example:
rf,vf = propagate_lagrangian(r0 = [1,0,0], v0 = [0,1,0], tof = pi/2, mu = 1)
pykep.
propagate_taylor
(*args)¶pykep.propagate_taylor(r0 = [1,0,0], v0 = [0,1,0], m0 = 100, thrust = [0,0,0], tof = pi/2, mu = 1, veff = 1, log10tol =-15, log10rtol = -15)
r: start position, x,y,z.
v: start velocity, vx,vy,vz.
m0: starting mass.
thrust: fixed inertial thrust, ux,uy,uz.
tof: propagation time.
mu: central body gravity constant.
veff: the product (Isp g0) defining the engine efficiency.
log10tol: the logarithm of the absolute tolerance passed to taylor propagator.
log10rtol: the logarithm of the relative tolerance passed to taylor propagator.
Returns a tuple (rf, vf, mf) containing the final position, velocity and mass after the propagation.
Example:
r,v,m = propagate_taylor(r0 = [1,0,0], v0 = [0,1,0], m0 = 100, thrust = [0,0,0], tof = pi/2, mu = 1, veff = 1, log10tol =-15, log10rtol = -15)
pykep.
fb_con
(*args)¶pykep.fb_con(vin,vout,pl)
vin: cartesian coordinates of the relative hyperbolic velocity before the fly-by
vout: vout, cartesian coordinates of the relative hyperbolic velocity after the fly-by
pl: fly-by planet
eq represents the violation of the equality constraint norm(vin)**2 = norm(vout)**2. ineq represents the violation of the inequality constraint on the hyperbola asymptote maximum deflection
Example:
v2_eq, delta_ineq = fb_con(vin, vout, planet_ss('earth'))
v2_eq = v2_eq / EARTH_VELOCITY**2
pykep.
fb_prop
(*args)¶pykep.fb_prop(v,v_pla,rp,beta,mu)
v: spacecarft velocity before the encounter (cartesian, absolute)
v-pla: planet inertial velocity at encounter (cartesian, absolute)
rp: fly-by radius
beta: fly-by plane orientation
mu: planet gravitational constant
Returns the cartesian components of the spacecarft velocity after the encounter. Example:
vout = fb_prop([1,0,0],[0,1,0],2,3.1415/2,1)
pykep.
fb_vel
(*args)¶pykep.fb_vel(vin,vout,pl)
vin: cartesian coordinates of the relative hyperbolic velocity before the fly-by
vout: vout, cartesian coordinates of the relative hyperbolic velocity after the fly-by
pl: fly-by planet
dV represents the norm of the delta-V needed to make a given fly-by possible. dV is necessary to fix the magnitude and the direction of vout.
Example:
dV = fb_vel(vin, vout, planet_ss('earth'))
pykep.
barker
(*args)¶pykep.barker(r1,r2,mu = 1.0)
r1: initial position (cartesian)
r2: final position (cartesian)
mu: gravity parameter
Returns the time of flight as evaluated by the Barker equation. Example:
t = barker([1,0,0],[0,1,0],1.0)
pykep.
ic2par
(*args)¶pykep.ic2par(r,v,mu = 1.0)
r: position (cartesian)
v: velocity (cartesian)
mu: gravity parameter
Returns the osculating keplerian elements a,e,i,W,w,E E is the eccentric anomaly for e<1, the Gudermannian for e>1 a is the semi-major axis always a positive quantity. NOTE: The routine is singular when the elements are not defined. Example:
el = ic2par([1,0,0],[0,1,0],1.0)
pykep.
par2ic
(*args)¶pykep.par2ic(E,mu = 1.0)
E: osculating keplerian elements a,e,i,W,w,E ( l, ND, rad, rad, rad, rad)
mu: gravity parameter (l^3/s^2)
Returns cartesian elements from Keplerian elements E is the eccentric anomaly for e<1, the Gudermannian for e>1 a is the semi-major axis always a positive quantity. Example:
r, v = pk.par2ic([1,0.3,0.1,0.1,0.2,0.2], 1)
pykep.
par2eq
(*args)¶pykep.par2eq(E, retrograde = False)
E: osculating keplerian elements a,e,i,W,w,E ( l, ND, rad, rad, rad, rad)
retrogade: uses the retrograde parameters. Default value is False.
Returns the modified equinoctial elements a(1-e^2),f,g,h,k,L L is the true mean longitude Example:
E = pk.par2eq(E = [1.0,0.1,0.2,0.3,0.4,0.5], retrograde = False)
pykep.
eq2par
(*args)¶pykep.eq2par(EQ, retrograde = False)
EQ: modified equinoctial elements a(1-e^2),f,g,h,k,L
retrogade: uses the retrograde parameters. Default value is False.
Returns the osculating Keplerian elements a,e,i,W,w,E E is the eccentric anomaly for e<1, the Gudermannian for e>1 a is the semi-major axis, always a positive quantity. L is the true longitude Example:
E = pk.eq2par(eq = [1,0.3,0.1,0.1,0.2,0.2], retrograde = False)
pykep.
eq2ic
(*args)¶pykep.eq2ic(EQ,mu = 1.0, retrograde = False)
EQ: modified equinoctial elements a(1-e^2),f,g,h,k,L
mu: gravity parameter (l^3/s^2)
Returns cartesian elements from Keplerian elements L is the true longitude Example:
r, v = pk.eq2ic(eq = [1,0.3,0.1,0.1,0.2,0.2], mu = 1, retrograde = False)
pykep.
ic2eq
(*args)¶pykep.ic2eq(r,v,mu = 1.0, retrogade = False)
r: position (cartesian)
v: velocity (cartesian)
mu: gravity parameter
retrogade: uses the retrograde parameters. Default value is False.
Returns the modified equinoctial elements a(1-e^2),f,g,h,k,L L is the true mean longitude Example:
E = ic2eq(r = [1,0,0], v = [0,1,0], mu =1.0)