02_jwst/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::{JUPITER_BARYCENTER, MOON, SUN},
10 frames::{EARTH_J2000, MOON_J2000},
11 },
12};
13use hifitime::{TimeUnits, Unit};
14use nyx::{
15 cosmic::{eclipse::EclipseLocator, Frame, Mass, MetaAlmanac, SRPData},
16 dynamics::{guidance::LocalFrame, OrbitalDynamics, SolarPressure, SpacecraftDynamics},
17 io::ExportCfg,
18 mc::MonteCarlo,
19 od::{msr::TrackingDataArc, prelude::KF, process::SpacecraftUncertainty, SpacecraftODProcess},
20 propagators::Propagator,
21 Spacecraft, State,
22};
23
24use std::{collections::BTreeMap, error::Error, sync::Arc};
25
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 OD process, and predicting for the analysis duration.
114
115 let ckf = KF::no_snc(jwst_estimate);
116
117 // Build the propagation instance for the OD process.
118 let prop = setup.with(jwst.with_stm(), almanac.clone());
119 let mut odp = SpacecraftODProcess::ckf(prop, ckf, BTreeMap::new(), None, almanac.clone());
120
121 // Define the prediction step, i.e. how often we want to know the covariance.
122 let step = 1_i64.minutes();
123 // Finally, predict, and export the trajectory with covariance to a parquet file.
124 odp.predict_for(step, prediction_duration)?;
125 odp.to_parquet(
126 &TrackingDataArc::default(),
127 "./02_jwst_covar_map.parquet",
128 ExportCfg::default(),
129 )?;
130
131 // === Monte Carlo framework ===
132 // Nyx comes with a complete multi-threaded Monte Carlo frame. It's blazing fast.
133
134 let my_mc = MonteCarlo::new(
135 jwst, // Nominal state
136 jwst_estimate.to_random_variable()?,
137 "02_jwst".to_string(), // Scenario name
138 None, // No specific seed specified, so one will be drawn from the computer's entropy.
139 );
140
141 let num_runs = 5_000;
142 let rslts = my_mc.run_until_epoch(
143 setup,
144 almanac.clone(),
145 jwst.epoch() + prediction_duration,
146 num_runs,
147 );
148
149 assert_eq!(rslts.runs.len(), num_runs);
150 // Finally, export these results, computing the eclipse percentage for all of these results.
151
152 // For all of the resulting trajectories, we'll want to compute the percentage of penumbra and umbra.
153 let eclipse_loc = EclipseLocator::cislunar(almanac.clone());
154 let umbra_event = eclipse_loc.to_umbra_event();
155 let penumbra_event = eclipse_loc.to_penumbra_event();
156
157 rslts.to_parquet(
158 "02_jwst_monte_carlo.parquet",
159 Some(vec![&umbra_event, &penumbra_event]),
160 ExportCfg::default(),
161 almanac,
162 )?;
163
164 Ok(())
165}