nyx_space/od/estimate/
kfestimate.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 super::{Estimate, State};
20use crate::cosmic::AstroError;
21use crate::linalg::allocator::Allocator;
22use crate::linalg::{DefaultAllocator, DimName, Matrix, OMatrix, OVector};
23use crate::mc::{MvnSpacecraft, StateDispersion};
24use crate::md::prelude::OrbitDual;
25use crate::md::StateParameter;
26use crate::Spacecraft;
27use nalgebra::Const;
28use nalgebra::SMatrix;
29use rand::SeedableRng;
30use rand_distr::Distribution;
31use rand_pcg::Pcg64Mcg;
32use std::cmp::PartialEq;
33use std::error::Error;
34use std::fmt;
35use std::ops::Mul;
36
37/// Kalman filter Estimate
38#[derive(Debug, Copy, Clone, PartialEq)]
39pub struct KfEstimate<T: State>
40where
41    DefaultAllocator: Allocator<<T as State>::Size>
42        + Allocator<<T as State>::Size, <T as State>::Size>
43        + Allocator<<T as State>::Size>
44        + Allocator<<T as State>::VecLength>,
45    <DefaultAllocator as Allocator<<T as State>::Size>>::Buffer<f64>: Copy,
46    <DefaultAllocator as Allocator<<T as State>::Size, <T as State>::Size>>::Buffer<f64>: Copy,
47{
48    /// The estimated state
49    pub nominal_state: T,
50    /// The state deviation
51    pub state_deviation: OVector<f64, <T as State>::Size>,
52    /// The Covariance of this estimate
53    pub covar: OMatrix<f64, <T as State>::Size, <T as State>::Size>,
54    /// The predicted covariance of this estimate
55    pub covar_bar: OMatrix<f64, <T as State>::Size, <T as State>::Size>,
56    /// Whether or not this is a predicted estimate from a time update, or an estimate from a measurement
57    pub predicted: bool,
58    /// The STM used to compute this Estimate
59    pub stm: OMatrix<f64, <T as State>::Size, <T as State>::Size>,
60}
61
62impl<T: State> KfEstimate<T>
63where
64    DefaultAllocator: Allocator<<T as State>::Size>
65        + Allocator<<T as State>::Size, <T as State>::Size>
66        + Allocator<<T as State>::Size>
67        + Allocator<<T as State>::VecLength>,
68    <DefaultAllocator as Allocator<<T as State>::Size>>::Buffer<f64>: Copy,
69    <DefaultAllocator as Allocator<<T as State>::Size, <T as State>::Size>>::Buffer<f64>: Copy,
70{
71    /// Initializes a new filter estimate from the nominal state (not dispersed) and the full covariance
72    pub fn from_covar(
73        nominal_state: T,
74        covar: OMatrix<f64, <T as State>::Size, <T as State>::Size>,
75    ) -> Self {
76        Self {
77            nominal_state,
78            state_deviation: OVector::<f64, <T as State>::Size>::zeros(),
79            covar,
80            covar_bar: covar,
81            predicted: true,
82            stm: OMatrix::<f64, <T as State>::Size, <T as State>::Size>::identity(),
83        }
84    }
85
86    /// Initializes a new filter estimate from the nominal state (not dispersed) and the diagonal of the covariance
87    pub fn from_diag(nominal_state: T, diag: OVector<f64, <T as State>::Size>) -> Self {
88        let covar = Matrix::from_diagonal(&diag);
89        Self {
90            nominal_state,
91            state_deviation: OVector::<f64, <T as State>::Size>::zeros(),
92            covar,
93            covar_bar: covar,
94            predicted: true,
95            stm: OMatrix::<f64, <T as State>::Size, <T as State>::Size>::identity(),
96        }
97    }
98}
99
100impl KfEstimate<Spacecraft> {
101    /// Generates an initial Kalman filter state estimate dispersed from the nominal state using the provided standard deviation parameters.
102    ///
103    /// The resulting estimate will have a diagonal covariance matrix constructed from the variances of each parameter.
104    /// *Limitation:* This method may not work correctly for all Keplerian orbital elements, refer to
105    /// <https://github.com/nyx-space/nyx/issues/339> for details.
106    pub fn disperse_from_diag(
107        nominal_state: Spacecraft,
108        dispersions: Vec<StateDispersion>,
109        seed: Option<u128>,
110    ) -> Result<Self, Box<dyn Error>> {
111        let generator = MvnSpacecraft::new(nominal_state, dispersions)?;
112
113        let mut rng = match seed {
114            Some(seed) => Pcg64Mcg::new(seed),
115            None => Pcg64Mcg::from_entropy(),
116        };
117        let dispersed_state = generator.sample(&mut rng);
118
119        // Compute the difference between both states
120        let delta_orbit = (nominal_state.orbit - dispersed_state.state.orbit).unwrap();
121
122        // Build the covariance as three times the absolute value of the error, squared.
123
124        let diag_data = [
125            (3.0 * delta_orbit.radius_km.x.abs()).powi(2),
126            (3.0 * delta_orbit.radius_km.y.abs()).powi(2),
127            (3.0 * delta_orbit.radius_km.z.abs()).powi(2),
128            (3.0 * delta_orbit.velocity_km_s.x.abs()).powi(2),
129            (3.0 * delta_orbit.velocity_km_s.y.abs()).powi(2),
130            (3.0 * delta_orbit.velocity_km_s.z.abs()).powi(2),
131            (3.0 * (nominal_state.srp.coeff_reflectivity
132                - dispersed_state.state.srp.coeff_reflectivity)
133                .abs())
134            .powi(2),
135            (3.0 * (nominal_state.drag.coeff_drag - dispersed_state.state.drag.coeff_drag).abs())
136                .powi(2),
137            (3.0 * (nominal_state.mass.prop_mass_kg - dispersed_state.state.mass.prop_mass_kg)
138                .abs())
139            .powi(2),
140        ];
141
142        let diag = OVector::<f64, Const<9>>::from_iterator(diag_data);
143
144        // Build the covar from the diagonal
145        let covar = Matrix::from_diagonal(&diag);
146
147        Ok(Self {
148            nominal_state: dispersed_state.state,
149            state_deviation: OVector::<f64, Const<9>>::zeros(),
150            covar,
151            covar_bar: covar,
152            predicted: true,
153            stm: OMatrix::<f64, Const<9>, Const<9>>::identity(),
154        })
155    }
156
157    /// Builds a multivariate random variable from this estimate's nominal state and covariance, zero mean.
158    pub fn to_random_variable(&self) -> Result<MvnSpacecraft, Box<dyn Error>> {
159        MvnSpacecraft::from_spacecraft_cov(self.nominal_state, self.covar, self.state_deviation)
160    }
161
162    /// Returns the 1-sigma uncertainty for a given parameter, in that parameter's unit
163    ///
164    /// This method uses the [OrbitDual] structure to compute the estimate in the hyperdual space
165    /// and rotate the nominal covariance into that space.
166    pub fn sigma_for(&self, param: StateParameter) -> Result<f64, AstroError> {
167        // Build the rotation matrix using Orbit Dual.
168        let mut rotmat = SMatrix::<f64, 1, 6>::zeros();
169        let orbit_dual = OrbitDual::from(self.nominal_state.orbit);
170
171        let xf_partial = orbit_dual.partial_for(param)?;
172        for (cno, val) in [
173            xf_partial.wtr_x(),
174            xf_partial.wtr_y(),
175            xf_partial.wtr_z(),
176            xf_partial.wtr_vx(),
177            xf_partial.wtr_vy(),
178            xf_partial.wtr_vz(),
179        ]
180        .iter()
181        .copied()
182        .enumerate()
183        {
184            rotmat[(0, cno)] = val;
185        }
186
187        Ok((rotmat * self.covar.fixed_view::<6, 6>(0, 0) * rotmat.transpose())[(0, 0)].sqrt())
188    }
189
190    /// Returns the 6x6 covariance (i.e. square of the sigma/uncertainty) of the SMA, ECC, INC, RAAN, AOP, and True Anomaly.
191    pub fn keplerian_covar(&self) -> SMatrix<f64, 6, 6> {
192        // Build the rotation matrix using Orbit Dual.
193        let mut rotmat = SMatrix::<f64, 6, 6>::zeros();
194        let orbit_dual = OrbitDual::from(self.nominal_state.orbit);
195        for (pno, param) in [
196            StateParameter::SMA,
197            StateParameter::Eccentricity,
198            StateParameter::Inclination,
199            StateParameter::RAAN,
200            StateParameter::AoP,
201            StateParameter::TrueAnomaly,
202        ]
203        .iter()
204        .copied()
205        .enumerate()
206        {
207            let xf_partial = orbit_dual.partial_for(param).unwrap();
208            for (cno, val) in [
209                xf_partial.wtr_x(),
210                xf_partial.wtr_y(),
211                xf_partial.wtr_z(),
212                xf_partial.wtr_vx(),
213                xf_partial.wtr_vy(),
214                xf_partial.wtr_vz(),
215            ]
216            .iter()
217            .copied()
218            .enumerate()
219            {
220                rotmat[(pno, cno)] = val;
221            }
222        }
223
224        rotmat * self.covar.fixed_view::<6, 6>(0, 0) * rotmat.transpose()
225    }
226}
227
228impl<T: State> Estimate<T> for KfEstimate<T>
229where
230    DefaultAllocator: Allocator<<T as State>::Size>
231        + Allocator<<T as State>::Size, <T as State>::Size>
232        + Allocator<<T as State>::Size>
233        + Allocator<<T as State>::VecLength>,
234    <DefaultAllocator as Allocator<<T as State>::Size>>::Buffer<f64>: Copy,
235    <DefaultAllocator as Allocator<<T as State>::Size, <T as State>::Size>>::Buffer<f64>: Copy,
236{
237    fn zeros(nominal_state: T) -> Self {
238        Self {
239            nominal_state,
240            state_deviation: OVector::<f64, <T as State>::Size>::zeros(),
241            covar: OMatrix::<f64, <T as State>::Size, <T as State>::Size>::zeros(),
242            covar_bar: OMatrix::<f64, <T as State>::Size, <T as State>::Size>::zeros(),
243            predicted: true,
244            stm: OMatrix::<f64, <T as State>::Size, <T as State>::Size>::identity(),
245        }
246    }
247
248    fn nominal_state(&self) -> T {
249        self.nominal_state
250    }
251
252    fn state_deviation(&self) -> OVector<f64, <T as State>::Size> {
253        self.state_deviation
254    }
255
256    fn covar(&self) -> OMatrix<f64, <T as State>::Size, <T as State>::Size> {
257        self.covar
258    }
259
260    fn predicted_covar(&self) -> OMatrix<f64, <T as State>::Size, <T as State>::Size> {
261        self.covar_bar
262    }
263
264    fn predicted(&self) -> bool {
265        self.predicted
266    }
267    fn stm(&self) -> &OMatrix<f64, <T as State>::Size, <T as State>::Size> {
268        &self.stm
269    }
270    fn set_state_deviation(&mut self, new_state: OVector<f64, <T as State>::Size>) {
271        self.state_deviation = new_state;
272    }
273    fn set_covar(&mut self, new_covar: OMatrix<f64, <T as State>::Size, <T as State>::Size>) {
274        self.covar = new_covar;
275    }
276}
277
278impl<T: State> fmt::Display for KfEstimate<T>
279where
280    DefaultAllocator: Allocator<<T as State>::Size>
281        + Allocator<<T as State>::Size, <T as State>::Size>
282        + Allocator<<T as State>::Size>
283        + Allocator<<T as State>::VecLength>,
284    <DefaultAllocator as Allocator<<T as State>::Size>>::Buffer<f64>: Copy,
285    <DefaultAllocator as Allocator<<T as State>::Size, <T as State>::Size>>::Buffer<f64>: Copy,
286{
287    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
288        let dim = <T as State>::Size::dim();
289        let word = if self.predicted {
290            "Prediction"
291        } else {
292            "Estimate"
293        };
294        let mut fmt_cov = Vec::with_capacity(dim);
295        for i in 0..dim {
296            let unit = if i < 3 {
297                "km"
298            } else if i < 6 {
299                "km/s"
300            } else {
301                ""
302            };
303            fmt_cov.push(format!("{:.6} {unit}", &self.covar[(i, i)]));
304        }
305        write!(
306            f,
307            "=== {} @ {} -- within 3 sigma: {} ===\nstate {}\nsigmas [{}]\n",
308            word,
309            &self.epoch(),
310            self.within_3sigma(),
311            &self.state(),
312            fmt_cov.join(", ")
313        )
314    }
315}
316
317impl<T: State> fmt::LowerExp for KfEstimate<T>
318where
319    DefaultAllocator: Allocator<<T as State>::Size>
320        + Allocator<<T as State>::Size, <T as State>::Size>
321        + Allocator<<T as State>::Size>
322        + Allocator<<T as State>::VecLength>,
323    <DefaultAllocator as Allocator<<T as State>::Size>>::Buffer<f64>: Copy,
324    <DefaultAllocator as Allocator<<T as State>::Size, <T as State>::Size>>::Buffer<f64>: Copy,
325{
326    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
327        write!(
328            f,
329            "=== PREDICTED: {} ===\nEstState {:e} Covariance {:e}\n=====================",
330            &self.predicted, &self.state_deviation, &self.covar
331        )
332    }
333}
334
335impl<T: State> Mul<f64> for KfEstimate<T>
336where
337    DefaultAllocator: Allocator<<T as State>::Size>
338        + Allocator<<T as State>::Size, <T as State>::Size>
339        + Allocator<<T as State>::Size>
340        + Allocator<<T as State>::VecLength>,
341    <DefaultAllocator as Allocator<<T as State>::Size>>::Buffer<f64>: Copy,
342    <DefaultAllocator as Allocator<<T as State>::Size, <T as State>::Size>>::Buffer<f64>: Copy,
343{
344    type Output = Self;
345
346    fn mul(mut self, rhs: f64) -> Self::Output {
347        self.covar *= rhs.powi(2);
348        self
349    }
350}
351
352#[cfg(test)]
353mod ut_kfest {
354    use crate::{
355        mc::StateDispersion, md::StateParameter, od::estimate::KfEstimate, Spacecraft,
356        GMAT_EARTH_GM,
357    };
358    use anise::{constants::frames::EARTH_J2000, prelude::Orbit};
359    use hifitime::Epoch;
360
361    #[test]
362    fn test_estimate_from_disp() {
363        let eme2k = EARTH_J2000.with_mu_km3_s2(GMAT_EARTH_GM);
364        let dt = Epoch::from_gregorian_tai_at_midnight(2020, 1, 1);
365        let initial_state = Spacecraft::builder()
366            .orbit(Orbit::keplerian(
367                22000.0, 0.01, 30.0, 80.0, 40.0, 0.0, dt, eme2k,
368            ))
369            .build();
370
371        let initial_estimate = KfEstimate::disperse_from_diag(
372            initial_state,
373            vec![
374                StateDispersion::builder()
375                    .param(StateParameter::SMA)
376                    .std_dev(1.1)
377                    .build(),
378                StateDispersion::zero_mean(StateParameter::Inclination, 0.2),
379                StateDispersion::zero_mean(StateParameter::RAAN, 0.2),
380                StateDispersion::zero_mean(StateParameter::AoP, 0.2),
381            ],
382            Some(0),
383        )
384        .unwrap();
385
386        let initial_state_dev = initial_estimate.nominal_state;
387
388        let (init_rss_pos_km, init_rss_vel_km_s, _) =
389            initial_state.rss(&initial_state_dev).unwrap();
390
391        let delta = (initial_state.orbit - initial_state_dev.orbit).unwrap();
392
393        println!("Truth initial state:\n{initial_state}\n{initial_state:x}");
394        println!("Filter initial state:\n{initial_state_dev}\n{initial_state_dev:x}");
395        println!(
396            "Initial state dev:\t{init_rss_pos_km:.6} km\t{init_rss_vel_km_s:.6} km/s\n{delta}",
397        );
398        println!("covariance: {:.6}", initial_estimate.covar);
399
400        // Check that the error is in the square root of the covariance
401        assert!(delta.radius_km.x < initial_estimate.covar[(0, 0)].sqrt());
402        assert!(delta.radius_km.y < initial_estimate.covar[(1, 1)].sqrt());
403        assert!(delta.radius_km.z < initial_estimate.covar[(2, 2)].sqrt());
404        assert!(delta.velocity_km_s.x < initial_estimate.covar[(3, 3)].sqrt());
405        assert!(delta.velocity_km_s.y < initial_estimate.covar[(4, 4)].sqrt());
406        assert!(delta.velocity_km_s.z < initial_estimate.covar[(5, 5)].sqrt());
407    }
408}