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