01_orbit_prop/
main.rs

1#![doc = include_str!("./README.md")]
2extern crate log;
3extern crate nyx_space as nyx;
4extern crate pretty_env_logger as pel;
5
6use anise::{
7    almanac::metaload::MetaFile,
8    constants::{
9        celestial_objects::{MOON, SUN},
10        frames::{EARTH_J2000, IAU_EARTH_FRAME},
11    },
12};
13use hifitime::{Epoch, Unit};
14use log::warn;
15use nyx::{
16    cosmic::{Mass, MetaAlmanac, Orbit, SRPData},
17    dynamics::{Harmonics, OrbitalDynamics, SolarPressure, SpacecraftDynamics},
18    io::{gravity::HarmonicsMem, ExportCfg},
19    od::GroundStation,
20    propagators::Propagator,
21    Spacecraft, State,
22};
23use polars::{
24    frame::column::ScalarColumn,
25    prelude::{df, AnyValue, ChunkCompareIneq, Column, DataType, Scalar},
26};
27
28use std::{error::Error, sync::Arc};
29
30fn main() -> Result<(), Box<dyn Error>> {
31    pel::init();
32    // Dynamics models require planetary constants and ephemerides to be defined.
33    // Let's start by grabbing those by using ANISE's latest MetaAlmanac.
34    // This will automatically download the DE440s planetary ephemeris,
35    // the daily-updated Earth Orientation Parameters, the high fidelity Moon orientation
36    // parameters (for the Moon Mean Earth and Moon Principal Axes frames), and the PCK11
37    // planetary constants kernels.
38    // For details, refer to https://github.com/nyx-space/anise/blob/master/data/latest.dhall.
39    // Note that we place the Almanac into an Arc so we can clone it cheaply and provide read-only
40    // references to many functions.
41    let almanac = Arc::new(MetaAlmanac::latest().map_err(Box::new)?);
42    // Define the orbit epoch
43    let epoch = Epoch::from_gregorian_utc_hms(2024, 2, 29, 12, 13, 14);
44
45    // Define the orbit.
46    // First we need to fetch the Earth J2000 from information from the Almanac.
47    // This allows the frame to include the gravitational parameters and the shape of the Earth,
48    // defined as a tri-axial ellipoid. Note that this shape can be changed manually or in the Almanac
49    // by loading a different set of planetary constants.
50    let earth_j2000 = almanac.frame_from_uid(EARTH_J2000)?;
51
52    let orbit =
53        Orbit::try_keplerian_altitude(300.0, 0.015, 68.5, 65.2, 75.0, 0.0, epoch, earth_j2000)?;
54    // Print in in Keplerian form.
55    println!("{orbit:x}");
56
57    // There are two ways to propagate an orbit. We can make a quick approximation assuming only two-body
58    // motion. This is a useful first order approximation but it isn't used in real-world applications.
59
60    // This approach is a feature of ANISE.
61    let future_orbit_tb = orbit.at_epoch(epoch + Unit::Day * 3)?;
62    println!("{future_orbit_tb:x}");
63
64    // Two body propagation relies solely on Kepler's laws, so only the true anomaly will change.
65    println!(
66        "SMA changed by {:.3e} km",
67        orbit.sma_km()? - future_orbit_tb.sma_km()?
68    );
69    println!(
70        "ECC changed by {:.3e}",
71        orbit.ecc()? - future_orbit_tb.ecc()?
72    );
73    println!(
74        "INC changed by {:.3e} deg",
75        orbit.inc_deg()? - future_orbit_tb.inc_deg()?
76    );
77    println!(
78        "RAAN changed by {:.3e} deg",
79        orbit.raan_deg()? - future_orbit_tb.raan_deg()?
80    );
81    println!(
82        "AOP changed by {:.3e} deg",
83        orbit.aop_deg()? - future_orbit_tb.aop_deg()?
84    );
85    println!(
86        "TA changed by {:.3} deg",
87        orbit.ta_deg()? - future_orbit_tb.ta_deg()?
88    );
89
90    // Nyx is used for high fidelity propagation, not Keplerian propagation as above.
91    // Nyx only propagates Spacecraft at the moment, which allows it to account for acceleration
92    // models such as solar radiation pressure.
93
94    // Let's build a cubesat sized spacecraft, with an SRP area of 10 cm^2 and a mass of 9.6 kg.
95    let sc = Spacecraft::builder()
96        .orbit(orbit)
97        .mass(Mass::from_dry_mass(9.60))
98        .srp(SRPData {
99            area_m2: 10e-4,
100            coeff_reflectivity: 1.1,
101        })
102        .build();
103    println!("{sc:x}");
104
105    // Set up the spacecraft dynamics.
106
107    // Specify that the orbital dynamics must account for the graviational pull of the Moon and the Sun.
108    // The gravity of the Earth will also be accounted for since the spaceraft in an Earth orbit.
109    let mut orbital_dyn = OrbitalDynamics::point_masses(vec![MOON, SUN]);
110
111    // We want to include the spherical harmonics, so let's download the gravitational data from the Nyx Cloud.
112    // We're using the JGM3 model here, which is the default in GMAT.
113    let mut jgm3_meta = MetaFile {
114        uri: "http://public-data.nyxspace.com/nyx/models/JGM3.cof.gz".to_string(),
115        crc32: Some(0xF446F027), // Specifying the CRC32 avoids redownloading it if it's cached.
116    };
117    // And let's download it if we don't have it yet.
118    jgm3_meta.process(true)?;
119
120    // Build the spherical harmonics.
121    // The harmonics must be computed in the body fixed frame.
122    // We're using the long term prediction of the Earth centered Earth fixed frame, IAU Earth.
123    let harmonics_21x21 = Harmonics::from_stor(
124        almanac.frame_from_uid(IAU_EARTH_FRAME)?,
125        HarmonicsMem::from_cof(&jgm3_meta.uri, 21, 21, true).unwrap(),
126    );
127
128    // Include the spherical harmonics into the orbital dynamics.
129    orbital_dyn.accel_models.push(harmonics_21x21);
130
131    // We define the solar radiation pressure, using the default solar flux and accounting only
132    // for the eclipsing caused by the Earth.
133    let srp_dyn = SolarPressure::default(EARTH_J2000, almanac.clone())?;
134
135    // Finalize setting up the dynamics, specifying the force models (orbital_dyn) separately from the
136    // acceleration models (SRP in this case). Use `from_models` to specify multiple accel models.
137    let dynamics = SpacecraftDynamics::from_model(orbital_dyn, srp_dyn);
138
139    println!("{dynamics}");
140
141    // Finally, let's propagate this orbit to the same epoch as above.
142    // The first returned value is the spacecraft state at the final epoch.
143    // The second value is the full trajectory where the step size is variable step used by the propagator.
144    let (future_sc, trajectory) = Propagator::default(dynamics)
145        .with(sc, almanac.clone())
146        .until_epoch_with_traj(future_orbit_tb.epoch)?;
147
148    println!("=== High fidelity propagation ===");
149    println!(
150        "SMA changed by {:.3} km",
151        orbit.sma_km()? - future_sc.orbit.sma_km()?
152    );
153    println!(
154        "ECC changed by {:.6}",
155        orbit.ecc()? - future_sc.orbit.ecc()?
156    );
157    println!(
158        "INC changed by {:.3e} deg",
159        orbit.inc_deg()? - future_sc.orbit.inc_deg()?
160    );
161    println!(
162        "RAAN changed by {:.3} deg",
163        orbit.raan_deg()? - future_sc.orbit.raan_deg()?
164    );
165    println!(
166        "AOP changed by {:.3} deg",
167        orbit.aop_deg()? - future_sc.orbit.aop_deg()?
168    );
169    println!(
170        "TA changed by {:.3} deg",
171        orbit.ta_deg()? - future_sc.orbit.ta_deg()?
172    );
173
174    // We also have access to the full trajectory throughout the propagation.
175    println!("{trajectory}");
176
177    // With the trajectory, let's build a few data products.
178
179    // 1. Export the trajectory as a CCSDS OEM version 2.0 file and as a parquet file, which includes the Keplerian orbital elements.
180
181    trajectory.to_oem_file(
182        "./01_cubesat_hf_prop.oem",
183        ExportCfg::builder().step(Unit::Minute * 2).build(),
184    )?;
185
186    trajectory.to_parquet_with_cfg(
187        "./01_cubesat_hf_prop.parquet",
188        ExportCfg::builder().step(Unit::Minute * 2).build(),
189        almanac.clone(),
190    )?;
191
192    // 2. Compare the difference in the radial-intrack-crosstrack frame between the high fidelity
193    // and Keplerian propagation. The RIC frame is commonly used to compute the difference in position
194    // and velocity of different spacecraft.
195    // 3. Compute the azimuth, elevation, range, and range-rate data of that spacecraft as seen from Boulder, CO, USA.
196
197    let boulder_station = GroundStation::from_point(
198        "Boulder, CO, USA".to_string(),
199        40.014984,   // latitude in degrees
200        -105.270546, // longitude in degrees
201        1.6550,      // altitude in kilometers
202        almanac.frame_from_uid(IAU_EARTH_FRAME)?,
203    );
204
205    // We iterate over the trajectory, grabbing a state every two minutes.
206    let mut offset_s = vec![];
207    let mut epoch_str = vec![];
208    let mut ric_x_km = vec![];
209    let mut ric_y_km = vec![];
210    let mut ric_z_km = vec![];
211    let mut ric_vx_km_s = vec![];
212    let mut ric_vy_km_s = vec![];
213    let mut ric_vz_km_s = vec![];
214
215    let mut azimuth_deg = vec![];
216    let mut elevation_deg = vec![];
217    let mut range_km = vec![];
218    let mut range_rate_km_s = vec![];
219    for state in trajectory.every(Unit::Minute * 2) {
220        // Try to compute the Keplerian/two body state just in time.
221        // This method occasionally fails to converge on an appropriate true anomaly
222        // from the mean anomaly. If that happens, we just skip this state.
223        // The high fidelity and Keplerian states diverge continuously, and we're curious
224        // about the divergence in this quick analysis.
225        let this_epoch = state.epoch();
226        match orbit.at_epoch(this_epoch) {
227            Ok(tb_then) => {
228                offset_s.push((this_epoch - orbit.epoch).to_seconds());
229                epoch_str.push(format!("{this_epoch}"));
230                // Compute the two body state just in time.
231                let ric = state.orbit.ric_difference(&tb_then)?;
232                ric_x_km.push(ric.radius_km.x);
233                ric_y_km.push(ric.radius_km.y);
234                ric_z_km.push(ric.radius_km.z);
235                ric_vx_km_s.push(ric.velocity_km_s.x);
236                ric_vy_km_s.push(ric.velocity_km_s.y);
237                ric_vz_km_s.push(ric.velocity_km_s.z);
238
239                // Compute the AER data for each state.
240                let aer = almanac.azimuth_elevation_range_sez(
241                    state.orbit,
242                    boulder_station.to_orbit(this_epoch, &almanac)?,
243                    None,
244                    None,
245                )?;
246                azimuth_deg.push(aer.azimuth_deg);
247                elevation_deg.push(aer.elevation_deg);
248                range_km.push(aer.range_km);
249                range_rate_km_s.push(aer.range_rate_km_s);
250            }
251            Err(e) => warn!("{} {e}", state.epoch()),
252        };
253    }
254
255    // Build the data frames.
256    let ric_df = df!(
257        "Offset (s)" => offset_s.clone(),
258        "Epoch" => epoch_str.clone(),
259        "RIC X (km)" => ric_x_km,
260        "RIC Y (km)" => ric_y_km,
261        "RIC Z (km)" => ric_z_km,
262        "RIC VX (km/s)" => ric_vx_km_s,
263        "RIC VY (km/s)" => ric_vy_km_s,
264        "RIC VZ (km/s)" => ric_vz_km_s,
265    )?;
266
267    println!("RIC difference at start\n{}", ric_df.head(Some(10)));
268    println!("RIC difference at end\n{}", ric_df.tail(Some(10)));
269
270    let aer_df = df!(
271        "Offset (s)" => offset_s.clone(),
272        "Epoch" => epoch_str.clone(),
273        "azimuth (deg)" => azimuth_deg,
274        "elevation (deg)" => elevation_deg,
275        "range (km)" => range_km,
276        "range rate (km/s)" => range_rate_km_s,
277    )?;
278
279    // Finally, let's see when the spacecraft is visible, assuming 15 degrees minimum elevation.
280    let mask = aer_df
281        .column("elevation (deg)")?
282        .gt(&Column::Scalar(ScalarColumn::new(
283            "elevation mask (deg)".into(),
284            Scalar::new(DataType::Float64, AnyValue::Float64(15.0)),
285            offset_s.len(),
286        )))?;
287    let cubesat_visible = aer_df.filter(&mask)?;
288
289    println!("{cubesat_visible}");
290
291    Ok(())
292}