Skip to main content

nyx_space/dynamics/
orbital.rs

1/*
2    Nyx, blazing fast astrodynamics
3    Copyright (C) 2018-onwards Christopher Rabotin <christopher.rabotin@gmail.com>
4
5    This program is free software: you can redistribute it and/or modify
6    it under the terms of the GNU Affero General Public License as published
7    by the Free Software Foundation, either version 3 of the License, or
8    (at your option) any later version.
9
10    This program is distributed in the hope that it will be useful,
11    but WITHOUT ANY WARRANTY; without even the implied warranty of
12    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
13    GNU Affero General Public License for more details.
14
15    You should have received a copy of the GNU Affero General Public License
16    along with this program.  If not, see <https://www.gnu.org/licenses/>.
17*/
18
19use super::{
20    AccelModel, DynamicsAlmanacSnafu, DynamicsAstroSnafu, DynamicsError, DynamicsPlanetarySnafu,
21};
22use crate::cosmic::{AstroPhysicsSnafu, Frame, Orbit};
23use crate::linalg::{Const, Matrix3, Matrix6, OVector, Vector3, Vector6};
24
25use anise::almanac::Almanac;
26use anise::astro::Aberration;
27use hyperdual::linalg::norm;
28use hyperdual::{extract_jacobian_and_result, hyperspace_from_vector, Float, OHyperdual};
29use serde::{Deserialize, Serialize};
30use serde_dhall::{SimpleType, StaticType};
31use snafu::ResultExt;
32use std::collections::HashMap;
33use std::f64;
34use std::fmt;
35use std::sync::Arc;
36
37pub use super::sph_harmonics::GravityField;
38
39/// `OrbitalDynamics` provides the equations of motion for any celestial dynamic, without state transition matrix computation.
40#[derive(Clone)]
41pub struct OrbitalDynamics {
42    pub accel_models: Vec<Arc<dyn AccelModel + Sync>>,
43}
44
45impl OrbitalDynamics {
46    /// Initializes the point masses gravities with the provided list of bodies
47    pub fn point_masses(celestial_objects: Vec<i32>) -> Self {
48        // Create the point masses
49        Self::new(vec![PointMasses::new(celestial_objects)])
50    }
51
52    /// Initializes a OrbitalDynamics which does not simulate the gravity pull of other celestial objects but the primary one.
53    pub fn two_body() -> Self {
54        Self::new(vec![])
55    }
56
57    /// Initialize orbital dynamics with a list of acceleration models
58    pub fn new(accel_models: Vec<Arc<dyn AccelModel + Sync>>) -> Self {
59        Self { accel_models }
60    }
61
62    /// Initialize new orbital mechanics with the provided model.
63    /// **Note:** Orbital dynamics _always_ include two body dynamics, these cannot be turned off.
64    pub fn from_model(accel_model: Arc<dyn AccelModel + Sync>) -> Self {
65        Self::new(vec![accel_model])
66    }
67}
68
69impl fmt::Display for OrbitalDynamics {
70    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
71        let models: Vec<String> = self.accel_models.iter().map(|x| format!("{x}")).collect();
72        write!(f, "Orbital dynamics: {}", models.join("; "))
73    }
74}
75
76impl OrbitalDynamics {
77    pub(crate) fn eom(
78        &self,
79        osc: &Orbit,
80        almanac: Arc<Almanac>,
81    ) -> Result<OVector<f64, Const<42>>, DynamicsError> {
82        // Still return something of size 42, but the STM will be zeros.
83        let body_acceleration = (-osc
84            .frame
85            .mu_km3_s2()
86            .context(AstroPhysicsSnafu)
87            .context(DynamicsAstroSnafu)?
88            / osc.rmag_km().powi(3))
89            * osc.radius_km;
90
91        let mut d_x = Vector6::from_iterator(
92            osc.velocity_km_s
93                .iter()
94                .chain(body_acceleration.iter())
95                .cloned(),
96        );
97
98        // Apply the acceleration models
99        for model in &self.accel_models {
100            let model_acc = model.eom(osc, almanac.clone())?;
101            for i in 0..3 {
102                d_x[i + 3] += model_acc[i];
103            }
104        }
105
106        Ok(OVector::<f64, Const<42>>::from_iterator(
107            d_x.iter()
108                .chain(OVector::<f64, Const<36>>::zeros().iter())
109                .cloned(),
110        ))
111    }
112
113    pub fn dual_eom(
114        &self,
115        _delta_t_s: f64,
116        osc: &Orbit,
117        almanac: Arc<Almanac>,
118    ) -> Result<(Vector6<f64>, Matrix6<f64>), DynamicsError> {
119        // Extract data from hyperspace
120        // Build full state vector with partials in the right position (hence building with all six components)
121        let state: Vector6<OHyperdual<f64, Const<7>>> =
122            hyperspace_from_vector(&osc.to_cartesian_pos_vel());
123
124        let radius = state.fixed_rows::<3>(0).into_owned();
125        let velocity = state.fixed_rows::<3>(3).into_owned();
126
127        // Code up math as usual
128        let rmag = norm(&radius);
129        let body_acceleration = radius
130            * (OHyperdual::<f64, Const<7>>::from_real(
131                -osc.frame
132                    .mu_km3_s2()
133                    .context(AstroPhysicsSnafu)
134                    .context(DynamicsAstroSnafu)?,
135            ) / rmag.powi(3));
136
137        // Extract result into Vector6 and Matrix6
138        let mut dx = Vector6::zeros();
139        let mut grad = Matrix6::zeros();
140        for i in 0..6 {
141            dx[i] = if i < 3 {
142                velocity[i].real()
143            } else {
144                body_acceleration[i - 3].real()
145            };
146            for j in 1..7 {
147                grad[(i, j - 1)] = if i < 3 {
148                    velocity[i][j]
149                } else {
150                    body_acceleration[i - 3][j]
151                };
152            }
153        }
154
155        // Apply the acceleration models
156        for model in &self.accel_models {
157            let (model_acc, model_grad) = model.dual_eom(osc, almanac.clone())?;
158            for i in 0..3 {
159                dx[i + 3] += model_acc[i];
160                for j in 1..4 {
161                    grad[(i + 3, j - 1)] += model_grad[(i, j - 1)];
162                }
163            }
164        }
165
166        // This function returns the time derivative of each function. The propagator will add this to the state vector (which has the previous STM).
167        // This is why we don't multiply the gradient (A matrix) with the previous STM
168        Ok((dx, grad))
169    }
170}
171
172/// PointMasses model
173#[derive(Clone, Debug, Serialize, Deserialize)]
174pub struct PointMasses {
175    pub celestial_objects: Vec<i32>,
176    /// Light-time correction computation if extra point masses are needed
177    pub correction: Option<Aberration>,
178}
179
180impl PointMasses {
181    /// Initializes the point masses gravities with the provided list of bodies
182    pub fn new(celestial_objects: Vec<i32>) -> Arc<Self> {
183        Arc::new(Self {
184            celestial_objects,
185            correction: None,
186        })
187    }
188
189    /// Initializes the point masses gravities with the provided list of bodies, and accounting for some light time correction
190    pub fn with_correction(celestial_objects: Vec<i32>, correction: Option<Aberration>) -> Self {
191        Self {
192            celestial_objects,
193            correction,
194        }
195    }
196}
197
198impl fmt::Display for PointMasses {
199    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
200        let masses: Vec<String> = self
201            .celestial_objects
202            .iter()
203            .map(|third_body| format!("{}", Frame::from_ephem_j2000(*third_body)))
204            .collect();
205        write!(f, "Point masses of {}", masses.join(", "))
206    }
207}
208
209impl AccelModel for PointMasses {
210    fn eom(&self, osc: &Orbit, almanac: Arc<Almanac>) -> Result<Vector3<f64>, DynamicsError> {
211        let mut d_x = Vector3::zeros();
212        // Get all of the position vectors between the center body and the third bodies
213        for third_body in self.celestial_objects.iter().copied() {
214            if osc.frame.ephem_origin_id_match(third_body) {
215                // Ignore the contribution of the integration frame, that's handled by OrbitalDynamics
216                continue;
217            }
218
219            let third_body_frame = almanac
220                .frame_info(osc.frame.with_ephem(third_body))
221                .context(DynamicsPlanetarySnafu {
222                    action: "planetary data from third body not loaded",
223                })?;
224
225            // Orbit of j-th body as seen from primary body
226            let st_ij = almanac
227                .transform(third_body_frame, osc.frame, osc.epoch, self.correction)
228                .context(DynamicsAlmanacSnafu {
229                    action: "computing third body gravitational pull",
230                })?;
231
232            let r_ij = st_ij.radius_km;
233            let r_ij3 = st_ij.rmag_km().powi(3);
234            let r_j = osc.radius_km - r_ij; // sc as seen from 3rd body
235            let r_j3 = r_j.norm().powi(3);
236            d_x += -third_body_frame
237                .mu_km3_s2()
238                .context(AstroPhysicsSnafu)
239                .context(DynamicsAstroSnafu)?
240                * (r_j / r_j3 + r_ij / r_ij3);
241        }
242        Ok(d_x)
243    }
244
245    fn dual_eom(
246        &self,
247        osc: &Orbit,
248        almanac: Arc<Almanac>,
249    ) -> Result<(Vector3<f64>, Matrix3<f64>), DynamicsError> {
250        // Build the hyperdual space of the radius vector
251        let radius: Vector3<OHyperdual<f64, Const<7>>> = hyperspace_from_vector(&osc.radius_km);
252        // Extract result into Vector6 and Matrix6
253        let mut fx = Vector3::zeros();
254        let mut grad = Matrix3::zeros();
255
256        // Get all of the position vectors between the center body and the third bodies
257        for third_body in self.celestial_objects.iter().copied() {
258            if osc.frame.ephem_origin_id_match(third_body) {
259                // Ignore the contribution of the integration frame, that's handled by OrbitalDynamics
260                continue;
261            }
262
263            let third_body_frame = almanac
264                .frame_info(osc.frame.with_ephem(third_body))
265                .context(DynamicsPlanetarySnafu {
266                    action: "planetary data from third body not loaded",
267                })?;
268
269            let gm_d = OHyperdual::<f64, Const<7>>::from_real(
270                -third_body_frame
271                    .mu_km3_s2()
272                    .context(AstroPhysicsSnafu)
273                    .context(DynamicsAstroSnafu)?,
274            );
275
276            // Orbit of j-th body as seen from primary body
277            let st_ij = almanac
278                .transform(third_body_frame, osc.frame, osc.epoch, self.correction)
279                .context(DynamicsAlmanacSnafu {
280                    action: "computing third body gravitational pull",
281                })?;
282
283            let r_ij: Vector3<OHyperdual<f64, Const<7>>> = hyperspace_from_vector(&st_ij.radius_km);
284            let r_ij3 = norm(&r_ij).powi(3);
285
286            // The difference leads to the dual parts nulling themselves out, so let's fix that.
287            let mut r_j = radius - r_ij; // sc as seen from 3rd body
288            r_j[0][1] = 1.0;
289            r_j[1][2] = 1.0;
290            r_j[2][3] = 1.0;
291
292            let r_j3 = norm(&r_j).powi(3);
293            let mut third_body_acc_d = r_j / r_j3 + r_ij / r_ij3;
294            third_body_acc_d[0] *= gm_d;
295            third_body_acc_d[1] *= gm_d;
296            third_body_acc_d[2] *= gm_d;
297
298            let (fxp, gradp) = extract_jacobian_and_result::<_, 3, 3, 7>(&third_body_acc_d);
299            fx += fxp;
300            grad += gradp;
301        }
302
303        Ok((fx, grad))
304    }
305}
306impl StaticType for PointMasses {
307    fn static_type() -> SimpleType {
308        let mut fields = HashMap::new();
309
310        fields.insert("celestial_objects".to_string(), Vec::<i32>::static_type());
311
312        // Manually define the record for Aberration right here
313        // instead of calling Aberration::static_type()
314        let aberration_fields = {
315            let mut f = HashMap::new();
316            f.insert("converged".to_string(), bool::static_type());
317            f.insert("stellar".to_string(), bool::static_type());
318            f.insert("transmit_mode".to_string(), bool::static_type());
319            SimpleType::Record(f)
320        };
321
322        fields.insert(
323            "correction".to_string(),
324            SimpleType::Optional(Box::new(aberration_fields)),
325        );
326
327        SimpleType::Record(fields)
328    }
329}