nyx_space/od/kalman/
filtering.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 crate::cosmic::AstroPhysicsSnafu;
20pub use crate::errors::NyxError;
21use crate::errors::StateAstroSnafu;
22use crate::linalg::allocator::Allocator;
23use crate::linalg::{DefaultAllocator, DimName, OMatrix, OVector};
24use crate::md::StateParameter;
25pub use crate::od::estimate::{Estimate, KfEstimate, Residual};
26use crate::od::prelude::KalmanVariant;
27use crate::od::process::ResidRejectCrit;
28pub use crate::od::snc::ProcessNoise;
29use crate::od::{ODDynamicsSnafu, ODError, ODStateSnafu, State};
30pub use crate::time::{Epoch, Unit};
31use snafu::prelude::*;
32
33use super::KalmanFilter;
34
35impl<T, A> KalmanFilter<T, A>
36where
37    A: DimName,
38    T: State,
39    DefaultAllocator: Allocator<<T as State>::Size>
40        + Allocator<<T as State>::VecLength>
41        + Allocator<A>
42        + Allocator<<T as State>::Size, <T as State>::Size>
43        + Allocator<A, A>
44        + Allocator<<T as State>::Size, A>
45        + Allocator<A, <T as State>::Size>,
46    <DefaultAllocator as Allocator<<T as State>::Size>>::Buffer<f64>: Copy,
47    <DefaultAllocator as Allocator<<T as State>::Size, <T as State>::Size>>::Buffer<f64>: Copy,
48{
49    /// Returns the previous estimate
50    pub fn previous_estimate(&self) -> &KfEstimate<T> {
51        &self.prev_estimate
52    }
53
54    pub fn set_previous_estimate(&mut self, est: &KfEstimate<T>) {
55        self.prev_estimate = *est;
56    }
57
58    /// Computes a time update/prediction (i.e. advances the filter estimate with the updated STM).
59    ///
60    /// May return a FilterError if the STM was not updated.
61    pub fn time_update(&mut self, nominal_state: T) -> Result<KfEstimate<T>, ODError> {
62        let stm = nominal_state.stm().context(ODDynamicsSnafu)?;
63        let mut covar_bar = stm * self.prev_estimate.covar * stm.transpose();
64
65        // Try to apply an SNC, if applicable
66        for (i, snc) in self.process_noise.iter().enumerate().rev() {
67            if let Some(mut snc_matrix) = snc.to_matrix(nominal_state.epoch()) {
68                if let Some(local_frame) = snc.local_frame {
69                    // Rotate the SNC from the definition frame into the state frame.
70                    let dcm = local_frame
71                        .dcm_to_inertial(nominal_state.orbit())
72                        .context(AstroPhysicsSnafu)
73                        .context(StateAstroSnafu {
74                            param: StateParameter::Epoch,
75                        })
76                        .context(ODStateSnafu {
77                            action: "rotating SNC from definition frame into state frame",
78                        })?;
79
80                    // Note: the SNC must be a diagonal matrix, so we only update the diagonals!
81                    match A::USIZE {
82                        3 => {
83                            let new_snc = dcm.rot_mat
84                                * snc_matrix.fixed_view::<3, 3>(0, 0)
85                                * dcm.rot_mat.transpose();
86                            for i in 0..A::USIZE {
87                                snc_matrix[(i, i)] = new_snc[(i, i)];
88                            }
89                        }
90                        6 => {
91                            let new_snc = dcm.state_dcm()
92                                * snc_matrix.fixed_view::<6, 6>(0, 0)
93                                * dcm.transpose().state_dcm();
94                            for i in 0..A::USIZE {
95                                snc_matrix[(i, i)] = new_snc[(i, i)];
96                            }
97                        }
98                        _ => {
99                            return Err(ODError::ODLimitation {
100                                action: "only process noises of size 3x3 or 6x6 are supported",
101                            })
102                        }
103                    }
104                }
105                // Check if we're using another SNC than the one before
106                if self.prev_used_snc != i {
107                    info!("Switched to {}-th {}", i, snc);
108                    self.prev_used_snc = i;
109                }
110
111                // Let's compute the Gamma matrix, an approximation of the time integral
112                // which assumes that the acceleration is constant between these two measurements.
113                let mut gamma = OMatrix::<f64, <T as State>::Size, A>::zeros();
114                let delta_t = (nominal_state.epoch() - self.prev_estimate.epoch()).to_seconds();
115                for blk in 0..A::dim() / 3 {
116                    for i in 0..3 {
117                        let idx_i = i + A::dim() * blk;
118                        let idx_j = i + 3 * blk;
119                        let idx_k = i + 3 + A::dim() * blk;
120                        // For first block
121                        // (0, 0) (1, 1) (2, 2) <=> \Delta t^2/2
122                        // (3, 0) (4, 1) (5, 2) <=> \Delta t
123                        // Second block
124                        // (6, 3) (7, 4) (8, 5) <=> \Delta t^2/2
125                        // (9, 3) (10, 4) (11, 5) <=> \Delta t
126                        // * \Delta t^2/2
127                        // (i, i) when blk = 0
128                        // (i + A::dim() * blk, i + 3) when blk = 1
129                        // (i + A::dim() * blk, i + 3 * blk)
130                        // * \Delta t
131                        // (i + 3, i) when blk = 0
132                        // (i + 3, i + 9) when blk = 1 (and I think i + 12 + 3)
133                        // (i + 3 + A::dim() * blk, i + 3 * blk)
134                        gamma[(idx_i, idx_j)] = delta_t.powi(2) / 2.0;
135                        gamma[(idx_k, idx_j)] = delta_t;
136                    }
137                }
138                // Let's add the process noise
139                covar_bar += &gamma * snc_matrix * &gamma.transpose();
140                // And break so we don't add any more process noise
141                break;
142            }
143        }
144
145        let state_bar = if matches!(self.variant, KalmanVariant::DeviationTracking) {
146            stm * self.prev_estimate.state_deviation
147        } else {
148            OVector::<f64, <T as State>::Size>::zeros()
149        };
150        let estimate = KfEstimate {
151            nominal_state,
152            state_deviation: state_bar,
153            covar: covar_bar,
154            covar_bar,
155            stm,
156            predicted: true,
157        };
158        self.prev_estimate = estimate;
159        // Update the prev epoch for all SNCs
160        for snc in &mut self.process_noise {
161            snc.prev_epoch = Some(self.prev_estimate.epoch());
162        }
163        Ok(estimate)
164    }
165
166    /// Computes the measurement update with a provided real observation and computed observation.
167    ///
168    /// May return a FilterError if the STM or sensitivity matrices were not updated.
169    pub fn measurement_update<M: DimName>(
170        &mut self,
171        nominal_state: T,
172        real_obs: OVector<f64, M>,
173        computed_obs: OVector<f64, M>,
174        r_k: OMatrix<f64, M, M>,
175        h_tilde: OMatrix<f64, M, <T as State>::Size>,
176        resid_rejection: Option<ResidRejectCrit>,
177    ) -> Result<
178        (
179            KfEstimate<T>,
180            Residual<M>,
181            Option<OMatrix<f64, <T as State>::Size, M>>,
182        ),
183        ODError,
184    >
185    where
186        DefaultAllocator: Allocator<M>
187            + Allocator<M, M>
188            + Allocator<M, <T as State>::Size>
189            + Allocator<<T as State>::Size, M>
190            + Allocator<nalgebra::Const<1>, M>,
191    {
192        let epoch = nominal_state.epoch();
193
194        // Grab the state transition matrix.
195        let stm = nominal_state.stm().context(ODDynamicsSnafu)?;
196
197        // Propagate the covariance.
198        let mut covar_bar = stm * self.prev_estimate.covar * stm.transpose();
199
200        // Apply any process noise as in a normal time update, if applicable
201        for (i, snc) in self.process_noise.iter().enumerate().rev() {
202            if let Some(mut snc_matrix) = snc.to_matrix(epoch) {
203                if let Some(local_frame) = snc.local_frame {
204                    // Rotate the SNC from the definition frame into the state frame.
205                    let dcm = local_frame
206                        .dcm_to_inertial(nominal_state.orbit())
207                        .context(AstroPhysicsSnafu)
208                        .context(StateAstroSnafu {
209                            param: StateParameter::Epoch,
210                        })
211                        .context(ODStateSnafu {
212                            action: "rotating SNC from definition frame into state frame",
213                        })?;
214
215                    // Note: the SNC must be a diagonal matrix, so we only update the diagonals!
216                    match A::USIZE {
217                        3 => {
218                            let new_snc = dcm.rot_mat
219                                * snc_matrix.fixed_view::<3, 3>(0, 0)
220                                * dcm.rot_mat.transpose();
221                            for i in 0..A::USIZE {
222                                snc_matrix[(i, i)] = new_snc[(i, i)];
223                            }
224                        }
225                        6 => {
226                            let new_snc = dcm.state_dcm()
227                                * snc_matrix.fixed_view::<6, 6>(0, 0)
228                                * dcm.transpose().state_dcm();
229                            for i in 0..A::USIZE {
230                                snc_matrix[(i, i)] = new_snc[(i, i)];
231                            }
232                        }
233                        _ => {
234                            return Err(ODError::ODLimitation {
235                                action: "only process noises of size 3x3 or 6x6 are supported",
236                            })
237                        }
238                    }
239                }
240                // Check if we're using another SNC than the one before
241                if self.prev_used_snc != i {
242                    info!("Switched to {}-th {}", i, snc);
243                    self.prev_used_snc = i;
244                }
245
246                // Let's compute the Gamma matrix, an approximation of the time integral
247                // which assumes that the acceleration is constant between these two measurements.
248                let mut gamma = OMatrix::<f64, <T as State>::Size, A>::zeros();
249                let delta_t = (nominal_state.epoch() - self.prev_estimate.epoch()).to_seconds();
250                for blk in 0..A::dim() / 3 {
251                    for i in 0..3 {
252                        let idx_i = i + A::dim() * blk;
253                        let idx_j = i + 3 * blk;
254                        let idx_k = i + 3 + A::dim() * blk;
255                        // For first block
256                        // (0, 0) (1, 1) (2, 2) <=> \Delta t^2/2
257                        // (3, 0) (4, 1) (5, 2) <=> \Delta t
258                        // Second block
259                        // (6, 3) (7, 4) (8, 5) <=> \Delta t^2/2
260                        // (9, 3) (10, 4) (11, 5) <=> \Delta t
261                        // * \Delta t^2/2
262                        // (i, i) when blk = 0
263                        // (i + A::dim() * blk, i + 3) when blk = 1
264                        // (i + A::dim() * blk, i + 3 * blk)
265                        // * \Delta t
266                        // (i + 3, i) when blk = 0
267                        // (i + 3, i + 9) when blk = 1 (and I think i + 12 + 3)
268                        // (i + 3 + A::dim() * blk, i + 3 * blk)
269                        gamma[(idx_i, idx_j)] = delta_t.powi(2) / 2.0;
270                        gamma[(idx_k, idx_j)] = delta_t;
271                    }
272                }
273                // Let's add the process noise
274                covar_bar += &gamma * snc_matrix * &gamma.transpose();
275                // And break so we don't add any more process noise
276                break;
277            }
278        }
279
280        // Project the propagated covariance into the measurement space.
281        let h_p_ht = &h_tilde * covar_bar * &h_tilde.transpose();
282
283        // Compute the innovation matrix (S_k).
284        let s_k = &h_p_ht + &r_k;
285
286        // Compute observation deviation/error (usually marked as y_i)
287        let prefit = real_obs.clone() - computed_obs.clone();
288
289        // Compute the prefit ratio for the automatic rejection.
290        // The measurement covariance is the square of the measurement itself.
291        // So we compute its Cholesky decomposition to return to the non squared values.
292        let r_k_chol = match s_k.clone().cholesky() {
293            Some(r_k_clone) => r_k_clone.l(),
294            None => {
295                // In very rare case, when there isn't enough noise in the measurements,
296                // the inverting of S_k fails. If so, we revert back to the nominal Kalman derivation.
297                r_k.clone().cholesky().ok_or(ODError::SingularNoiseRk)?.l()
298            }
299        };
300
301        // Compute the ratio as the average of each component of the prefit over the square root of the measurement
302        // matrix r_k. Refer to ODTK MathSpec equation 4.10.
303        let ratio = s_k
304            .diagonal()
305            .iter()
306            .copied()
307            .enumerate()
308            .map(|(idx, r)| prefit[idx] / r.sqrt())
309            .sum::<f64>()
310            / (M::USIZE as f64);
311
312        if let Some(resid_reject) = resid_rejection {
313            if ratio.abs() > resid_reject.num_sigmas {
314                // Reject this whole measurement and perform only a time update
315                let pred_est = self.time_update(nominal_state)?;
316                return Ok((
317                    pred_est,
318                    Residual::rejected(
319                        epoch,
320                        prefit,
321                        ratio,
322                        r_k_chol.diagonal(),
323                        real_obs,
324                        computed_obs,
325                    ),
326                    None,
327                ));
328            }
329        }
330
331        // Invert the innovation covariance.
332        let s_k_inv = match s_k.try_inverse() {
333            Some(s_k_inv) => s_k_inv,
334            None => return Err(ODError::SingularKalmanGain),
335        };
336
337        let gain = covar_bar * &h_tilde.transpose() * &s_k_inv;
338
339        // Compute the state estimate, depends on the variant.
340        let (state_hat, res) = match self.variant {
341            KalmanVariant::ReferenceUpdate => {
342                // In EKF, the state hat is actually the state deviation. We trust the gain to be correct,
343                // so we just apply it directly to the prefit residual.
344                let state_hat = &gain * &prefit;
345                let postfit = &prefit - (&h_tilde * state_hat);
346                (
347                    state_hat,
348                    Residual::accepted(
349                        epoch,
350                        prefit,
351                        postfit,
352                        ratio,
353                        r_k_chol.diagonal(),
354                        real_obs,
355                        computed_obs,
356                    ),
357                )
358            }
359            KalmanVariant::DeviationTracking => {
360                // Time update
361                let state_bar = stm * self.prev_estimate.state_deviation;
362                let postfit = &prefit - (&h_tilde * state_bar);
363                (
364                    state_bar + &gain * &postfit,
365                    Residual::accepted(
366                        epoch,
367                        prefit,
368                        postfit,
369                        ratio,
370                        r_k_chol.diagonal(),
371                        real_obs,
372                        computed_obs,
373                    ),
374                )
375            }
376        };
377
378        // Compute covariance (Joseph update)
379        let first_term =
380            OMatrix::<f64, <T as State>::Size, <T as State>::Size>::identity() - &gain * &h_tilde;
381        let covar =
382            first_term * covar_bar * first_term.transpose() + &gain * &r_k * &gain.transpose();
383
384        // And wrap up
385        let estimate = KfEstimate {
386            nominal_state,
387            state_deviation: state_hat,
388            covar,
389            covar_bar,
390            stm,
391            predicted: false,
392        };
393
394        self.prev_estimate = estimate;
395        // Update the prev epoch for all SNCs
396        for snc in &mut self.process_noise {
397            snc.prev_epoch = Some(self.prev_estimate.epoch());
398        }
399
400        Ok((estimate, res, Some(gain)))
401    }
402
403    pub fn replace_state(&self) -> bool {
404        matches!(self.variant, KalmanVariant::ReferenceUpdate)
405    }
406
407    /// Overwrites all of the process noises to the one provided
408    pub fn set_process_noise(&mut self, snc: ProcessNoise<A>) {
409        self.process_noise = vec![snc];
410    }
411}