nyx_space/mc/
multivariate.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 std::error::Error;
20
21use super::{DispersedState, StateDispersion};
22use crate::errors::StateError;
23use crate::md::prelude::{BPlane, OrbitDual};
24use crate::md::{AstroSnafu, StateParameter};
25use crate::{pseudo_inverse, NyxError, Spacecraft, State};
26use nalgebra::{DMatrix, DVector, SMatrix, SVector};
27use rand_distr::{Distribution, Normal};
28use snafu::ResultExt;
29
30/// A multivariate spacecraft state generator for Monte Carlo analyses. Ensures that the covariance is properly applied on all provided state variables.
31pub struct MvnSpacecraft {
32    /// The template state
33    pub template: Spacecraft,
34    pub dispersions: Vec<StateDispersion>,
35    /// The mean of the multivariate normal distribution
36    pub mean: SVector<f64, 9>,
37    /// The dot product \sqrt{\vec s} \cdot \vec v, where S is the singular values and V the V matrix from the SVD decomp of the covariance of multivariate normal distribution
38    pub sqrt_s_v: SMatrix<f64, 9, 9>,
39    /// The standard normal distribution used to seed the multivariate normal distribution
40    pub std_norm_distr: Normal<f64>,
41}
42
43impl MvnSpacecraft {
44    /// Creates a new mulivariate state generator from a mean and covariance on the set of state parameters.
45    /// The covariance must be positive semi definite.
46    ///
47    /// # Algorithm
48    /// This function will build the rotation matrix to rotate the requested dispersions into the Spacecraft state space using [OrbitDual].
49    /// If there are any dispersions on the Cr and Cd, then these are dispersed independently (because they are iid).
50    pub fn new(
51        template: Spacecraft,
52        dispersions: Vec<StateDispersion>,
53    ) -> Result<Self, Box<dyn Error>> {
54        let mut cov = SMatrix::<f64, 9, 9>::zeros();
55        let mut mean = SVector::<f64, 9>::zeros();
56
57        let orbit_dual = OrbitDual::from(template.orbit);
58        let mut b_plane = None;
59        for obj in &dispersions {
60            if obj.param.is_b_plane() {
61                b_plane = Some(
62                    BPlane::from_dual(orbit_dual)
63                        .context(AstroSnafu)
64                        .map_err(Box::new)?,
65                );
66                break;
67            }
68        }
69
70        let num_orbital = dispersions
71            .iter()
72            .filter(|disp| disp.param.is_orbital())
73            .count();
74
75        if num_orbital > 0 {
76            // Build the rotation matrix from the orbital dispersions to the Cartesian state.
77            let mut jac = DMatrix::from_element(num_orbital, 6, 0.0);
78            let mut covar = DMatrix::from_element(num_orbital, num_orbital, 0.0);
79            let mut means = DVector::from_element(num_orbital, 0.0);
80            let orbit_dual = OrbitDual::from(template.orbit);
81            let mut rno = 0;
82            for disp in &dispersions {
83                if disp.param.is_orbital() {
84                    let partial = if disp.param.is_b_plane() {
85                        match disp.param {
86                            StateParameter::BdotR => b_plane.unwrap().b_r,
87                            StateParameter::BdotT => b_plane.unwrap().b_t,
88                            StateParameter::BLTOF => b_plane.unwrap().ltof_s,
89                            _ => unreachable!(),
90                        }
91                    } else {
92                        orbit_dual.partial_for(disp.param).context(AstroSnafu)?
93                    };
94
95                    for (cno, val) in [
96                        partial.wtr_x(),
97                        partial.wtr_y(),
98                        partial.wtr_z(),
99                        partial.wtr_vx(),
100                        partial.wtr_vy(),
101                        partial.wtr_vz(),
102                    ]
103                    .iter()
104                    .copied()
105                    .enumerate()
106                    {
107                        jac[(rno, cno)] = val;
108                    }
109                    covar[(rno, rno)] = disp.std_dev.unwrap_or(0.0).powi(2);
110                    means[rno] = disp.mean.unwrap_or(0.0);
111                    rno += 1;
112                }
113            }
114
115            // Now that we have the Jacobian that rotates from the Cartesian elements to the dispersions parameters,
116            // let's compute the inverse of this Jacobian to rotate from the dispersion params into the Cartesian elements.
117            let jac_inv = pseudo_inverse!(&jac)?;
118
119            // Rotate the orbital covariance back into the Cartesian state space, making this a 6x6.
120            let orbit_cov = &jac_inv * &covar * jac_inv.transpose();
121            // Rotate the means into the Cartesian space
122            let cartesian_mean = jac_inv * means;
123
124            for ii in 0..6 {
125                for jj in 0..6 {
126                    cov[(ii, jj)] = orbit_cov[(ii, jj)];
127                }
128                mean[ii] = cartesian_mean[ii];
129            }
130        };
131
132        if dispersions.len() > num_orbital {
133            for disp in &dispersions {
134                if disp.param.is_orbital() {
135                    continue;
136                } else {
137                    match disp.param {
138                        StateParameter::Cr => {
139                            cov[(7, 7)] = disp.mean.unwrap_or(0.0).powi(2);
140                        }
141                        StateParameter::Cd => {
142                            cov[(8, 8)] = disp.mean.unwrap_or(0.0).powi(2);
143                        }
144                        StateParameter::DryMass | StateParameter::PropMass => {
145                            cov[(9, 9)] = disp.mean.unwrap_or(0.0).powi(2);
146                        }
147                        _ => return Err(Box::new(StateError::ReadOnly { param: disp.param })),
148                    }
149                }
150            }
151        }
152
153        // At this point, the cov matrix is a 9x9 with all dispersions transformed into the Cartesian state space.
154
155        let svd = cov.svd(false, true);
156        if svd.v_t.is_none() {
157            return Err(Box::new(NyxError::CovarianceMatrixNotPsd));
158        }
159
160        let sqrt_s = svd.singular_values.map(|x| x.sqrt());
161        let mut sqrt_s_v_t = svd.v_t.unwrap().transpose();
162
163        for (i, mut col) in sqrt_s_v_t.column_iter_mut().enumerate() {
164            col *= sqrt_s[i];
165        }
166
167        Ok(Self {
168            template,
169            dispersions,
170            mean,
171            sqrt_s_v: sqrt_s_v_t,
172            std_norm_distr: Normal::new(0.0, 1.0).unwrap(),
173        })
174    }
175
176    /// Same as `new` but with a zero mean
177    pub fn zero_mean(
178        template: Spacecraft,
179        mut dispersions: Vec<StateDispersion>,
180    ) -> Result<Self, Box<dyn Error>> {
181        for disp in &mut dispersions {
182            disp.mean = Some(0.0);
183        }
184
185        Self::new(template, dispersions)
186    }
187
188    /// Initializes a new multivariate distribution using the state data in the spacecraft state space.
189    pub fn from_spacecraft_cov(
190        template: Spacecraft,
191        cov: SMatrix<f64, 9, 9>,
192        mean: SVector<f64, 9>,
193    ) -> Result<Self, Box<dyn Error>> {
194        // Check that covariance is PSD by ensuring that all the eigenvalues are positive or nil
195        match cov.eigenvalues() {
196            None => return Err(Box::new(NyxError::CovarianceMatrixNotPsd)),
197            Some(evals) => {
198                for eigenval in &evals {
199                    if *eigenval < 0.0 {
200                        return Err(Box::new(NyxError::CovarianceMatrixNotPsd));
201                    }
202                }
203            }
204        };
205
206        let svd = cov.svd(false, true);
207        if svd.v_t.is_none() {
208            return Err(Box::new(NyxError::CovarianceMatrixNotPsd));
209        }
210
211        let s = svd.singular_values;
212        // Item by item multiplication
213        let mut sqrt_s_v = svd.v_t.unwrap().transpose();
214        for (i, mut col) in sqrt_s_v.column_iter_mut().enumerate() {
215            col *= s[i].sqrt();
216        }
217
218        let dispersions = vec![
219            StateDispersion::builder()
220                .param(StateParameter::X)
221                .std_dev(cov[(0, 0)])
222                .build(),
223            StateDispersion::builder()
224                .param(StateParameter::Y)
225                .std_dev(cov[(1, 1)])
226                .build(),
227            StateDispersion::builder()
228                .param(StateParameter::Z)
229                .std_dev(cov[(2, 2)])
230                .build(),
231            StateDispersion::builder()
232                .param(StateParameter::VX)
233                .std_dev(cov[(3, 3)])
234                .build(),
235            StateDispersion::builder()
236                .param(StateParameter::VY)
237                .std_dev(cov[(4, 4)])
238                .build(),
239            StateDispersion::builder()
240                .param(StateParameter::VZ)
241                .std_dev(cov[(5, 5)])
242                .build(),
243            StateDispersion::builder()
244                .param(StateParameter::Cr)
245                .std_dev(cov[(6, 6)])
246                .build(),
247            StateDispersion::builder()
248                .param(StateParameter::Cd)
249                .std_dev(cov[(7, 7)])
250                .build(),
251            StateDispersion::builder()
252                .param(StateParameter::PropMass)
253                .std_dev(cov[(8, 8)])
254                .build(),
255        ];
256
257        Ok(Self {
258            template,
259            dispersions,
260            mean,
261            sqrt_s_v,
262            std_norm_distr: Normal::new(0.0, 1.0).unwrap(),
263        })
264    }
265}
266
267impl Distribution<DispersedState<Spacecraft>> for MvnSpacecraft {
268    fn sample<R: rand::Rng + ?Sized>(&self, rng: &mut R) -> DispersedState<Spacecraft> {
269        // Generate the vector representing the state
270        let x_rng = SVector::<f64, 9>::from_fn(|_, _| self.std_norm_distr.sample(rng));
271        let x = self.sqrt_s_v * x_rng + self.mean;
272        let mut state = self.template;
273
274        // Set the new state data
275        for (coord, val) in x.iter().copied().enumerate() {
276            if coord < 3 {
277                state.orbit.radius_km[coord] += val;
278            } else if coord < 6 {
279                state.orbit.velocity_km_s[coord % 3] += val;
280            } else if coord == 6 {
281                state.srp.coeff_reflectivity += val;
282            } else if coord == 7 {
283                state.drag.coeff_drag += val;
284            } else if coord == 8 {
285                state.mass.prop_mass_kg += val;
286            }
287        }
288
289        let mut actual_dispersions = Vec::new();
290        for disp in &self.dispersions {
291            // Compute the delta
292            let delta = self.template.value(disp.param).unwrap() - state.value(disp.param).unwrap();
293            actual_dispersions.push((disp.param, delta));
294        }
295
296        DispersedState {
297            state,
298            actual_dispersions,
299        }
300    }
301}
302
303#[cfg(test)]
304mod multivariate_ut {
305    use super::*;
306    use crate::Spacecraft;
307    use crate::GMAT_EARTH_GM;
308
309    #[test]
310    fn disperse_r_mag() {
311        use anise::constants::frames::EARTH_J2000;
312        use anise::prelude::Orbit;
313
314        use crate::time::Epoch;
315        use rand_pcg::Pcg64Mcg;
316
317        // Ensure that this worked: a 3 sigma deviation around 1 km means we shouldn't have 99.7% of samples within those bounds.
318        // Create a reproducible fast seed
319        let seed = 0;
320
321        let eme2k = EARTH_J2000.with_mu_km3_s2(GMAT_EARTH_GM);
322
323        let dt = Epoch::from_gregorian_utc_at_midnight(2021, 1, 31);
324        let state = Spacecraft::builder()
325            .orbit(Orbit::keplerian(
326                8_191.93, 1e-6, 12.85, 306.614, 314.19, 99.887_7, dt, eme2k,
327            ))
328            .build();
329
330        // Check that we can modify the radius magnitude
331        let std_dev = 1.0;
332        let generator = MvnSpacecraft::new(
333            state,
334            vec![StateDispersion::builder()
335                .param(StateParameter::Rmag)
336                .std_dev(std_dev)
337                .build()],
338        )
339        .unwrap();
340
341        let rng = Pcg64Mcg::new(seed);
342        let init_rmag = state.orbit.rmag_km();
343        let cnt_too_far: u16 = generator
344            .sample_iter(rng)
345            .take(1000)
346            .map(|dispersed_state| {
347                if (init_rmag - dispersed_state.state.orbit.rmag_km()).abs() >= 3.0 * std_dev {
348                    1
349                } else {
350                    0
351                }
352            })
353            .sum::<u16>();
354
355        // We specified a seed so we know exactly what to expect and we've reset the seed to 0.
356        assert_eq!(
357            cnt_too_far,
358            6, // Mathematically, this should be 3!
359            "Should have about 3% of samples being more than 3 sigma away, got {cnt_too_far}"
360        );
361    }
362
363    #[test]
364    fn disperse_full_cartesian() {
365        use anise::constants::frames::EARTH_J2000;
366        use anise::prelude::Orbit;
367
368        use crate::time::Epoch;
369        use crate::Spacecraft;
370        use crate::GMAT_EARTH_GM;
371
372        use rand_pcg::Pcg64Mcg;
373
374        let eme2k = EARTH_J2000.with_mu_km3_s2(GMAT_EARTH_GM);
375
376        let dt = Epoch::from_gregorian_utc_at_midnight(2021, 1, 31);
377        let state = Orbit::keplerian(8_191.93, 1e-6, 12.85, 306.614, 314.19, 99.887_7, dt, eme2k);
378
379        let std_dev = [10.0, 10.0, 10.0, 0.2, 0.2, 0.2, 0.0, 0.0, 0.0];
380
381        let generator = MvnSpacecraft::new(
382            Spacecraft {
383                orbit: state,
384                ..Default::default()
385            },
386            vec![
387                StateDispersion::builder()
388                    .param(StateParameter::X)
389                    .std_dev(10.0)
390                    .build(),
391                StateDispersion::builder()
392                    .param(StateParameter::Y)
393                    .std_dev(10.0)
394                    .build(),
395                StateDispersion::builder()
396                    .param(StateParameter::Z)
397                    .std_dev(10.0)
398                    .build(),
399                StateDispersion::builder()
400                    .param(StateParameter::VX)
401                    .std_dev(0.2)
402                    .build(),
403                StateDispersion::builder()
404                    .param(StateParameter::VY)
405                    .std_dev(0.2)
406                    .build(),
407                StateDispersion::builder()
408                    .param(StateParameter::VZ)
409                    .std_dev(0.2)
410                    .build(),
411            ],
412        )
413        .unwrap();
414
415        // Ensure that this worked: a 3 sigma deviation around 1 km means we shouldn't have 99.7% of samples within those bounds.
416        // Create a reproducible fast seed
417        let seed = 0;
418        let rng = Pcg64Mcg::new(seed);
419
420        let cnt_too_far: u16 = generator
421            .sample_iter(rng)
422            .take(1000)
423            .map(|dispersed_state| {
424                let mut cnt = 0;
425                for (idx, val_std_dev) in std_dev.iter().take(6).enumerate() {
426                    let cur_val = dispersed_state.state.to_vector()[idx];
427                    let nom_val = state.to_cartesian_pos_vel()[idx];
428                    if (cur_val - nom_val).abs() > *val_std_dev {
429                        cnt += 1;
430                    }
431                }
432                cnt
433            })
434            .sum::<u16>();
435
436        // We specified a seed so we know exactly what to expect
437        assert_eq!(
438            cnt_too_far / 6,
439            312,
440            "Should have about 3% of samples being more than 3 sigma away, got {cnt_too_far}"
441        );
442    }
443
444    #[test]
445    fn disperse_raan_only() {
446        use anise::constants::frames::EARTH_J2000;
447        use anise::prelude::Orbit;
448
449        use crate::time::Epoch;
450        use rand_pcg::Pcg64Mcg;
451
452        let eme2k = EARTH_J2000.with_mu_km3_s2(GMAT_EARTH_GM);
453
454        let dt = Epoch::from_gregorian_utc_at_midnight(2021, 1, 31);
455        let state = Spacecraft::builder()
456            .orbit(Orbit::keplerian(
457                8_100.0, 1e-6, 12.85, 356.614, 14.19, 199.887_7, dt, eme2k,
458            ))
459            .build();
460
461        let angle_sigma_deg = 0.2;
462
463        let generator = MvnSpacecraft::new(
464            state,
465            vec![StateDispersion::zero_mean(
466                StateParameter::RAAN,
467                angle_sigma_deg,
468            )],
469        )
470        .unwrap();
471
472        // Ensure that this worked: a 3 sigma deviation around 1 km means we shouldn't have 99.7% of samples within those bounds.
473        // Create a reproducible fast seed
474        let seed = 0;
475        let rng = Pcg64Mcg::new(seed);
476
477        let cnt_too_far: u16 = generator
478            .sample_iter(rng)
479            .take(1000)
480            .map(|dispersed_state| {
481                // For all other orbital parameters, make sure that we have not changed things dramatically.
482                for param in [StateParameter::SMA, StateParameter::Inclination] {
483                    let orig = state.value(param).unwrap();
484                    let new = dispersed_state.state.value(param).unwrap();
485                    let prct_change = 100.0 * (orig - new).abs() / orig;
486                    assert!(
487                        prct_change < 5.0,
488                        "{param} changed by {prct_change:.3} % ({orig:.3e} -> {new:.3e})"
489                    );
490                }
491
492                if (dispersed_state.actual_dispersions[0].1).abs() > 3.0 * angle_sigma_deg {
493                    1
494                } else {
495                    0
496                }
497            })
498            .sum::<u16>();
499
500        // We specified a seed so we know exactly what to expect
501        // Consider: https://github.com/nyx-space/nyx/issues/339
502        assert_eq!(
503            cnt_too_far,
504            7, // This is about twice too high
505            "Should have about 3% of samples being more than 3 sigma away, got {cnt_too_far}"
506        );
507    }
508
509    #[ignore = "https://github.com/nyx-space/nyx/issues/339"]
510    #[test]
511    fn disperse_keplerian() {
512        use anise::constants::frames::EARTH_J2000;
513        use anise::prelude::Orbit;
514
515        use crate::time::Epoch;
516        use rand_pcg::Pcg64Mcg;
517
518        let eme2k = EARTH_J2000.with_mu_km3_s2(GMAT_EARTH_GM);
519
520        let dt = Epoch::from_gregorian_utc_at_midnight(2021, 1, 31);
521        let state = Spacecraft::builder()
522            .orbit(Orbit::keplerian(
523                8_100.0, 1e-6, 12.85, 356.614, 14.19, 199.887_7, dt, eme2k,
524            ))
525            .build();
526
527        let sma_sigma_km = 10.0;
528        let inc_sigma_deg = 0.15;
529        let angle_sigma_deg = 0.02;
530
531        let generator = MvnSpacecraft::new(
532            state,
533            vec![
534                StateDispersion::zero_mean(StateParameter::SMA, sma_sigma_km),
535                StateDispersion::zero_mean(StateParameter::Inclination, inc_sigma_deg),
536                StateDispersion::zero_mean(StateParameter::RAAN, angle_sigma_deg),
537                StateDispersion::zero_mean(StateParameter::AoP, angle_sigma_deg),
538            ],
539        )
540        .unwrap();
541
542        // Ensure that this worked: a 3 sigma deviation around 1 km means we shouldn't have 99.7% of samples within those bounds.
543        // Create a reproducible fast seed
544        let seed = 0;
545        let rng = Pcg64Mcg::new(seed);
546
547        let cnt_too_far: u16 = generator
548            .sample_iter(rng)
549            .take(1000)
550            .map(|dispersed_state| {
551                if (dispersed_state.actual_dispersions[0].1).abs() > 3.0 * sma_sigma_km
552                    || (dispersed_state.actual_dispersions[1].1).abs() > 3.0 * inc_sigma_deg
553                    || (dispersed_state.actual_dispersions[2].1).abs() > 3.0 * angle_sigma_deg
554                    || (dispersed_state.actual_dispersions[3].1).abs() > 3.0 * angle_sigma_deg
555                {
556                    1
557                } else {
558                    0
559                }
560            })
561            .sum::<u16>();
562
563        // We specified a seed so we know exactly what to expect
564        assert_eq!(
565            dbg!(cnt_too_far) / 3,
566            3,
567            "Should have about 3% of samples being more than 3 sigma away, got {cnt_too_far}"
568        );
569    }
570}