Skip to content

High-Precision Orbit Propagator

propagate(state, begin, end=..., *, duration=..., duration_secs=..., duration_days=..., output_phi=False, propsettings=..., satproperties=...)

High-precision orbit propagator

Propagate orbits with high-precision force modeling via adaptive Runge-Kutta methods (default is order 9/8).

Parameters:

Name Type Description Default
state NDArray[float64]

6-element numpy array representing satellite GCRF position and velocity in meters and meters/second

required
begin time

Time at which satellite is at input state

required
end time

Time at which new position and velocity will be computed

...
duration duration

Duration from "begin" at which new position & velocity will be computed

...
duration_secs float

Duration in seconds from "begin" at which new position and velocity will be computed

...
duration_days float

Duration in days from "begin" at which new position and velocity will be computed

...
output_phi bool

Output 6x6 state transition matrix between begin and end times

False
propsettings propsettings

Settings for the propagation; if omitted, defaults are used

...
satproperties satproperties_static

Drag and radiation pressure susceptibility of satellite

...

Returns:

Type Description
propresult

Propagation result object holding state outputs, statistics, and dense output if requested

Notes

Propagates satellite ephemeris (position, velocity in GCRF & time) to new time and outputs new position and velocity via Runge-Kutta integration. Inputs and outputs are all in the Geocentric Celestial Reference Frame (GCRF).

Included forces:

  • Earth gravity with higher-order zonal terms
  • Sun, Moon gravity
  • Radiation pressure
  • Atmospheric drag: NRL-MSISE 2000 density model, with option to include space weather effects

End time must be set by keyword argument, either explicitly or by duration. Solid Earth tides are not (yet) included in the model.

Example
import numpy as np

# Define initial state in GCRF (position in meters, velocity in m/s)
state = np.array([6.781e6, 0, 0, 0, 7.5e3, 0])
t0 = satkit.time(2024, 1, 1)

# Propagate forward by 1 day
result = satkit.propagate(state, t0, duration_days=1.0)
print(f"End position: {result.pos} m")
print(f"End velocity: {result.vel} m/s")

# Interpolate at intermediate time
t_mid = t0 + satkit.duration(hours=12)
mid_state = result.interp(t_mid)

propresult

Results of a satellite propagation

This class lets the user access results of the satellite propagation

Notes:

  • If "enable_interp" is set to True in the propagation settings, the propresult object can be used to interpolate solutions at any time between the begin and end times of the propagation via the "interp" method

pos property

GCRF position of satellite, meters

Returns:

Type Description
NDArray[float64]

npt.ArrayLike[float]: 3-element numpy array representing GCRF position (meters) at end of propagation

vel property

GCRF velocity of satellite, meters/second

Returns:

Type Description
NDArray[float64]

npt.ArrayLike[float]: 3-element numpy array representing GCRF velocity in meters/second at end of propagation

state property

6-element end state (pos + vel) of satellite in meters & meters/second

Returns:

Type Description
NDArray[float64]

npt.ArrayLike[float]: 6-element numpy array representing state of satellite in meters & meters/second

state_end property

6-element state (pos + vel) of satellite in meters & meters/second at end of propagation

Notes: * This is the same as the "state" property

Returns:

Type Description
NDArray[float64]

npt.ArrayLike[float]: 6-element numpy array representing state of satellite in meters & meters/second

state_begin property

6-element state (pos + vel) of satellite in meters & meters/second at begin of propagation Returns: npt.NDArray[np.float64]: 6-element numpy array representing state of satellite in meters & meters/second at begin of propagation

time property

Time at which state is valid

Returns:

Type Description
time

satkit.time: Time at which state is valid

time_end property

Time at which state is valid

Notes: * This is identical to "time" property

Returns:

Type Description
time

satkit.time: Time at which state is valid

time_begin property

Time at which state_begin is valid

Returns:

Type Description
time

satkit.time: Time at which state_begin is valid

stats property

Statistics of propagation

Returns:

Name Type Description
propstats propstats

Object containing statistics of propagation

phi property

State transition matrix

Returns:

Type Description
NDArray[float64] | None

npt.ArrayLike[np.float64] | None: 6x6 numpy array representing state transition matrix or None if not computed

interp(time, output_phi=False)

Interpolate state at given time

Parameters:

Name Type Description Default
time time

Time at which to interpolate state

required
output_phi bool

Output 6x6 state transition matrix at the interpolated time Default is False

False

Returns:

Type Description
NDArray[float64] | tuple[NDArray[float64], NDArray[float64]]

