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
use super::hyperdual::linalg::norm;
use super::hyperdual::{extract_jacobian_and_result, hyperspace_from_vector, Float, Hyperdual};
use super::{AccelModel, Dynamics, NyxError};
use crate::celestia::{Bodies, Cosm, Frame, LTCorr, Orbit};
use crate::dimensions::{
DimName, Matrix3, Matrix6, Vector3, Vector6, VectorN, U3, U36, U4, U42, U6, U7,
};
use crate::State;
use std::f64;
use std::sync::Arc;
pub use super::sph_harmonics::Harmonics;
#[derive(Clone)]
pub struct OrbitalDynamics<'a> {
pub accel_models: Vec<Arc<dyn AccelModel + Sync + 'a>>,
}
impl<'a> OrbitalDynamics<'a> {
pub fn point_masses(bodies: &[Bodies], cosm: Arc<Cosm>) -> Arc<Self> {
Self::new(vec![PointMasses::new(bodies, cosm)])
}
pub fn two_body() -> Arc<Self> {
Self::new(vec![])
}
pub fn new(accel_models: Vec<Arc<dyn AccelModel + Sync + 'a>>) -> Arc<Self> {
Arc::new(Self::new_raw(accel_models))
}
pub fn new_raw(accel_models: Vec<Arc<dyn AccelModel + Sync + 'a>>) -> Self {
Self { accel_models }
}
pub fn from_model(accel_model: Arc<dyn AccelModel + Sync + 'a>) -> Arc<Self> {
Self::new(vec![accel_model])
}
pub fn add_model(&mut self, accel_model: Arc<dyn AccelModel + Sync + 'a>) {
self.accel_models.push(accel_model);
}
pub fn with_model(self, accel_model: Arc<dyn AccelModel + Sync + 'a>) -> Arc<Self> {
let mut me = self.clone();
me.add_model(accel_model);
Arc::new(me)
}
}
impl<'a> Dynamics for OrbitalDynamics<'a> {
type HyperdualSize = U7;
type StateType = Orbit;
fn eom(
&self,
delta_t_s: f64,
state: &VectorN<f64, U42>,
ctx: &Orbit,
) -> Result<VectorN<f64, U42>, NyxError> {
let (new_state, new_stm) = if ctx.stm.is_some() {
let pos_vel = state.fixed_rows::<U6>(0).into_owned();
let (state, grad) = self.eom_grad(delta_t_s, &pos_vel, ctx)?;
let stm_dt = ctx.stm() * grad;
let mut stm_as_vec = VectorN::<f64, U36>::zeros();
let mut stm_idx = 0;
for i in 0..U6::dim() {
for j in 0..U6::dim() {
stm_as_vec[(stm_idx, 0)] = stm_dt[(i, j)];
stm_idx += 1;
}
}
(state, stm_as_vec)
} else {
let osc = ctx.ctor_from(delta_t_s, state);
let body_acceleration = (-osc.frame.gm() / osc.rmag().powi(3)) * osc.radius();
let mut d_x = Vector6::from_iterator(
osc.velocity()
.iter()
.chain(body_acceleration.iter())
.cloned(),
);
for model in &self.accel_models {
let model_acc = model.eom(&osc)?;
for i in 0..3 {
d_x[i + 3] += model_acc[i];
}
}
(d_x, VectorN::<f64, U36>::zeros())
};
Ok(VectorN::<f64, U42>::from_iterator(
new_state.iter().chain(new_stm.iter()).cloned(),
))
}
fn dual_eom(
&self,
_delta_t_s: f64,
state: &VectorN<Hyperdual<f64, U7>, U6>,
ctx: &Orbit,
) -> Result<(Vector6<f64>, Matrix6<f64>), NyxError> {
let radius = state.fixed_rows::<U3>(0).into_owned();
let velocity = state.fixed_rows::<U3>(3).into_owned();
let rmag = norm(&radius);
let body_acceleration =
radius * (Hyperdual::<f64, U7>::from_real(-ctx.frame.gm()) / rmag.powi(3));
let mut fx = Vector6::zeros();
let mut grad = Matrix6::zeros();
for i in 0..U6::dim() {
fx[i] = if i < 3 {
velocity[i].real()
} else {
body_acceleration[i - 3].real()
};
for j in 1..U7::dim() {
grad[(i, j - 1)] = if i < 3 {
velocity[i][j]
} else {
body_acceleration[i - 3][j]
};
}
}
for model in &self.accel_models {
let (model_acc, model_grad) = model.dual_eom(&radius, ctx)?;
for i in 0..U3::dim() {
fx[i + 3] += model_acc[i];
for j in 1..U4::dim() {
grad[(i + 3, j - 1)] += model_grad[(i, j - 1)];
}
}
}
Ok((fx, grad))
}
}
pub struct PointMasses {
pub bodies: Vec<Frame>,
pub cosm: Arc<Cosm>,
pub correction: LTCorr,
}
impl PointMasses {
pub fn new(bodies: &[Bodies], cosm: Arc<Cosm>) -> Arc<Self> {
Arc::new(Self::with_correction(bodies, cosm, LTCorr::None))
}
pub fn with_correction(bodies: &[Bodies], cosm: Arc<Cosm>, correction: LTCorr) -> Self {
let mut refs = Vec::new();
for body in bodies {
refs.push(cosm.frame_from_ephem_path(&body.ephem_path()));
}
Self {
bodies: refs,
cosm,
correction,
}
}
pub fn specific(body_names: &[String], cosm: Arc<Cosm>, correction: LTCorr) -> Self {
let mut refs = Vec::new();
for body in body_names {
refs.push(cosm.frame(body));
}
Self {
bodies: refs,
cosm,
correction,
}
}
}
impl AccelModel for PointMasses {
fn eom(&self, osc: &Orbit) -> Result<Vector3<f64>, NyxError> {
let mut d_x = Vector3::zeros();
for third_body in &self.bodies {
if third_body == &osc.frame {
continue;
}
let st_ij = self.cosm.celestial_state(
&third_body.ephem_path(),
osc.dt,
osc.frame,
self.correction,
);
let r_ij = st_ij.radius();
let r_ij3 = st_ij.rmag().powi(3);
let r_j = osc.radius() - r_ij;
let r_j3 = r_j.norm().powi(3);
d_x += -third_body.gm() * (r_j / r_j3 + r_ij / r_ij3);
}
Ok(d_x)
}
fn dual_eom(
&self,
state: &VectorN<Hyperdual<f64, U7>, U3>,
osc: &Orbit,
) -> Result<(Vector3<f64>, Matrix3<f64>), NyxError> {
let radius = state.fixed_rows::<U3>(0).into_owned();
let mut fx = Vector3::zeros();
let mut grad = Matrix3::zeros();
for third_body in &self.bodies {
if third_body == &osc.frame {
continue;
}
let gm_d = Hyperdual::<f64, U7>::from_real(-third_body.gm());
let st_ij = self.cosm.celestial_state(
&third_body.ephem_path(),
osc.dt,
osc.frame,
self.correction,
);
let r_ij: Vector3<Hyperdual<f64, U7>> = hyperspace_from_vector(&st_ij.radius());
let r_ij3 = norm(&r_ij).powi(3) / gm_d;
let mut r_j = radius - r_ij;
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) / gm_d;
let third_body_acc_d = r_j / r_j3 + r_ij / r_ij3;
let (fxp, gradp) = extract_jacobian_and_result::<_, U3, U3, _>(&third_body_acc_d);
fx += fxp;
grad += gradp;
}
Ok((fx, grad))
}
}