Coordinate Frame Transforms¶
The satkit.frametransform module provides functions for transforming between various coordinate
frames used in satellite tracking and orbit determination. These include multiple variations of "inertial"
coordinate frames, and multiple versions of "Earth-fixed" coordinate frames.
Some notes:
- Most of the algorithms in this module are from the book "Fundamentals of Astrodynamics and Applications" by David Vallado.
- The frame transforms are defined as arbitrary rotations in a 3-dimensional space. The rotations are a function of time, and are represented as quaternions.
- The rotation from the Geocentric Celestial Reference Frame (GCRF) to the Earth-Centered Inertial (ECI) frame is defined by the International Astronomical Union (IAU), available at https://www.iers.org/. See IERS Technical Note 36 for the latest values.
frametransform
¶
Transformations between coordinate frames, and associated utility functions
Coordinate frame transforms are mostly pulled from Vallado: https://www.google.com/books/edition/Fundamentals_of_Astrodynamics_and_Applic/PJLlWzMBKjkC?hl=en&gbpv=0
or the IERS: https://www.iers.org/
gmst(*args, **kwargs)
¶
Greenwich Mean Sidereal Time
Notes
- GMST is the angle between the vernal equinox and the Greenwich meridian
- Vallado algorithm 15
- GMST = 67310.5481 + (876600h + 8640184.812866) * tᵤₜ₁ * (0.983104 + tᵤₜ₁ * −6.2e−6)
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
tm
|
time | datetime
|
scalar time at which to calculate output |
required |
Returns:
| Type | Description |
|---|---|
|
Greenwich Mean Sidereal Time, radians, at input time |
gast(*args, **kwargs)
¶
Greenwich Apparent Sidereal Time
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
tm
|
time
|
scalar, list, or numpy array of astro.time or datetime.datetime representing time at which to calculate output |
required |
Returns:
| Type | Description |
|---|---|
|
Greenwich apparent sidereal time, radians, at input time(s) |
earth_rotation_angle(*args, **kwargs)
¶
Earth Rotation Angle
Notes
- See: IERS Technical Note 36, Chapter 5, Equation 5.15
- Calculation Details:
- Let t be UT1 Julian date
- let f be fractional component of t (fraction of day)
- ERA = 2𝜋 ((0.7790572732640 + f + 0.00273781191135448 * (t - 2451545.0))
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
tm (satkit.time|datetime.datetime
|
Time[s] at which to calculate Earth Rotation Angle |
required |
Returns:
| Type | Description |
|---|---|
|
Earth Rotation Angle at input time[s] in radians |
qitrf2tirs(*args, **kwargs)
¶
Rotation from Terrestrial Intermediate Reference System to Celestial Intermediate Reference Systems
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
tm
|
time | ArrayLike[time] | datetime | ArrayLike[datetime]
|
Time[s] at which to calculate the quaternion |
required |
Returns:
| Type | Description |
|---|---|
|
Quaternion representing rotation from ITRF to TIRS at input time(s) |
qteme2gcrf(*args, **kwargs)
¶
Rotation from True Equator Mean Equinox (TEME) to Geocentric Celestial Reference Frame (GCRF)
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
tm
|
time | datetime
|
Time[s] at which to calculate the quaternion |
required |
Returns:
| Type | Description |
|---|---|
|
Quaternion representing rotation from TEME to GCRF at input time(s) |
qcirs2gcrf(*args, **kwargs)
¶
Rotation from Celestial Intermediate Reference System to Geocentric Celestial Reference Frame
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
tm
|
time | ArrayLike[time] | datetime | ArrayLike[datetime]
|
Time[s] at which to calculate the quaternion |
required |
Returns:
| Type | Description |
|---|---|
|
Quaternion representing rotation from CIRS to GCRF at input time(s) |
qtirs2cirs(*args, **kwargs)
¶
Rotation from Terrestrial Intermediate Reference System (TIRS) to the Celestial Intermediate Reference System (CIRS)
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
tm
|
time | datetime
|
Time[s] at which to calculate the quaternion |
required |
Returns:
| Type | Description |
|---|---|
|
Quaternion representing rotation from TIRS to CIRS at input time(s) |
qgcrf2itrf_approx(*args, **kwargs)
¶
Quaternion representing approximate rotation from the Geocentric Celestial Reference Frame (GCRF) to the International Terrestrial Reference Frame (ITRF)
Notes
- Accurate to approx. 1 arcsec
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
tm
|
time | datetime
|
Time[s] at which to calculate the quaternion |
required |
Returns:
| Type | Description |
|---|---|
|
Quaternion representing rotation from GCRF to ITRF at input time(s) |
qitrf2gcrf_approx(*args, **kwargs)
¶
Quaternion representing approximate rotation from the International Terrestrial Reference Frame (ITRF) to the Geocentric Celestial Reference Frame (GCRF)
Notes
- Accurate to approx. 1 arcsec
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
tm
|
time | datetime
|
Time[s] at which to calculate the quaternion |
required |
Returns:
| Type | Description |
|---|---|
|
Quaternion representing rotation from ITRF to GCRF at input time(s) |
qgcrf2itrf(*args, **kwargs)
¶
Quaternion representing rotation from the Geocentric Celestial Reference Frame (GCRF) to the International Terrestrial Reference Frame (ITRF)
Notes
- Uses full IAU2010 Reduction
- See IERS Technical Note 36, Chapter 5
- Does not include solid tides, ocean tides
- Very computationally expensive
- Velocity transforms: this quaternion rotates position vectors
between ITRF and GCRF but is not sufficient for velocity on
its own. ITRF is a rotating frame, so the velocity transform
picks up an extra
omega_earth x rterm (~470 m/s at LEO). Use :func:itrf_to_gcrf_state/ :func:gcrf_to_itrf_statefor full state (position + velocity) transforms.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
tm
|
time | datetime
|
Time[s] at which to calculate the quaternion |
required |
Returns:
| Type | Description |
|---|---|
|
Quaternion representing rotation from GCRF to ITRF at input time(s) |
qitrf2gcrf(*args, **kwargs)
¶
Quaternion representing rotation from the International Terrestrial Reference Frame (ITRF) to the Geocentric Celestial Reference Frame (GCRF)
Notes
- Uses full IAU2010 Reduction
- See IERS Technical Note 36, Chapter 5
- Does not include solid tides, ocean tides
- Very computationally expensive
- Velocity transforms: this quaternion rotates position vectors
between ITRF and GCRF but is not sufficient for velocity on
its own. ITRF is a rotating frame, so the velocity transform
picks up an extra
omega_earth x rterm (~470 m/s at LEO). Use :func:itrf_to_gcrf_state/ :func:gcrf_to_itrf_statefor full state (position + velocity) transforms.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
tm
|
satkit.time datetime.datetime
|
Time[s] at which to calculate the quaternion |
required |
Returns: Quaternion representing rotation from ITRF to GCRF at input time(s)
qteme2itrf(*args, **kwargs)
¶
Quaternion representing rotation from the True Equator Mean Equinox (TEME) frame to the International Terrestrial Reference Frame (ITRF)
Notes
- This is equation 3-90 in Vallado
- TEME is the output frame of the SGP4 propagator used to compute position from two-line element sets.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
tm
|
time | datetime
|
Time[s] at which to calculate the quaternion |
required |
Returns:
| Type | Description |
|---|---|
|
Quaternion representing rotation from TEME to ITRF at input time(s) |
earth_orientation_params(time)
¶
Get Earth Orientation Parameters at given instant
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
time
|
time
|
Instant at which to query parameters |
required |
Returns:
| Type | Description |
|---|---|
tuple[float, float, float, float, float, float]
|
(float, float, float, float, float, float) | None: Tuple with following elements: 0 : (UT1 - UTC) in seconds 1 : X polar motion in arcsecs 2 : Y polar motion in arcsecs 3 : LOD: instantaneous rate of change in (UT1-UTC), msec/day 4 : dX wrt IAU-2000A nutation, milli-arcsecs 5 : dY wrt IAU-2000A nutation, milli-arcsecs |
Notes
- Returns None if the time is before the range of available EOP data
- For times after the last available EOP data, the last entry's values are returned (constant extrapolation)
- EOP data is available from 1962 to current, with predictions ~4 months ahead
- See: https://www.iers.org/IERS/EN/DataProducts/EarthOrientationData/eop.html
qmod2gcrf(*args, **kwargs)
¶
Quaternion rotating Mean-of-Date (MOD) → GCRF at the given time(s).
Mean-of-Date accounts for precession but not nutation.
qtod2mod_approx(*args, **kwargs)
¶
Approximate True-of-Date (TOD) → Mean-of-Date (MOD) rotation.
Accounts for nutation only.
to_gcrf(frame, pos, vel)
¶
Return the 3x3 DCM that transforms a vector from a satellite-local orbital frame into GCRF at the current state.
This is the unified dispatch for satellite-local orbital frames. Supported values:
frame.GCRF— returns the 3x3 identity matrix (trivial case)frame.LVLH— Local Vertical / Local Horizontalframe.RTN— Radial / In-track / Cross-track (= RSW = RTN)frame.NTW— Normal-to-velocity / Tangent / Cross-track
For an arbitrary frame-to-frame rotation, compose with
:func:from_gcrf::
# NTW -> RIC
dcm = sk.frametransform.from_gcrf(sk.frame.RTN, pos, vel) @ \
sk.frametransform.to_gcrf(sk.frame.NTW, pos, vel)
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
frame
|
frame
|
Source satellite-local frame |
required |
pos
|
ArrayLike
|
3-element position vector in GCRF [m] |
required |
vel
|
ArrayLike
|
3-element velocity vector in GCRF [m/s] |
required |
Returns:
| Type | Description |
|---|---|
NDArray[float64]
|
numpy.ndarray: 3x3 rotation matrix (frame → GCRF) |
Raises:
| Type | Description |
|---|---|
RuntimeError
|
if |
itrf_to_gcrf_state(pos_itrf, vel_itrf, time)
¶
Transform a satellite state (position + velocity) from ITRF to GCRF.
Unlike the raw :func:qitrf2gcrf quaternion, this function correctly
handles the Earth-rotation contribution to velocity. A point at rest
on Earth's surface has zero velocity in ITRF but ~465 m/s in GCRF at
the equator, and this function accounts for that term.
The IAU 2010 ITRF → GCRF reduction decomposes into three stages:
polar motion (ITRF → TIRS), Earth rotation about the CIO polar axis
(TIRS → CIRS), and precession-nutation (CIRS → GCRF). The
Earth-rotation sweep term omega_earth x r is computed in
TIRS — not ITRF or GCRF — because TIRS is defined such that
Earth's rotation axis is exactly along its +z axis. Computing the
sweep anywhere else would introduce either a polar-motion-sized
error (~0.3 arcsec in ITRF) or a precession-sized error (tens of
degrees in GCRF).
Implementation:
- Rotate
pos_itrfandvel_itrfinto TIRS via polar motion. - Add
omega_earth x r_tirsto the velocity in TIRS, whereomega_earth = (0, 0, OMEGA_EARTH)exactly. - Rotate TIRS → CIRS → GCRF via the full IAU 2010 chain.
Uses the full IAU 2010 reduction (polar motion + Earth rotation + precession-nutation with dX/dY corrections from Earth orientation parameters).
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
pos_itrf
|
ArrayLike
|
3-element position vector in ITRF [m] |
required |
vel_itrf
|
ArrayLike
|
3-element velocity vector as observed in ITRF [m/s] (zero for a point at rest on Earth's surface) |
required |
time
|
time
|
Epoch of the state |
required |
Returns:
| Type | Description |
|---|---|
NDArray[float64]
|
A 2-tuple |
NDArray[float64]
|
state expressed in GCRF. |
Example
import satkit as sk
import numpy as np
# Geostationary satellite, stationary in ITRF
t = sk.time(2024, 1, 1)
pos_itrf = np.array([42164.17e3, 0.0, 0.0])
vel_itrf = np.array([0.0, 0.0, 0.0])
pos_gcrf, vel_gcrf = sk.frametransform.itrf_to_gcrf_state(
pos_itrf, vel_itrf, t)
# |vel_gcrf| ≈ 3075 m/s (the GEO orbital speed)
gcrf_to_itrf_state(pos_gcrf, vel_gcrf, time)
¶
Transform a satellite state (position + velocity) from GCRF to ITRF.
Inverse of :func:itrf_to_gcrf_state. Rotates the state through
GCRF → CIRS → TIRS, subtracts the Earth-rotation omega_earth x r
term in TIRS (where Earth's rotation axis is exactly along +z),
then applies inverse polar motion to reach ITRF. A geostationary
satellite (whose GCRF velocity is pure orbital motion) produces
zero velocity in ITRF. Uses the full IAU 2010 reduction.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
pos_gcrf
|
ArrayLike
|
3-element position vector in GCRF [m] |
required |
vel_gcrf
|
ArrayLike
|
3-element velocity vector in GCRF [m/s] |
required |
time
|
time
|
Epoch of the state |
required |
Returns:
| Type | Description |
|---|---|
NDArray[float64]
|
A 2-tuple |
NDArray[float64]
|
velocity as observed in ITRF. |
from_gcrf(frame, pos, vel)
¶
Return the 3x3 DCM that transforms a vector from GCRF into a satellite-local orbital frame at the current state.
Transpose of :func:to_gcrf. See that function for the list of
supported frames, composition examples, and error conditions.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
frame
|
frame
|
Destination satellite-local frame |
required |
pos
|
ArrayLike
|
3-element position vector in GCRF [m] |
required |
vel
|
ArrayLike
|
3-element velocity vector in GCRF [m/s] |
required |
Returns:
| Type | Description |
|---|---|
NDArray[float64]
|
numpy.ndarray: 3x3 rotation matrix (GCRF → frame) |
Raises:
| Type | Description |
|---|---|
RuntimeError
|
if |
disable_eop_time_warning()
¶
Disable the warning printed to stderr when Earth Orientation Parameters (EOP) are not available for a given time.
Notes
- This function is used to disable the warning printed when EOP are not available for a given time.
- If not disabled, warning will be shown only once per library load,