Name |
Type |
Description |
---|---|---|
class |
The base class for all planets (cannot be instantiated) |
|
class |
A simple planet with Keplerian ephemerides |
|
class |
A solar system planet using jpl low-precision ephemerides |
|
class |
An Earth artificial satellite from its TLE (ephemerides are computed via the SGP4 propagator) |
|
class |
A planet with ephemerides computed using the JPL SPICE Toolbox (requires BUILD_SPICE option active when building from cmake) |
|
class |
A planet from the MPCORB database (keplerian ephemerides) |
|
class |
An asteroid from the GTOC2 competition (keplerian ephemerides) |
|
class |
An asteroid from the GTOC5 competition (keplerian ephemerides) |
|
class |
A Jupiter moon from the GTOC6 competition (keplerian ephemerides) |
|
class |
An asteroid from the GTOC7 competition (keplerian ephemerides) |
pykep.planet.
_base
(*args)¶The base class for all planets, it cannot be instantiated as it contains pure virtual methods
__init__
(*args)¶pykep.planet._base(mu_central body, mu_self, radius, self_radius, name)
mu_central_body: Gravity parameter of the central body (this is not used to compute the ephemerides)
mu_self: Gravity parameter of the target
radius: Radius of target body
self_radius: Safe radius of target body
name: Body name
eph
(*args)¶pykep.planet._base.eph(when)
when: a pykep.epoch
indicating the epoch at which the ephemerides are needed, it can also be a double in which case its interpreted as a mjd2000
Retuns a tuple containing the planet position and velocity in SI units
Note
This is a pure virtual method and must be reimplemented in the derived class
Example:
r,v = earth.eph(epoch(5433), 'mjd2000')
r,v = earth.eph(5433)
osculating_elements
(*args)¶Retruns a tuple containing the six osculating keplerian elements a,e,i,W,w,M at the reference epoch (SI units used). If no epoch is passed 0. is assumed
Example:
elem = earth.osculating_elements()
elem = earth.osculating_elements(epoch(2345.3, 'mjd2000'))
compute_period
(*args)¶The planet orbital period in [sec]
Example:
T = earth.compute_period()
T = earth.compute_period(epoch(2345.3, 'mjd2000'))
human_readable_extra
(*args)¶pykep.planet._base.human_readable_extra()
Retuns a string with extra information on the problem.
Note
This is a virtual method and can be reimplemented in the derived class. In which case __repr__ will append its returned value
mu_central_body
¶The central body gravity parameter in [m^3/s^2]
Example:
mu = earth.mu_central_body
mu_self
¶The body gravity parameter in [m^3/s^2]
Example:
mu_pla = earth.mu_self
radius
¶The planet radius in [m]
Example:
R = earth.radius
safe_radius
¶The body safe radius in [m](distance at which it is safe for spacecraft to fly-by)
Example:
Rs = earth.safe_radius
earth.safe_radius = 1.05
name
¶The body Name
Example:
name = earth.name
pykep.planet.
keplerian
(*args)¶A planet with Keplerian ephemerides, derives from pykep.planet._base
__init__
(*args)¶pykep.planet.keplerian(when,orbital_elements, mu_central_body, mu_self,radius, safe_radius [, name = ‘unknown’])
pykep.planet.keplerian(when,r,v, mu_central_body, mu_self,radius, safe_radius [, name = ‘unknown’])
when: a pykep.epoch
indicating the orbital elements reference epoch
orbital_elements: a sequence of six containing a,e,i,W,w,M (SI units, i.e. meters and radiants)
r,v: position and velocity of an object at when (SI units)
mu_central_body: gravity parameter of the central body (SI units, i.e. m^2/s^3)
mu_self: gravity parameter of the planet (SI units, i.e. m^2/s^3)
radius: body radius (SI units, i.e. meters)
safe_radius: mimimual radius that is safe during a fly-by of the planet (SI units, i.e. m)
name: body name
Note
All classes having Keplerian ephemerides as pykep.planet.mpcorb
inherit from this (c++) class
Example:
earth = planet(epoch(54000,"mjd"),(9.99e-01 * AU, 1.67e-02, 8.85e-04 * DEG2RAD, 1.75e+02 * DEG2RAD, 2.87e+02 * DEG2RAD, 2.57e+02 * DEG2RAD), MU_SUN, 398600e9, 6378000, 6900000, 'Earth')"
pykep.planet.
jpl_lp
(*args)¶A solar system planet that uses the JPL low-precision ephemerides, derives from pykep.planet._base
__init__
(*args)¶pykep.planet.jpl_lp(name)
name: string containing the common planet name (e.g. ‘earth’, ‘venus’ etc.)
Example:
earth = planet.jpl_lp('earth')
pykep.planet.
mpcorb
(*args)¶A planet from the MPCORB database, derives from pykep.planet.keplerian
__init__
(*args)¶pykep.planet.mpcorb(line)
line: a line from the MPCORB database file
Example:
apophis = planet.mpcorb('99942 19.2 0.15 K107N 202.49545 126.41859 204.43202 3.33173 0.1911104 1.11267324 0.9223398 1 MPO164109 1397 2 2004-2008 0.40 M-v 3Eh MPCAPO C802 (99942) Apophis 20080109')
H
¶The asteroid absolute magnitude. This is assuming an albedo of 0.25 and using the formula at www.physics.sfasu.edu/astro/asteroids/sizemagnitude.html Example:
H = apophis.H
n_observations
¶Number of observations made on the asteroid Example:
R = apophis.n_observations
n_oppositions
¶Numper of oppositions the asteroid has been observed in. Example:
R = apophis.n_oppositions
year_of_discovery
¶The year the asteroid was first discovered. In case the asteroid has been observed only once (n_observations), this number is, instead, the Arc Length in days Example:
R = apophis.year_of_discovery
pykep.planet.
tle
(*args)¶An Earth satellite defined from the TLE format, derives from pykep.planet._base
__init__
(*args)¶pykep.planet.tle(line1, line2)
line1: string containing the first line of a TLE (69 well formatted chars)
line2: string containing the second line of a TLE (69 well formatted chars)
Example:
line1 = '1 23177U 94040C 06175.45752052 .00000386 00000-0 76590-3 0 95'
line2 = '2 23177 7.0496 179.8238 7258491 296.0482 8.3061 2.25906668 97438'
arianne = planet.tle(line1, line2)
pykep.planet.
spice
(*args)¶A planet to use the ephemerides from a loaded SPICE kernel
__init__
(*args)¶pykep.planet.spice(target, observer, ref_frame, aberrations, mu_central_body, mu_self, radius, self_radius)
target (string
): Target body (see NAIF IDs)
observer (string
): Observer body (see NAIF IDs
ref_frame (string
): The reference frame (see SPICE supported frames)
aberrations (string
): Aberration correction type (see spkezr_c docs <http://naif.jpl.nasa.gov/pub/naif/toolkit_docs/C/cspice/spkezr_c.html>`_)
mu_central_body (`float`
): Gravity parameter of the central body (this is not necessary to compute the ephemerides).
mu_self (`float`
): Gravity parameter of the target (this is not necessary to compute the ephemerides).
radius (`float`
): Radius of target body (this is not necessary to compute the ephemerides).
self_radius (`float`
): Safe radius of target body (this is not necessary to compute the ephemerides).
Example:
import pykep as pk
pk.load_spice_kernel("de405.bsp")
earth = pk.planet.spice('EARTH', 'SUN', 'ECLIPJ2000', 'NONE', MU_SUN, MU_EARTH, EARTH_R, EARTH_R * 1.05)
r,v = earth.eph(pk.epoch(9500))
Note
The presence in memory of the needed SPICE kernel files are only checked upon call to the ephemerides method.
As a consequence this object can still be constructed with invalid arguments. Only later the ephemerides call
will fail throwing an excpetion. To load SPICE kernels use the pykep.utils.load_spice_kernel()
utility.
Note
mu_central_body must be set to compute the period of the orbital elements.
pykep.planet.
gtoc2
(*args)¶An asteroid from gtoc2, derives from pykep.planet.keplerian
__init__
(*args)¶pykep.planet.gtoc2(ast_id)
ast_id: Construct from a consecutive id from 0 to 910 (Earth).The order is that of the original data file from JPL
Group 1: 0 - 95
Group 2: 96 - 271
Group 3: 272 - 571
Group 4: 572 - 909
Earth: 910
Example:
earth_gtoc2 = planet.gtoc2(910)
pykep.planet.
gtoc5
(*args)¶An asteroid from gtoc5, derives from pykep.planet.keplerian
__init__
(*args)¶pykep.planet.gtoc5(ast_id)
ast_id: a consecutive id from 1 to 7076 (Earth). The order is that of the originaldata file distributed by the Russian, see the (gtoc portal). Earth is 7076
Example:
russian_ast = planet.gtoc5(1)
pykep.planet.
gtoc6
(*args)¶A Jupiter moon from gtoc6, derives from pykep.planet.keplerian
__init__
(*args)¶pykep.planet.gtoc6(name)
name: string containing the common planet name (e.g. ‘io’, ‘europa’, ‘callisto’ or ‘ganymede’)
Example:
io = planet.gtoc6('io')
pykep.planet.
gtoc7
(*args)¶An asteroid from gtoc7, derives from pykep.planet.keplerian
__init__
(*args)¶pykep.planet.gtoc7(ast_id)
ast_id: a consecutive id from 0 (Earth) to 16256. The order is that of the original file , see the (gtoc portal).
Example:
earth = planet.gtoc7(0)