1use 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::{Float, OHyperdual, extract_jacobian_and_result, hyperspace_from_vector};
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
37#[cfg(feature = "python")]
38use pyo3::prelude::*;
39
40pub use super::gravity_field::GravityField;
41
42#[derive(Clone)]
44pub struct OrbitalDynamics {
45 pub accel_models: Vec<Arc<dyn AccelModel + Sync>>,
46}
47
48impl OrbitalDynamics {
49 pub fn point_masses(celestial_objects: Vec<i32>) -> Self {
51 Self::new(vec![PointMasses::new(celestial_objects)])
53 }
54
55 pub fn two_body() -> Self {
57 Self::new(vec![])
58 }
59
60 pub fn new(accel_models: Vec<Arc<dyn AccelModel + Sync>>) -> Self {
62 Self { accel_models }
63 }
64
65 pub fn from_model(accel_model: Arc<dyn AccelModel + Sync>) -> Self {
68 Self::new(vec![accel_model])
69 }
70}
71
72impl fmt::Display for OrbitalDynamics {
73 fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
74 let models: Vec<String> = self.accel_models.iter().map(|x| format!("{x}")).collect();
75 write!(f, "Orbital dynamics: {}", models.join("; "))
76 }
77}
78
79impl OrbitalDynamics {
80 pub(crate) fn eom(
81 &self,
82 osc: &Orbit,
83 almanac: Arc<Almanac>,
84 ) -> Result<OVector<f64, Const<42>>, DynamicsError> {
85 let body_acceleration = (-osc
87 .frame
88 .mu_km3_s2()
89 .context(AstroPhysicsSnafu)
90 .context(DynamicsAstroSnafu)?
91 / osc.rmag_km().powi(3))
92 * osc.radius_km;
93
94 let mut d_x = Vector6::from_iterator(
95 osc.velocity_km_s
96 .iter()
97 .chain(body_acceleration.iter())
98 .cloned(),
99 );
100
101 for model in &self.accel_models {
103 let model_acc = model.eom(osc, almanac.clone())?;
104 for i in 0..3 {
105 d_x[i + 3] += model_acc[i];
106 }
107 }
108
109 Ok(OVector::<f64, Const<42>>::from_iterator(
110 d_x.iter()
111 .chain(OVector::<f64, Const<36>>::zeros().iter())
112 .cloned(),
113 ))
114 }
115
116 pub fn dual_eom(
117 &self,
118 _delta_t_s: f64,
119 osc: &Orbit,
120 almanac: Arc<Almanac>,
121 ) -> Result<(Vector6<f64>, Matrix6<f64>), DynamicsError> {
122 let state: Vector6<OHyperdual<f64, Const<7>>> =
125 hyperspace_from_vector(&osc.to_cartesian_pos_vel());
126
127 let radius = state.fixed_rows::<3>(0).into_owned();
128 let velocity = state.fixed_rows::<3>(3).into_owned();
129
130 let rmag = norm(&radius);
132 let body_acceleration = radius
133 * (OHyperdual::<f64, Const<7>>::from_real(
134 -osc.frame
135 .mu_km3_s2()
136 .context(AstroPhysicsSnafu)
137 .context(DynamicsAstroSnafu)?,
138 ) / rmag.powi(3));
139
140 let mut dx = Vector6::zeros();
142 let mut grad = Matrix6::zeros();
143 for i in 0..6 {
144 dx[i] = if i < 3 {
145 velocity[i].real()
146 } else {
147 body_acceleration[i - 3].real()
148 };
149 for j in 1..7 {
150 grad[(i, j - 1)] = if i < 3 {
151 velocity[i][j]
152 } else {
153 body_acceleration[i - 3][j]
154 };
155 }
156 }
157
158 for model in &self.accel_models {
160 let (model_acc, model_grad) = model.gradient(osc, almanac.clone())?;
161 for i in 0..3 {
162 dx[i + 3] += model_acc[i];
163 for j in 1..4 {
164 grad[(i + 3, j - 1)] += model_grad[(i, j - 1)];
165 }
166 }
167 }
168
169 Ok((dx, grad))
172 }
173}
174
175#[derive(Clone, Debug, Serialize, Deserialize)]
177#[cfg_attr(feature = "python", pyclass(from_py_object, get_all, set_all))]
178pub struct PointMasses {
179 pub celestial_objects: Vec<i32>,
180 pub correction: Option<Aberration>,
182}
183
184impl PointMasses {
185 pub fn new(celestial_objects: Vec<i32>) -> Arc<Self> {
187 Arc::new(Self {
188 celestial_objects,
189 correction: None,
190 })
191 }
192
193 pub fn with_correction(celestial_objects: Vec<i32>, correction: Option<Aberration>) -> Self {
195 Self {
196 celestial_objects,
197 correction,
198 }
199 }
200}
201
202impl fmt::Display for PointMasses {
203 fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
204 let masses: Vec<String> = self
205 .celestial_objects
206 .iter()
207 .map(|third_body| format!("{}", Frame::from_ephem_j2000(*third_body)))
208 .collect();
209 write!(f, "Point masses of {}", masses.join(", "))
210 }
211}
212
213impl AccelModel for PointMasses {
214 fn eom(&self, osc: &Orbit, almanac: Arc<Almanac>) -> Result<Vector3<f64>, DynamicsError> {
215 let mut d_x = Vector3::zeros();
216 for third_body in self.celestial_objects.iter().copied() {
218 if osc.frame.ephem_origin_id_match(third_body) {
219 continue;
221 }
222
223 let third_body_frame = almanac
224 .frame_info(osc.frame.with_ephem(third_body))
225 .context(DynamicsPlanetarySnafu {
226 action: "planetary data from third body not loaded",
227 })?;
228
229 let st_ij = almanac
231 .transform(third_body_frame, osc.frame, osc.epoch, self.correction)
232 .context(DynamicsAlmanacSnafu {
233 action: "computing third body gravitational pull",
234 })?;
235
236 let r_ij = st_ij.radius_km;
237 let r_ij3 = st_ij.rmag_km().powi(3);
238 let r_j = osc.radius_km - r_ij; let r_j3 = r_j.norm().powi(3);
240 d_x += -third_body_frame
241 .mu_km3_s2()
242 .context(AstroPhysicsSnafu)
243 .context(DynamicsAstroSnafu)?
244 * (r_j / r_j3 + r_ij / r_ij3);
245 }
246 Ok(d_x)
247 }
248
249 fn gradient(
250 &self,
251 osc: &Orbit,
252 almanac: Arc<Almanac>,
253 ) -> Result<(Vector3<f64>, Matrix3<f64>), DynamicsError> {
254 let radius: Vector3<OHyperdual<f64, Const<7>>> = hyperspace_from_vector(&osc.radius_km);
256 let mut fx = Vector3::zeros();
258 let mut grad = Matrix3::zeros();
259
260 for third_body in self.celestial_objects.iter().copied() {
262 if osc.frame.ephem_origin_id_match(third_body) {
263 continue;
265 }
266
267 let third_body_frame = almanac
268 .frame_info(osc.frame.with_ephem(third_body))
269 .context(DynamicsPlanetarySnafu {
270 action: "planetary data from third body not loaded",
271 })?;
272
273 let gm_d = OHyperdual::<f64, Const<7>>::from_real(
274 -third_body_frame
275 .mu_km3_s2()
276 .context(AstroPhysicsSnafu)
277 .context(DynamicsAstroSnafu)?,
278 );
279
280 let st_ij = almanac
282 .transform(third_body_frame, osc.frame, osc.epoch, self.correction)
283 .context(DynamicsAlmanacSnafu {
284 action: "computing third body gravitational pull",
285 })?;
286
287 let r_ij: Vector3<OHyperdual<f64, Const<7>>> = hyperspace_from_vector(&st_ij.radius_km);
288 let r_ij3 = norm(&r_ij).powi(3);
289
290 let mut r_j = radius - r_ij; r_j[0][1] = 1.0;
293 r_j[1][2] = 1.0;
294 r_j[2][3] = 1.0;
295
296 let r_j3 = norm(&r_j).powi(3);
297 let mut third_body_acc_d = r_j / r_j3 + r_ij / r_ij3;
298 third_body_acc_d[0] *= gm_d;
299 third_body_acc_d[1] *= gm_d;
300 third_body_acc_d[2] *= gm_d;
301
302 let (fxp, gradp) = extract_jacobian_and_result::<_, 3, 3, 7>(&third_body_acc_d);
303 fx += fxp;
304 grad += gradp;
305 }
306
307 Ok((fx, grad))
308 }
309}
310impl StaticType for PointMasses {
311 fn static_type() -> SimpleType {
312 let mut fields = HashMap::new();
313
314 fields.insert("celestial_objects".to_string(), Vec::<i32>::static_type());
315
316 let aberration_fields = {
319 let mut f = HashMap::new();
320 f.insert("converged".to_string(), bool::static_type());
321 f.insert("stellar".to_string(), bool::static_type());
322 f.insert("transmit_mode".to_string(), bool::static_type());
323 SimpleType::Record(f)
324 };
325
326 fields.insert(
327 "correction".to_string(),
328 SimpleType::Optional(Box::new(aberration_fields)),
329 );
330
331 SimpleType::Record(fields)
332 }
333}