Struct OrbitalDynamics

Source
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

Source

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 65)
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_from_uid(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(StateParameter::SMA, 42_164.0, 5.0), // 5 km
58        Objective::within_tolerance(StateParameter::Eccentricity, 0.001, 5e-5),
59        Objective::within_tolerance(StateParameter::Inclination, 0.05, 1e-2),
60    ];
61
62    let ruggiero_ctrl = Ruggiero::from_max_eclipse(objectives, sc, 0.2)?;
63    println!("{ruggiero_ctrl}");
64
65    let mut orbital_dyn = OrbitalDynamics::point_masses(vec![MOON, SUN]);
66
67    let mut jgm3_meta = MetaFile {
68        uri: "http://public-data.nyxspace.com/nyx/models/JGM3.cof.gz".to_string(),
69        crc32: Some(0xF446F027), // Specifying the CRC32 avoids redownloading it if it's cached.
70    };
71    jgm3_meta.process(true)?;
72
73    let harmonics = Harmonics::from_stor(
74        almanac.frame_from_uid(IAU_EARTH_FRAME)?,
75        HarmonicsMem::from_cof(&jgm3_meta.uri, 8, 8, true)?,
76    );
77    orbital_dyn.accel_models.push(harmonics);
78
79    let srp_dyn = SolarPressure::default(EARTH_J2000, almanac.clone())?;
80    let sc_dynamics = SpacecraftDynamics::from_model(orbital_dyn, srp_dyn)
81        .with_guidance_law(ruggiero_ctrl.clone());
82
83    println!("{sc_dynamics}");
84
85    // Finally, let's use the Monte Carlo framework built into Nyx to propagate spacecraft.
86
87    // Let's start by defining the dispersion.
88    // 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.
89    // Note that additional validation on the MVN is in progress -- https://github.com/nyx-space/nyx/issues/339.
90    let mc_rv = MvnSpacecraft::new(
91        sc,
92        vec![StateDispersion::zero_mean(StateParameter::SMA, 3.0)],
93    )?;
94
95    let my_mc = MonteCarlo::new(
96        sc, // Nominal state
97        mc_rv,
98        "03_geo_sk".to_string(), // Scenario name
99        None, // No specific seed specified, so one will be drawn from the computer's entropy.
100    );
101
102    // Build the propagator setup.
103    let setup = Propagator::rk89(
104        sc_dynamics.clone(),
105        IntegratorOptions::builder()
106            .min_step(10.0_f64.seconds())
107            .error_ctrl(ErrorControl::RSSCartesianStep)
108            .build(),
109    );
110
111    let num_runs = 25;
112    let rslts = my_mc.run_until_epoch(setup, almanac.clone(), sc.epoch() + prop_time, num_runs);
113
114    assert_eq!(rslts.runs.len(), num_runs);
115
116    // For all of the resulting trajectories, we'll want to compute the percentage of penumbra and umbra.
117
118    rslts.to_parquet(
119        "03_geo_sk.parquet",
120        Some(vec![
121            &EclipseLocator::cislunar(almanac.clone()).to_penumbra_event()
122        ]),
123        ExportCfg::default(),
124        almanac,
125    )?;
126
127    Ok(())
128}
More examples
Hide additional examples
examples/03_geo_analysis/raise.rs (line 83)
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_from_uid(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(StateParameter::SMA, 42_165.0, 20.0),
69        Objective::within_tolerance(StateParameter::Eccentricity, 0.001, 5e-5),
70        Objective::within_tolerance(StateParameter::Inclination, 0.05, 1e-2),
71    ];
72
73    // Ensure that we only thrust if we have more than 20% illumination.
74    let ruggiero_ctrl = Ruggiero::from_max_eclipse(objectives, sc, 0.2).unwrap();
75    println!("{ruggiero_ctrl}");
76
77    // Define the high fidelity dynamics
78
79    // Set up the spacecraft dynamics.
80
81    // Specify that the orbital dynamics must account for the graviational pull of the Moon and the Sun.
82    // The gravity of the Earth will also be accounted for since the spaceraft in an Earth orbit.
83    let mut orbital_dyn = OrbitalDynamics::point_masses(vec![MOON, SUN]);
84
85    // We want to include the spherical harmonics, so let's download the gravitational data from the Nyx Cloud.
86    // We're using the JGM3 model here, which is the default in GMAT.
87    let mut jgm3_meta = MetaFile {
88        uri: "http://public-data.nyxspace.com/nyx/models/JGM3.cof.gz".to_string(),
89        crc32: Some(0xF446F027), // Specifying the CRC32 avoids redownloading it if it's cached.
90    };
91    // And let's download it if we don't have it yet.
92    jgm3_meta.process(true)?;
93
94    // Build the spherical harmonics.
95    // The harmonics must be computed in the body fixed frame.
96    // We're using the long term prediction of the Earth centered Earth fixed frame, IAU Earth.
97    let harmonics = Harmonics::from_stor(
98        almanac.frame_from_uid(IAU_EARTH_FRAME)?,
99        HarmonicsMem::from_cof(&jgm3_meta.uri, 8, 8, true).unwrap(),
100    );
101
102    // Include the spherical harmonics into the orbital dynamics.
103    orbital_dyn.accel_models.push(harmonics);
104
105    // We define the solar radiation pressure, using the default solar flux and accounting only
106    // for the eclipsing caused by the Earth.
107    let srp_dyn = SolarPressure::default(EARTH_J2000, almanac.clone())?;
108
109    // Finalize setting up the dynamics, specifying the force models (orbital_dyn) separately from the
110    // acceleration models (SRP in this case). Use `from_models` to specify multiple accel models.
111    let sc_dynamics = SpacecraftDynamics::from_model(orbital_dyn, srp_dyn)
112        .with_guidance_law(ruggiero_ctrl.clone());
113
114    println!("{:x}", orbit);
115
116    // We specify a minimum step in the propagator because the Ruggiero control would otherwise drive this step very low.
117    let (final_state, traj) = Propagator::rk89(
118        sc_dynamics.clone(),
119        IntegratorOptions::builder()
120            .min_step(10.0_f64.seconds())
121            .error_ctrl(ErrorControl::RSSCartesianStep)
122            .build(),
123    )
124    .with(sc, almanac.clone())
125    .for_duration_with_traj(prop_time)?;
126
127    let prop_usage = sc.mass.prop_mass_kg - final_state.mass.prop_mass_kg;
128    println!("{:x}", final_state.orbit);
129    println!("prop usage: {:.3} kg", prop_usage);
130
131    // Finally, export the results for analysis, including the penumbra percentage throughout the orbit raise.
132    traj.to_parquet(
133        "./03_geo_raise.parquet",
134        Some(vec![
135            &EclipseLocator::cislunar(almanac.clone()).to_penumbra_event()
136        ]),
137        ExportCfg::default(),
138        almanac,
139    )?;
140
141    for status_line in ruggiero_ctrl.status(&final_state) {
142        println!("{status_line}");
143    }
144
145    ruggiero_ctrl
146        .achieved(&final_state)
147        .expect("objective not achieved");
148
149    Ok(())
150}
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 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}
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_from_uid(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 = Harmonics::from_stor(
91        almanac.frame_from_uid(IAU_EARTH_FRAME)?,
92        HarmonicsMem::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        Some(vec![
155            &EclipseLocator::cislunar(almanac.clone()).to_penumbra_event()
156        ]),
157        ExportCfg::builder().step(analysis_step).build(),
158        almanac.clone(),
159    )?;
160
161    // 2. Compute the latitude, longitude, and altitude throughout the trajectory by rotating the spacecraft position into the Earth body fixed frame.
162
163    // We iterate over the trajectory, grabbing a state every two minutes.
164    let mut offset_s = vec![];
165    let mut epoch_str = vec![];
166    let mut longitude_deg = vec![];
167    let mut latitude_deg = vec![];
168    let mut altitude_km = vec![];
169
170    for state in trajectory.every(analysis_step) {
171        // Convert the GEO bird state into the body fixed frame, and keep track of its latitude, longitude, and altitude.
172        // These define the GEO stationkeeping box.
173
174        let this_epoch = state.epoch();
175
176        offset_s.push((this_epoch - orbit.epoch).to_seconds());
177        epoch_str.push(this_epoch.to_isoformat());
178
179        let state_bf = almanac.transform_to(state.orbit, IAU_EARTH_FRAME, None)?;
180        let (lat_deg, long_deg, alt_km) = state_bf.latlongalt()?;
181        longitude_deg.push(long_deg);
182        latitude_deg.push(lat_deg);
183        altitude_km.push(alt_km);
184    }
185
186    println!(
187        "Longitude changed by {:.3} deg -- Box is 0.1 deg E-W",
188        orig_long_deg - longitude_deg.last().unwrap()
189    );
190
191    println!(
192        "Latitude changed by {:.3} deg -- Box is 0.05 deg N-S",
193        orig_lat_deg - latitude_deg.last().unwrap()
194    );
195
196    println!(
197        "Altitude changed by {:.3} km -- Box is 30 km",
198        orig_alt_km - altitude_km.last().unwrap()
199    );
200
201    // Build the station keeping data frame.
202    let mut sk_df = df!(
203        "Offset (s)" => offset_s.clone(),
204        "Epoch (UTC)" => epoch_str.clone(),
205        "Longitude E-W (deg)" => longitude_deg,
206        "Latitude N-S (deg)" => latitude_deg,
207        "Altitude (km)" => altitude_km,
208
209    )?;
210
211    // Create a file to write the Parquet to
212    let file = File::create("./03_geo_lla.parquet").expect("Could not create file");
213
214    // Create a ParquetWriter and write the DataFrame to the file
215    ParquetWriter::new(file).finish(&mut sk_df)?;
216
217    Ok(())
218}
examples/04_lro_od/main.rs (line 112)
33fn main() -> Result<(), Box<dyn Error>> {
34    pel::init();
35
36    // ====================== //
37    // === ALMANAC SET UP === //
38    // ====================== //
39
40    // Dynamics models require planetary constants and ephemerides to be defined.
41    // Let's start by grabbing those by using ANISE's MetaAlmanac.
42
43    let data_folder: PathBuf = [env!("CARGO_MANIFEST_DIR"), "examples", "04_lro_od"]
44        .iter()
45        .collect();
46
47    let meta = data_folder.join("lro-dynamics.dhall");
48
49    // Load this ephem in the general Almanac we're using for this analysis.
50    let mut almanac = MetaAlmanac::new(meta.to_string_lossy().to_string())
51        .map_err(Box::new)?
52        .process(true)
53        .map_err(Box::new)?;
54
55    let mut moon_pc = almanac.planetary_data.get_by_id(MOON)?;
56    moon_pc.mu_km3_s2 = 4902.74987;
57    almanac.planetary_data.set_by_id(MOON, moon_pc)?;
58
59    let mut earth_pc = almanac.planetary_data.get_by_id(EARTH)?;
60    earth_pc.mu_km3_s2 = 398600.436;
61    almanac.planetary_data.set_by_id(EARTH, earth_pc)?;
62
63    // Save this new kernel for reuse.
64    // In an operational context, this would be part of the "Lock" process, and should not change throughout the mission.
65    almanac
66        .planetary_data
67        .save_as(&data_folder.join("lro-specific.pca"), true)?;
68
69    // Lock the almanac (an Arc is a read only structure).
70    let almanac = Arc::new(almanac);
71
72    // Orbit determination requires a Trajectory structure, which can be saved as parquet file.
73    // In our case, the trajectory comes from the BSP file, so we need to build a Trajectory from the almanac directly.
74    // To query the Almanac, we need to build the LRO frame in the J2000 orientation in our case.
75    // Inspecting the LRO BSP in the ANISE GUI shows us that NASA has assigned ID -85 to LRO.
76    let lro_frame = Frame::from_ephem_j2000(-85);
77
78    // To build the trajectory we need to provide a spacecraft template.
79    let sc_template = Spacecraft::builder()
80        .mass(Mass::from_dry_and_prop_masses(1018.0, 900.0)) // Launch masses
81        .srp(SRPData {
82            // SRP configuration is arbitrary, but we will be estimating it anyway.
83            area_m2: 3.9 * 2.7,
84            coeff_reflectivity: 0.96,
85        })
86        .orbit(Orbit::zero(MOON_J2000)) // Setting a zero orbit here because it's just a template
87        .build();
88    // Now we can build the trajectory from the BSP file.
89    // We'll arbitrarily set the tracking arc to 24 hours with a five second time step.
90    let traj_as_flown = Traj::from_bsp(
91        lro_frame,
92        MOON_J2000,
93        almanac.clone(),
94        sc_template,
95        5.seconds(),
96        Some(Epoch::from_str("2024-01-01 00:00:00 UTC")?),
97        Some(Epoch::from_str("2024-01-02 00:00:00 UTC")?),
98        Aberration::LT,
99        Some("LRO".to_string()),
100    )?;
101
102    println!("{traj_as_flown}");
103
104    // ====================== //
105    // === MODEL MATCHING === //
106    // ====================== //
107
108    // Set up the spacecraft dynamics.
109
110    // Specify that the orbital dynamics must account for the graviational pull of the Earth and the Sun.
111    // The gravity of the Moon will also be accounted for since the spaceraft in a lunar orbit.
112    let mut orbital_dyn = OrbitalDynamics::point_masses(vec![EARTH, SUN, JUPITER_BARYCENTER]);
113
114    // We want to include the spherical harmonics, so let's download the gravitational data from the Nyx Cloud.
115    // We're using the GRAIL JGGRX model.
116    let mut jggrx_meta = MetaFile {
117        uri: "http://public-data.nyxspace.com/nyx/models/Luna_jggrx_1500e_sha.tab.gz".to_string(),
118        crc32: Some(0x6bcacda8), // Specifying the CRC32 avoids redownloading it if it's cached.
119    };
120    // And let's download it if we don't have it yet.
121    jggrx_meta.process(true)?;
122
123    // Build the spherical harmonics.
124    // The harmonics must be computed in the body fixed frame.
125    // We're using the long term prediction of the Moon principal axes frame.
126    let moon_pa_frame = MOON_PA_FRAME.with_orient(31008);
127    let sph_harmonics = Harmonics::from_stor(
128        almanac.frame_from_uid(moon_pa_frame)?,
129        HarmonicsMem::from_shadr(&jggrx_meta.uri, 80, 80, true)?,
130    );
131
132    // Include the spherical harmonics into the orbital dynamics.
133    orbital_dyn.accel_models.push(sph_harmonics);
134
135    // We define the solar radiation pressure, using the default solar flux and accounting only
136    // for the eclipsing caused by the Earth and Moon.
137    // Note that by default, enabling the SolarPressure model will also enable the estimation of the coefficient of reflectivity.
138    let srp_dyn = SolarPressure::new(vec![EARTH_J2000, MOON_J2000], almanac.clone())?;
139
140    // Finalize setting up the dynamics, specifying the force models (orbital_dyn) separately from the
141    // acceleration models (SRP in this case). Use `from_models` to specify multiple accel models.
142    let dynamics = SpacecraftDynamics::from_model(orbital_dyn, srp_dyn);
143
144    println!("{dynamics}");
145
146    // Now we can build the propagator.
147    let setup = Propagator::default_dp78(dynamics.clone());
148
149    // For reference, let's build the trajectory with Nyx's models from that LRO state.
150    let (sim_final, traj_as_sim) = setup
151        .with(*traj_as_flown.first(), almanac.clone())
152        .until_epoch_with_traj(traj_as_flown.last().epoch())?;
153
154    println!("SIM INIT:  {:x}", traj_as_flown.first());
155    println!("SIM FINAL: {sim_final:x}");
156    // Compute RIC difference between SIM and LRO ephem
157    let sim_lro_delta = sim_final
158        .orbit
159        .ric_difference(&traj_as_flown.last().orbit)?;
160    println!("{traj_as_sim}");
161    println!(
162        "SIM v LRO - RIC Position (m): {:.3}",
163        sim_lro_delta.radius_km * 1e3
164    );
165    println!(
166        "SIM v LRO - RIC Velocity (m/s): {:.3}",
167        sim_lro_delta.velocity_km_s * 1e3
168    );
169
170    traj_as_sim.ric_diff_to_parquet(
171        &traj_as_flown,
172        "./04_lro_sim_truth_error.parquet",
173        ExportCfg::default(),
174    )?;
175
176    // ==================== //
177    // === OD SIMULATOR === //
178    // ==================== //
179
180    // After quite some time trying to exactly match the model, we still end up with an oscillatory difference on the order of 150 meters between the propagated state
181    // and the truth LRO state.
182
183    // Therefore, we will actually run an estimation from a dispersed LRO state.
184    // The sc_seed is the true LRO state from the BSP.
185    let sc_seed = *traj_as_flown.first();
186
187    // Load the Deep Space Network ground stations.
188    // Nyx allows you to build these at runtime but it's pretty static so we can just load them from YAML.
189    let ground_station_file: PathBuf = [
190        env!("CARGO_MANIFEST_DIR"),
191        "examples",
192        "04_lro_od",
193        "dsn-network.yaml",
194    ]
195    .iter()
196    .collect();
197
198    let devices = GroundStation::load_named(ground_station_file)?;
199
200    // Typical OD software requires that you specify your own tracking schedule or you'll have overlapping measurements.
201    // Nyx can build a tracking schedule for you based on the first station with access.
202    let trkconfg_yaml: PathBuf = [
203        env!("CARGO_MANIFEST_DIR"),
204        "examples",
205        "04_lro_od",
206        "tracking-cfg.yaml",
207    ]
208    .iter()
209    .collect();
210
211    let configs: BTreeMap<String, TrkConfig> = TrkConfig::load_named(trkconfg_yaml)?;
212
213    // Build the tracking arc simulation to generate a "standard measurement".
214    let mut trk = TrackingArcSim::<Spacecraft, GroundStation>::new(
215        devices.clone(),
216        traj_as_flown.clone(),
217        configs,
218    )?;
219
220    trk.build_schedule(almanac.clone())?;
221    let arc = trk.generate_measurements(almanac.clone())?;
222    // Save the simulated tracking data
223    arc.to_parquet_simple("./04_lro_simulated_tracking.parquet")?;
224
225    // We'll note that in our case, we have continuous coverage of LRO when the vehicle is not behind the Moon.
226    println!("{arc}");
227
228    // Now that we have simulated measurements, we'll run the orbit determination.
229
230    // ===================== //
231    // === OD ESTIMATION === //
232    // ===================== //
233
234    let sc = SpacecraftUncertainty::builder()
235        .nominal(sc_seed)
236        .frame(LocalFrame::RIC)
237        .x_km(0.5)
238        .y_km(0.5)
239        .z_km(0.5)
240        .vx_km_s(5e-3)
241        .vy_km_s(5e-3)
242        .vz_km_s(5e-3)
243        .build();
244
245    // Build the filter initial estimate, which we will reuse in the filter.
246    let initial_estimate = sc.to_estimate()?;
247
248    println!("== FILTER STATE ==\n{sc_seed:x}\n{initial_estimate}");
249
250    let kf = KF::new(
251        // Increase the initial covariance to account for larger deviation.
252        initial_estimate,
253        // Until https://github.com/nyx-space/nyx/issues/351, we need to specify the SNC in the acceleration of the Moon J2000 frame.
254        SNC3::from_diagonal(10 * Unit::Minute, &[1e-12, 1e-12, 1e-12]),
255    );
256
257    // We'll set up the OD process to reject measurements whose residuals are move than 3 sigmas away from what we expect.
258    let mut odp = SpacecraftODProcess::ckf(
259        setup.with(initial_estimate.state().with_stm(), almanac.clone()),
260        kf,
261        devices,
262        Some(ResidRejectCrit::default()),
263        almanac.clone(),
264    );
265
266    odp.process_arc(&arc)?;
267
268    let ric_err = traj_as_flown
269        .at(odp.estimates.last().unwrap().epoch())?
270        .orbit
271        .ric_difference(&odp.estimates.last().unwrap().orbital_state())?;
272    println!("== RIC at end ==");
273    println!("RIC Position (m): {}", ric_err.radius_km * 1e3);
274    println!("RIC Velocity (m/s): {}", ric_err.velocity_km_s * 1e3);
275
276    odp.to_parquet(&arc, "./04_lro_od_results.parquet", ExportCfg::default())?;
277
278    // In our case, we have the truth trajectory from NASA.
279    // So we can compute the RIC state difference between the real LRO ephem and what we've just estimated.
280    // Export the OD trajectory first.
281    let od_trajectory = odp.to_traj()?;
282    // Build the RIC difference.
283    od_trajectory.ric_diff_to_parquet(
284        &traj_as_flown,
285        "./04_lro_od_truth_error.parquet",
286        ExportCfg::default(),
287    )?;
288
289    Ok(())
290}
examples/01_orbit_prop/main.rs (line 109)
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}
Source

