nyx_space/dynamics/
orbital.rs

1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
/*
    Nyx, blazing fast astrodynamics
    Copyright (C) 2018-onwards Christopher Rabotin <christopher.rabotin@gmail.com>

    This program is free software: you can redistribute it and/or modify
    it under the terms of the GNU Affero General Public License as published
    by the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    This program is distributed in the hope that it will be useful,
    but WITHOUT ANY WARRANTY; without even the implied warranty of
    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
    GNU Affero General Public License for more details.

    You should have received a copy of the GNU Affero General Public License
    along with this program.  If not, see <https://www.gnu.org/licenses/>.
*/

use super::{
    AccelModel, DynamicsAlmanacSnafu, DynamicsAstroSnafu, DynamicsError, DynamicsPlanetarySnafu,
};
use crate::cosmic::{AstroPhysicsSnafu, Frame, Orbit};
use crate::linalg::{Const, Matrix3, Matrix6, OVector, Vector3, Vector6};

use anise::almanac::Almanac;
use anise::astro::Aberration;
use hyperdual::linalg::norm;
use hyperdual::{extract_jacobian_and_result, hyperspace_from_vector, Float, OHyperdual};
use snafu::ResultExt;
use std::f64;
use std::fmt;
use std::sync::Arc;

pub use super::sph_harmonics::Harmonics;

/// `OrbitalDynamics` provides the equations of motion for any celestial dynamic, without state transition matrix computation.
#[derive(Clone)]

pub struct OrbitalDynamics {
    pub accel_models: Vec<Arc<dyn AccelModel + Sync>>,
}

impl OrbitalDynamics {
    /// Initializes the point masses gravities with the provided list of bodies
    pub fn point_masses(celestial_objects: Vec<i32>) -> Self {
        // Create the point masses
        Self::new(vec![PointMasses::new(celestial_objects)])
    }

    /// Initializes a OrbitalDynamics which does not simulate the gravity pull of other celestial objects but the primary one.
    pub fn two_body() -> Self {
        Self::new(vec![])
    }

    /// Initialize orbital dynamics with a list of acceleration models
    pub fn new(accel_models: Vec<Arc<dyn AccelModel + Sync>>) -> Self {
        Self { accel_models }
    }

    /// Initialize new orbital mechanics with the provided model.
    /// **Note:** Orbital dynamics _always_ include two body dynamics, these cannot be turned off.
    pub fn from_model(accel_model: Arc<dyn AccelModel + Sync>) -> Self {
        Self::new(vec![accel_model])
    }
}

impl fmt::Display for OrbitalDynamics {
    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
        let models: Vec<String> = self.accel_models.iter().map(|x| format!("{x}")).collect();
        write!(f, "Orbital dynamics: {}", models.join("; "))
    }
}

impl OrbitalDynamics {
    pub(crate) fn eom(
        &self,
        osc: &Orbit,
        almanac: Arc<Almanac>,
    ) -> Result<OVector<f64, Const<42>>, DynamicsError> {
        // Still return something of size 42, but the STM will be zeros.
        let body_acceleration = (-osc
            .frame
            .mu_km3_s2()
            .context(AstroPhysicsSnafu)
            .context(DynamicsAstroSnafu)?
            / osc.rmag_km().powi(3))
            * osc.radius_km;

        let mut d_x = Vector6::from_iterator(
            osc.velocity_km_s
                .iter()
                .chain(body_acceleration.iter())
                .cloned(),
        );

        // Apply the acceleration models
        for model in &self.accel_models {
            let model_acc = model.eom(osc, almanac.clone())?;
            for i in 0..3 {
                d_x[i + 3] += model_acc[i];
            }
        }

        Ok(OVector::<f64, Const<42>>::from_iterator(
            d_x.iter()
                .chain(OVector::<f64, Const<36>>::zeros().iter())
                .cloned(),
        ))
    }

    pub fn dual_eom(
        &self,
        _delta_t_s: f64,
        osc: &Orbit,
        almanac: Arc<Almanac>,
    ) -> Result<(Vector6<f64>, Matrix6<f64>), DynamicsError> {
        // Extract data from hyperspace
        // Build full state vector with partials in the right position (hence building with all six components)
        let state: Vector6<OHyperdual<f64, Const<7>>> =
            hyperspace_from_vector(&osc.to_cartesian_pos_vel());

        let radius = state.fixed_rows::<3>(0).into_owned();
        let velocity = state.fixed_rows::<3>(3).into_owned();

        // Code up math as usual
        let rmag = norm(&radius);
        let body_acceleration = radius
            * (OHyperdual::<f64, Const<7>>::from_real(
                -osc.frame
                    .mu_km3_s2()
                    .context(AstroPhysicsSnafu)
                    .context(DynamicsAstroSnafu)?,
            ) / rmag.powi(3));

        // Extract result into Vector6 and Matrix6
        let mut dx = Vector6::zeros();
        let mut grad = Matrix6::zeros();
        for i in 0..6 {
            dx[i] = if i < 3 {
                velocity[i].real()
            } else {
                body_acceleration[i - 3].real()
            };
            for j in 1..7 {
                grad[(i, j - 1)] = if i < 3 {
                    velocity[i][j]
                } else {
                    body_acceleration[i - 3][j]
                };
            }
        }

        // Apply the acceleration models
        for model in &self.accel_models {
            // let (model_acc, model_grad) = model.dual_eom(&radius, osc)?;
            let (model_acc, model_grad) = model.dual_eom(osc, almanac.clone())?;
            for i in 0..3 {
                dx[i + 3] += model_acc[i];
                for j in 1..4 {
                    grad[(i + 3, j - 1)] += model_grad[(i, j - 1)];
                }
            }
        }

        // This function returns the time derivative of each function. The propagator will add this to the state vector (which has the previous STM).
        // This is why we don't multiply the gradient (A matrix) with the previous STM
        Ok((dx, grad))
    }
}

