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