pub fn two_body() -> Self

Initializes a OrbitalDynamics which does not simulate the gravity pull of other celestial objects but the primary one.

Source

pub fn new(accel_models: Vec<Arc<dyn AccelModel + Sync>>) -> Self

Initialize orbital dynamics with a list of acceleration models

Source

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.

Source§

impl OrbitalDynamics

Source

pub fn dual_eom( &self, _delta_t_s: f64, osc: &Orbit, almanac: Arc<Almanac>, ) -> Result<(Vector6<f64>, Matrix6<f64>), DynamicsError>

Trait Implementations§

Source§

impl Clone for OrbitalDynamics

Source§

fn clone(&self) -> OrbitalDynamics

Returns a copy of the value. Read more
1.0.0 · Source§

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
Source§

impl Display for OrbitalDynamics

Source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more

Auto Trait Implementations§

Blanket Implementations§

Source§

impl<T> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<T> CloneToUninit for T
where T: Clone,

Source§

unsafe fn clone_to_uninit(&self, dst: *mut u8)

🔬This is a nightly-only experimental API. (clone_to_uninit)
Performs copy-assignment from self to dst. Read more
Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

§

impl<T> Instrument for T

§

fn instrument(self, span: Span) -> Instrumented<Self>

Instruments this type with the provided [Span], returning an Instrumented wrapper. Read more
§

