Skip to main content

sgp4/
lib.rs

1//! This crate implements the SGP4 algorithm for satellite propagation
2//!
3//! It also provides methods to parse Two-Line Element Set (TLE) and Orbit Mean-Elements Message (OMM) data.
4//!
5//! A UTF-8 transcript of the mathematical expressions that compose SGP4
6//! can be found in the repository [https://github.com/neuromorphicsystems/sgp4](https://github.com/neuromorphicsystems/sgp4).
7//!
8//! # Example
9//!
10//! The following standalone program downloads the lastest Galileo OMMs from Celestrak
11//! and predicts the satellites' positions and velocities after 12 h and 24 h.
12//!
13//! ```
14//! fn main() -> anyhow::Result<()> {
15//!     let mut response = ureq::get("https://celestrak.com/NORAD/elements/gp.php")
16//!         .query("GROUP", "galileo")
17//!         .query("FORMAT", "json")
18//!         .call()?;
19//!     let elements_vec: Vec<sgp4::Elements> = response.body_mut().read_json()?;
20//!     for elements in &elements_vec {
21//!         println!("{}", elements.object_name.as_ref().unwrap());
22//!         let constants = sgp4::Constants::from_elements(elements)?;
23//!         for hours in &[12, 24] {
24//!             println!("    t = {} min", hours * 60);
25//!             let prediction = constants.propagate(sgp4::MinutesSinceEpoch((hours * 60) as f64))?;
26//!             println!("        r = {:?} km", prediction.position);
27//!             println!("        ṙ = {:?} km.s⁻¹", prediction.velocity);
28//!         }
29//!     }
30//!     Ok(())
31//! }
32//! ```
33//! More examples can be found in the repository [https://github.com/neuromorphicsystems/sgp4/tree/master/examples](https://github.com/neuromorphicsystems/sgp4/tree/master/examples).
34//!
35
36#![cfg_attr(docsrs, feature(doc_cfg))]
37#![cfg_attr(not(feature = "std"), no_std)]
38
39#[cfg(not(any(feature = "std", feature = "libm")))]
40compile_error!("either feature \"std\" or feature \"libm\" must be enabled");
41
42#[cfg(all(feature = "std", feature = "libm"))]
43compile_error!("feature \"std\" and feature \"libm\" cannot be enabled at the same time");
44
45#[cfg(feature = "alloc")]
46extern crate alloc;
47
48#[cfg(not(feature = "std"))]
49use num_traits::Float;
50
51mod deep_space;
52mod gp;
53mod model;
54mod near_earth;
55mod propagator;
56mod third_body;
57mod tle;
58pub use chrono;
59pub use deep_space::ResonanceState;
60pub use gp::Error;
61pub use model::afspc_epoch_to_sidereal_time;
62pub use model::iau_epoch_to_sidereal_time;
63pub use model::Geopotential;
64pub use model::WGS72;
65pub use model::WGS84;
66pub use propagator::Constants;
67pub use propagator::Orbit;
68pub use propagator::Prediction;
69pub use tle::julian_years_since_j2000;
70pub use tle::julian_years_since_j2000_afspc_compatibility_mode;
71pub use tle::Classification;
72pub use tle::DatetimeToMinutesSinceEpochError;
73pub use tle::Elements;
74pub use tle::Error as TleError;
75pub use tle::ErrorLine as TleErrorLine;
76pub use tle::ErrorWhat as TleErrorWhat;
77pub use tle::MinutesSinceEpoch;
78pub use tle::MinutesSinceEpochToDatetimeError;
79
80#[cfg(feature = "alloc")]
81#[cfg_attr(docsrs, doc(cfg(feature = "alloc")))]
82pub use tle::parse_2les;
83
84#[cfg(feature = "alloc")]
85#[cfg_attr(docsrs, doc(cfg(feature = "alloc")))]
86pub use tle::parse_3les;
87
88/// Represents a propagation error caused by orbital elements divergence
89#[derive(Debug, Clone, PartialEq)]
90pub enum KozaiElementsError {
91    NegativeKozaiMeanMotion,
92    NegativeBrouwerMeanMotion,
93}
94
95impl core::fmt::Display for KozaiElementsError {
96    fn fmt(&self, formatter: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
97        match self {
98            KozaiElementsError::NegativeKozaiMeanMotion => formatter
99                .write_str("The Kozai mean motion calculated from epoch elements is negative"),
100            KozaiElementsError::NegativeBrouwerMeanMotion => formatter
101                .write_str("The Brouwer mean motion calculated from epoch elements is negative"),
102        }
103    }
104}
105
106#[cfg(feature = "std")]
107impl std::error::Error for KozaiElementsError {}
108
109/// The orbit used to generate epoch constants has an invalid eccentricity
110#[derive(Debug, Clone, PartialEq)]
111pub struct OutOfRangeEpochEccentricity(pub f64);
112
113impl core::fmt::Display for OutOfRangeEpochEccentricity {
114    fn fmt(&self, formatter: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
115        formatter.write_fmt(core::format_args!(
116            "The epoch eccentricity ({}) is outside the range [0, 1[",
117            self.0
118        ))
119    }
120}
121
122#[cfg(feature = "std")]
123impl std::error::Error for OutOfRangeEpochEccentricity {}
124
125/// Errors returned when creating epoch contants from elements
126#[derive(Debug, Clone, PartialEq)]
127pub enum ElementsError {
128    KozaiElementsError(KozaiElementsError),
129    OutOfRangeEpochEccentricity(OutOfRangeEpochEccentricity),
130}
131
132impl core::fmt::Display for ElementsError {
133    fn fmt(&self, formatter: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
134        match self {
135            ElementsError::KozaiElementsError(error) => error.fmt(formatter),
136            ElementsError::OutOfRangeEpochEccentricity(error) => error.fmt(formatter),
137        }
138    }
139}
140
141impl From<KozaiElementsError> for ElementsError {
142    fn from(value: KozaiElementsError) -> Self {
143        Self::KozaiElementsError(value)
144    }
145}
146
147impl From<OutOfRangeEpochEccentricity> for ElementsError {
148    fn from(value: OutOfRangeEpochEccentricity) -> Self {
149        Self::OutOfRangeEpochEccentricity(value)
150    }
151}
152
153#[cfg(feature = "std")]
154impl std::error::Error for ElementsError {}
155
156impl Orbit {
157    /// Creates a new Brouwer orbit representation from Kozai elements
158    ///
159    /// If the Kozai orbital elements are obtained from a TLE or OMM,
160    /// the convenience function [sgp4::Constants::from_elements](struct.Constants.html#method.from_elements)
161    /// can be used instead of manually mapping the `Elements` fields to the `Constants::new` parameters.
162    ///
163    /// # Arguments
164    ///
165    /// * `geopotential` - The model of Earth gravity to use in the conversion
166    /// * `inclination` - Angle between the equator and the orbit plane in rad
167    /// * `right_ascension` - Angle between vernal equinox and the point where the orbit crosses the equatorial plane in rad
168    /// * `eccentricity` - The shape of the orbit
169    /// * `argument_of_perigee` - Angle between the ascending node and the orbit's point of closest approach to the earth in rad
170    /// * `mean_anomaly` - Angle of the satellite location measured from perigee in rad
171    /// * `kozai_mean_motion` - Mean orbital angular velocity in rad.min⁻¹ (Kozai convention)
172    ///
173    /// # Example
174    ///
175    /// ```
176    /// # fn main() -> anyhow::Result<()> {
177    /// let elements = sgp4::Elements::from_tle(
178    ///     Some("ISS (ZARYA)".to_owned()),
179    ///     "1 25544U 98067A   20194.88612269 -.00002218  00000-0 -31515-4 0  9992".as_bytes(),
180    ///     "2 25544  51.6461 221.2784 0001413  89.1723 280.4612 15.49507896236008".as_bytes(),
181    /// )?;
182    /// let orbit_0 = sgp4::Orbit::from_kozai_elements(
183    ///     &sgp4::WGS84,
184    ///     elements.inclination * (core::f64::consts::PI / 180.0),
185    ///     elements.right_ascension * (core::f64::consts::PI / 180.0),
186    ///     elements.eccentricity,
187    ///     elements.argument_of_perigee * (core::f64::consts::PI / 180.0),
188    ///     elements.mean_anomaly * (core::f64::consts::PI / 180.0),
189    ///     elements.mean_motion * (core::f64::consts::PI / 720.0),
190    /// )?;
191    /// #     Ok(())
192    /// # }
193    /// ```
194    pub fn from_kozai_elements(
195        geopotential: &Geopotential,
196        inclination: f64,
197        right_ascension: f64,
198        eccentricity: f64,
199        argument_of_perigee: f64,
200        mean_anomaly: f64,
201        kozai_mean_motion: f64,
202    ) -> core::result::Result<Self, KozaiElementsError> {
203        if kozai_mean_motion <= 0.0 {
204            Err(KozaiElementsError::NegativeKozaiMeanMotion)
205        } else {
206            let mean_motion = {
207                // a₁ = (kₑ / n₀)²ᐟ³
208                let a1 = (geopotential.ke / kozai_mean_motion).powf(2.0 / 3.0);
209
210                //      3      3 cos²I₀ - 1
211                // p₀ = - J₂ ---------------
212                //      4       (1 − e₀²)³ᐟ²
213                let p0 = 0.75 * geopotential.j2 * (3.0 * inclination.cos().powi(2) - 1.0)
214                    / (1.0 - eccentricity.powi(2)).powf(3.0 / 2.0);
215
216                // 𝛿₁ = p₀ / a₁²
217                let d1 = p0 / a1.powi(2);
218
219                // 𝛿₀ = p₀ / (a₁ (1 - ¹/₃ 𝛿₁ - 𝛿₁² - ¹³⁴/₈₁ 𝛿₁³))²
220                let d0 = p0
221                    / (a1 * (1.0 - d1.powi(2) - d1 * (1.0 / 3.0 + 134.0 * d1.powi(2) / 81.0)))
222                        .powi(2);
223
224                //         n₀
225                // n₀" = ------
226                //       1 + 𝛿₀
227                kozai_mean_motion / (1.0 + d0)
228            };
229            if mean_motion <= 0.0 {
230                Err(KozaiElementsError::NegativeBrouwerMeanMotion)
231            } else {
232                Ok(propagator::Orbit {
233                    inclination,
234                    right_ascension,
235                    eccentricity,
236                    argument_of_perigee,
237                    mean_anomaly,
238                    mean_motion,
239                })
240            }
241        }
242    }
243}
244
245impl Constants {
246    /// Initializes a new propagator from epoch quantities
247    ///
248    /// If the orbital elements are obtained from a TLE or OMM,
249    /// the convenience function [sgp4::Constants::from_elements](struct.Constants.html#method.from_elements)
250    /// can be used instead of manually mapping the `Elements` fields to the `Constants::new` parameters.
251    ///
252    /// # Arguments
253    ///
254    /// * `geopotential` - The model of Earth gravity to use in the conversion
255    /// * `epoch_to_sidereal_time` - The function to use to convert the J2000 epoch to sidereal time
256    /// * `epoch` - The number of years since UTC 1 January 2000 12h00 (J2000)
257    /// * `drag_term` - The radiation pressure coefficient in earth radii⁻¹ (B*)
258    /// * `orbit_0` - The Brouwer orbital elements at epoch
259    ///
260    /// # Example
261    ///
262    /// ```
263    /// # fn main() -> anyhow::Result<()> {
264    /// let elements = sgp4::Elements::from_tle(
265    ///     Some("ISS (ZARYA)".to_owned()),
266    ///     "1 25544U 98067A   20194.88612269 -.00002218  00000-0 -31515-4 0  9992".as_bytes(),
267    ///     "2 25544  51.6461 221.2784 0001413  89.1723 280.4612 15.49507896236008".as_bytes(),
268    /// )?;
269    /// let constants = sgp4::Constants::new(
270    ///     sgp4::WGS84,
271    ///     sgp4::iau_epoch_to_sidereal_time,
272    ///     elements.epoch(),
273    ///     elements.drag_term,
274    ///     sgp4::Orbit::from_kozai_elements(
275    ///         &sgp4::WGS84,
276    ///         elements.inclination * (core::f64::consts::PI / 180.0),
277    ///         elements.right_ascension * (core::f64::consts::PI / 180.0),
278    ///         elements.eccentricity,
279    ///         elements.argument_of_perigee * (core::f64::consts::PI / 180.0),
280    ///         elements.mean_anomaly * (core::f64::consts::PI / 180.0),
281    ///         elements.mean_motion * (core::f64::consts::PI / 720.0),
282    ///     )?,
283    /// )?;
284    /// #     Ok(())
285    /// # }
286    /// ```
287    pub fn new(
288        geopotential: Geopotential,
289        epoch_to_sidereal_time: impl Fn(f64) -> f64,
290        epoch: f64,
291        drag_term: f64,
292        orbit_0: propagator::Orbit,
293    ) -> core::result::Result<Self, OutOfRangeEpochEccentricity> {
294        if orbit_0.eccentricity < 0.0 || orbit_0.eccentricity >= 1.0 {
295            Err(OutOfRangeEpochEccentricity(orbit_0.eccentricity))
296        } else {
297            // p₁ = cos I₀
298            let p1 = orbit_0.inclination.cos();
299
300            // p₂ = 1 − e₀²
301            let p2 = 1.0 - orbit_0.eccentricity.powi(2);
302
303            // k₆ = 3 p₁² - 1
304            let k6 = 3.0 * p1.powi(2) - 1.0;
305
306            // a₀" = (kₑ / n₀")²ᐟ³
307            let a0 = (geopotential.ke / orbit_0.mean_motion).powf(2.0 / 3.0);
308
309            // p₃ = a₀" (1 - e₀)
310            let p3 = a0 * (1.0 - orbit_0.eccentricity);
311            let (s, p6) = {
312                // p₄ = aₑ (p₃ - 1)
313                let p4 = geopotential.ae * (p3 - 1.0);
314
315                // p₅ = │ 20      if p₄ < 98
316                //      │ p₄ - 78 if 98 ≤ p₄ < 156
317                //      │ 78      otherwise
318                let p5 = if p4 < 98.0 {
319                    20.0
320                } else if p4 < 156.0 {
321                    p4 - 78.0
322                } else {
323                    78.0
324                };
325                (
326                    // s = p₅ / aₑ + 1
327                    p5 / geopotential.ae + 1.0,
328                    // p₆ = ((120 - p₅) / aₑ)⁴
329                    ((120.0 - p5) / geopotential.ae).powi(4),
330                )
331            };
332
333            // ξ = 1 / (a₀" - s)
334            let xi = 1.0 / (a0 - s);
335
336            // p₇ = p₆ ξ⁴
337            let p7 = p6 * xi.powi(4);
338
339            // η = a₀" e₀ ξ
340            let eta = a0 * orbit_0.eccentricity * xi;
341
342            // p₈ = |1 - η²|
343            let p8 = (1.0 - eta.powi(2)).abs();
344
345            // p₉ = p₇ / p₈⁷ᐟ²
346            let p9 = p7 / p8.powf(3.5);
347
348            // C₁ = B* p₉ n₀" (a₀" (1 + ³/₂ η² + e₀ η (4 + η²))
349            //      + ³/₈ J₂ ξ k₆ (8 + 3 η² (8 + η²)) / p₈)
350            let c1 = drag_term
351                * (p9
352                    * orbit_0.mean_motion
353                    * (a0
354                        * (1.0
355                            + 1.5 * eta.powi(2)
356                            + orbit_0.eccentricity * eta * (4.0 + eta.powi(2)))
357                        + 0.375 * geopotential.j2 * xi / p8
358                            * k6
359                            * (8.0 + 3.0 * eta.powi(2) * (8.0 + eta.powi(2)))));
360
361            // p₁₀ = (a₀" p₂)⁻²
362            let p10 = 1.0 / (a0 * p2).powi(2);
363
364            // β₀ = p₂¹ᐟ²
365            let b0 = p2.sqrt();
366
367            // p₁₁ = ³/₂ J₂ p₁₀ n₀"
368            let p11 = 1.5 * geopotential.j2 * p10 * orbit_0.mean_motion;
369
370            // p₁₂ = ¹/₂ p₁₁ J₂ p₁₀
371            let p12 = 0.5 * p11 * geopotential.j2 * p10;
372
373            // p₁₃ = - ¹⁵/₃₂ J₄ p₁₀² n₀"
374            let p13 = -0.46875 * geopotential.j4 * p10.powi(2) * orbit_0.mean_motion;
375
376            // p₁₄ = - p₁₁ p₁ + (¹/₂ p₁₂ (4 - 19 p₁²) + 2 p₁₃ (3 - 7 p₁²)) p₁
377            let p14 = -p11 * p1
378                + (0.5 * p12 * (4.0 - 19.0 * p1.powi(2)) + 2.0 * p13 * (3.0 - 7.0 * p1.powi(2)))
379                    * p1;
380
381            // k₁₄ = - ¹/₂ p₁₁ (1 - 5 p₁²) + ¹/₁₆ p₁₂ (7 - 114 p₁² + 395 p₁⁴)
382            let k14 = -0.5 * p11 * (1.0 - 5.0 * p1.powi(2))
383                + 0.0625 * p12 * (7.0 - 114.0 * p1.powi(2) + 395.0 * p1.powi(4))
384                + p13 * (3.0 - 36.0 * p1.powi(2) + 49.0 * p1.powi(4));
385
386            // p₁₅ = n₀" + ¹/₂ p₁₁ β₀ k₆ + ¹/₁₆ p₁₂ β₀ (13 - 78 p₁² + 137 p₁⁴)
387            let p15 = orbit_0.mean_motion
388                + 0.5 * p11 * b0 * k6
389                + 0.0625 * p12 * b0 * (13.0 - 78.0 * p1.powi(2) + 137.0 * p1.powi(4));
390
391            // C₄ = 2 B* n₀" p₉ a₀" p₂ (
392            //      η (2 + ¹/₂ η²)
393            //      + e₀ (¹/₂ + 2 η²)
394            //      - J₂ ξ / (a p₈) (-3 k₆ (1 - 2 e₀ η + η² (³/₂ - ¹/₂ e₀ η))
395            //      + ³/₄ (1 - p₁²) (2 η² - e₀ η (1 + η²)) cos 2 ω₀)
396            let c4 = drag_term
397                * (2.0
398                    * orbit_0.mean_motion
399                    * p9
400                    * a0
401                    * p2
402                    * (eta * (2.0 + 0.5 * eta.powi(2))
403                        + orbit_0.eccentricity * (0.5 + 2.0 * eta.powi(2))
404                        - geopotential.j2 * xi / (a0 * p8)
405                            * (-3.0
406                                * k6
407                                * (1.0 - 2.0 * orbit_0.eccentricity * eta
408                                    + eta.powi(2) * (1.5 - 0.5 * orbit_0.eccentricity * eta))
409                                + 0.75
410                                    * (1.0 - p1.powi(2))
411                                    * (2.0 * eta.powi(2)
412                                        - orbit_0.eccentricity * eta * (1.0 + eta.powi(2)))
413                                    * (2.0 * orbit_0.argument_of_perigee).cos())));
414
415            // k₀ = - ⁷/₂ p₂ p₁₁ p₁ C₁
416            let k0 = 3.5 * p2 * (-p11 * p1) * c1;
417
418            // k₁ = ³/₂ C₁
419            let k1 = 1.5 * c1;
420
421            if orbit_0.mean_motion > 2.0 * core::f64::consts::PI / 225.0 {
422                Ok(near_earth::constants(
423                    geopotential,
424                    drag_term,
425                    orbit_0,
426                    p1,
427                    a0,
428                    s,
429                    xi,
430                    eta,
431                    c1,
432                    c4,
433                    k0,
434                    k1,
435                    k6,
436                    k14,
437                    p2,
438                    p3,
439                    p7,
440                    p9,
441                    p14,
442                    p15,
443                ))
444            } else {
445                Ok(deep_space::constants(
446                    geopotential,
447                    epoch_to_sidereal_time,
448                    epoch,
449                    orbit_0,
450                    p1,
451                    a0,
452                    c1,
453                    b0,
454                    c4,
455                    k0,
456                    k1,
457                    k14,
458                    p2,
459                    p14,
460                    p15,
461                ))
462            }
463        }
464    }
465
466    /// Initializes a new propagator from an `Elements` object
467    ///
468    /// This is the recommended method to initialize a propagator from a TLE or OMM.
469    /// The WGS84 model, the IAU sidereal time expression and the accurate UTC to J2000 expression are used.
470    ///
471    /// # Arguments
472    ///
473    /// * `elements` - Orbital elements and drag term parsed from a TLE or OMM
474    ///
475    /// # Example
476    ///
477    /// ```
478    /// # fn main() -> anyhow::Result<()> {
479    /// let constants = sgp4::Constants::from_elements(
480    ///     &sgp4::Elements::from_tle(
481    ///         Some("ISS (ZARYA)".to_owned()),
482    ///         "1 25544U 98067A   20194.88612269 -.00002218  00000-0 -31515-4 0  9992".as_bytes(),
483    ///         "2 25544  51.6461 221.2784 0001413  89.1723 280.4612 15.49507896236008".as_bytes(),
484    ///     )?,
485    /// )?;
486    /// #     Ok(())
487    /// # }
488    /// ```
489    pub fn from_elements(elements: &Elements) -> core::result::Result<Self, ElementsError> {
490        Ok(Constants::new(
491            WGS84,
492            iau_epoch_to_sidereal_time,
493            elements.epoch(),
494            elements.drag_term,
495            Orbit::from_kozai_elements(
496                &WGS84,
497                elements.inclination * (core::f64::consts::PI / 180.0),
498                elements.right_ascension * (core::f64::consts::PI / 180.0),
499                elements.eccentricity,
500                elements.argument_of_perigee * (core::f64::consts::PI / 180.0),
501                elements.mean_anomaly * (core::f64::consts::PI / 180.0),
502                elements.mean_motion * (core::f64::consts::PI / 720.0),
503            )?,
504        )?)
505    }
506
507    /// Initializes a new propagator from an `Elements` object
508    ///
509    /// This method should be used if compatibility with the AFSPC implementation is needed.
510    /// The WGS72 model, the AFSPC sidereal time expression and the AFSPC UTC to J2000 expression are used.
511    ///
512    /// # Arguments
513    ///
514    /// * `elements` - Orbital elements and drag term parsed from a TLE or OMM
515    ///
516    /// # Example
517    ///
518    /// ```
519    /// # fn main() -> anyhow::Result<()> {
520    /// let constants = sgp4::Constants::from_elements_afspc_compatibility_mode(
521    ///     &sgp4::Elements::from_tle(
522    ///         Some("ISS (ZARYA)".to_owned()),
523    ///         "1 25544U 98067A   20194.88612269 -.00002218  00000-0 -31515-4 0  9992".as_bytes(),
524    ///         "2 25544  51.6461 221.2784 0001413  89.1723 280.4612 15.49507896236008".as_bytes(),
525    ///     )?,
526    /// )?;
527    /// #     Ok(())
528    /// # }
529    /// ```
530    pub fn from_elements_afspc_compatibility_mode(
531        elements: &Elements,
532    ) -> core::result::Result<Self, ElementsError> {
533        Ok(Constants::new(
534            WGS72,
535            afspc_epoch_to_sidereal_time,
536            elements.epoch_afspc_compatibility_mode(),
537            elements.drag_term,
538            Orbit::from_kozai_elements(
539                &WGS72,
540                elements.inclination * (core::f64::consts::PI / 180.0),
541                elements.right_ascension * (core::f64::consts::PI / 180.0),
542                elements.eccentricity,
543                elements.argument_of_perigee * (core::f64::consts::PI / 180.0),
544                elements.mean_anomaly * (core::f64::consts::PI / 180.0),
545                elements.mean_motion * (core::f64::consts::PI / 720.0),
546            )?,
547        )?)
548    }
549
550    /// Returns the initial deep space resonance integrator state
551    ///
552    /// For most orbits, SGP4 propagation is stateless.
553    /// That is, predictions at different times are always calculated from the epoch quantities.
554    /// No calculations are saved by propagating to successive times sequentially.
555    ///
556    /// However, resonant deep space orbits (geosynchronous or Molniya) use an integrator
557    /// to estimate the resonance effects of Earth gravity, with a 720 min time step. If the propagation
558    /// times are monotonic, a few operations per prediction can be saved by re-using the integrator state.
559    ///
560    /// The high-level API `Constants::propagate` re-initializes the state with each propagation for simplicity.
561    /// `Constants::initial_state` and `Constants::propagate_from_state` can be used together
562    /// to speed up resonant deep space satellite propagation.
563    /// For non-deep space or non-resonant orbits, their behavior is identical to `Constants::propagate`.
564    ///
565    /// See `Constants::propagate_from_state` for an example.
566    pub fn initial_state(&self) -> Option<ResonanceState> {
567        match &self.method {
568            propagator::Method::NearEarth { .. } => None,
569            propagator::Method::DeepSpace { resonant, .. } => match resonant {
570                propagator::Resonant::No { .. } => None,
571                propagator::Resonant::Yes { lambda_0, .. } => {
572                    Some(ResonanceState::new(self.orbit_0.mean_motion, *lambda_0))
573                }
574            },
575        }
576    }
577
578    /// Calculates the SGP4 position and velocity predictions
579    ///
580    /// This is an advanced API which results in marginally faster propagation than `Constants::propagate` in some cases
581    /// (see `Constants::initial_state` for details), at the cost of added complexity for the user.
582    ///
583    /// The propagation times must be monotonic if the same resonance state is used repeatedly.
584    /// The `afspc_compatibility_mode` makes a difference only if the satellite is on a Lyddane deep space orbit
585    /// (period greater than 225 min and inclination smaller than 0.2 rad).
586    ///
587    /// # Arguments
588    ///
589    /// * `t` - The number of minutes since epoch (can be positive, negative or zero)
590    /// * `state` - The deep space propagator state returned by `Constants::initial_state`
591    /// * `afspc_compatibility_mode` - Set to true if compatibility with the AFSPC implementation is needed
592    ///
593    /// # Example
594    ///
595    /// ```
596    /// # fn main() -> anyhow::Result<()> {
597    /// let elements = sgp4::Elements::from_tle(
598    ///     Some("MOLNIYA 1-36".to_owned()),
599    ///     "1 08195U 75081A   06176.33215444  .00000099  00000-0  11873-3 0   813".as_bytes(),
600    ///     "2 08195  64.1586 279.0717 6877146 264.7651  20.2257  2.00491383225656".as_bytes(),
601    /// )?;
602    /// let constants = sgp4::Constants::from_elements(&elements)?;
603    /// let mut state = constants.initial_state();
604    /// for days in 0..7 {
605    ///     println!("t = {} min", days * 60 * 24);
606    ///     let prediction =
607    ///         constants.propagate_from_state(sgp4::MinutesSinceEpoch((days * 60 * 24) as f64), state.as_mut(), false)?;
608    ///     println!("    r = {:?} km", prediction.position);
609    ///     println!("    ṙ = {:?} km.s⁻¹", prediction.velocity);
610    /// }
611    /// #     Ok(())
612    /// # }
613    /// ```
614    #[allow(clippy::many_single_char_names)]
615    pub fn propagate_from_state(
616        &self,
617        t: MinutesSinceEpoch,
618        state: Option<&mut ResonanceState>,
619        afspc_compatibility_mode: bool,
620    ) -> core::result::Result<Prediction, gp::Error> {
621        // p₂₂ = Ω₀ + Ω̇ t + k₀ t²
622        let p22 =
623            self.orbit_0.right_ascension + self.right_ascension_dot * t.0 + self.k0 * t.0.powi(2);
624
625        // p₂₃ = ω₀ + ω̇ t
626        let p23 = self.orbit_0.argument_of_perigee + self.argument_of_perigee_dot * t.0;
627        let (orbit, a, p32, p33, p34, p35, p36) = match &self.method {
628            propagator::Method::NearEarth {
629                a0,
630                k2,
631                k3,
632                k4,
633                k5,
634                k6,
635                high_altitude,
636            } => {
637                assert!(
638                    state.is_none(),
639                    "state must be None with a near-earth propagator",
640                );
641                self.near_earth_orbital_elements(
642                    *a0,
643                    *k2,
644                    *k3,
645                    *k4,
646                    *k5,
647                    *k6,
648                    high_altitude,
649                    t.0,
650                    p22,
651                    p23,
652                )
653            }
654            propagator::Method::DeepSpace {
655                eccentricity_dot,
656                inclination_dot,
657                solar_perturbations,
658                lunar_perturbations,
659                resonant,
660            } => self.deep_space_orbital_elements(
661                *eccentricity_dot,
662                *inclination_dot,
663                solar_perturbations,
664                lunar_perturbations,
665                resonant,
666                state,
667                t.0,
668                p22,
669                p23,
670                afspc_compatibility_mode,
671            ),
672        }?;
673
674        // p₃₇ = 1 / (a (1 - e²))
675        let p37 = 1.0 / (a * (1.0 - orbit.eccentricity.powi(2)));
676
677        // aₓₙ = e cos ω
678        let axn = orbit.eccentricity * orbit.argument_of_perigee.cos();
679
680        // aᵧₙ = e sin ω + p₃₇ p₃₂
681        let ayn = orbit.eccentricity * orbit.argument_of_perigee.sin() + p37 * p32;
682
683        // p₃₈ = M + ω + p₃₇ p₃₅ aₓₙ rem 2π
684        let p38 = (orbit.mean_anomaly + orbit.argument_of_perigee + p37 * p35 * axn)
685            % (2.0 * core::f64::consts::PI);
686
687        // (E + ω)₀ = p₃₈
688        let mut ew = p38;
689        for _ in 0..10 {
690            //             p₃₈ - aᵧₙ cos (E + ω)ᵢ + aₓₙ sin (E + ω)ᵢ - (E + ω)ᵢ
691            // Δ(E + ω)ᵢ = ---------------------------------------------------
692            //                   1 - cos (E + ω)ᵢ aₓₙ - sin (E + ω)ᵢ aᵧₙ
693            let delta = (p38 - ayn * ew.cos() + axn * ew.sin() - ew)
694                / (1.0 - ew.cos() * axn - ew.sin() * ayn);
695
696            if delta.abs() < 1.0e-12 {
697                break;
698            }
699
700            // (E + ω)ᵢ₊₁ = (E + ω)ᵢ + Δ(E + ω)ᵢ|[-0.95, 0.95]
701            ew += delta.clamp(-0.95, 0.95);
702        }
703
704        // p₃₉ = aₓₙ² + aᵧₙ²
705        let p39 = axn.powi(2) + ayn.powi(2);
706
707        // pₗ = a (1 - p₃₉)
708        let pl = a * (1.0 - p39);
709        if pl < 0.0 {
710            Err(gp::Error::NegativeSemiLatusRectum { t: t.0 })
711        } else {
712            // p₄₀ = aₓₙ sin(E + ω) - aᵧₙ cos(E + ω)
713            let p40 = axn * ew.sin() - ayn * ew.cos();
714
715            // r = a (1 - aₓₙ cos(E + ω) + aᵧₙ sin(E + ω))
716            let r = a * (1.0 - (axn * ew.cos() + ayn * ew.sin()));
717
718            // ṙ = a¹ᐟ² p₄₀ / r
719            let r_dot = a.sqrt() * p40 / r;
720
721            // β = (1 - p₃₉)¹ᐟ²
722            let b = (1.0 - p39).sqrt();
723
724            // p₄₁ = p₄₀ / (1 + β)
725            let p41 = p40 / (1.0 + b);
726
727            // p₄₂ = a / r (sin(E + ω) - aᵧₙ - aₓₙ p₄₁)
728            let p42 = a / r * (ew.sin() - ayn - axn * p41);
729
730            // p₄₃ = a / r (cos(E + ω) - aₓₙ + aᵧₙ p₄₁)
731            let p43 = a / r * (ew.cos() - axn + ayn * p41);
732
733            //           p₄₂
734            // u = tan⁻¹ ---
735            //           p₄₃
736            let u = p42.atan2(p43);
737
738            // p₄₄ = 2 p₄₃ p₄₂
739            let p44 = 2.0 * p43 * p42;
740
741            // p₄₅ = 1 - 2 p₄₂²
742            let p45 = 1.0 - 2.0 * p42.powi(2);
743
744            // p₄₆ = (¹/₂ J₂ / pₗ) / pₗ
745            let p46 = 0.5 * self.geopotential.j2 / pl / pl;
746
747            // rₖ = r (1 - ³/₂ p₄₆ β p₃₆) + ¹/₂ (¹/₂ J₂ / pₗ) p₃₃ p₄₅
748            let rk = r * (1.0 - 1.5 * p46 * b * p36)
749                + 0.5 * (0.5 * self.geopotential.j2 / pl) * p33 * p45;
750
751            // uₖ = u - ¹/₄ p₄₆ p₃₄ p₄₄
752            let uk = u - 0.25 * p46 * p34 * p44;
753
754            // Iₖ = I + ³/₂ p₄₆ cos I sin I p₄₅
755            let inclination_k = orbit.inclination
756                + 1.5 * p46 * orbit.inclination.cos() * orbit.inclination.sin() * p45;
757
758            // Ωₖ = Ω + ³/₂ p₄₆ cos I p₄₄
759            let right_ascension_k =
760                orbit.right_ascension + 1.5 * p46 * orbit.inclination.cos() * p44;
761
762            // ṙₖ = ṙ + n (¹/₂ J₂ / pₗ) p₃₃ / kₑ
763            let rk_dot = r_dot
764                - orbit.mean_motion * (0.5 * self.geopotential.j2 / pl) * p33 * p44
765                    / self.geopotential.ke;
766
767            // rḟₖ = pₗ¹ᐟ² / r + n (¹/₂ J₂ / pₗ) (p₃₃ p₄₅ + ³/₂ p₃₆) / kₑ
768            let rfk_dot = pl.sqrt() / r
769                + orbit.mean_motion * (0.5 * self.geopotential.j2 / pl) * (p33 * p45 + 1.5 * p36)
770                    / self.geopotential.ke;
771
772            // u₀ = - sin Ωₖ cos Iₖ sin uₖ + cos Ωₖ cos uₖ
773            let u0 = -right_ascension_k.sin() * inclination_k.cos() * uk.sin()
774                + right_ascension_k.cos() * uk.cos();
775            // u₁ = cos Ωₖ cos Iₖ sin uₖ + sin Ωₖ cos uₖ
776            let u1 = right_ascension_k.cos() * inclination_k.cos() * uk.sin()
777                + right_ascension_k.sin() * uk.cos();
778            // u₂ = sin Iₖ sin uₖ
779            let u2 = inclination_k.sin() * uk.sin();
780            Ok(Prediction {
781                position: [
782                    // r₀ = rₖ u₀ aₑ
783                    rk * u0 * self.geopotential.ae,
784                    // r₁ = rₖ u₁ aₑ
785                    rk * u1 * self.geopotential.ae,
786                    // r₂ = rₖ u₂ aₑ
787                    rk * u2 * self.geopotential.ae,
788                ],
789                velocity: [
790                    // ṙ₀ = (ṙₖ u₀ + rḟₖ (- sin Ωₖ cos Iₖ cos uₖ - cos Ωₖ sin uₖ)) aₑ kₑ / 60
791                    (rk_dot * u0
792                        + rfk_dot
793                            * (-right_ascension_k.sin() * inclination_k.cos() * uk.cos()
794                                - right_ascension_k.cos() * uk.sin()))
795                        * (self.geopotential.ae * self.geopotential.ke / 60.0),
796                    // ṙ₁ = (ṙₖ u₁ + rḟₖ (cos Ωₖ cos Iₖ cos uₖ - sin Ωₖ sin uₖ)) aₑ kₑ / 60
797                    (rk_dot * u1
798                        + rfk_dot
799                            * (right_ascension_k.cos() * inclination_k.cos() * uk.cos()
800                                - right_ascension_k.sin() * uk.sin()))
801                        * (self.geopotential.ae * self.geopotential.ke / 60.0),
802                    // ṙ₂ = (ṙₖ u₂ + rḟₖ (sin Iₖ cos uₖ)) aₑ kₑ / 60
803                    (rk_dot * u2 + rfk_dot * (inclination_k.sin() * uk.cos()))
804                        * (self.geopotential.ae * self.geopotential.ke / 60.0),
805                ],
806            })
807        }
808    }
809
810    /// Calculates the SGP4 position and velocity predictions
811    ///
812    /// This is the recommended method to propagate epoch orbital elements.
813    ///
814    /// # Arguments
815    /// `t` - The number of minutes since epoch (can be positive, negative or zero)
816    ///
817    /// # Example
818    ///
819    /// ```
820    /// # fn main() -> anyhow::Result<()> {
821    /// let constants = sgp4::Constants::from_elements_afspc_compatibility_mode(
822    ///     &sgp4::Elements::from_tle(
823    ///         Some("ISS (ZARYA)".to_owned()),
824    ///         "1 25544U 98067A   20194.88612269 -.00002218  00000-0 -31515-4 0  9992".as_bytes(),
825    ///         "2 25544  51.6461 221.2784 0001413  89.1723 280.4612 15.49507896236008".as_bytes(),
826    ///     )?,
827    /// )?;
828    /// let prediction = constants.propagate(sgp4::MinutesSinceEpoch(60.0 * 24.0));
829    /// #     Ok(())
830    /// # }
831    /// ```
832    pub fn propagate(&self, t: MinutesSinceEpoch) -> core::result::Result<Prediction, gp::Error> {
833        self.propagate_from_state(t, self.initial_state().as_mut(), false)
834    }
835
836    /// Calculates the SGP4 position and velocity predictions
837    ///
838    /// This method should be used if compatibility with the AFSPC implementation is needed.
839    /// Its behavior is different from `Constants::propagate`
840    /// only if the satellite is on a Lyddane deep space orbit
841    /// (period greater than 225 min and inclination smaller than 0.2 rad).
842    ///
843    /// # Arguments
844    /// `t` - The number of minutes since epoch (can be positive, negative or zero)
845    ///
846    /// # Example
847    ///
848    /// ```
849    /// # fn main() -> anyhow::Result<()> {
850    /// let constants = sgp4::Constants::from_elements_afspc_compatibility_mode(
851    ///     &sgp4::Elements::from_tle(
852    ///         Some("ISS (ZARYA)".to_owned()),
853    ///         "1 25544U 98067A   20194.88612269 -.00002218  00000-0 -31515-4 0  9992".as_bytes(),
854    ///         "2 25544  51.6461 221.2784 0001413  89.1723 280.4612 15.49507896236008".as_bytes(),
855    ///     )?,
856    /// )?;
857    /// let prediction = constants.propagate_afspc_compatibility_mode(sgp4::MinutesSinceEpoch(60.0 * 24.0));
858    /// #     Ok(())
859    /// # }
860    /// ```
861    pub fn propagate_afspc_compatibility_mode(
862        &self,
863        t: MinutesSinceEpoch,
864    ) -> core::result::Result<Prediction, gp::Error> {
865        self.propagate_from_state(t, self.initial_state().as_mut(), true)
866    }
867}