API Reference¶
Core Functions¶
- sgp4jax.tle_to_satrec(line1, line2, gravity=None)¶
Parse TLE and initialize a SatRec ready for propagation.
- Parameters:
line1 (
str) – First TLE line.line2 (
str) – Second TLE line.gravity (
GravityConstants|None) – Gravity constants. Default WGS72.
- Returns:
Initialized satellite record ready for propagation.
- Return type:
- sgp4jax.tles_to_satrec(tles, gravity=None)¶
Parse an array of TLEs and return a batched SatRec.
- Parameters:
tles (
list[list[str]]) – TLE lines with shape(n_sat, 2). Each row contains[line1, line2]as strings.gravity (
GravityConstants|None) – Gravity constants. Default WGS72.
- Returns:
A single
SatRecwhose fields are stacked arrays with a leading dimension ofn_sat, ready for use withjax.vmap.- Return type:
- sgp4jax.satrec_from_elements(inclo, nodeo, ecco, argpo, mo, no_kozai, *, bstar=0.0, ndot=0.0, nddot=0.0, epoch_jd=2451545.0, gravity=None)¶
Initialize a SatRec directly from orbital elements.
A user-friendly alternative to parsing a TLE string when the orbital elements are already known. Calls the internal
sgp4initroutine to compute all derived SGP4 quantities, so the returned record is ready for use withpropagate(),propagate_jd(), and all GCRF/covariance helpers.- Parameters:
inclo (
float) – Inclination, rad.nodeo (
float) – Right ascension of ascending node, rad.ecco (
float) – Eccentricity (0 ≤ e < 1).argpo (
float) – Argument of perigee, rad.mo (
float) – Mean anomaly at epoch, rad.no_kozai (
float) – Mean motion at epoch, rad/min.bstar (
float) – SGP4 drag term B*, km⁻¹. Default 0.ndot (
float) – First derivative of mean motion, rad/min². Default 0.nddot (
float) – Second derivative of mean motion, rad/min³. Default 0.epoch_jd (
float) – Epoch as a full Julian date (e.g. 2451545.0 = J2000.0). Internally split into whole + fractional parts to preserve float64 precision. Default J2000.0.gravity (
GravityConstants|None) – Gravity model constants. DefaultWGS72.
- Returns:
Fully initialised satellite record ready for propagation.
- Return type:
Examples
>>> import sgp4jax, jax.numpy as jnp >>> sat = sgp4jax.satrec_from_elements( ... inclo=0.9006, ... nodeo=1.2217, ... ecco=0.0004, ... argpo=4.6194, ... mo=3.6160, ... no_kozai=0.0672, ... bstar=2.5e-4, ... epoch_jd=2458924.686, ... ) >>> r, v, err = sgp4jax.propagate(sat, jnp.array(60.0))
- sgp4jax.propagate(satrec, tsince)¶
Propagate satellite to time tsince (minutes from epoch).
- Parameters:
- Return type:
- Returns:
r (jax.Array, shape (3,)) – Position in TEME frame, km.
v (jax.Array, shape (3,)) – Velocity in TEME frame, km/s.
error (jax.Array) – Error code (0 = success).
- sgp4jax.propagate_jd(satrec, jd, fr)¶
Propagate satellite to Julian Date (jd + fr) in TEME.
Converts the Julian Date to minutes-since-epoch and calls
propagate().- Parameters:
- Return type:
- Returns:
r (jax.Array, shape (3,)) – Position in TEME frame, km.
v (jax.Array, shape (3,)) – Velocity in TEME frame, km/s.
error (jax.Array) – Error code (0 = success).
Specialized Propagators¶
These propagators are optimised for specific orbit families and are
significantly faster than propagate() for homogeneous
batches. Each eliminates the dead-branch computation present in the
general propagator:
propagate_leo()— near-earth orbits (method=0), removes the deep-space integrator entirely.propagate_sdp4_nr()— deep-space no-resonance orbits (method=1,irez=0, e.g. GPS/GLONASS/Galileo MEO), replaces the 64-step resonance scan with five scalar multiplications.propagate_mixed()— heterogeneous batches of any orbit type; groups satellites by type internally and dispatches each group to the appropriate specialised propagator.
- sgp4jax.propagate_leo(satrec, tsince)¶
Propagate a near-earth (LEO) satellite to time tsince (minutes from epoch).
This is a stripped-down version of sgp4 for near-earth satellites only. The deep-space resonance integrator (dspace) and lunar-solar periodic perturbations (dpper) are removed entirely, reducing XLA graph size and eliminating the dominant runtime cost for LEO orbits.
Warning: Do not use for deep-space satellites (orbital period > 225 min, mean motion < ~6.4 rev/day). Results will be incorrect for GEO/HEO/MEO objects because the deep-space secular and periodic corrections are skipped.
- Parameters:
- Return type:
- Returns:
r (jax.Array, shape (3,)) – Position in TEME frame, km.
v (jax.Array, shape (3,)) – Velocity in TEME frame, km/s.
error (jax.Array) – Error code (0 = success).
- sgp4jax.propagate_jd_leo(satrec, jd, fr)¶
Propagate a near-earth satellite to Julian Date (jd + fr) in TEME.
LEO/near-earth only variant — deep-space code paths are absent. See
propagate_leo()for the performance trade-offs and limitations.- Parameters:
satrec (
SatRec) – Initialized SatRec fromtle_to_satrec(). Near-earth only.jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date (UTC), integer/whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date (UTC), fractional part (scalar).
- Return type:
- Returns:
r (jax.Array, shape (3,)) – Position in TEME frame, km.
v (jax.Array, shape (3,)) – Velocity in TEME frame, km/s.
error (jax.Array) – Error code (0 = success).
- sgp4jax.propagate_sdp4_nr(satrec, tsince)¶
Propagate a deep-space, no-resonance (irez=0) satellite to time tsince.
This is a stripped-down version of sgp4 for deep-space satellites with irez=0 (no resonance). This captures MEO navigation satellites (GPS, GLONASS, Galileo, BeiDou MEO) and other deep-space objects that do not fall in the 24-hour synchronous (irez=1) or 12-hour Molniya (irez=2) resonance bands.
Simplifications over the full sgp4:
The resonance integrator (dspace / lax.scan) is removed entirely. For irez=0, dspace reduces to five scalar secular drift terms (dedt, didt, domdt, dnodt, dmdt), which are applied directly.
Deep-space satellites always have isimp=1, so the full-drag computation branch is omitted.
All
jnp.where(is_deep, ...)guards collapse to unconditional deep-space assignments.dpper (lunar-solar periodics) is retained unchanged.
Warning
Do not use for near-earth satellites (period < 225 min), irez=1 (GEO/synchronous), or irez=2 (Molniya/12h resonant) satellites — resonance and drag terms will be missing or zero.
- Parameters:
- Return type:
- Returns:
r (jax.Array, shape (3,)) – Position in TEME frame, km.
v (jax.Array, shape (3,)) – Velocity in TEME frame, km/s.
error (jax.Array) – Error code (0 = success).
- sgp4jax.propagate_jd_sdp4_nr(satrec, jd, fr)¶
Propagate a deep-space no-resonance (irez=0) satellite to Julian Date (jd + fr) in TEME.
Deep-space irez=0 only variant — resonance integrator is absent. See
propagate_sdp4_nr()for the performance trade-offs and limitations.- Parameters:
satrec (
SatRec) – Initialized SatRec fromtle_to_satrec(). Deep-space,irez=0only.jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date (UTC), integer/whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date (UTC), fractional part (scalar).
- Return type:
- Returns:
r (jax.Array, shape (3,)) – Position in TEME frame, km.
v (jax.Array, shape (3,)) – Velocity in TEME frame, km/s.
error (jax.Array) – Error code (0 = success).
- sgp4jax.propagate_mixed(satrec_batch, times)¶
Propagate a heterogeneous batch of satellites over a shared array of times.
Groups satellites by type (near-earth / deep-space-no-resonance / deep-space-resonant) and dispatches each group to the appropriate specialized propagator, avoiding dead-branch computation for every group. Results are reassembled into the original satellite ordering.
Note
This function is not JIT-compilable as a whole and does not compose with
jax.gradorjax.vmap. For JIT / AD / vmap compatibility, group satellites by type and call the specialized propagators directly (propagate_leo(),propagate_sdp4_nr(),propagate()).- Parameters:
- Return type:
- Returns:
r (jax.Array, shape (N, M, 3)) – Positions in TEME frame, km.
v (jax.Array, shape (N, M, 3)) – Velocities in TEME frame, km/s.
error (jax.Array, shape (N, M)) – Error codes (0 = success).
Frame Transformations¶
- sgp4jax.teme_to_gcrf(r_teme, v_teme, jd, fr)¶
Transform position/velocity from TEME to GCRF (≈ICRS).
Replicates Skyfield’s
EarthSatellite._at()transformation chain:M = N @ P @ B, thenR = rot_z(theta_GMST1982 - GAST) @ Mwith the final GCRF vectors obtained viaR^T @ r_TEME.The input Julian date
jd + fris treated as UT1 time (≈UTC). TT and TDB are derived internally using an approximate delta_T model.- Parameters:
r_teme (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Position in TEME frame, km.v_teme (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Velocity in TEME frame, km/s.jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, integer/whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, fractional part (scalar).
- Return type:
- Returns:
r_gcrf (jax.Array, shape (3,)) – Position in GCRF frame, km.
v_gcrf (jax.Array, shape (3,)) – Velocity in GCRF frame, km/s.
- sgp4jax.itrf_to_gcrf(r_itrf, jd, fr)¶
Transform position from ITRF (Earth-fixed) to GCRF (≈ICRS).
The rotation is
R = M^T @ rot_z(GAST)whereM = N @ P @ Bis the combined nutation-precession-bias matrix and GAST is the Greenwich Apparent Sidereal Time. Matches Skyfield’sITRS.rotation_at(t).Tto sub-millimetre precision.The input Julian date
jd + fris treated as UT1 (≈UTC).- Parameters:
r_itrf (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Position in ITRF frame, km.jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date (UT1), integer/whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date (UT1), fractional part (scalar).
- Returns:
r_gcrf – Position in GCRF frame, km.
- Return type:
- sgp4jax.gcrf_to_itrf(r_gcrf, jd, fr)¶
Transform position from GCRF (≈ICRS) to ITRF (Earth-fixed).
Inverse of
itrf_to_gcrf(). The rotation isR = rot_z(-GAST) @ MwhereM = N @ P @ B. Matches Skyfield’sITRS.rotation_at(t)to sub-millimetre precision.The input Julian date
jd + fris treated as UT1 (≈UTC).- Parameters:
r_gcrf (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Position in GCRF frame, km.jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date (UT1), integer/whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date (UT1), fractional part (scalar).
- Returns:
r_itrf – Position in ITRF frame, km.
- Return type:
- sgp4jax.propagate_gcrf(satrec, tsince)¶
Propagate satellite and return GCRF position/velocity.
- Parameters:
- Return type:
- Returns:
r_gcrf (jax.Array, shape (3,)) – Position in GCRF frame, km.
v_gcrf (jax.Array, shape (3,)) – Velocity in GCRF frame, km/s.
error (jax.Array) – Error code (0 = success).
- sgp4jax.propagate_jd_gcrf(satrec, jd, fr)¶
Propagate satellite to Julian Date (jd + fr) and return GCRF.
The input time
jd + fris interpreted as UTC, matching the time scale of TLE epochs and Measurement Set timestamps. The same UTC value is used for SGP4 propagation and, as an approximation to UT1, for the TEME → GCRF frame transformation (the UT1 − UTC difference is at most 0.9 s, introducing < 1 m frame error).- Parameters:
- Return type:
- Returns:
r_gcrf (jax.Array, shape (3,)) – Position in GCRF frame, km.
v_gcrf (jax.Array, shape (3,)) – Velocity in GCRF frame, km/s.
error (jax.Array) – Error code (0 = success).
Covariance Transformations¶
Transform covariance matrices between the RIC (Radial, In-track, Cross-track)
frame, TEME Cartesian, and orbital element space. All Jacobians are computed
via jax.jacobian(), so all transforms are fully differentiable.
Two element models are provided:
6-element —
(inclo, nodeo, ecco, argpo, mo, no_kozai)— Keplerian Jacobian, no drag.7-element —
(inclo, nodeo, ecco, argpo, mo, no_kozai, bstar)— SGP4 Jacobian through the fullsgp4init → sgp4pipeline; bstar’s influence on the trajectory is captured.
Empirical Covariance from TLE Age¶
- sgp4jax.tle_ric_covariance(satrec, jd, fr, *, sigma_r0=0.05, sigma_t0=0.3, sigma_n0=0.05)¶
Empirical RIC position-velocity covariance estimate based on TLE age.
Implements a Vallado-style error-growth model. At epoch the dominant uncertainty is in-track (along-track timing). Both in-track and (to a lesser degree) radial errors grow with time as tracking information becomes stale; in-track growth is further scaled by the drag coefficient
bstarrelative to the LEO population median.The velocity block is derived from the position block via the Keplerian approximation
σ_ṽ ≈ n · σ_r, where n is the mean motion. Position-velocity cross-terms are set to zero (diagonal model); they are small relative to the diagonal terms over typical TLE-age timescales and the uncertainty in the growth-rate parameters themselves dominates.Error-growth model (all in km):
σ_R(Δt) = σ_R0 + 0.05 · |Δt| radial σ_T(Δt) = σ_T0 + γ_T · |Δt| in-track σ_N(Δt) = σ_N0 + 0.05 · |Δt| cross-track γ_T = 0.5 + 2.0 · (|bstar| / 3.6×10⁻⁴) km/day
The in-track growth rate
γ_Thas two contributions:Base term (0.5 km/day) — from mean-motion uncertainty, present even at zero drag.
Drag term (2.0 km/day at median LEO
bstar) — from atmospheric-density uncertainty amplified by the drag coefficient. The referencebstarof 3.6×10⁻⁴ km⁻¹ is the empirical median over 13 901 active LEO objects (March 2026).
Typical 1-σ in-track error at median
bstar:Fresh TLE (Δt = 0): 0.3 km
1 day old: 2.8 km
3 days old: 8.1 km
7 days old: 18 km (treat as effectively unusable)
Note
This model is calibrated for well-tracked LEO objects (radar cross-section ≳ 10 cm, altitude 200–2000 km). High area-to-mass objects (debris, balloons) degrade 5–10× faster. When a CDM covariance is available from space-track it should be preferred over this estimate.
- Parameters:
satrec (
SatRec) – Initialized SatRec fromtle_to_satrec().jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Target Julian date, whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Target Julian date, fractional part (scalar).sigma_r0 (
float) – 1-σ radial position uncertainty at TLE epoch, km. Default 0.050 (50 m), reflecting typical OD residuals.sigma_t0 (
float) – 1-σ in-track position uncertainty at TLE epoch, km. Default 0.300.sigma_n0 (
float) – 1-σ cross-track position uncertainty at TLE epoch, km. Default 0.050.
- Returns:
cov_ric – Covariance matrix in the RIC frame, block ordered
[R, T, N, Ṙ, Ṫ, Ṅ]. Units: km² (position block), km²/s² (velocity block). Cross-terms are zero.- Return type:
- sgp4jax.tle_bstar_sigma(satrec, jd, fr, *, bstar_frac0=0.3, bstar_floor=1e-05, bstar_growth_per_day=0.1)¶
Empirical 1-σ uncertainty on
bstarbased on TLE age.bstaris estimated by fitting TLE residuals to tracking data. Its uncertainty has two components:Epoch uncertainty — from OD fit residuals and short-term atmospheric variability (~10 % diurnal cycle). Modelled as a fraction of the fitted |bstar| value with a minimum floor for low-drag objects.
Age growth — atmospheric density varies with solar activity (F10.7 index, 27-day solar-rotation cycle, solar-cycle envelope). A TLE’s fitted bstar does not track these changes, so uncertainty grows at roughly 10 % of |bstar| per day. By 7 days the bstar estimate is effectively unconstrained; by 3 days it has roughly doubled from its epoch value.
Model (units: km⁻¹):
σ_bstar₀ = max(bstar_frac0 · |bstar|, bstar_floor) σ_bstar(Δt) = σ_bstar₀ + bstar_growth_per_day · |bstar| · |Δt|
Typical values at median LEO
bstar(3.6×10⁻⁴ km⁻¹):Fresh TLE (Δt = 0): 1.1×10⁻⁴ km⁻¹ (≈ 30 %)
3 days old: 2.2×10⁻⁴ km⁻¹ (≈ 60 %)
7 days old: 3.6×10⁻⁴ km⁻¹ (≈ 100 %)
This σ is suitable as the standard deviation of a Gaussian or half-normal prior on bstar in Bayesian TLE fitting. For a log-normal prior on |bstar|, use this value divided by |bstar| as the log-scale σ.
- Parameters:
satrec (
SatRec) – Initialized SatRec fromtle_to_satrec().jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Target Julian date, whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Target Julian date, fractional part (scalar).bstar_frac0 (
float) – Fractional uncertainty in bstar at TLE epoch. Default 0.30 (30 %).bstar_floor (
float) – Minimum 1-σ at epoch, km⁻¹. Prevents near-zero bstar from producing an unrealistically tight prior. Default 1×10⁻⁵ km⁻¹.bstar_growth_per_day (
float) – Fractional growth in bstar uncertainty per day, relative to|bstar|. Default 0.10 (10 %/day), reflecting typical atmospheric-density variability from solar activity.
- Returns:
sigma_bstar – Scalar 1-σ uncertainty on
bstar, km⁻¹.- Return type:
Frame utilities¶
- sgp4jax.ric_rotation(r, v)¶
Build the 3×3 rotation matrix from the RIC frame to Cartesian (TEME).
The columns of the returned matrix are the RIC unit vectors expressed in the Cartesian frame:
T[:, 0] = R̂ (radial) T[:, 1] = Î (in-track) T[:, 2] = Ĉ (cross-track)
So a position expressed in RIC coordinates transforms to Cartesian as
x_cart = T @ x_ric.
- sgp4jax.cov_ric_to_teme(cov_ric, r, v)¶
Transform a 6×6 state covariance from the RIC frame to TEME Cartesian.
- Parameters:
cov_ric (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Covariance matrix in RIC frame. Block structure:[[Σ_pos, Σ_pv], [Σ_vp, Σ_vel]]where each block is(3, 3).r (
Array) – Position in TEME at the epoch of the covariance, km.v (
Array) – Velocity in TEME at the epoch of the covariance, km/s.
- Returns:
Covariance in TEME frame.
- Return type:
- sgp4jax.cov_teme_to_ric(cov_teme, r, v)¶
Transform a 6×6 TEME Cartesian covariance to the RIC frame.
6-element model¶
- sgp4jax.elements_jacobian(satrec, jd, fr)¶
Jacobian of the Keplerian state map with respect to orbital elements.
Computes J = ∂(r, v) / ∂(inclo, nodeo, ecco, argpo, mo, no_kozai) evaluated at the given Julian date using the pure two-body Keplerian map.
- Parameters:
- Returns:
J – Jacobian matrix. Rows correspond to
(r_x, r_y, r_z, v_x, v_y, v_z); columns to(inclo, nodeo, ecco, argpo, mo, no_kozai).- Return type:
- sgp4jax.cov_elements_to_teme(cov_elements, satrec, jd, fr)¶
Transform a 6×6 orbital element covariance to TEME Cartesian covariance.
Element ordering:
(inclo, nodeo, ecco, argpo, mo, no_kozai).Uses the linear propagation rule Σ_rv = J Σ_el J^T where J is the Keplerian Jacobian ∂(r,v)/∂(elements) at the given time.
- Parameters:
cov_elements (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Covariance in element space.satrec (
SatRec) – Scalar SatRec fromtle_to_satrec().jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, fractional part (scalar).
- Returns:
Covariance in TEME Cartesian frame.
- Return type:
- sgp4jax.cov_teme_to_elements(cov_teme, satrec, jd, fr)¶
Transform a 6×6 TEME Cartesian covariance to orbital element covariance.
Element ordering:
(inclo, nodeo, ecco, argpo, mo, no_kozai).Inverts the Keplerian Jacobian: Σ_el = J⁻¹ Σ_rv J⁻ᵀ.
- Parameters:
cov_teme (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Covariance in TEME Cartesian frame.satrec (
SatRec) – Scalar SatRec fromtle_to_satrec().jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, fractional part (scalar).
- Returns:
Covariance in element space. Element ordering:
(inclo, nodeo, ecco, argpo, mo, no_kozai).- Return type:
- sgp4jax.cov_ric_to_elements(cov_ric, satrec, jd, fr)¶
Transform a 6×6 RIC covariance to orbital element covariance.
Full chain: RIC → TEME → elements. The Keplerian state at
(jd, fr)is used both to define the RIC frame and to evaluate the Jacobian.Element ordering:
(inclo, nodeo, ecco, argpo, mo, no_kozai).- Parameters:
cov_ric (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Covariance in RIC frame.satrec (
SatRec) – Scalar SatRec fromtle_to_satrec().jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, fractional part (scalar).
- Returns:
Covariance in element space. Element ordering:
(inclo, nodeo, ecco, argpo, mo, no_kozai).- Return type:
- sgp4jax.cov_elements_to_ric(cov_elements, satrec, jd, fr)¶
Transform a 6×6 orbital element covariance to the RIC frame.
Full chain: elements → TEME → RIC. Inverse of
cov_ric_to_elements().Element ordering:
(inclo, nodeo, ecco, argpo, mo, no_kozai).- Parameters:
cov_elements (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Covariance in element space. Element ordering:(inclo, nodeo, ecco, argpo, mo, no_kozai).satrec (
SatRec) – Scalar SatRec fromtle_to_satrec().jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, fractional part (scalar).
- Returns:
Covariance in RIC frame.
- Return type:
7-element model (with drag)¶
- sgp4jax.elements7_jacobian(satrec, jd, fr)¶
Jacobian of the SGP4 state map with respect to the 7-element model.
Computes J = ∂(r, v) / ∂(inclo, nodeo, ecco, argpo, mo, no_kozai, bstar) by differentiating through the full
sgp4init → sgp4pipeline.Unlike the 6-element
elements_jacobian(), this captures bstar’s influence on the trajectory via the drag terms in SGP4.- Parameters:
- Returns:
J – Jacobian matrix. Rows correspond to
(r_x, r_y, r_z, v_x, v_y, v_z); columns to(inclo, nodeo, ecco, argpo, mo, no_kozai, bstar).- Return type:
- sgp4jax.cov_elements7_to_teme(cov_elements7, satrec, jd, fr)¶
Transform a 7×7 element covariance to TEME Cartesian covariance.
Element ordering:
(inclo, nodeo, ecco, argpo, mo, no_kozai, bstar).Uses the linear propagation rule Σ_rv = J Σ_el7 Jᵀ where J is the SGP4 Jacobian ∂(r,v)/∂(elements7), shape
(6, 7).- Parameters:
cov_elements7 (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Covariance in 7-element space.satrec (
SatRec) – Scalar SatRec fromtle_to_satrec().jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, fractional part (scalar).
- Returns:
Covariance in TEME Cartesian frame.
- Return type:
- sgp4jax.cov_teme_to_elements7(cov_teme, satrec, jd, fr)¶
Transform a 6×6 TEME Cartesian covariance to 7-element covariance.
Element ordering:
(inclo, nodeo, ecco, argpo, mo, no_kozai, bstar).Uses the right pseudo-inverse J† = Jᵀ(JJᵀ)⁻¹ to invert the 6×7 Jacobian: Σ_el7 = J† Σ_rv (J†)ᵀ.
Warning
The result is rank-deficient (rank ≤ 6). A 6-dimensional Cartesian state cannot uniquely constrain all 7 element dimensions — bstar is only weakly observable from a single-epoch state snapshot. For a full bstar estimate, propagate and fit over multiple epochs.
- Parameters:
cov_teme (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Covariance in TEME Cartesian frame.satrec (
SatRec) – Scalar SatRec fromtle_to_satrec().jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, fractional part (scalar).
- Returns:
Covariance in 7-element space. Rank ≤ 6.
- Return type:
- sgp4jax.cov_elements7_to_ric(cov_elements7, satrec, jd, fr)¶
Transform a 7×7 element covariance to the RIC frame.
Full chain: elements7 → TEME → RIC. The SGP4 state at
(jd, fr)defines the RIC frame.Element ordering:
(inclo, nodeo, ecco, argpo, mo, no_kozai, bstar).- Parameters:
cov_elements7 (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Covariance in 7-element space.satrec (
SatRec) – Scalar SatRec fromtle_to_satrec().jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, fractional part (scalar).
- Returns:
Covariance in RIC frame.
- Return type:
- sgp4jax.cov_ric_to_elements7(cov_ric, satrec, jd, fr)¶
Transform a 6×6 RIC covariance to 7-element covariance.
Full chain: RIC → TEME → elements7. The SGP4 state at
(jd, fr)defines the RIC frame.Element ordering:
(inclo, nodeo, ecco, argpo, mo, no_kozai, bstar).Warning
The result is rank-deficient (rank ≤ 6). See
cov_teme_to_elements7()for details.- Parameters:
cov_ric (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Covariance in RIC frame.satrec (
SatRec) – Scalar SatRec fromtle_to_satrec().jd (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, whole part (scalar).fr (
Union[Array,ndarray,bool,number,bool,int,float,complex]) – Julian date, fractional part (scalar).
- Returns:
Covariance in 7-element space. Rank ≤ 6.
- Return type:
Keplerian Propagation (no perturbations)¶
Pure two-body Keplerian motion using only the six classical TLE elements (inclination, RAAN, eccentricity, argument of perigee, mean anomaly, mean motion). All perturbations are ignored. Fully differentiable.
Note
TLE elements are Brouwer mean elements. Because SGP4 adds short-period corrections when converting to osculating elements, Keplerian positions will differ from SGP4 by O(1–50 km). This is expected behaviour, not a bug.
- sgp4jax.kepler_gcrf_positions(satrec, times_jd)¶
Propagate one satellite to multiple times using pure Keplerian motion.
No perturbations are applied. Only the six TLE orbital elements (
inclo,nodeo,ecco,argpo,mo,no_kozai) and the epoch (jdsatepoch,jdsatepochF) are used.The function is JIT-compiled and differentiable with respect to all fields of satrec and with respect to times_jd.
- Parameters:
- Return type:
- Returns:
r_gcrf (jax.Array, shape (n_times, 3)) – GCRF positions, km.
v_gcrf (jax.Array, shape (n_times, 3)) – GCRF velocities, km/s.
Notes
The TEME → GCRF rotation uses the input Julian dates as UT1 (UTC is an approximation accurate to < 1 s). For sub-arcsecond frame accuracy call
utc_to_ut1()first and pass the UT1 dates.
- sgp4jax.kepler_gcrf_positions_multi(satrec, times_jd)¶
Propagate multiple satellites to multiple times using pure Keplerian motion.
No perturbations are applied. The outer
vmapiterates over satellites, the inner over times.The function is JIT-compiled and differentiable with respect to all fields of satrec and with respect to times_jd.
- Parameters:
- Return type:
- Returns:
r_gcrf (jax.Array, shape (n_sat, n_times, 3)) – GCRF positions, km.
v_gcrf (jax.Array, shape (n_sat, n_times, 3)) – GCRF velocities, km/s.
Notes
The TEME → GCRF rotation uses the input Julian dates as UT1 (UTC is an approximation accurate to < 1 s). For sub-arcsecond frame accuracy call
utc_to_ut1()first and pass the UT1 dates.
Batch GCRF Convenience Functions¶
These functions propagate one or more satellites to an array of UTC Julian dates and return positions and velocities in the GCRF frame. Choose the variant that matches the orbit type of your satellite batch for best performance:
- sgp4jax.gcrf_positions(satrec, times_jd)¶
Propagate a single satellite to multiple UTC Julian dates.
- Parameters:
- Return type:
- Returns:
r_gcrf (jax.Array, shape (n_times, 3)) – Positions in GCRF frame, km.
v_gcrf (jax.Array, shape (n_times, 3)) – Velocities in GCRF frame, km/s.
- sgp4jax.gcrf_positions_multi(satrec, times_jd)¶
Propagate multiple satellites to multiple UTC Julian dates.
- Parameters:
- Return type:
- Returns:
r_gcrf (jax.Array, shape (n_sat, n_times, 3)) – Positions in GCRF frame, km.
v_gcrf (jax.Array, shape (n_sat, n_times, 3)) – Velocities in GCRF frame, km/s.
- sgp4jax.gcrf_positions_multi_leo(satrec, times_jd)¶
Propagate a near-earth (LEO) satellite batch to multiple UTC Julian dates in GCRF.
Use this for homogeneous batches of near-earth satellites (
method=0). For heterogeneous batches, usegcrf_positions_mixed().- Parameters:
- Return type:
- Returns:
r_gcrf (jax.Array, shape (n_sat, n_times, 3)) – Positions in GCRF frame, km.
v_gcrf (jax.Array, shape (n_sat, n_times, 3)) – Velocities in GCRF frame, km/s.
- sgp4jax.gcrf_positions_multi_sdp4_nr(satrec, times_jd)¶
Propagate a deep-space no-resonance satellite batch to multiple UTC Julian dates in GCRF.
Use this for homogeneous batches of deep-space irez=0 satellites (GPS/GLONASS/Galileo/BeiDou MEO constellations). For heterogeneous batches, use
gcrf_positions_mixed().- Parameters:
- Return type:
- Returns:
r_gcrf (jax.Array, shape (n_sat, n_times, 3)) – Positions in GCRF frame, km.
v_gcrf (jax.Array, shape (n_sat, n_times, 3)) – Velocities in GCRF frame, km/s.
- sgp4jax.gcrf_positions_mixed(satrec_batch, times_jd)¶
Propagate a heterogeneous satellite batch to multiple UTC Julian dates in GCRF.
Groups satellites by type and dispatches each group to its specialized propagator (near-earth / deep-space-no-resonance / deep-space-resonant), then rotates all results to the GCRF frame. Results are reassembled in the original satellite ordering.
Note
Like
propagate_mixed(), this function is not JIT-compilable as a whole.- Parameters:
- Return type:
- Returns:
r_gcrf (jax.Array, shape (N, M, 3)) – Positions in GCRF frame, km.
v_gcrf (jax.Array, shape (N, M, 3)) – Velocities in GCRF frame, km/s.
Earth Orientation (IERS)¶
These functions manage the IERS Bulletin A table used for UTC → UT1 conversion, which is needed for accurate ITRF ↔ GCRF frame transformations.
- sgp4jax.update_iers_table(cache_path=None, url=None)¶
Download the latest IERS finals2000A table and refresh the cache.
Fetches the file from the IERS data centre (falling back to USNO), parses the Bulletin A UT1-UTC values, writes a compact
.npzcache, and immediately loads the result into the module-level interpolation table used byutc_to_ut1().
- sgp4jax.load_iers_table(cache_path=None)¶
Load the IERS table from the local
.npzcache.Called automatically at module import if the default cache exists. Re-call with an explicit cache_path to load from a non-default location or to reload after
update_iers_table().- Parameters:
cache_path (
str|Path|None) – Path to the.npzcache file. Default:~/.cache/sgp4jax/finals2000A.npz.- Raises:
FileNotFoundError – Cache file not found. Call
update_iers_table()to download it.- Return type:
- sgp4jax.utc_to_ut1(jd_utc, fr_utc)¶
Convert a UTC Julian date to UT1 using the IERS finals2000A table.
Interpolates DUT1 = UT1 − UTC (seconds) from the Bulletin A column of the loaded IERS table and adjusts the fractional Julian date. Matches Skyfield’s
t.dut1values to better than 1 µs for all dates within the table’s coverage (~1972 to ~1 year ahead).The function is JAX-traceable and works inside
jax.jitandjax.vmapafter the table has been loaded.- Parameters:
- Return type:
- Returns:
jd_ut1 (jax.Array) – UT1 Julian date, whole part (same value as
jd_utc).fr_ut1 (jax.Array) – UT1 Julian date, fractional part.
- Raises:
RuntimeError – IERS table not loaded — call
update_iers_table()orload_iers_table()first.
Data Structures¶
- class sgp4jax.SatRec(bstar: Array, ecco: Array, argpo: Array, inclo: Array, mo: Array, no_kozai: Array, nodeo: Array, ndot: Array, nddot: Array, j2: Array, j3: Array, j4: Array, j3oj2: Array, xke: Array, mu: Array, radiusearthkm: Array, tumin: Array, jdsatepoch: Array, jdsatepochF: Array, no_unkozai: Array, a: Array, alta: Array, altp: Array, con41: Array, cc1: Array, cc4: Array, cc5: Array, d2: Array, d3: Array, d4: Array, delmo: Array, eta: Array, argpdot: Array, omgcof: Array, sinmao: Array, t2cof: Array, t3cof: Array, t4cof: Array, t5cof: Array, x1mth2: Array, x7thm1: Array, mdot: Array, nodedot: Array, xlcof: Array, xmcof: Array, nodecf: Array, aycof: Array, gsto: Array, d2201: Array, d2211: Array, d3210: Array, d3222: Array, d4410: Array, d4422: Array, d5220: Array, d5232: Array, d5421: Array, d5433: Array, dedt: Array, del1: Array, del2: Array, del3: Array, didt: Array, dmdt: Array, dnodt: Array, domdt: Array, e3: Array, ee2: Array, peo: Array, pgho: Array, pho: Array, pinco: Array, plo: Array, se2: Array, se3: Array, sgh2: Array, sgh3: Array, sgh4: Array, sh2: Array, sh3: Array, si2: Array, si3: Array, sl2: Array, sl3: Array, sl4: Array, xfact: Array, xgh2: Array, xgh3: Array, xgh4: Array, xh2: Array, xh3: Array, xi2: Array, xi3: Array, xl2: Array, xl3: Array, xl4: Array, xlamo: Array, xli: Array, xni: Array, zmol: Array, zmos: Array, atime: Array, method: Array, isimp: Array, irez: Array)¶
JAX-compatible satellite record.
All fields are jnp.ndarray scalars (float64). String flags are encoded as floats:
method: 0.0 = near-earth, 1.0 = deep-space
isimp: 0.0 = full perturbations, 1.0 = simplified
irez: 0.0 = no resonance, 1.0 = synchronous, 2.0 = half-day
Gravity Constants¶
- class sgp4jax._constants.GravityConstants(tumin: float, mu: float, radiusearthkm: float, xke: float, j2: float, j3: float, j4: float, j3oj2: float)¶
Earth gravity model constants for use with SGP4.
- sgp4jax.WGS84 = (13.446851082044981, 398600.5, 6378.137, 0.07436685316871385, 0.00108262998905, -2.53215306e-06, -1.61098761e-06, -0.0023388905587420003)¶
Earth gravity model constants for use with SGP4.
- sgp4jax.WGS72¶
WGS-72 Earth gravity constants (default for SGP4). Field names identical to
GravityConstants.
- sgp4jax.WGS72OLD¶
WGS-72 “old” Earth gravity constants (Vallado original). Field names identical to
GravityConstants.