fn in_current_span(self) -> Instrumented<Self>

Instruments this type with the current Span, returning an Instrumented wrapper. Read more
Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

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

Source§

impl<T> IntoEither for T

Source§

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 more
Source§

fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
where F: FnOnce(&Self) -> bool,

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

§

const ALIGN: usize

The alignment of pointer.
§

type Init = T

The type for initializers.
§

unsafe fn init(init: <T as Pointable>::Init) -> usize

Initializes a with the given initializer. Read more
§

unsafe fn deref<'a>(ptr: usize) -> &'a T

Dereferences the given pointer. Read more
§

unsafe fn deref_mut<'a>(ptr: usize) -> &'a mut T

Mutably dereferences the given pointer. Read more
§

unsafe fn drop(ptr: usize)

Drops the object pointed to by the given pointer. Read more
Source§

impl<T> Same for T

Source§

type Output = T

Should always be Self
§

impl<SS, SP> SupersetOf<SS> for SP
where SS: SubsetOf<SP>,

§

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

Checks if self is actually part of its subset T (and can be converted to it).
§

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

The inclusion map: converts self to the equivalent element of its superset.
Source§

impl<T> ToOwned for T
where T: Clone,

Source§

type Owned = T

The resulting type after obtaining ownership.
Source§

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
Source§

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
Source§

impl<T> ToString for T
where T: Display + ?Sized,

Source§

fn to_string(&self) -> String

Converts the given value to a String. Read more
Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.
§

impl<V, T> VZip<V> for T
where V: MultiLane<T>,

§

fn vzip(self) -> V

§

impl<T> WithSubscriber for T

§

fn with_subscriber<S>(self, subscriber: S) -> WithDispatch<Self>
where S: Into<Dispatch>,

Attaches the provided Subscriber to this type, returning a [WithDispatch] wrapper. Read more
§

fn with_current_subscriber(self) -> WithDispatch<Self>

Attaches the current default Subscriber to this type, returning a [WithDispatch] wrapper. Read more
§

impl<T> ErasedDestructor for T
where T: 'static,

§

impl<T> MaybeSendSync for T