The utility module

Name

Type

Description

pykep.util.read_satcat()

function

Reads the SATCAT database into a dictionary

pykep.util.read_tle()

function

Reads a TLE file. These can be downloaded from Celestrack or, if you have an account, Space Track

pykep.util.load_spice_kernel()

function

Loads in memory a kernel from the JPL SPICE toolbox (requires BUILD_SPICE option active when building from cmake)

pykep.util.load_gravity_model()

function

Loads a spherical harmonics gravity model

pykep.util.gravity_spherical_harmonic()

function

Calculates the gravitational acceleration due to the provided spherical harmonic gravity model

Detailed Documentation

pykep.util.read_satcat(*args)

pykep.read_satcat(“my_dir/satcat.txt”)

This function reads the satcat catalogue, as can be downloaded from http://www.celestrak.com/NORAD/elements/ and returns a dictionary keyed with the satellite international designator (e.g. “1958-002B”). Each entry is a named tuple with the following fields:

noradn, multnameflag, payloadflag, operationstatus, name, ownership, launchdate, launchsite, decay, period, incl, apogee, perigee, radarA, orbitstatus

see http://celestrak.com/satcat/satcat-format.asp for the meaning of all entries.

Example:

satcat = pykep.read_satcat("my_dir/satcat.txt")
cross_section = satcat["1958-002B"].radarA

pykep.util.read_tle(*args)

planet_list = pykep.read_tle(tle_file, verbose=False, with_name=True)

  • tle_file: A string containin the file name (assumed to be in the working directory)

  • verbose: Activates some screen output to show the progress.

  • with_name: When True the TLE files contain three lines, the first one containing the satellite common name

  • [out] planet_list: a list of pykep tle_planets having as name their international designator (e.g. “1958-002B”).

This function reads a Two-Line-Element file as taken from the NORAD database http://www.celestrak.com/NORAD/elements/ or equivalent and returns a list of pykep planet_tle objects

Note

The name of each of the instantiated planets will be its international designator and thus a valid key to the satcat dictionary

Example:

planet_list = pykep.read_tle("my_dir/cosmos.tle")
satcat = pykep.read_satcat("my_dir/satcat.txt")
cross_sections = [satcat[pl.name()].radarA for pl in  planet_list]

pykep.util.load_spice_kernel(*args)

pykep.util.load_spice_kernel(file_name)

  • file_name: string containing the kernel file to load

Loads the SPICE kernel specified by the filename into memory.

Note

The kernel must be in memory before its used, for example, when computing a planets.spice ephemerides

Example:

util.load_spice_kernel('de432s.bsp')

pykep.util.load_gravity_model(*args)

Load spherical harmonics model.

Args:
  • fp (str): path to model file.

Returns:
  • planet_radius (float): Radius of central body (m)

  • gravity_param (float): Gravitational parameter of central body (m3/s2)

  • c_coefficients (array-like): Array of C_(n, m) coefficients.

  • s_coefficients (array-like): Array of S_(n, m) coefficients.

  • degree (int): Maximum degree of model.

  • order (int): Maximum order of model.

Example:

r, mu, c, s, n, m = pykep.util.load_gravity_model('gravity_models/Earth/egm96.txt')
File format:

The input file must follow a specific format. The first row contains planetary radius in km and gravitational parameter in km3/s2. The second row until the end contains normalised C and S coefficients for every degree and order. The file is comma-delimited. Gravity models can be downloaded from http://pds-geosciences.wustl.edu/dataserv/gravity_models.htm (NASA) and https://earth.esa.int/web/guest/data_access/browse-data-products (ESA).

File:

1| radius,mu,degree,order
2| 0,0,C(0,0),S(0,0)
3| 1,0,C(1,0),S(1,0)
4| ...        

pykep.util.gravity_spherical_harmonic(*args)

Calculate the gravitational acceleration due to the spherical harmonics gravity model supplied.

Args:
  • x (array-like): (N x 3) array of Cartesian coordinates of satellite in frame defined by gravity model.

  • r_planet (float): Equatorial radius of central body.

  • mu (float): Gravitational parameter of central body.

  • c (array-like): Two-dimensional normalised C coefficient array: C[degree, order] = C_(n, m).

  • s (array-like): Two-dimensional normalised S coefficient array: S[degree, order] = S_(n, m).

  • n_max (int): Degree up to which to calculate the gravitational acceleration.

  • m_max (int): Order up to which to calculate the gravitational acceleration. Cannot be higher than n_max.

Returns:
  • acc (array-like): (N x 3) array of the gravitational acceleration.

Example:

r, mu, c, s, n, m = pykep.util.load_gravity_model('gravity_models/Moon/glgm3150.txt')
x = numpy.array([[459.5, 795.8773461, 1591.754692], [-459.5, -795.8773461, -1591.754692]])

acc_zonal = pykep.util.gravity_spherical_harmonic(x, r, mu, c, s, 20, 0)
acc_sqr = pykep.util.gravity_spherical_harmonic(x, r, mu, c, s, 20, 20)
acc_high = pykep.util.gravity_spherical_harmonic(x, r, mu, c, s, 900, 900)

Note

This model was taken from a report by NASA: https://ntrs.nasa.gov/archive/nasa/casi.ntrs.nasa.gov/20160011252.pdf This is the normalised gottlieb algorithm, as coded in MATLAB in the report and transferred to Python.