pub struct OrbitalDynamics {
pub accel_models: Vec<Arc<dyn AccelModel + Sync>>,
}Expand description
OrbitalDynamics provides the equations of motion for any celestial dynamic, without state transition matrix computation.
Fields§
§accel_models: Vec<Arc<dyn AccelModel + Sync>>Implementations§
Source§impl OrbitalDynamics
impl OrbitalDynamics
Sourcepub fn point_masses(celestial_objects: Vec<i32>) -> Self
pub fn point_masses(celestial_objects: Vec<i32>) -> Self
Initializes the point masses gravities with the provided list of bodies
Examples found in repository?
examples/03_geo_analysis/stationkeeping.rs (line 77)
28fn main() -> Result<(), Box<dyn Error>> {
29 pel::init();
30 // Set up the dynamics like in the orbit raise.
31 let almanac = Arc::new(MetaAlmanac::latest().map_err(Box::new)?);
32 let epoch = Epoch::from_gregorian_utc_hms(2024, 2, 29, 12, 13, 14);
33
34 // Define the GEO orbit, and we're just going to maintain it very tightly.
35 let earth_j2000 = almanac.frame_info(EARTH_J2000)?;
36 let orbit = Orbit::try_keplerian(42164.0, 1e-5, 0., 163.0, 75.0, 0.0, epoch, earth_j2000)?;
37 println!("{orbit:x}");
38
39 let sc = Spacecraft::builder()
40 .orbit(orbit)
41 .mass(Mass::from_dry_and_prop_masses(1000.0, 1000.0)) // 1000 kg of dry mass and prop, totalling 2.0 tons
42 .srp(SRPData::from_area(3.0 * 6.0)) // Assuming 1 kW/m^2 or 18 kW, giving a margin of 4.35 kW for on-propulsion consumption
43 .thruster(Thruster {
44 // "NEXT-STEP" row in Table 2
45 isp_s: 4435.0,
46 thrust_N: 0.472,
47 })
48 .mode(GuidanceMode::Thrust) // Start thrusting immediately.
49 .build();
50
51 // Set up the spacecraft dynamics like in the orbit raise example.
52
53 let prop_time = 30.0 * Unit::Day;
54
55 // Define the guidance law -- we're just using a Ruggiero controller as demonstrated in AAS-2004-5089.
56 let objectives = &[
57 Objective::within_tolerance(
58 StateParameter::Element(OrbitalElement::SemiMajorAxis),
59 42_165.0,
60 20.0,
61 ),
62 Objective::within_tolerance(
63 StateParameter::Element(OrbitalElement::Eccentricity),
64 0.001,
65 5e-5,
66 ),
67 Objective::within_tolerance(
68 StateParameter::Element(OrbitalElement::Inclination),
69 0.05,
70 1e-2,
71 ),
72 ];
73
74 let ruggiero_ctrl = Ruggiero::from_max_eclipse(objectives, sc, 0.2)?;
75 println!("{ruggiero_ctrl}");
76
77 let mut orbital_dyn = OrbitalDynamics::point_masses(vec![MOON, SUN]);
78
79 let mut jgm3_meta = MetaFile {
80 uri: "http://public-data.nyxspace.com/nyx/models/JGM3.cof.gz".to_string(),
81 crc32: Some(0xF446F027), // Specifying the CRC32 avoids redownloading it if it's cached.
82 };
83 jgm3_meta.process(true)?;
84
85 let harmonics = GravityField::from_stor(
86 almanac.frame_info(IAU_EARTH_FRAME)?,
87 GravityFieldData::from_cof(&jgm3_meta.uri, 8, 8, true)?,
88 );
89 orbital_dyn.accel_models.push(harmonics);
90
91 let srp_dyn = SolarPressure::default_flux(EARTH_J2000, almanac.clone())?;
92 let sc_dynamics = SpacecraftDynamics::from_model(orbital_dyn, srp_dyn)
93 .with_guidance_law(ruggiero_ctrl.clone());
94
95 println!("{sc_dynamics}");
96
97 // Finally, let's use the Monte Carlo framework built into Nyx to propagate spacecraft.
98
99 // Let's start by defining the dispersion.
100 // The MultivariateNormal structure allows us to define the dispersions in any of the orbital parameters, but these are applied directly in the Cartesian state space.
101 // Note that additional validation on the MVN is in progress -- https://github.com/nyx-space/nyx/issues/339.
102 let mc_rv = MvnSpacecraft::new(
103 sc,
104 vec![StateDispersion::zero_mean(
105 StateParameter::Element(OrbitalElement::SemiMajorAxis),
106 3.0,
107 )],
108 )?;
109
110 let my_mc = MonteCarlo::new(
111 sc, // Nominal state
112 mc_rv,
113 "03_geo_sk".to_string(), // Scenario name
114 None, // No specific seed specified, so one will be drawn from the computer's entropy.
115 );
116
117 // Build the propagator setup.
118 let setup = Propagator::rk89(
119 sc_dynamics.clone(),
120 IntegratorOptions::builder()
121 .min_step(10.0_f64.seconds())
122 .error_ctrl(ErrorControl::RSSCartesianStep)
123 .build(),
124 );
125
126 let num_runs = 25;
127 let rslts = my_mc.run_until_epoch(setup, almanac.clone(), sc.epoch() + prop_time, num_runs);
128
129 assert_eq!(rslts.runs.len(), num_runs);
130
131 rslts.to_parquet("03_geo_sk.parquet", ExportCfg::default())?;
132
133 Ok(())
134}More examples
examples/03_geo_analysis/raise.rs (line 95)
27fn main() -> Result<(), Box<dyn Error>> {
28 pel::init();
29
30 // Dynamics models require planetary constants and ephemerides to be defined.
31 // Let's start by grabbing those by using ANISE's latest MetaAlmanac.
32 // This will automatically download the DE440s planetary ephemeris,
33 // the daily-updated Earth Orientation Parameters, the high fidelity Moon orientation
34 // parameters (for the Moon Mean Earth and Moon Principal Axes frames), and the PCK11
35 // planetary constants kernels.
36 // For details, refer to https://github.com/nyx-space/anise/blob/master/data/latest.dhall.
37 // Note that we place the Almanac into an Arc so we can clone it cheaply and provide read-only
38 // references to many functions.
39 let almanac = Arc::new(MetaAlmanac::latest().map_err(Box::new)?);
40 // Fetch the EME2000 frame from the Almabac
41 let eme2k = almanac.frame_info(EARTH_J2000).unwrap();
42 // Define the orbit epoch
43 let epoch = Epoch::from_gregorian_utc_hms(2024, 2, 29, 12, 13, 14);
44
45 // Build the spacecraft itself.
46 // Using slide 6 of https://aerospace.org/sites/default/files/2018-11/Davis-Mayberry_HPSEP_11212018.pdf
47 // for the "next gen" SEP characteristics.
48
49 // GTO start
50 let orbit = Orbit::keplerian(24505.9, 0.725, 7.05, 0.0, 0.0, 0.0, epoch, eme2k);
51
52 let sc = Spacecraft::builder()
53 .orbit(orbit)
54 .mass(Mass::from_dry_and_prop_masses(1000.0, 1000.0)) // 1000 kg of dry mass and prop, totalling 2.0 tons
55 .srp(SRPData::from_area(3.0 * 6.0)) // Assuming 1 kW/m^2 or 18 kW, giving a margin of 4.35 kW for on-propulsion consumption
56 .thruster(Thruster {
57 // "NEXT-STEP" row in Table 2
58 isp_s: 4435.0,
59 thrust_N: 0.472,
60 })
61 .mode(GuidanceMode::Thrust) // Start thrusting immediately.
62 .build();
63
64 let prop_time = 180.0 * Unit::Day;
65
66 // Define the guidance law -- we're just using a Ruggiero controller as demonstrated in AAS-2004-5089.
67 let objectives = &[
68 Objective::within_tolerance(
69 StateParameter::Element(OrbitalElement::SemiMajorAxis),
70 42_165.0,
71 20.0,
72 ),
73 Objective::within_tolerance(
74 StateParameter::Element(OrbitalElement::Eccentricity),
75 0.001,
76 5e-5,
77 ),
78 Objective::within_tolerance(
79 StateParameter::Element(OrbitalElement::Inclination),
80 0.05,
81 1e-2,
82 ),
83 ];
84
85 // Ensure that we only thrust if we have more than 20% illumination.
86 let ruggiero_ctrl = Ruggiero::from_max_eclipse(objectives, sc, 0.2).unwrap();
87 println!("{ruggiero_ctrl}");
88
89 // Define the high fidelity dynamics
90
91 // Set up the spacecraft dynamics.
92
93 // Specify that the orbital dynamics must account for the graviational pull of the Moon and the Sun.
94 // The gravity of the Earth will also be accounted for since the spaceraft in an Earth orbit.
95 let mut orbital_dyn = OrbitalDynamics::point_masses(vec![MOON, SUN]);
96
97 // We want to include the spherical harmonics, so let's download the gravitational data from the Nyx Cloud.
98 // We're using the JGM3 model here, which is the default in GMAT.
99 let mut jgm3_meta = MetaFile {
100 uri: "http://public-data.nyxspace.com/nyx/models/JGM3.cof.gz".to_string(),
101 crc32: Some(0xF446F027), // Specifying the CRC32 avoids redownloading it if it's cached.
102 };
103 // And let's download it if we don't have it yet.
104 jgm3_meta.process(true)?;
105
106 // Build the spherical harmonics.
107 // The harmonics must be computed in the body fixed frame.
108 // We're using the long term prediction of the Earth centered Earth fixed frame, IAU Earth.
109 let harmonics = GravityField::from_stor(
110 almanac.frame_info(IAU_EARTH_FRAME)?,
111 GravityFieldData::from_cof(&jgm3_meta.uri, 8, 8, true).unwrap(),
112 );
113
114 // Include the spherical harmonics into the orbital dynamics.
115 orbital_dyn.accel_models.push(harmonics);
116
117 // We define the solar radiation pressure, using the default solar flux and accounting only
118 // for the eclipsing caused by the Earth.
119 let srp_dyn = SolarPressure::default_flux(EARTH_J2000, almanac.clone())?;
120
121 // Finalize setting up the dynamics, specifying the force models (orbital_dyn) separately from the
122 // acceleration models (SRP in this case). Use `from_models` to specify multiple accel models.
123 let sc_dynamics = SpacecraftDynamics::from_model(orbital_dyn, srp_dyn)
124 .with_guidance_law(ruggiero_ctrl.clone());
125
126 println!("{orbit:x}");
127
128 // We specify a minimum step in the propagator because the Ruggiero control would otherwise drive this step very low.
129 let (final_state, traj) = Propagator::rk89(
130 sc_dynamics.clone(),
131 IntegratorOptions::builder()
132 .min_step(10.0_f64.seconds())
133 .error_ctrl(ErrorControl::RSSCartesianStep)
134 .build(),
135 )
136 .with(sc, almanac.clone())
137 .for_duration_with_traj(prop_time)?;
138
139 let prop_usage = sc.mass.prop_mass_kg - final_state.mass.prop_mass_kg;
140 println!("{:x}", final_state.orbit);
141 println!("prop usage: {prop_usage:.3} kg");
142
143 // Finally, export the results for analysis, including the penumbra percentage throughout the orbit raise.
144 traj.to_parquet("./03_geo_raise.parquet", ExportCfg::default())?;
145
146 for status_line in ruggiero_ctrl.status(&final_state) {
147 println!("{status_line}");
148 }
149
150 ruggiero_ctrl
151 .achieved(&final_state)
152 .expect("objective not achieved");
153
154 Ok(())
155}examples/02_jwst_covar_monte_carlo/main.rs (line 99)
26fn main() -> Result<(), Box<dyn Error>> {
27 pel::init();
28 // Dynamics models require planetary constants and ephemerides to be defined.
29 // Let's start by grabbing those by using ANISE's latest MetaAlmanac.
30 // For details, refer to https://github.com/nyx-space/anise/blob/master/data/latest.dhall.
31
32 // Download the regularly update of the James Webb Space Telescope reconstucted (or definitive) ephemeris.
33 // Refer to https://naif.jpl.nasa.gov/pub/naif/JWST/kernels/spk/aareadme.txt for details.
34 let mut latest_jwst_ephem = MetaFile {
35 uri: "https://naif.jpl.nasa.gov/pub/naif/JWST/kernels/spk/jwst_rec.bsp".to_string(),
36 crc32: None,
37 };
38 latest_jwst_ephem.process(true)?;
39
40 // Load this ephem in the general Almanac we're using for this analysis.
41 let almanac = Arc::new(
42 MetaAlmanac::latest()
43 .map_err(Box::new)?
44 .load_from_metafile(latest_jwst_ephem, true)?,
45 );
46
47 // By loading this ephemeris file in the ANISE GUI or ANISE CLI, we can find the NAIF ID of the JWST
48 // in the BSP. We need this ID in order to query the ephemeris.
49 const JWST_NAIF_ID: i32 = -170;
50 // Let's build a frame in the J2000 orientation centered on the JWST.
51 const JWST_J2000: Frame = Frame::from_ephem_j2000(JWST_NAIF_ID);
52
53 // Since the ephemeris file is updated regularly, we'll just grab the latest state in the ephem.
54 let (earliest_epoch, latest_epoch) = almanac.spk_domain(JWST_NAIF_ID)?;
55 println!("JWST defined from {earliest_epoch} to {latest_epoch}");
56 // Fetch the state, printing it in the Earth J2000 frame.
57 let jwst_orbit = almanac.transform(JWST_J2000, EARTH_J2000, latest_epoch, None)?;
58 println!("{jwst_orbit:x}");
59
60 // Build the spacecraft
61 // SRP area assumed to be the full sunshield and mass if 6200.0 kg, c.f. https://webb.nasa.gov/content/about/faqs/facts.html
62 // SRP Coefficient of reflectivity assumed to be that of Kapton, i.e. 2 - 0.44 = 1.56, table 1 from https://amostech.com/TechnicalPapers/2018/Poster/Bengtson.pdf
63 let jwst = Spacecraft::builder()
64 .orbit(jwst_orbit)
65 .srp(SRPData {
66 area_m2: 21.197 * 14.162,
67 coeff_reflectivity: 1.56,
68 })
69 .mass(Mass::from_dry_mass(6200.0))
70 .build();
71
72 // Build up the spacecraft uncertainty builder.
73 // We can use the spacecraft uncertainty structure to build this up.
74 // We start by specifying the nominal state (as defined above), then the uncertainty in position and velocity
75 // in the RIC frame. We could also specify the Cr, Cd, and mass uncertainties, but these aren't accounted for until
76 // Nyx can also estimate the deviation of the spacecraft parameters.
77 let jwst_uncertainty = SpacecraftUncertainty::builder()
78 .nominal(jwst)
79 .frame(LocalFrame::RIC)
80 .x_km(0.5)
81 .y_km(0.3)
82 .z_km(1.5)
83 .vx_km_s(1e-4)
84 .vy_km_s(0.6e-3)
85 .vz_km_s(3e-3)
86 .build();
87
88 println!("{jwst_uncertainty}");
89
90 // Build the Kalman filter estimate.
91 // Note that we could have used the KfEstimate structure directly (as seen throughout the OD integration tests)
92 // but this approach requires quite a bit more boilerplate code.
93 let jwst_estimate = jwst_uncertainty.to_estimate()?;
94
95 // Set up the spacecraft dynamics.
96 // We'll use the point masses of the Earth, Sun, Jupiter (barycenter, because it's in the DE440), and the Moon.
97 // We'll also enable solar radiation pressure since the James Webb has a huge and highly reflective sun shield.
98
99 let orbital_dyn = OrbitalDynamics::point_masses(vec![MOON, SUN, JUPITER_BARYCENTER]);
100 let srp_dyn = SolarPressure::new(vec![EARTH_J2000, MOON_J2000], almanac.clone())?;
101
102 // Finalize setting up the dynamics.
103 let dynamics = SpacecraftDynamics::from_model(orbital_dyn, srp_dyn);
104
105 // Build the propagator set up to use for the whole analysis.
106 let setup = Propagator::default(dynamics);
107
108 // All of the analysis will use this duration.
109 let prediction_duration = 6.5 * Unit::Day;
110
111 // === Covariance mapping ===
112 // For the covariance mapping / prediction, we'll use the common orbit determination approach.
113 // This is done by setting up a spacecraft Kalman filter OD process, and predicting for the analysis duration.
114
115 // Build the propagation instance for the OD process.
116 let odp = SpacecraftKalmanOD::new(
117 setup.clone(),
118 KalmanVariant::DeviationTracking,
119 None,
120 BTreeMap::new(),
121 almanac.clone(),
122 );
123
124 // The prediction step is 1 minute by default, configured in the OD process, i.e. how often we want to know the covariance.
125 assert_eq!(odp.max_step, 1_i64.minutes());
126 // Finally, predict, and export the trajectory with covariance to a parquet file.
127 let od_sol = odp.predict_for(jwst_estimate, prediction_duration)?;
128 od_sol.to_parquet("./02_jwst_covar_map.parquet", ExportCfg::default())?;
129
130 // === Monte Carlo framework ===
131 // Nyx comes with a complete multi-threaded Monte Carlo frame. It's blazing fast.
132
133 let my_mc = MonteCarlo::new(
134 jwst, // Nominal state
135 jwst_estimate.to_random_variable()?,
136 "02_jwst".to_string(), // Scenario name
137 None, // No specific seed specified, so one will be drawn from the computer's entropy.
138 );
139
140 let num_runs = 5_000;
141 let rslts = my_mc.run_until_epoch(
142 setup,
143 almanac.clone(),
144 jwst.epoch() + prediction_duration,
145 num_runs,
146 );
147
148 assert_eq!(rslts.runs.len(), num_runs);
149 // Finally, export these results, computing the eclipse percentage for all of these results.
150
151 rslts.to_parquet("02_jwst_monte_carlo.parquet", ExportCfg::default())?;
152
153 Ok(())
154}examples/05_cislunar_spacecraft_link_od/main.rs (line 82)
34fn main() -> Result<(), Box<dyn Error>> {
35 pel::init();
36
37 // ====================== //
38 // === ALMANAC SET UP === //
39 // ====================== //
40
41 let manifest_dir =
42 PathBuf::from(std::env::var("CARGO_MANIFEST_DIR").unwrap_or(".".to_string()));
43
44 let out = manifest_dir.join("data/04_output/");
45
46 let almanac = Arc::new(
47 Almanac::new(
48 &manifest_dir
49 .join("data/01_planetary/pck08.pca")
50 .to_string_lossy(),
51 )
52 .unwrap()
53 .load(
54 &manifest_dir
55 .join("data/01_planetary/de440s.bsp")
56 .to_string_lossy(),
57 )
58 .unwrap(),
59 );
60
61 let eme2k = almanac.frame_info(EARTH_J2000).unwrap();
62 let moon_iau = almanac.frame_info(IAU_MOON_FRAME).unwrap();
63
64 let epoch = Epoch::from_gregorian_tai(2021, 5, 29, 19, 51, 16, 852_000);
65 let nrho = Orbit::cartesian(
66 166_473.631_302_239_7,
67 -274_715.487_253_382_7,
68 -211_233.210_176_686_7,
69 0.933_451_604_520_018_4,
70 0.436_775_046_841_900_9,
71 -0.082_211_021_250_348_95,
72 epoch,
73 eme2k,
74 );
75
76 let tx_nrho_sc = Spacecraft::from(nrho);
77
78 let state_luna = almanac.transform_to(nrho, MOON_J2000, None).unwrap();
79 println!("Start state (dynamics: Earth, Moon, Sun gravity):\n{state_luna}");
80
81 let bodies = vec![EARTH, SUN];
82 let dynamics = SpacecraftDynamics::new(OrbitalDynamics::point_masses(bodies));
83
84 let setup = Propagator::rk89(
85 dynamics,
86 IntegratorOptions::builder().max_step(0.5.minutes()).build(),
87 );
88
89 /* == Propagate the NRHO vehicle == */
90 let prop_time = 1.1 * state_luna.period().unwrap();
91
92 let (nrho_final, mut tx_traj) = setup
93 .with(tx_nrho_sc, almanac.clone())
94 .for_duration_with_traj(prop_time)
95 .unwrap();
96
97 tx_traj.name = Some("NRHO Tx SC".to_string());
98
99 println!("{tx_traj}");
100
101 /* == Propagate an LLO vehicle == */
102 let llo_orbit =
103 Orbit::try_keplerian_altitude(110.0, 1e-4, 90.0, 0.0, 0.0, 0.0, epoch, moon_iau).unwrap();
104
105 let llo_sc = Spacecraft::builder().orbit(llo_orbit).build();
106
107 let (_, llo_traj) = setup
108 .with(llo_sc, almanac.clone())
109 .until_epoch_with_traj(nrho_final.epoch())
110 .unwrap();
111
112 // Export the subset of the first two hours.
113 llo_traj
114 .clone()
115 .filter_by_offset(..2.hours())
116 .to_parquet_simple(out.join("05_caps_llo_truth.pq"))?;
117
118 /* == Setup the interlink == */
119
120 let mut measurement_types = IndexSet::new();
121 measurement_types.insert(MeasurementType::Range);
122 measurement_types.insert(MeasurementType::Doppler);
123
124 let mut stochastics = IndexMap::new();
125
126 let sa45_csac_allan_dev = 1e-11;
127
128 stochastics.insert(
129 MeasurementType::Range,
130 StochasticNoise::from_hardware_range_km(
131 sa45_csac_allan_dev,
132 10.0.seconds(),
133 link_specific::ChipRate::StandardT4B,
134 link_specific::SN0::Average,
135 ),
136 );
137
138 stochastics.insert(
139 MeasurementType::Doppler,
140 StochasticNoise::from_hardware_doppler_km_s(
141 sa45_csac_allan_dev,
142 10.0.seconds(),
143 link_specific::CarrierFreq::SBand,
144 link_specific::CN0::Average,
145 ),
146 );
147
148 let interlink = InterlinkTxSpacecraft {
149 traj: tx_traj,
150 measurement_types,
151 integration_time: None,
152 timestamp_noise_s: None,
153 ab_corr: Aberration::LT,
154 stochastic_noises: Some(stochastics),
155 };
156
157 // Devices are the transmitter, which is our NRHO vehicle.
158 let mut devices = BTreeMap::new();
159 devices.insert("NRHO Tx SC".to_string(), interlink);
160
161 let mut configs = BTreeMap::new();
162 configs.insert(
163 "NRHO Tx SC".to_string(),
164 TrkConfig::builder()
165 .strands(vec![Strand {
166 start: epoch,
167 end: nrho_final.epoch(),
168 }])
169 .build(),
170 );
171
172 let mut trk_sim =
173 TrackingArcSim::with_seed(devices.clone(), llo_traj.clone(), configs, 0).unwrap();
174 println!("{trk_sim}");
175
176 let trk_data = trk_sim.generate_measurements(almanac.clone()).unwrap();
177 println!("{trk_data}");
178
179 trk_data
180 .to_parquet_simple(out.clone().join("nrho_interlink_msr.pq"))
181 .unwrap();
182
183 // Run a truth OD where we estimate the LLO position
184 let llo_uncertainty = SpacecraftUncertainty::builder()
185 .nominal(llo_sc)
186 .x_km(1.0)
187 .y_km(1.0)
188 .z_km(1.0)
189 .vx_km_s(1e-3)
190 .vy_km_s(1e-3)
191 .vz_km_s(1e-3)
192 .build();
193
194 let mut proc_devices = devices.clone();
195
196 // Define the initial estimate, randomized, seed for reproducibility
197 let mut initial_estimate = llo_uncertainty.to_estimate_randomized(Some(0)).unwrap();
198 // Inflate the covariance -- https://github.com/nyx-space/nyx/issues/339
199 initial_estimate.covar *= 2.5;
200
201 // Increase the noise in the devices to accept more measurements.
202
203 for link in proc_devices.values_mut() {
204 for noise in &mut link.stochastic_noises.as_mut().unwrap().values_mut() {
205 *noise.white_noise.as_mut().unwrap() *= 3.0;
206 }
207 }
208
209 let init_err = initial_estimate
210 .orbital_state()
211 .ric_difference(&llo_orbit)
212 .unwrap();
213
214 println!("initial estimate:\n{initial_estimate}");
215 println!("RIC errors = {init_err}",);
216
217 let odp = InterlinkKalmanOD::new(
218 setup.clone(),
219 KalmanVariant::ReferenceUpdate,
220 Some(ResidRejectCrit::default()),
221 proc_devices,
222 almanac.clone(),
223 );
224
225 // Shrink the data to process.
226 let arc = trk_data.filter_by_offset(..2.hours());
227
228 let od_sol = odp.process_arc(initial_estimate, &arc).unwrap();
229
230 println!("{od_sol}");
231
232 od_sol
233 .to_parquet(
234 out.join("05_caps_interlink_od_sol.pq"),
235 ExportCfg::default(),
236 )
237 .unwrap();
238
239 let od_traj = od_sol.to_traj().unwrap();
240
241 od_traj
242 .ric_diff_to_parquet(
243 &llo_traj,
244 out.join("05_caps_interlink_llo_est_error.pq"),
245 ExportCfg::default(),
246 )
247 .unwrap();
248
249 let final_est = od_sol.estimates.last().unwrap();
250 assert!(final_est.within_3sigma(), "should be within 3 sigma");
251
252 println!("ESTIMATE\n{final_est:x}\n");
253 let truth = llo_traj.at(final_est.epoch()).unwrap();
254 println!("TRUTH\n{truth:x}");
255
256 let final_err = truth
257 .orbit
258 .ric_difference(&final_est.orbital_state())
259 .unwrap();
260 println!("ERROR {final_err}");
261
262 // Build the residuals versus reference plot.
263 let rvr_sol = odp
264 .process_arc(initial_estimate, &arc.resid_vs_ref_check())
265 .unwrap();
266
267 rvr_sol
268 .to_parquet(
269 out.join("05_caps_interlink_resid_v_ref.pq"),
270 ExportCfg::default(),
271 )
272 .unwrap();
273
274 let final_rvr = rvr_sol.estimates.last().unwrap();
275
276 println!("RMAG error {:.3} m", final_err.rmag_km() * 1e3);
277 println!(
278 "Pure prop error {:.3} m",
279 final_rvr
280 .orbital_state()
281 .ric_difference(&final_est.orbital_state())
282 .unwrap()
283 .rmag_km()
284 * 1e3
285 );
286
287 Ok(())
288}examples/06_lunar_orbit_determination/main.rs (line 92)
37fn main() -> Result<(), Box<dyn Error>> {
38 pel::init();
39
40 // ====================== //
41 // === ALMANAC SET UP === //
42 // ====================== //
43
44 // Dynamics models require planetary constants and ephemerides to be defined.
45 // Let's start by grabbing those by using ANISE's MetaAlmanac.
46
47 let data_folder: PathBuf = [
48 env!("CARGO_MANIFEST_DIR"),
49 "examples",
50 "06_lunar_orbit_determination",
51 ]
52 .iter()
53 .collect();
54
55 let meta = data_folder.join("metaalmanac.dhall");
56
57 // Load this ephem in the general Almanac we're using for this analysis.
58 let almanac = MetaAlmanac::new(meta.to_string_lossy().as_ref())
59 .map_err(Box::new)?
60 .process(true)
61 .map_err(Box::new)?;
62
63 // Lock the almanac (an Arc is a read only structure).
64 let almanac = Arc::new(almanac);
65
66 // Build a nominal trajectory
67 // TODO: Switch this to a sequence once the OD over a spacecraft sequence is implemented.
68
69 let epoch = Epoch::from_gregorian_utc_at_noon(2024, 2, 29);
70 let moon_j2000 = almanac.frame_info(MOON_J2000)?;
71
72 // To build the trajectory we need to provide a spacecraft template.
73 let orbiter = Spacecraft::builder()
74 .mass(Mass::from_dry_and_prop_masses(1018.0, 900.0))
75 .srp(SRPData {
76 area_m2: 3.9 * 2.7,
77 coeff_reflectivity: 0.96,
78 })
79 .orbit(Orbit::try_keplerian_altitude(
80 150.0, 0.00212, 33.6, 45.0, 45.0, 0.0, epoch, moon_j2000,
81 )?) // Setting a zero orbit here because it's just a template
82 .build();
83
84 // ========================== //
85 // === BUILD NOMINAL TRAJ === //
86 // ========================== //
87
88 // Set up the spacecraft dynamics.
89
90 // Specify that the orbital dynamics must account for the graviational pull of the Earth and the Sun.
91 // The gravity of the Moon will also be accounted for since the spaceraft in a lunar orbit.
92 let mut orbital_dyn = OrbitalDynamics::point_masses(vec![EARTH, SUN, JUPITER_BARYCENTER]);
93
94 // We want to include the spherical harmonics, so let's download the gravitational data from the Nyx Cloud.
95 // We're using the GRAIL JGGRX model.
96 let mut jggrx_meta = MetaFile {
97 uri: "http://public-data.nyxspace.com/nyx/models/Luna_jggrx_1500e_sha.tab.gz".to_string(),
98 crc32: Some(0x6bcacda8), // Specifying the CRC32 avoids redownloading it if it's cached.
99 };
100 // And let's download it if we don't have it yet.
101 jggrx_meta.process(true)?;
102
103 // Build the spherical harmonics.
104 // The harmonics must be computed in the body fixed frame.
105 // We're using the long term prediction of the Moon principal axes frame.
106 let moon_pa_frame = MOON_PA_FRAME.with_orient(31008);
107 let sph_harmonics = GravityField::from_stor(
108 almanac.frame_info(moon_pa_frame)?,
109 GravityFieldData::from_shadr(&jggrx_meta.uri, 80, 80, true)?,
110 );
111
112 // Include the spherical harmonics into the orbital dynamics.
113 orbital_dyn.accel_models.push(sph_harmonics);
114
115 // We define the solar radiation pressure, using the default solar flux and accounting only
116 // for the eclipsing caused by the Earth and Moon.
117 // Note that by default, enabling the SolarPressure model will also enable the estimation of the coefficient of reflectivity.
118 let srp_dyn = SolarPressure::new(vec![MOON_J2000], almanac.clone())?;
119
120 // Finalize setting up the dynamics, specifying the force models (orbital_dyn) separately from the
121 // acceleration models (SRP in this case). Use `from_models` to specify multiple accel models.
122 let dynamics = SpacecraftDynamics::from_model(orbital_dyn, srp_dyn);
123
124 println!("{dynamics}");
125
126 let setup = Propagator::rk89(dynamics.clone(), IntegratorOptions::default());
127
128 let truth_traj = setup
129 .with(orbiter, almanac.clone())
130 .for_duration_with_traj(Unit::Day * 2)?
131 .1;
132
133 // ==================== //
134 // === OD SIMULATOR === //
135 // ==================== //
136
137 // Load the Deep Space Network ground stations.
138 // Nyx allows you to build these at runtime but it's pretty static so we can just load them from YAML.
139 let ground_station_file = data_folder.join("dsn-network.yaml");
140 let devices = GroundStation::load_named(ground_station_file)?;
141
142 let proc_devices = devices.clone();
143
144 // Typical OD software requires that you specify your own tracking schedule or you'll have overlapping measurements.
145 // Nyx can build a tracking schedule for you based on the first station with access.
146 let configs: BTreeMap<String, TrkConfig> =
147 TrkConfig::load_named(data_folder.join("tracking-cfg.yaml"))?;
148
149 // Build the tracking arc simulation to generate a "standard measurement".
150 let mut trk = TrackingArcSim::<Spacecraft, GroundStation>::with_seed(
151 devices.clone(),
152 truth_traj.clone(),
153 configs,
154 123, // Set a seed for reproducibility
155 )?;
156
157 trk.build_schedule(almanac.clone())?;
158 let arc = trk.generate_measurements(almanac.clone())?;
159 // Save the simulated tracking data
160 arc.to_parquet_simple("./data/04_output/06_lunar_simulated_tracking.parquet")?;
161
162 // We'll note that in our case, we have continuous coverage of LRO when the vehicle is not behind the Moon.
163 println!("{arc}");
164
165 // Now that we have simulated measurements, we'll run the orbit determination.
166
167 // ===================== //
168 // === OD ESTIMATION === //
169 // ===================== //
170
171 let sc = SpacecraftUncertainty::builder()
172 .nominal(orbiter)
173 .frame(LocalFrame::RIC)
174 .x_km(0.5)
175 .y_km(0.5)
176 .z_km(0.5)
177 .vx_km_s(5e-3)
178 .vy_km_s(5e-3)
179 .vz_km_s(5e-3)
180 .build();
181
182 // Build the filter initial estimate, which we will reuse in the filter.
183 let initial_estimate = sc.to_estimate()?;
184
185 println!("== FILTER STATE ==\n{orbiter:x}\n{initial_estimate}");
186
187 // Build the SNC in the Moon J2000 frame, specified as a velocity noise over time.
188 let process_noise = ProcessNoise3D::from_velocity_km_s(
189 &[1e-13, 1e-13, 1e-13],
190 1 * Unit::Hour,
191 10 * Unit::Minute,
192 None,
193 );
194
195 println!("{process_noise}");
196
197 // We'll set up the OD process to reject measurements whose residuals are move than 3 sigmas away from what we expect.
198 let odp = SpacecraftKalmanScalarOD::new(
199 setup,
200 KalmanVariant::ReferenceUpdate,
201 Some(ResidRejectCrit::default()),
202 proc_devices,
203 almanac.clone(),
204 )
205 .with_process_noise(process_noise);
206
207 let od_sol = odp.process_arc(initial_estimate, &arc)?;
208
209 let final_est = od_sol.estimates.last().unwrap();
210
211 println!("{final_est}");
212
213 let ric_err = truth_traj
214 .at(final_est.epoch())?
215 .orbit
216 .ric_difference(&final_est.orbital_state())?;
217 println!("== RIC at end ==");
218 println!("RIC Position (m): {:.3}", ric_err.radius_km * 1e3);
219 println!("RIC Velocity (m/s): {:.3}", ric_err.velocity_km_s * 1e3);
220
221 println!(
222 "Num residuals rejected: #{}",
223 od_sol.rejected_residuals().len()
224 );
225 println!(
226 "Percentage within +/-3: {}",
227 od_sol.residual_ratio_within_threshold(3.0).unwrap()
228 );
229 println!("Ratios normal? {}", od_sol.is_normal(None).unwrap());
230
231 od_sol.to_parquet(
232 "./data/04_output/06_lunar_od_results.parquet",
233 ExportCfg::default(),
234 )?;
235
236 let od_trajectory = od_sol.to_traj()?;
237 // Build the RIC difference.
238 od_trajectory.ric_diff_to_parquet(
239 &truth_traj,
240 "./data/04_output/06_lunar_od_truth_error.parquet",
241 ExportCfg::default(),
242 )?;
243
244 Ok(())
245}examples/03_geo_analysis/drift.rs (line 76)
26fn main() -> Result<(), Box<dyn Error>> {
27 pel::init();
28 // Dynamics models require planetary constants and ephemerides to be defined.
29 // Let's start by grabbing those by using ANISE's latest MetaAlmanac.
30 // This will automatically download the DE440s planetary ephemeris,
31 // the daily-updated Earth Orientation Parameters, the high fidelity Moon orientation
32 // parameters (for the Moon Mean Earth and Moon Principal Axes frames), and the PCK11
33 // planetary constants kernels.
34 // For details, refer to https://github.com/nyx-space/anise/blob/master/data/latest.dhall.
35 // Note that we place the Almanac into an Arc so we can clone it cheaply and provide read-only
36 // references to many functions.
37 let almanac = Arc::new(MetaAlmanac::latest().map_err(Box::new)?);
38 // Define the orbit epoch
39 let epoch = Epoch::from_gregorian_utc_hms(2024, 2, 29, 12, 13, 14);
40
41 // Define the orbit.
42 // First we need to fetch the Earth J2000 from information from the Almanac.
43 // This allows the frame to include the gravitational parameters and the shape of the Earth,
44 // defined as a tri-axial ellipoid. Note that this shape can be changed manually or in the Almanac
45 // by loading a different set of planetary constants.
46 let earth_j2000 = almanac.frame_info(EARTH_J2000)?;
47
48 // Placing this GEO bird just above Colorado.
49 // In theory, the eccentricity is zero, but in practice, it's about 1e-5 to 1e-6 at best.
50 let orbit = Orbit::try_keplerian(42164.0, 1e-5, 0., 163.0, 75.0, 0.0, epoch, earth_j2000)?;
51 // Print in in Keplerian form.
52 println!("{orbit:x}");
53
54 let state_bf = almanac.transform_to(orbit, IAU_EARTH_FRAME, None)?;
55 let (orig_lat_deg, orig_long_deg, orig_alt_km) = state_bf.latlongalt()?;
56
57 // Nyx is used for high fidelity propagation, not Keplerian propagation as above.
58 // Nyx only propagates Spacecraft at the moment, which allows it to account for acceleration
59 // models such as solar radiation pressure.
60
61 // Let's build a cubesat sized spacecraft, with an SRP area of 10 cm^2 and a mass of 9.6 kg.
62 let sc = Spacecraft::builder()
63 .orbit(orbit)
64 .mass(Mass::from_dry_mass(9.60))
65 .srp(SRPData {
66 area_m2: 10e-4,
67 coeff_reflectivity: 1.1,
68 })
69 .build();
70 println!("{sc:x}");
71
72 // Set up the spacecraft dynamics.
73
74 // Specify that the orbital dynamics must account for the graviational pull of the Moon and the Sun.
75 // The gravity of the Earth will also be accounted for since the spaceraft in an Earth orbit.
76 let mut orbital_dyn = OrbitalDynamics::point_masses(vec![MOON, SUN]);
77
78 // We want to include the spherical harmonics, so let's download the gravitational data from the Nyx Cloud.
79 // We're using the JGM3 model here, which is the default in GMAT.
80 let mut jgm3_meta = MetaFile {
81 uri: "http://public-data.nyxspace.com/nyx/models/JGM3.cof.gz".to_string(),
82 crc32: Some(0xF446F027), // Specifying the CRC32 avoids redownloading it if it's cached.
83 };
84 // And let's download it if we don't have it yet.
85 jgm3_meta.process(true)?;
86
87 // Build the spherical harmonics.
88 // The harmonics must be computed in the body fixed frame.
89 // We're using the long term prediction of the Earth centered Earth fixed frame, IAU Earth.
90 let harmonics_21x21 = GravityField::from_stor(
91 almanac.frame_info(IAU_EARTH_FRAME)?,
92 GravityFieldData::from_cof(&jgm3_meta.uri, 21, 21, true).unwrap(),
93 );
94
95 // Include the spherical harmonics into the orbital dynamics.
96 orbital_dyn.accel_models.push(harmonics_21x21);
97
98 // We define the solar radiation pressure, using the default solar flux and accounting only
99 // for the eclipsing caused by the Earth and Moon.
100 let srp_dyn = SolarPressure::new(vec![EARTH_J2000, MOON_J2000], almanac.clone())?;
101
102 // Finalize setting up the dynamics, specifying the force models (orbital_dyn) separately from the
103 // acceleration models (SRP in this case). Use `from_models` to specify multiple accel models.
104 let dynamics = SpacecraftDynamics::from_model(orbital_dyn, srp_dyn);
105
106 println!("{dynamics}");
107
108 // Finally, let's propagate this orbit to the same epoch as above.
109 // The first returned value is the spacecraft state at the final epoch.
110 // The second value is the full trajectory where the step size is variable step used by the propagator.
111 let (future_sc, trajectory) = Propagator::default(dynamics)
112 .with(sc, almanac.clone())
113 .until_epoch_with_traj(epoch + Unit::Century * 0.03)?;
114
115 println!("=== High fidelity propagation ===");
116 println!(
117 "SMA changed by {:.3} km",
118 orbit.sma_km()? - future_sc.orbit.sma_km()?
119 );
120 println!(
121 "ECC changed by {:.6}",
122 orbit.ecc()? - future_sc.orbit.ecc()?
123 );
124 println!(
125 "INC changed by {:.3e} deg",
126 orbit.inc_deg()? - future_sc.orbit.inc_deg()?
127 );
128 println!(
129 "RAAN changed by {:.3} deg",
130 orbit.raan_deg()? - future_sc.orbit.raan_deg()?
131 );
132 println!(
133 "AOP changed by {:.3} deg",
134 orbit.aop_deg()? - future_sc.orbit.aop_deg()?
135 );
136 println!(
137 "TA changed by {:.3} deg",
138 orbit.ta_deg()? - future_sc.orbit.ta_deg()?
139 );
140
141 // We also have access to the full trajectory throughout the propagation.
142 println!("{trajectory}");
143
144 println!("Spacecraft params after 3 years without active control:\n{future_sc:x}");
145
146 // With the trajectory, let's build a few data products.
147
148 // 1. Export the trajectory as a parquet file, which includes the Keplerian orbital elements.
149
150 let analysis_step = Unit::Minute * 5;
151
152 trajectory.to_parquet(
153 "./03_geo_hf_prop.parquet",
154 ExportCfg::builder().step(analysis_step).build(),
155 )?;
156
157 // 2. Compute the latitude, longitude, and altitude throughout the trajectory by rotating the spacecraft position into the Earth body fixed frame.
158
159 // We iterate over the trajectory, grabbing a state every two minutes.
160 let mut offset_s = vec![];
161 let mut epoch_str = vec![];
162 let mut longitude_deg = vec![];
163 let mut latitude_deg = vec![];
164 let mut altitude_km = vec![];
165
166 for state in trajectory.every(analysis_step) {
167 // Convert the GEO bird state into the body fixed frame, and keep track of its latitude, longitude, and altitude.
168 // These define the GEO stationkeeping box.
169
170 let this_epoch = state.epoch();
171
172 offset_s.push((this_epoch - orbit.epoch).to_seconds());
173 epoch_str.push(this_epoch.to_isoformat());
174
175 let state_bf = almanac.transform_to(state.orbit, IAU_EARTH_FRAME, None)?;
176 let (lat_deg, long_deg, alt_km) = state_bf.latlongalt()?;
177 longitude_deg.push(long_deg);
178 latitude_deg.push(lat_deg);
179 altitude_km.push(alt_km);
180 }
181
182 println!(
183 "Longitude changed by {:.3} deg -- Box is 0.1 deg E-W",
184 orig_long_deg - longitude_deg.last().unwrap()
185 );
186
187 println!(
188 "Latitude changed by {:.3} deg -- Box is 0.05 deg N-S",
189 orig_lat_deg - latitude_deg.last().unwrap()
190 );
191
192 println!(
193 "Altitude changed by {:.3} km -- Box is 30 km",
194 orig_alt_km - altitude_km.last().unwrap()
195 );
196
197 // Build the station keeping data frame.
198 let mut sk_df = df!(
199 "Offset (s)" => offset_s.clone(),
200 "Epoch (UTC)" => epoch_str.clone(),
201 "Longitude E-W (deg)" => longitude_deg,
202 "Latitude N-S (deg)" => latitude_deg,
203 "Altitude (km)" => altitude_km,
204
205 )?;
206
207 // Create a file to write the Parquet to
208 let file = File::create("./03_geo_lla.parquet").expect("Could not create file");
209
210 // Create a ParquetWriter and write the DataFrame to the file
211 ParquetWriter::new(file).finish(&mut sk_df)?;
212
213 Ok(())
214}Additional examples can be found in:
Sourcepub fn two_body() -> Self
pub fn two_body() -> Self
Initializes a OrbitalDynamics which does not simulate the gravity pull of other celestial objects but the primary one.
Sourcepub fn new(accel_models: Vec<Arc<dyn AccelModel + Sync>>) -> Self
pub fn new(accel_models: Vec<Arc<dyn AccelModel + Sync>>) -> Self
Initialize orbital dynamics with a list of acceleration models
Sourcepub fn from_model(accel_model: Arc<dyn AccelModel + Sync>) -> Self
pub fn from_model(accel_model: Arc<dyn AccelModel + Sync>) -> Self
Initialize new orbital mechanics with the provided model. Note: Orbital dynamics always include two body dynamics, these cannot be turned off.
Trait Implementations§
Source§impl Clone for OrbitalDynamics
impl Clone for OrbitalDynamics
Source§fn clone(&self) -> OrbitalDynamics
fn clone(&self) -> OrbitalDynamics
Returns a duplicate of the value. Read more
1.0.0 · Source§fn clone_from(&mut self, source: &Self)
fn clone_from(&mut self, source: &Self)
Performs copy-assignment from
source. Read moreAuto Trait Implementations§
impl Freeze for OrbitalDynamics
impl !RefUnwindSafe for OrbitalDynamics
impl Send for OrbitalDynamics
impl Sync for OrbitalDynamics
impl Unpin for OrbitalDynamics
impl UnsafeUnpin for OrbitalDynamics
impl !UnwindSafe for OrbitalDynamics
Blanket Implementations§
Source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
Source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Mutably borrows from an owned value. Read more
Source§impl<T> CloneToUninit for Twhere
T: Clone,
impl<T> CloneToUninit for Twhere
T: Clone,
§impl<T> Instrument for T
impl<T> Instrument for T
§fn instrument(self, span: Span) -> Instrumented<Self>
fn instrument(self, span: Span) -> Instrumented<Self>
§fn in_current_span(self) -> Instrumented<Self>
fn in_current_span(self) -> Instrumented<Self>
Source§impl<T> IntoEither for T
impl<T> IntoEither for T
Source§fn into_either(self, into_left: bool) -> Either<Self, Self>
fn into_either(self, into_left: bool) -> Either<Self, Self>
Converts
self into a Left variant of Either<Self, Self>
if into_left is true.
Converts self into a Right variant of Either<Self, Self>
otherwise. Read moreSource§fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
Converts
self into a Left variant of Either<Self, Self>
if into_left(&self) returns true.
Converts self into a Right variant of Either<Self, Self>
otherwise. Read more§impl<T> Pointable for T
impl<T> Pointable for T
§impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
§fn to_subset(&self) -> Option<SS>
fn to_subset(&self) -> Option<SS>
The inverse inclusion map: attempts to construct
self from the equivalent element of its
superset. Read more§fn is_in_subset(&self) -> bool
fn is_in_subset(&self) -> bool
Checks if
self is actually part of its subset T (and can be converted to it).§fn to_subset_unchecked(&self) -> SS
fn to_subset_unchecked(&self) -> SS
Use with care! Same as
self.to_subset but without any property checks. Always succeeds.§fn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
The inclusion map: converts
self to the equivalent element of its superset.