Struct sgp4::Constants

source ·
pub struct Constants { /* private fields */ }
Expand description

Propagator variables calculated from epoch quantities and used during propagation

Constants can be initialized from general perturbation elements. They are not mutated during propagation, which means they can be used by different threads in parallel (for example to generate predictions at different times).

Implementations

Initializes a new propagator from epoch quantities

If the orbital elements are obtained from a TLE or OMM, the convenience function sgp4::Constants::from_elements can be used instead of manually mapping the Elements fields to the Constants::new parameters.

Arguments
  • geopotential - The model of Earth gravity to use in the conversion
  • epoch_to_sidereal_time - The function to use to convert the J2000 epoch to sidereal time
  • epoch - The number of years since UTC 1 January 2000 12h00 (J2000)
  • drag_term - The radiation pressure coefficient in earth radii⁻¹ (B*)
  • orbit_0 - The Brouwer orbital elements at epoch
Example
let elements = sgp4::Elements::from_tle(
    Some("ISS (ZARYA)".to_owned()),
    "1 25544U 98067A   20194.88612269 -.00002218  00000-0 -31515-4 0  9992".as_bytes(),
    "2 25544  51.6461 221.2784 0001413  89.1723 280.4612 15.49507896236008".as_bytes(),
)?;
let constants = sgp4::Constants::new(
    sgp4::WGS84,
    sgp4::iau_epoch_to_sidereal_time,
    elements.epoch(),
    elements.drag_term,
    sgp4::Orbit::from_kozai_elements(
        &sgp4::WGS84,
        elements.inclination * (core::f64::consts::PI / 180.0),
        elements.right_ascension * (core::f64::consts::PI / 180.0),
        elements.eccentricity,
        elements.argument_of_perigee * (core::f64::consts::PI / 180.0),
        elements.mean_anomaly * (core::f64::consts::PI / 180.0),
        elements.mean_motion * (core::f64::consts::PI / 720.0),
    )?,
)?;

Initializes a new propagator from an Elements object

This is the recommended method to initialize a propagator from a TLE or OMM. The WGS84 model, the IAU sidereal time expression and the accurate UTC to J2000 expression are used.

Arguments
  • elements - Orbital elements and drag term parsed from a TLE or OMM
Example
let constants = sgp4::Constants::from_elements(
    &sgp4::Elements::from_tle(
        Some("ISS (ZARYA)".to_owned()),
        "1 25544U 98067A   20194.88612269 -.00002218  00000-0 -31515-4 0  9992".as_bytes(),
        "2 25544  51.6461 221.2784 0001413  89.1723 280.4612 15.49507896236008".as_bytes(),
    )?,
)?;

Initializes a new propagator from an Elements object

This method should be used if compatibility with the AFSPC implementation is needed. The WGS72 model, the AFSPC sidereal time expression and the AFSPC UTC to J2000 expression are used.

Arguments
  • elements - Orbital elements and drag term parsed from a TLE or OMM
Example
let constants = sgp4::Constants::from_elements_afspc_compatibility_mode(
    &sgp4::Elements::from_tle(
        Some("ISS (ZARYA)".to_owned()),
        "1 25544U 98067A   20194.88612269 -.00002218  00000-0 -31515-4 0  9992".as_bytes(),
        "2 25544  51.6461 221.2784 0001413  89.1723 280.4612 15.49507896236008".as_bytes(),
    )?,
)?;

Returns the initial deep space resonance integrator state

For most orbits, SGP4 propagation is stateless. That is, predictions at different times are always calculated from the epoch quantities. No calculations are saved by propagating to successive times sequentially.

However, resonant deep space orbits (geosynchronous or Molniya) use an integrator to estimate the resonance effects of Earth gravity, with a 720 min time step. If the propagation times are monotonic, a few operations per prediction can be saved by re-using the integrator state.

The high-level API Constants::propagate re-initializes the state with each propagation for simplicity. Constants::initial_state and Constants::propagate_from_state can be used together to speed up resonant deep space satellite propagation. For non-deep space or non-resonant orbits, their behavior is identical to Constants::propagate.

See Constants::propagate_from_state for an example.

Calculates the SGP4 position and velocity predictions

This is an advanced API which results in marginally faster propagation than Constants::propagate in some cases (see Constants::initial_state for details), at the cost of added complexity for the user.

The propagation times must be monotonic if the same resonance state is used repeatedly. The afspc_compatibility_mode makes a difference only if the satellite is on a Lyddane deep space orbit (period greater than 225 min and inclination smaller than 0.2 rad).

Arguments
  • t - The number of minutes since epoch (can be positive, negative or zero)
  • state - The deep space propagator state returned by Constants::initial_state
  • afspc_compatibility_mode - Set to true if compatibility with the AFSPC implementation is needed
Example
let elements = sgp4::Elements::from_tle(
    Some("MOLNIYA 1-36".to_owned()),
    "1 08195U 75081A   06176.33215444  .00000099  00000-0  11873-3 0   813".as_bytes(),
    "2 08195  64.1586 279.0717 6877146 264.7651  20.2257  2.00491383225656".as_bytes(),
)?;
let constants = sgp4::Constants::from_elements(&elements)?;
let mut state = constants.initial_state();
for days in 0..7 {
    println!("t = {} min", days * 60 * 24);
    let prediction =
        constants.propagate_from_state((days * 60 * 24) as f64, state.as_mut(), false)?;
    println!("    r = {:?} km", prediction.position);
    println!("    ṙ = {:?} km.s⁻¹", prediction.velocity);
}

Calculates the SGP4 position and velocity predictions

This is the recommended method to propagate epoch orbital elements.

Arguments

t - The number of minutes since epoch (can be positive, negative or zero)

Example
let constants = sgp4::Constants::from_elements_afspc_compatibility_mode(
    &sgp4::Elements::from_tle(
        Some("ISS (ZARYA)".to_owned()),
        "1 25544U 98067A   20194.88612269 -.00002218  00000-0 -31515-4 0  9992".as_bytes(),
        "2 25544  51.6461 221.2784 0001413  89.1723 280.4612 15.49507896236008".as_bytes(),
    )?,
)?;
let prediction = constants.propagate(60.0 * 24.0);

Calculates the SGP4 position and velocity predictions

This method should be used if compatibility with the AFSPC implementation is needed. Its behavior is different from Constants::propagate only if the satellite is on a Lyddane deep space orbit (period greater than 225 min and inclination smaller than 0.2 rad).

Arguments

t - The number of minutes since epoch (can be positive, negative or zero)

Example
let constants = sgp4::Constants::from_elements_afspc_compatibility_mode(
    &sgp4::Elements::from_tle(
        Some("ISS (ZARYA)".to_owned()),
        "1 25544U 98067A   20194.88612269 -.00002218  00000-0 -31515-4 0  9992".as_bytes(),
        "2 25544  51.6461 221.2784 0001413  89.1723 280.4612 15.49507896236008".as_bytes(),
    )?,
)?;
let prediction = constants.propagate_afspc_compatibility_mode(60.0 * 24.0);

Trait Implementations

Deserialize this value from the given Serde deserializer. Read more
Serialize this value into the given Serde serializer. Read more

Auto Trait Implementations

Blanket Implementations

Gets the TypeId of self. Read more
Immutably borrows from an owned value. Read more
Mutably borrows from an owned value. Read more

Returns the argument unchanged.

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

The type returned in the event of a conversion error.
Performs the conversion.
The type returned in the event of a conversion error.
Performs the conversion.