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 snafu::ResultExt;
30use std::f64;
31use std::fmt;
32use std::sync::Arc;
33
34pub use super::sph_harmonics::Harmonics;
35
36/// `OrbitalDynamics` provides the equations of motion for any celestial dynamic, without state transition matrix computation.
37#[derive(Clone)]
38
39pub struct OrbitalDynamics {
40    pub accel_models: Vec<Arc<dyn AccelModel + Sync>>,
41}
42
43impl OrbitalDynamics {
44    /// Initializes the point masses gravities with the provided list of bodies
45    pub fn point_masses(celestial_objects: Vec<i32>) -> Self {
46        // Create the point masses
47        Self::new(vec![PointMasses::new(celestial_objects)])
48    }
49
50    /// Initializes a OrbitalDynamics which does not simulate the gravity pull of other celestial objects but the primary one.
51    pub fn two_body() -> Self {
52        Self::new(vec![])
53    }
54
55    /// Initialize orbital dynamics with a list of acceleration models
56    pub fn new(accel_models: Vec<Arc<dyn AccelModel + Sync>>) -> Self {
57        Self { accel_models }
58    }
59
60    /// Initialize new orbital mechanics with the provided model.
61    /// **Note:** Orbital dynamics _always_ include two body dynamics, these cannot be turned off.
62    pub fn from_model(accel_model: Arc<dyn AccelModel + Sync>) -> Self {
63        Self::new(vec![accel_model])
64    }
65}
66
67impl fmt::Display for OrbitalDynamics {
68    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
69        let models: Vec<String> = self.accel_models.iter().map(|x| format!("{x}")).collect();
70        write!(f, "Orbital dynamics: {}", models.join("; "))
71    }
72}
73
74impl OrbitalDynamics {
75    pub(crate) fn eom(
76        &self,
77        osc: &Orbit,
78        almanac: Arc<Almanac>,
79    ) -> Result<OVector<f64, Const<42>>, DynamicsError> {
80        // Still return something of size 42, but the STM will be zeros.
81        let body_acceleration = (-osc
82            .frame
83            .mu_km3_s2()
84            .context(AstroPhysicsSnafu)
85            .context(DynamicsAstroSnafu)?
86            / osc.rmag_km().powi(3))
87            * osc.radius_km;
88
89        let mut d_x = Vector6::from_iterator(
90            osc.velocity_km_s
91                .iter()
92                .chain(body_acceleration.iter())
93                .cloned(),
94        );
95
96        // Apply the acceleration models
97        for model in &self.accel_models {
98            let model_acc = model.eom(osc, almanac.clone())?;
99            for i in 0..3 {
100                d_x[i + 3] += model_acc[i];
101            }
102        }
103
104        Ok(OVector::<f64, Const<42>>::from_iterator(
105            d_x.iter()
106                .chain(OVector::<f64, Const<36>>::zeros().iter())
107                .cloned(),
108        ))
109    }
110
111    pub fn dual_eom(
112        &self,
113        _delta_t_s: f64,
114        osc: &Orbit,
115        almanac: Arc<Almanac>,
116    ) -> Result<(Vector6<f64>, Matrix6<f64>), DynamicsError> {
117        // Extract data from hyperspace
118        // Build full state vector with partials in the right position (hence building with all six components)
119        let state: Vector6<OHyperdual<f64, Const<7>>> =
120            hyperspace_from_vector(&osc.to_cartesian_pos_vel());
121
122        let radius = state.fixed_rows::<3>(0).into_owned();
123        let velocity = state.fixed_rows::<3>(3).into_owned();
124
125        // Code up math as usual
126        let rmag = norm(&radius);
127        let body_acceleration = radius
128            * (OHyperdual::<f64, Const<7>>::from_real(
129                -osc.frame
130                    .mu_km3_s2()
131                    .context(AstroPhysicsSnafu)
132                    .context(DynamicsAstroSnafu)?,
133            ) / rmag.powi(3));
134
135        // Extract result into Vector6 and Matrix6
136        let mut dx = Vector6::zeros();
137        let mut grad = Matrix6::zeros();
138        for i in 0..6 {
139            dx[i] = if i < 3 {
140                velocity[i].real()
141            } else {
142                body_acceleration[i - 3].real()
143            };
144            for j in 1..7 {
145                grad[(i, j - 1)] = if i < 3 {
146                    velocity[i][j]
147                } else {
148                    body_acceleration[i - 3][j]
149                };
150            }
151        }
152
153        // Apply the acceleration models
154        for model in &self.accel_models {
155            let (model_acc, model_grad) = model.dual_eom(osc, almanac.clone())?;
156            for i in 0..3 {
157                dx[i + 3] += model_acc[i];
158                for j in 1..4 {
159                    grad[(i + 3, j - 1)] += model_grad[(i, j - 1)];
160                }
161            }
162        }
163
164        // This function returns the time derivative of each function. The propagator will add this to the state vector (which has the previous STM).
165        // This is why we don't multiply the gradient (A matrix) with the previous STM
166        Ok((dx, grad))
167    }
168}
169
170/// PointMasses model
171pub struct PointMasses {
172    pub celestial_objects: Vec<i32>,
173    /// Light-time correction computation if extra point masses are needed
174    pub correction: Option<Aberration>,
175}
176
177impl PointMasses {
178    /// Initializes the point masses gravities with the provided list of bodies
179    pub fn new(celestial_objects: Vec<i32>) -> Arc<Self> {
180        Arc::new(Self {
181            celestial_objects,
182            correction: None,
183        })
184    }
185
186    /// Initializes the point masses gravities with the provided list of bodies, and accounting for some light time correction
187    pub fn with_correction(celestial_objects: Vec<i32>, correction: Aberration) -> Self {
188        Self {
189            celestial_objects,
190            correction: Some(correction),
191        }
192    }
193}
194
195impl fmt::Display for PointMasses {
196    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
197        let masses: Vec<String> = self
198            .celestial_objects
199            .iter()
200            .map(|third_body| format!("{}", Frame::from_ephem_j2000(*third_body)))
201            .collect();
202        write!(f, "Point masses of {}", masses.join(", "))
203    }
204}
205
206impl AccelModel for PointMasses {
207    fn eom(&self, osc: &Orbit, almanac: Arc<Almanac>) -> Result<Vector3<f64>, DynamicsError> {
208        let mut d_x = Vector3::zeros();
209        // Get all of the position vectors between the center body and the third bodies
210        for third_body in self.celestial_objects.iter().copied() {
211            if osc.frame.ephem_origin_id_match(third_body) {
212                // Ignore the contribution of the integration frame, that's handled by OrbitalDynamics
213                continue;
214            }
215
216            let third_body_frame = almanac
217                .frame_from_uid(osc.frame.with_ephem(third_body))
218                .context(DynamicsPlanetarySnafu {
219                    action: "planetary data from third body not loaded",
220                })?;
221
222            // Orbit of j-th body as seen from primary body
223            let st_ij = almanac
224                .transform(third_body_frame, osc.frame, osc.epoch, self.correction)
225                .context(DynamicsAlmanacSnafu {
226                    action: "computing third body gravitational pull",
227                })?;
228
229            let r_ij = st_ij.radius_km;
230            let r_ij3 = st_ij.rmag_km().powi(3);
231            let r_j = osc.radius_km - r_ij; // sc as seen from 3rd body
232            let r_j3 = r_j.norm().powi(3);
233            d_x += -third_body_frame
234                .mu_km3_s2()
235                .context(AstroPhysicsSnafu)
236                .context(DynamicsAstroSnafu)?
237                * (r_j / r_j3 + r_ij / r_ij3);
238        }
239        Ok(d_x)
240    }
241
242    fn dual_eom(
243        &self,
244        osc: &Orbit,
245        almanac: Arc<Almanac>,
246    ) -> Result<(Vector3<f64>, Matrix3<f64>), DynamicsError> {
247        // Build the hyperdual space of the radius vector
248        let radius: Vector3<OHyperdual<f64, Const<7>>> = hyperspace_from_vector(&osc.radius_km);
249        // Extract result into Vector6 and Matrix6
250        let mut fx = Vector3::zeros();
251        let mut grad = Matrix3::zeros();
252
253        // Get all of the position vectors between the center body and the third bodies
254        for third_body in &self.celestial_objects {
255            let third_body_frame = almanac
256                .frame_from_uid(Frame::from_ephem_j2000(*third_body))
257                .context(DynamicsPlanetarySnafu {
258                    action: "planetary data from third body not loaded",
259                })?;
260
261            if osc.frame.ephem_origin_match(third_body_frame) {
262                // Ignore the contribution of the integration frame, that's handled by OrbitalDynamics
263                continue;
264            }
265
266            let gm_d = OHyperdual::<f64, Const<7>>::from_real(
267                -third_body_frame
268                    .mu_km3_s2()
269                    .context(AstroPhysicsSnafu)
270                    .context(DynamicsAstroSnafu)?,
271            );
272
273            // Orbit of j-th body as seen from primary body
274            let st_ij = almanac
275                .transform(third_body_frame, osc.frame, osc.epoch, self.correction)
276                .context(DynamicsAlmanacSnafu {
277                    action: "computing third body gravitational pull",
278                })?;
279
280            let r_ij: Vector3<OHyperdual<f64, Const<7>>> = hyperspace_from_vector(&st_ij.radius_km);
281            let r_ij3 = norm(&r_ij).powi(3);
282
283            // The difference leads to the dual parts nulling themselves out, so let's fix that.
284            let mut r_j = radius - r_ij; // sc as seen from 3rd body
285            r_j[0][1] = 1.0;
286            r_j[1][2] = 1.0;
287            r_j[2][3] = 1.0;
288
289            let r_j3 = norm(&r_j).powi(3);
290            let mut third_body_acc_d = r_j / r_j3 + r_ij / r_ij3;
291            third_body_acc_d[0] *= gm_d;
292            third_body_acc_d[1] *= gm_d;
293            third_body_acc_d[2] *= gm_d;
294
295            let (fxp, gradp) = extract_jacobian_and_result::<_, 3, 3, 7>(&third_body_acc_d);
296            fx += fxp;
297            grad += gradp;
298        }
299
300        Ok((fx, grad))
301    }
302}