1pub use crate::errors::NyxError;
20use crate::linalg::allocator::Allocator;
21use crate::linalg::{DefaultAllocator, DimName, OMatrix, OVector, U3};
22pub use crate::od::estimate::{Estimate, KfEstimate, Residual};
23use crate::od::process::ResidRejectCrit;
24pub use crate::od::snc::SNC;
25use crate::od::{Filter, ODDynamicsSnafu, ODError, State};
26pub use crate::time::{Epoch, Unit};
27use snafu::prelude::*;
28
29#[derive(Debug, Clone)]
34#[allow(clippy::upper_case_acronyms)]
35pub struct KF<T, A, M>
36where
37 A: DimName,
38 M: DimName,
39 T: State,
40 DefaultAllocator: Allocator<M>
41 + Allocator<<T as State>::Size>
42 + Allocator<<T as State>::VecLength>
43 + Allocator<A>
44 + Allocator<M, M>
45 + Allocator<M, <T as State>::Size>
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 <DefaultAllocator as Allocator<<T as State>::Size>>::Buffer<f64>: Copy,
51 <DefaultAllocator as Allocator<<T as State>::Size, <T as State>::Size>>::Buffer<f64>: Copy,
52{
53 pub prev_estimate: KfEstimate<T>,
55 pub process_noise: Vec<SNC<A>>,
57 pub ekf: bool,
60 h_tilde: OMatrix<f64, M, <T as State>::Size>,
61 h_tilde_updated: bool,
62 prev_used_snc: usize,
63}
64
65impl<T, A, M> KF<T, A, M>
66where
67 A: DimName,
68 M: DimName,
69 T: State,
70 DefaultAllocator: Allocator<M>
71 + Allocator<<T as State>::Size>
72 + Allocator<<T as State>::VecLength>
73 + Allocator<A>
74 + Allocator<M, M>
75 + Allocator<M, <T as State>::Size>
76 + Allocator<<T as State>::Size, M>
77 + Allocator<<T as State>::Size, <T as State>::Size>
78 + Allocator<A, A>
79 + Allocator<<T as State>::Size, A>
80 + Allocator<A, <T as State>::Size>,
81 <DefaultAllocator as Allocator<<T as State>::Size>>::Buffer<f64>: Copy,
82 <DefaultAllocator as Allocator<<T as State>::Size, <T as State>::Size>>::Buffer<f64>: Copy,
83{
84 pub fn new(initial_estimate: KfEstimate<T>, process_noise: SNC<A>) -> Self {
86 assert_eq!(
87 A::dim() % 3,
88 0,
89 "SNC can only be applied to accelerations multiple of 3"
90 );
91
92 let mut process_noise = process_noise;
94 process_noise.init_epoch = Some(initial_estimate.epoch());
95
96 Self {
97 prev_estimate: initial_estimate,
98 process_noise: vec![process_noise],
99 ekf: false,
100 h_tilde: OMatrix::<f64, M, <T as State>::Size>::zeros(),
101 h_tilde_updated: false,
102 prev_used_snc: 0,
103 }
104 }
105
106 pub fn with_sncs(initial_estimate: KfEstimate<T>, process_noises: Vec<SNC<A>>) -> Self {
110 assert_eq!(
111 A::dim() % 3,
112 0,
113 "SNC can only be applied to accelerations multiple of 3"
114 );
115 let mut process_noises = process_noises;
116 for snc in &mut process_noises {
118 snc.init_epoch = Some(initial_estimate.epoch());
119 }
120
121 Self {
122 prev_estimate: initial_estimate,
123 process_noise: process_noises,
124 ekf: false,
125 h_tilde: OMatrix::<f64, M, <T as State>::Size>::zeros(),
126 h_tilde_updated: false,
127 prev_used_snc: 0,
128 }
129 }
130}
131
132impl<T, M> KF<T, U3, M>
133where
134 M: DimName,
135 T: State,
136 DefaultAllocator: Allocator<M>
137 + Allocator<<T as State>::Size>
138 + Allocator<<T as State>::VecLength>
139 + Allocator<M, M>
140 + Allocator<M, <T as State>::Size>
141 + Allocator<<T as State>::Size, M>
142 + Allocator<<T as State>::Size, <T as State>::Size>
143 + Allocator<U3, U3>
144 + Allocator<<T as State>::Size, U3>
145 + Allocator<U3, <T as State>::Size>,
146 <DefaultAllocator as Allocator<<T as State>::Size>>::Buffer<f64>: Copy,
147 <DefaultAllocator as Allocator<<T as State>::Size, <T as State>::Size>>::Buffer<f64>: Copy,
148{
149 pub fn no_snc(initial_estimate: KfEstimate<T>) -> Self {
151 Self {
152 prev_estimate: initial_estimate,
153 process_noise: Vec::new(),
154 ekf: false,
155 h_tilde: OMatrix::<f64, M, <T as State>::Size>::zeros(),
156 h_tilde_updated: false,
157 prev_used_snc: 0,
158 }
159 }
160}
161
162impl<T, A, M> Filter<T, A, M> for KF<T, A, M>
163where
164 A: DimName,
165 M: DimName,
166 T: State,
167 DefaultAllocator: Allocator<M>
168 + Allocator<<T as State>::Size>
169 + Allocator<<T as State>::VecLength>
170 + Allocator<A>
171 + Allocator<M, M>
172 + Allocator<M, <T as State>::Size>
173 + Allocator<<T as State>::Size, M>
174 + Allocator<<T as State>::Size, <T as State>::Size>
175 + Allocator<A, A>
176 + Allocator<<T as State>::Size, A>
177 + Allocator<A, <T as State>::Size>
178 + Allocator<nalgebra::Const<1>, M>,
179 <DefaultAllocator as Allocator<<T as State>::Size>>::Buffer<f64>: Copy,
180 <DefaultAllocator as Allocator<<T as State>::Size, <T as State>::Size>>::Buffer<f64>: Copy,
181{
182 type Estimate = KfEstimate<T>;
183
184 fn previous_estimate(&self) -> &Self::Estimate {
186 &self.prev_estimate
187 }
188
189 fn set_previous_estimate(&mut self, est: &Self::Estimate) {
190 self.prev_estimate = *est;
191 }
192
193 fn update_h_tilde(&mut self, h_tilde: OMatrix<f64, M, <T as State>::Size>) {
196 self.h_tilde = h_tilde;
197 self.h_tilde_updated = true;
198 }
199
200 fn time_update(&mut self, nominal_state: T) -> Result<Self::Estimate, ODError> {
204 let stm = nominal_state.stm().context(ODDynamicsSnafu)?;
205 let mut covar_bar = stm * self.prev_estimate.covar * stm.transpose();
206
207 for (i, snc) in self.process_noise.iter().enumerate().rev() {
209 if let Some(snc_matrix) = snc.to_matrix(nominal_state.epoch()) {
210 if self.prev_used_snc != i {
212 info!("Switched to {}-th {}", i, snc);
213 self.prev_used_snc = i;
214 }
215
216 let mut gamma = OMatrix::<f64, <T as State>::Size, A>::zeros();
219 let delta_t = (nominal_state.epoch() - self.prev_estimate.epoch()).to_seconds();
220 for blk in 0..A::dim() / 3 {
221 for i in 0..3 {
222 let idx_i = i + A::dim() * blk;
223 let idx_j = i + 3 * blk;
224 let idx_k = i + 3 + A::dim() * blk;
225 gamma[(idx_i, idx_j)] = delta_t.powi(2) / 2.0;
240 gamma[(idx_k, idx_j)] = delta_t;
241 }
242 }
243 covar_bar += &gamma * snc_matrix * &gamma.transpose();
245 break;
247 }
248 }
249
250 let state_bar = if self.ekf {
251 OVector::<f64, <T as State>::Size>::zeros()
252 } else {
253 stm * self.prev_estimate.state_deviation
254 };
255 let estimate = KfEstimate {
256 nominal_state,
257 state_deviation: state_bar,
258 covar: covar_bar,
259 covar_bar,
260 stm,
261 predicted: true,
262 };
263 self.prev_estimate = estimate;
264 for snc in &mut self.process_noise {
266 snc.prev_epoch = Some(self.prev_estimate.epoch());
267 }
268 Ok(estimate)
269 }
270
271 fn measurement_update(
275 &mut self,
276 nominal_state: T,
277 real_obs: &OVector<f64, M>,
278 computed_obs: &OVector<f64, M>,
279 r_k: OMatrix<f64, M, M>,
280 resid_rejection: Option<ResidRejectCrit>,
281 ) -> Result<(Self::Estimate, Residual<M>), ODError> {
282 if !self.h_tilde_updated {
283 return Err(ODError::SensitivityNotUpdated);
284 }
285
286 let stm = nominal_state.stm().context(ODDynamicsSnafu)?;
287
288 let epoch = nominal_state.epoch();
289
290 let covar_bar = stm * self.prev_estimate.covar * stm.transpose();
291
292 let h_tilde_t = &self.h_tilde.transpose();
293 let h_p_ht = &self.h_tilde * covar_bar * h_tilde_t;
294
295 let s_k = &h_p_ht + &r_k;
296
297 let prefit = real_obs - computed_obs;
299
300 let r_k_chol = match s_k.clone().cholesky() {
304 Some(r_k_clone) => r_k_clone.l(),
305 None => {
306 r_k.clone().cholesky().ok_or(ODError::SingularNoiseRk)?.l()
309 }
310 };
311
312 let ratio = s_k
315 .diagonal()
316 .iter()
317 .copied()
318 .enumerate()
319 .map(|(idx, r)| prefit[idx] / r.sqrt())
320 .sum::<f64>()
321 / (M::USIZE as f64);
322
323 if let Some(resid_reject) = resid_rejection {
324 if ratio.abs() > resid_reject.num_sigmas {
325 let pred_est = self.time_update(nominal_state)?;
327 return Ok((
328 pred_est,
329 Residual::rejected(epoch, prefit, ratio, r_k_chol.diagonal()),
330 ));
331 }
332 }
333
334 let mut innovation_covar = h_p_ht + &s_k;
340 if !innovation_covar.try_inverse_mut() {
341 return Err(ODError::SingularKalmanGain);
342 }
343
344 let gain = covar_bar * h_tilde_t * &innovation_covar;
345
346 let (state_hat, res) = if self.ekf {
348 let state_hat = &gain * &prefit;
349 let postfit = &prefit - (&self.h_tilde * state_hat);
350 (
351 state_hat,
352 Residual::accepted(epoch, prefit, postfit, ratio, r_k_chol.diagonal()),
353 )
354 } else {
355 let state_bar = stm * self.prev_estimate.state_deviation;
357 let postfit = &prefit - (&self.h_tilde * state_bar);
358 (
359 state_bar + &gain * &postfit,
360 Residual::accepted(epoch, prefit, postfit, ratio, r_k_chol.diagonal()),
361 )
362 };
363
364 let first_term = OMatrix::<f64, <T as State>::Size, <T as State>::Size>::identity()
366 - &gain * &self.h_tilde;
367 let covar =
368 first_term * covar_bar * first_term.transpose() + &gain * &s_k * &gain.transpose();
369
370 let estimate = KfEstimate {
372 nominal_state,
373 state_deviation: state_hat,
374 covar,
375 covar_bar,
376 stm,
377 predicted: false,
378 };
379
380 self.h_tilde_updated = false;
381 self.prev_estimate = estimate;
382 for snc in &mut self.process_noise {
384 snc.prev_epoch = Some(self.prev_estimate.epoch());
385 }
386 Ok((estimate, res))
387 }
388
389 fn is_extended(&self) -> bool {
390 self.ekf
391 }
392
393 fn set_extended(&mut self, status: bool) {
394 self.ekf = status;
395 }
396
397 fn set_process_noise(&mut self, snc: SNC<A>) {
399 self.process_noise = vec![snc];
400 }
401}