npt.ArrayLike[np.float64] | tuple[npt.ArrayLike[np.float64], npt.ArrayLike[np.float64]]: 6-element vector representing state at given time. if output_phi, also output 6x6 state transition matrix at given time

Example
# After propagation with enable_interp=True
result = satkit.propagate(state, t0, duration_days=1.0)
t_mid = t0 + satkit.duration(hours=12)
mid_state = result.interp(t_mid)
print(f"Position at 12h: {mid_state[0:3]} m")

propsettings

This class contains settings used in the high-precision orbit propagator part of the "satkit" python toolbox

Notes:

  • Default settings:

    • abs_error: 1e-8
    • rel_error: 1e-8
    • gravity_order: 4
    • use_spaceweather: True
    • enable_interp: True

    • enable_interp enables high-precision interpolation of state between begin and end times via the returned function, it is enabled by default. There is a small increase in computational efficiency if set to false

abs_error property writable

Maximum absolute value of error for any element in propagated state following ODE integration

Returns:

Name Type Description
float float

Maximum absolute value of error for any element in propagated state following ODE integration, default is 1e-8

rel_error property writable

Maximum relative error of any element in propagated state following ODE integration

Returns:

Name Type Description
float float

Maximum relative error of any element in propagated state following ODE integration, default is 1e-8

gravity_order property writable

Earth gravity order to use in ODE integration

Returns:

Name Type Description
int int

Earth gravity order to use in ODE integration, default is 4

use_spaceweather property writable

Use space weather data when computing atmospheric density for drag forces

Notes:

  • Space weather data can have a large effect on the density of the atmosphere
  • This can be important for accurate drag force calculations
  • Space weather data is updated every 3 hours. Most-recent data can be downloaded with satkit.utils.update_datafiles()
  • Default value is True

Returns:

Name Type Description
bool bool

Indicate whether or not space weather data should be used when computing atmospheric density for drag forces

enable_interp property writable

Store intermediate data that allows for fast high-precision interpolation of state between begin and end times If not needed, there is a small computational advantage if set to False

__init__(*, abs_error=1e-08, rel_error=1e-08, gravity_order=4, use_spaceweather=True, enable_interp=True)

Create propagation settings object used to configure high-precision orbit propagator

Parameters:

Name Type Description Default
abs_error float

Maximum absolute value of error for any element in propagated state following ODE integration. Default is 1e-8

1e-08
rel_error float

Maximum relative error of any element in propagated state following ODE integration. Default is 1e-8

1e-08
gravity_order int

Earth gravity order to use in ODE integration. Default is 4

4
use_spaceweather bool

Use space weather data when computing atmospheric density for drag forces. Default is True

True
enable_interp bool

Store intermediate data that allows for fast high-precision interpolation of state between begin and end times. Default is True

True

Returns:

Type Description
None

New propsettings object with default settings

Example
settings = satkit.propsettings(
    gravity_order=16,
    abs_error=1e-10,
    rel_error=1e-10,
)

precompute_terms(begin, end, step=None)

Precompute terms for fast interpolation of state between begin and end times

This can be used, for example, to compute sun and moon positions only once if propagating many satellites over the same time period

Parameters:

Name Type Description Default
begin time

Begin time of propagation

required
end time

End time of propagation

required
step duration

Step size for interpolation. Default = 60 seconds

None

propstats

Statistics of a satellite propagation

num_eval property

Number of function evaluations

num_accept property

Number of accepted steps in adaptive RK integrator

num_reject property

Number of rejected steps in adaptive RK integrator

satproperties_static

Satellite properties relevant for drag and radiation pressure

This class lets the satellite radiation pressure and drag parameters be set to static values for duration of propagation

Attributes:

Name Type Description
cdaoverm float

Coefficient of drag times area over mass in m^2/kg

craoverm float

Coefficient of radiation pressure times area over mass in m^2/kg

cdaoverm property writable

Coefficient of drag times area over mass. Units are m^2/kg

craoverm property writable

Coefficient of radiation pressure times area over mass. Units are m^2/kg

__init__(cdaoverm=0, craoverm=0)

Create a satproperties_static object with given craoverm and cdaoverm in m^2/kg

Parameters:

Name Type Description Default
cdaoverm float

Coefficient of drag times area over mass in m^2/kg

0
craoverm float

Coefficient of radiation pressure times area over mass in m^2/kg

0

Notes:

  • The two arguments can be passed as positional arguments or as keyword arguments

Example:

properties = satproperties_static(craoverm = 0.5, cdaoverm = 0.4)

# or with same output

properties = satproperties_static(0.5, 0.4)