/// PointMasses model
pub struct PointMasses {
    pub celestial_objects: Vec<i32>,
    /// Light-time correction computation if extra point masses are needed
    pub correction: Option<Aberration>,
}

impl PointMasses {
    /// Initializes the point masses gravities with the provided list of bodies
    pub fn new(celestial_objects: Vec<i32>) -> Arc<Self> {
        Arc::new(Self {
            celestial_objects,
            correction: None,
        })
    }

    /// Initializes the point masses gravities with the provided list of bodies, and accounting for some light time correction
    pub fn with_correction(celestial_objects: Vec<i32>, correction: Aberration) -> Self {
        Self {
            celestial_objects,
            correction: Some(correction),
        }
    }
}

impl fmt::Display for PointMasses {
    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
        let masses: Vec<String> = self
            .celestial_objects
            .iter()
            .map(|third_body| format!("{}", Frame::from_ephem_j2000(*third_body)))
            .collect();
        write!(f, "Point masses of {}", masses.join(", "))
    }
}

impl AccelModel for PointMasses {
    fn eom(&self, osc: &Orbit, almanac: Arc<Almanac>) -> Result<Vector3<f64>, DynamicsError> {
        let mut d_x = Vector3::zeros();
        // Get all of the position vectors between the center body and the third bodies
        for third_body in self.celestial_objects.iter().copied() {
            if osc.frame.ephem_origin_id_match(third_body) {
                // Ignore the contribution of the integration frame, that's handled by OrbitalDynamics
                continue;
            }

            let third_body_frame = almanac
                .frame_from_uid(osc.frame.with_ephem(third_body))
                .context(DynamicsPlanetarySnafu {
                    action: "planetary data from third body not loaded",
                })?;

            // Orbit of j-th body as seen from primary body
            let st_ij = almanac
                .transform(third_body_frame, osc.frame, osc.epoch, self.correction)
                .context(DynamicsAlmanacSnafu {
                    action: "computing third body gravitational pull",
                })?;

            let r_ij = st_ij.radius_km;
            let r_ij3 = st_ij.rmag_km().powi(3);
            let r_j = osc.radius_km - r_ij; // sc as seen from 3rd body
            let r_j3 = r_j.norm().powi(3);
            d_x += -third_body_frame
                .mu_km3_s2()
                .context(AstroPhysicsSnafu)
                .context(DynamicsAstroSnafu)?
                * (r_j / r_j3 + r_ij / r_ij3);
        }
        Ok(d_x)
    }

    fn dual_eom(
        &self,
        osc: &Orbit,
        almanac: Arc<Almanac>,
    ) -> Result<(Vector3<f64>, Matrix3<f64>), DynamicsError> {
        // Build the hyperdual space of the radius vector
        let radius: Vector3<OHyperdual<f64, Const<7>>> = hyperspace_from_vector(&osc.radius_km);
        // Extract result into Vector6 and Matrix6
        let mut fx = Vector3::zeros();
        let mut grad = Matrix3::zeros();

        // Get all of the position vectors between the center body and the third bodies
        for third_body in &self.celestial_objects {
            let third_body_frame = almanac
                .frame_from_uid(Frame::from_ephem_j2000(*third_body))
                .context(DynamicsPlanetarySnafu {
                    action: "planetary data from third body not loaded",
                })?;

            if osc.frame.ephem_origin_match(third_body_frame) {
                // Ignore the contribution of the integration frame, that's handled by OrbitalDynamics
                continue;
            }

            let gm_d = OHyperdual::<f64, Const<7>>::from_real(
                -third_body_frame
                    .mu_km3_s2()
                    .context(AstroPhysicsSnafu)
                    .context(DynamicsAstroSnafu)?,
            );

            // Orbit of j-th body as seen from primary body
            let st_ij = almanac
                .transform(third_body_frame, osc.frame, osc.epoch, self.correction)
                .context(DynamicsAlmanacSnafu {
                    action: "computing third body gravitational pull",
                })?;

            let r_ij: Vector3<OHyperdual<f64, Const<7>>> = hyperspace_from_vector(&st_ij.radius_km);
            let r_ij3 = norm(&r_ij).powi(3);

            // The difference leads to the dual parts nulling themselves out, so let's fix that.
            let mut r_j = radius - r_ij; // sc as seen from 3rd body
            r_j[0][1] = 1.0;
            r_j[1][2] = 1.0;
            r_j[2][3] = 1.0;

            let r_j3 = norm(&r_j).powi(3);
            let mut third_body_acc_d = r_j / r_j3 + r_ij / r_ij3;
            third_body_acc_d[0] *= gm_d;
            third_body_acc_d[1] *= gm_d;
            third_body_acc_d[2] *= gm_d;

            let (fxp, gradp) = extract_jacobian_and_result::<_, 3, 3, 7>(&third_body_acc_d);
            fx += fxp;
            grad += gradp;
        }

        Ok((fx, grad))
    }
}