brahe 1.4.0

Brahe is a modern satellite dynamics library for research and engineering applications designed to be easy-to-learn, high-performance, and quick-to-deploy. The north-star of the development is enabling users to solve meaningful problems and answer questions quickly, easily, and correctly.
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
/*!
 * Sensitivity matrix computation for parameter estimation and consider covariance analysis.
 *
 * This module provides trait-based interfaces for computing sensitivity matrices (∂f/∂p),
 * with implementations for both analytical and numerical (finite difference) methods.
 *
 * Sensitivity matrices describe how a function's output changes with respect to consider
 * parameters, which is essential for orbit determination with consider parameters and
 * covariance analysis.
 */

use crate::math::jacobian::{DifferenceMethod, PerturbationStrategy};
use nalgebra::{DMatrix, DVector, SMatrix, SVector};

// ============================================================================
// Type Aliases for Sensitivity Functions
// ============================================================================

/// Sensitivity function type for static-sized systems.
type SSensitivityFn<const S: usize, const P: usize> =
    Box<dyn Fn(f64, &SVector<f64, S>, &SVector<f64, P>) -> SMatrix<f64, S, P> + Send + Sync>;

/// Sensitivity function type for dynamic-sized systems.
type DSensitivityFn = Box<dyn Fn(f64, &DVector<f64>, &DVector<f64>) -> DMatrix<f64> + Send + Sync>;

/// Dynamics function type for static-sized sensitivity computation.
type SDynamicsWithParams<const S: usize, const P: usize> =
    Box<dyn Fn(f64, &SVector<f64, S>, &SVector<f64, P>) -> SVector<f64, S> + Send + Sync>;

/// Dynamics function type for dynamic-sized sensitivity computation.
type DDynamicsWithParams =
    Box<dyn Fn(f64, &DVector<f64>, &DVector<f64>) -> DVector<f64> + Send + Sync>;

/// Trait for static-sized sensitivity providers.
///
/// Computes the sensitivity matrix ∂f/∂p where f is the dynamics function
/// and p are the consider parameters.
///
/// # Thread Safety
///
/// Requires `Send + Sync` for thread-safe integrator usage.
pub trait SSensitivityProvider<const S: usize, const P: usize>: Send + Sync {
    /// Compute the sensitivity matrix at the given time, state, and parameters.
    ///
    /// # Arguments
    /// - `t`: Current time
    /// - `state`: State vector at time t
    /// - `params`: Consider parameters
    ///
    /// # Returns
    /// Sensitivity matrix ∂f/∂p (S×P matrix)
    fn compute(
        &self,
        t: f64,
        state: &SVector<f64, S>,
        params: &SVector<f64, P>,
    ) -> SMatrix<f64, S, P>;
}

/// Trait for dynamic-sized sensitivity providers.
///
/// Computes the sensitivity matrix ∂f/∂p where f is the dynamics function
/// and p are the consider parameters.
pub trait DSensitivityProvider: Send + Sync {
    /// Compute the sensitivity matrix at the given time, state, and parameters.
    ///
    /// # Arguments
    /// - `t`: Current time
    /// - `state`: State vector at time t
    /// - `params`: Consider parameters
    ///
    /// # Returns
    /// Sensitivity matrix ∂f/∂p (S×P matrix)
    fn compute(&self, t: f64, state: &DVector<f64>, params: &DVector<f64>) -> DMatrix<f64>;
}

/// Analytical sensitivity provider for static-sized systems.
///
/// Uses a user-provided function that directly computes the analytical sensitivity matrix.
pub struct SAnalyticSensitivity<const S: usize, const P: usize> {
    sensitivity_fn: SSensitivityFn<S, P>,
}

impl<const S: usize, const P: usize> SAnalyticSensitivity<S, P> {
    /// Create a new analytical sensitivity provider.
    pub fn new(sensitivity_fn: SSensitivityFn<S, P>) -> Self {
        Self { sensitivity_fn }
    }
}

impl<const S: usize, const P: usize> SSensitivityProvider<S, P> for SAnalyticSensitivity<S, P> {
    fn compute(
        &self,
        t: f64,
        state: &SVector<f64, S>,
        params: &SVector<f64, P>,
    ) -> SMatrix<f64, S, P> {
        (self.sensitivity_fn)(t, state, params)
    }
}

/// Analytical sensitivity provider for dynamic-sized systems.
///
/// Uses a user-provided function that directly computes the analytical sensitivity matrix.
pub struct DAnalyticSensitivity {
    sensitivity_fn: DSensitivityFn,
}

impl DAnalyticSensitivity {
    /// Create a new analytical sensitivity provider.
    pub fn new(sensitivity_fn: DSensitivityFn) -> Self {
        Self { sensitivity_fn }
    }
}

impl DSensitivityProvider for DAnalyticSensitivity {
    fn compute(&self, t: f64, state: &DVector<f64>, params: &DVector<f64>) -> DMatrix<f64> {
        (self.sensitivity_fn)(t, state, params)
    }
}

/// Numerical sensitivity provider for static-sized systems using finite differences.
///
/// Computes the sensitivity matrix numerically by perturbing the parameters
/// and evaluating the dynamics function.
pub struct SNumericalSensitivity<const S: usize, const P: usize> {
    dynamics_fn: SDynamicsWithParams<S, P>,
    method: DifferenceMethod,
    strategy: PerturbationStrategy,
}

impl<const S: usize, const P: usize> SNumericalSensitivity<S, P> {
    /// Create a new numerical sensitivity provider with default settings.
    ///
    /// Uses central differences with adaptive perturbation strategy.
    pub fn new(dynamics_fn: SDynamicsWithParams<S, P>) -> Self {
        Self {
            dynamics_fn,
            method: DifferenceMethod::Central,
            strategy: PerturbationStrategy::Adaptive {
                scale_factor: 1.0,
                min_value: 1.0,
            },
        }
    }

    /// Create with central differences (default, most accurate).
    pub fn central(dynamics_fn: SDynamicsWithParams<S, P>) -> Self {
        Self::new(dynamics_fn)
    }

    /// Create with forward differences (faster but less accurate).
    pub fn forward(dynamics_fn: SDynamicsWithParams<S, P>) -> Self {
        Self {
            dynamics_fn,
            method: DifferenceMethod::Forward,
            strategy: PerturbationStrategy::Adaptive {
                scale_factor: 1.0,
                min_value: 1.0,
            },
        }
    }

    /// Create with backward differences.
    pub fn backward(dynamics_fn: SDynamicsWithParams<S, P>) -> Self {
        Self {
            dynamics_fn,
            method: DifferenceMethod::Backward,
            strategy: PerturbationStrategy::Adaptive {
                scale_factor: 1.0,
                min_value: 1.0,
            },
        }
    }

    /// Set the perturbation strategy.
    pub fn with_strategy(mut self, strategy: PerturbationStrategy) -> Self {
        self.strategy = strategy;
        self
    }

    fn compute_perturbation(&self, value: f64) -> f64 {
        match self.strategy {
            PerturbationStrategy::Adaptive {
                scale_factor,
                min_value,
            } => {
                let eps = f64::EPSILON;
                scale_factor * eps.sqrt() * value.abs().max(min_value)
            }
            PerturbationStrategy::Fixed(h) => h,
            PerturbationStrategy::Percentage(pct) => value.abs() * pct,
        }
    }
}

impl<const S: usize, const P: usize> SSensitivityProvider<S, P> for SNumericalSensitivity<S, P> {
    fn compute(
        &self,
        t: f64,
        state: &SVector<f64, S>,
        params: &SVector<f64, P>,
    ) -> SMatrix<f64, S, P> {
        let mut sensitivity = SMatrix::<f64, S, P>::zeros();

        for j in 0..P {
            let h = self.compute_perturbation(params[j]);

            let column = match self.method {
                DifferenceMethod::Forward => {
                    let mut params_plus = *params;
                    params_plus[j] += h;
                    let f_plus = (self.dynamics_fn)(t, state, &params_plus);
                    let f_0 = (self.dynamics_fn)(t, state, params);
                    (f_plus - f_0) / h
                }
                DifferenceMethod::Central => {
                    let mut params_plus = *params;
                    let mut params_minus = *params;
                    params_plus[j] += h;
                    params_minus[j] -= h;
                    let f_plus = (self.dynamics_fn)(t, state, &params_plus);
                    let f_minus = (self.dynamics_fn)(t, state, &params_minus);
                    (f_plus - f_minus) / (2.0 * h)
                }
                DifferenceMethod::Backward => {
                    let mut params_minus = *params;
                    params_minus[j] -= h;
                    let f_0 = (self.dynamics_fn)(t, state, params);
                    let f_minus = (self.dynamics_fn)(t, state, &params_minus);
                    (f_0 - f_minus) / h
                }
            };

            sensitivity.set_column(j, &column);
        }

        sensitivity
    }
}

/// Numerical sensitivity provider for dynamic-sized systems using finite differences.
///
/// Computes the sensitivity matrix numerically by perturbing the parameters
/// and evaluating the dynamics function.
pub struct DNumericalSensitivity {
    dynamics_fn: DDynamicsWithParams,
    method: DifferenceMethod,
    strategy: PerturbationStrategy,
}

impl DNumericalSensitivity {
    /// Create a new numerical sensitivity provider with default settings.
    ///
    /// Uses central differences with adaptive perturbation strategy.
    pub fn new(dynamics_fn: DDynamicsWithParams) -> Self {
        Self {
            dynamics_fn,
            method: DifferenceMethod::Central,
            strategy: PerturbationStrategy::Adaptive {
                scale_factor: 1.0,
                min_value: 1.0,
            },
        }
    }

    /// Create with central differences (default, most accurate).
    pub fn central(dynamics_fn: DDynamicsWithParams) -> Self {
        Self::new(dynamics_fn)
    }

    /// Create with forward differences (faster but less accurate).
    pub fn forward(dynamics_fn: DDynamicsWithParams) -> Self {
        Self {
            dynamics_fn,
            method: DifferenceMethod::Forward,
            strategy: PerturbationStrategy::Adaptive {
                scale_factor: 1.0,
                min_value: 1.0,
            },
        }
    }

    /// Create with backward differences.
    pub fn backward(dynamics_fn: DDynamicsWithParams) -> Self {
        Self {
            dynamics_fn,
            method: DifferenceMethod::Backward,
            strategy: PerturbationStrategy::Adaptive {
                scale_factor: 1.0,
                min_value: 1.0,
            },
        }
    }

    /// Set the perturbation strategy.
    pub fn with_strategy(mut self, strategy: PerturbationStrategy) -> Self {
        self.strategy = strategy;
        self
    }

    fn compute_perturbation(&self, value: f64) -> f64 {
        match self.strategy {
            PerturbationStrategy::Adaptive {
                scale_factor,
                min_value,
            } => {
                let eps = f64::EPSILON;
                scale_factor * eps.sqrt() * value.abs().max(min_value)
            }
            PerturbationStrategy::Fixed(h) => h,
            PerturbationStrategy::Percentage(pct) => value.abs() * pct,
        }
    }
}

impl DSensitivityProvider for DNumericalSensitivity {
    fn compute(&self, t: f64, state: &DVector<f64>, params: &DVector<f64>) -> DMatrix<f64> {
        let s = state.len();
        let p = params.len();
        let mut sensitivity = DMatrix::<f64>::zeros(s, p);

        for j in 0..p {
            let h = self.compute_perturbation(params[j]);

            let column = match self.method {
                DifferenceMethod::Forward => {
                    let mut params_plus = params.clone();
                    params_plus[j] += h;
                    let f_plus = (self.dynamics_fn)(t, state, &params_plus);
                    let f_0 = (self.dynamics_fn)(t, state, params);
                    (f_plus - f_0) / h
                }
                DifferenceMethod::Central => {
                    let mut params_plus = params.clone();
                    let mut params_minus = params.clone();
                    params_plus[j] += h;
                    params_minus[j] -= h;
                    let f_plus = (self.dynamics_fn)(t, state, &params_plus);
                    let f_minus = (self.dynamics_fn)(t, state, &params_minus);
                    (f_plus - f_minus) / (2.0 * h)
                }
                DifferenceMethod::Backward => {
                    let mut params_minus = params.clone();
                    params_minus[j] -= h;
                    let f_0 = (self.dynamics_fn)(t, state, params);
                    let f_minus = (self.dynamics_fn)(t, state, &params_minus);
                    (f_0 - f_minus) / h
                }
            };

            sensitivity.set_column(j, &column);
        }

        sensitivity
    }
}

#[cfg(test)]
mod tests {
    use super::*;
    use approx::assert_abs_diff_eq;

    #[test]
    fn test_dynamic_numerical_sensitivity() {
        // Dynamics: f(t, x, p) = p[0] * x
        let dynamics = |_t: f64, state: &DVector<f64>, params: &DVector<f64>| -> DVector<f64> {
            params[0] * state
        };

        let provider = DNumericalSensitivity::central(Box::new(dynamics));

        let state = DVector::from_vec(vec![1.0, 2.0]);
        let params = DVector::from_vec(vec![3.0]);

        let sens = provider.compute(0.0, &state, &params);

        // ∂f/∂p = x, so sensitivity should be [[1], [2]]
        assert_eq!(sens.nrows(), 2);
        assert_eq!(sens.ncols(), 1);
        assert_abs_diff_eq!(sens[(0, 0)], 1.0, epsilon = 1e-6);
        assert_abs_diff_eq!(sens[(1, 0)], 2.0, epsilon = 1e-6);
    }

    #[test]
    fn test_dynamic_analytical_sensitivity() {
        // Analytical sensitivity: ∂f/∂p = x
        let sensitivity_fn =
            |_t: f64, state: &DVector<f64>, _params: &DVector<f64>| -> DMatrix<f64> {
                DMatrix::from_column_slice(state.len(), 1, state.as_slice())
            };

        let provider = DAnalyticSensitivity::new(Box::new(sensitivity_fn));

        let state = DVector::from_vec(vec![1.0, 2.0]);
        let params = DVector::from_vec(vec![3.0]);

        let sens = provider.compute(0.0, &state, &params);

        assert_eq!(sens.nrows(), 2);
        assert_eq!(sens.ncols(), 1);
        assert_abs_diff_eq!(sens[(0, 0)], 1.0, epsilon = 1e-10);
        assert_abs_diff_eq!(sens[(1, 0)], 2.0, epsilon = 1e-10);
    }

    #[test]
    fn test_static_numerical_sensitivity() {
        // Dynamics: f(t, x, p) = [p[0] * x[0], p[1] * x[1]]
        let dynamics =
            |_t: f64, state: &SVector<f64, 2>, params: &SVector<f64, 2>| -> SVector<f64, 2> {
                SVector::<f64, 2>::new(params[0] * state[0], params[1] * state[1])
            };

        let provider = SNumericalSensitivity::central(Box::new(dynamics));

        let state = SVector::<f64, 2>::new(1.0, 2.0);
        let params = SVector::<f64, 2>::new(3.0, 4.0);

        let sens = provider.compute(0.0, &state, &params);

        // ∂f/∂p = [[x[0], 0], [0, x[1]]]
        assert_abs_diff_eq!(sens[(0, 0)], 1.0, epsilon = 1e-6);
        assert_abs_diff_eq!(sens[(0, 1)], 0.0, epsilon = 1e-6);
        assert_abs_diff_eq!(sens[(1, 0)], 0.0, epsilon = 1e-6);
        assert_abs_diff_eq!(sens[(1, 1)], 2.0, epsilon = 1e-6);
    }

    #[test]
    fn test_static_analytical_sensitivity() {
        // Analytical sensitivity
        let sensitivity_fn =
            |_t: f64, state: &SVector<f64, 2>, _params: &SVector<f64, 2>| -> SMatrix<f64, 2, 2> {
                SMatrix::<f64, 2, 2>::new(state[0], 0.0, 0.0, state[1])
            };

        let provider = SAnalyticSensitivity::new(Box::new(sensitivity_fn));

        let state = SVector::<f64, 2>::new(1.0, 2.0);
        let params = SVector::<f64, 2>::new(3.0, 4.0);

        let sens = provider.compute(0.0, &state, &params);

        assert_abs_diff_eq!(sens[(0, 0)], 1.0, epsilon = 1e-10);
        assert_abs_diff_eq!(sens[(0, 1)], 0.0, epsilon = 1e-10);
        assert_abs_diff_eq!(sens[(1, 0)], 0.0, epsilon = 1e-10);
        assert_abs_diff_eq!(sens[(1, 1)], 2.0, epsilon = 1e-10);
    }

    #[test]
    fn test_forward_vs_central_difference() {
        // Simple quadratic dynamics to test difference methods
        let dynamics = |_t: f64, state: &DVector<f64>, params: &DVector<f64>| -> DVector<f64> {
            DVector::from_vec(vec![params[0] * params[0] * state[0]])
        };

        let state = DVector::from_vec(vec![1.0]);
        let params = DVector::from_vec(vec![2.0]);

        let forward = DNumericalSensitivity::forward(Box::new(dynamics));
        let central = DNumericalSensitivity::central(Box::new(dynamics));

        let sens_forward = forward.compute(0.0, &state, &params);
        let sens_central = central.compute(0.0, &state, &params);

        // Analytical: ∂f/∂p = 2*p*x = 4
        // Central should be more accurate
        assert_abs_diff_eq!(sens_central[(0, 0)], 4.0, epsilon = 1e-8);
        assert_abs_diff_eq!(sens_forward[(0, 0)], 4.0, epsilon = 1e-4);
    }

    #[test]
    fn test_static_numerical_sensitivity_forward() {
        // Dynamics: f(t, x, p) = [p[0] * x[0], p[1] * x[1]]
        let dynamics =
            |_t: f64, state: &SVector<f64, 2>, params: &SVector<f64, 2>| -> SVector<f64, 2> {
                SVector::<f64, 2>::new(params[0] * state[0], params[1] * state[1])
            };

        let provider = SNumericalSensitivity::forward(Box::new(dynamics));

        let state = SVector::<f64, 2>::new(1.0, 2.0);
        let params = SVector::<f64, 2>::new(3.0, 4.0);

        let sens = provider.compute(0.0, &state, &params);

        // ∂f/∂p = [[x[0], 0], [0, x[1]]]
        assert_abs_diff_eq!(sens[(0, 0)], 1.0, epsilon = 1e-6);
        assert_abs_diff_eq!(sens[(0, 1)], 0.0, epsilon = 1e-6);
        assert_abs_diff_eq!(sens[(1, 0)], 0.0, epsilon = 1e-6);
        assert_abs_diff_eq!(sens[(1, 1)], 2.0, epsilon = 1e-6);
    }

    #[test]
    fn test_static_numerical_sensitivity_backward() {
        // Dynamics: f(t, x, p) = [p[0] * x[0], p[1] * x[1]]
        let dynamics =
            |_t: f64, state: &SVector<f64, 2>, params: &SVector<f64, 2>| -> SVector<f64, 2> {
                SVector::<f64, 2>::new(params[0] * state[0], params[1] * state[1])
            };

        let provider = SNumericalSensitivity::backward(Box::new(dynamics));

        let state = SVector::<f64, 2>::new(1.0, 2.0);
        let params = SVector::<f64, 2>::new(3.0, 4.0);

        let sens = provider.compute(0.0, &state, &params);

        // ∂f/∂p = [[x[0], 0], [0, x[1]]]
        assert_abs_diff_eq!(sens[(0, 0)], 1.0, epsilon = 1e-6);
        assert_abs_diff_eq!(sens[(0, 1)], 0.0, epsilon = 1e-6);
        assert_abs_diff_eq!(sens[(1, 0)], 0.0, epsilon = 1e-6);
        assert_abs_diff_eq!(sens[(1, 1)], 2.0, epsilon = 1e-6);
    }

    #[test]
    fn test_static_numerical_sensitivity_with_strategy_fixed() {
        // Dynamics: f(t, x, p) = [p[0] * x[0], p[1] * x[1]]
        let dynamics =
            |_t: f64, state: &SVector<f64, 2>, params: &SVector<f64, 2>| -> SVector<f64, 2> {
                SVector::<f64, 2>::new(params[0] * state[0], params[1] * state[1])
            };

        let provider = SNumericalSensitivity::new(Box::new(dynamics))
            .with_strategy(PerturbationStrategy::Fixed(1e-6));

        let state = SVector::<f64, 2>::new(1.0, 2.0);
        let params = SVector::<f64, 2>::new(3.0, 4.0);

        let sens = provider.compute(0.0, &state, &params);

        // ∂f/∂p = [[x[0], 0], [0, x[1]]]
        assert_abs_diff_eq!(sens[(0, 0)], 1.0, epsilon = 1e-5);
        assert_abs_diff_eq!(sens[(0, 1)], 0.0, epsilon = 1e-5);
        assert_abs_diff_eq!(sens[(1, 0)], 0.0, epsilon = 1e-5);
        assert_abs_diff_eq!(sens[(1, 1)], 2.0, epsilon = 1e-5);
    }

    #[test]
    fn test_static_numerical_sensitivity_with_strategy_percentage() {
        // Dynamics: f(t, x, p) = [p[0] * x[0], p[1] * x[1]]
        let dynamics =
            |_t: f64, state: &SVector<f64, 2>, params: &SVector<f64, 2>| -> SVector<f64, 2> {
                SVector::<f64, 2>::new(params[0] * state[0], params[1] * state[1])
            };

        let provider = SNumericalSensitivity::new(Box::new(dynamics))
            .with_strategy(PerturbationStrategy::Percentage(1e-6));

        let state = SVector::<f64, 2>::new(1.0, 2.0);
        let params = SVector::<f64, 2>::new(3.0, 4.0);

        let sens = provider.compute(0.0, &state, &params);

        // ∂f/∂p = [[x[0], 0], [0, x[1]]]
        assert_abs_diff_eq!(sens[(0, 0)], 1.0, epsilon = 1e-5);
        assert_abs_diff_eq!(sens[(0, 1)], 0.0, epsilon = 1e-5);
        assert_abs_diff_eq!(sens[(1, 0)], 0.0, epsilon = 1e-5);
        assert_abs_diff_eq!(sens[(1, 1)], 2.0, epsilon = 1e-5);
    }

    #[test]
    fn test_dynamic_numerical_sensitivity_backward() {
        // Dynamics: f(t, x, p) = p[0] * x
        let dynamics = |_t: f64, state: &DVector<f64>, params: &DVector<f64>| -> DVector<f64> {
            params[0] * state
        };

        let provider = DNumericalSensitivity::backward(Box::new(dynamics));

        let state = DVector::from_vec(vec![1.0, 2.0]);
        let params = DVector::from_vec(vec![3.0]);

        let sens = provider.compute(0.0, &state, &params);

        // ∂f/∂p = x, so sensitivity should be [[1], [2]]
        assert_eq!(sens.nrows(), 2);
        assert_eq!(sens.ncols(), 1);
        assert_abs_diff_eq!(sens[(0, 0)], 1.0, epsilon = 1e-6);
        assert_abs_diff_eq!(sens[(1, 0)], 2.0, epsilon = 1e-6);
    }

    #[test]
    fn test_dynamic_numerical_sensitivity_with_strategy() {
        // Dynamics: f(t, x, p) = p[0] * x
        let dynamics = |_t: f64, state: &DVector<f64>, params: &DVector<f64>| -> DVector<f64> {
            params[0] * state
        };

        let provider = DNumericalSensitivity::new(Box::new(dynamics))
            .with_strategy(PerturbationStrategy::Fixed(1e-6));

        let state = DVector::from_vec(vec![1.0, 2.0]);
        let params = DVector::from_vec(vec![3.0]);

        let sens = provider.compute(0.0, &state, &params);

        // ∂f/∂p = x, so sensitivity should be [[1], [2]]
        assert_eq!(sens.nrows(), 2);
        assert_eq!(sens.ncols(), 1);
        assert_abs_diff_eq!(sens[(0, 0)], 1.0, epsilon = 1e-5);
        assert_abs_diff_eq!(sens[(1, 0)], 2.0, epsilon = 1e-5);
    }

    #[test]
    fn test_dynamic_numerical_sensitivity_fixed_perturbation() {
        // Test that Fixed perturbation strategy is used correctly
        // Dynamics: f(t, x, p) = p[0]^2 * x
        let dynamics = |_t: f64, state: &DVector<f64>, params: &DVector<f64>| -> DVector<f64> {
            params[0] * params[0] * state
        };

        let provider = DNumericalSensitivity::central(Box::new(dynamics))
            .with_strategy(PerturbationStrategy::Fixed(1e-5));

        let state = DVector::from_vec(vec![1.0]);
        let params = DVector::from_vec(vec![2.0]);

        let sens = provider.compute(0.0, &state, &params);

        // ∂f/∂p = 2*p*x = 4
        assert_abs_diff_eq!(sens[(0, 0)], 4.0, epsilon = 1e-4);
    }

    #[test]
    fn test_dynamic_numerical_sensitivity_percentage_perturbation() {
        // Test that Percentage perturbation strategy is used correctly
        // Dynamics: f(t, x, p) = p[0]^2 * x
        let dynamics = |_t: f64, state: &DVector<f64>, params: &DVector<f64>| -> DVector<f64> {
            params[0] * params[0] * state
        };

        let provider = DNumericalSensitivity::central(Box::new(dynamics))
            .with_strategy(PerturbationStrategy::Percentage(1e-6));

        let state = DVector::from_vec(vec![1.0]);
        let params = DVector::from_vec(vec![2.0]);

        let sens = provider.compute(0.0, &state, &params);

        // ∂f/∂p = 2*p*x = 4
        assert_abs_diff_eq!(sens[(0, 0)], 4.0, epsilon = 1e-4);
    }

    #[test]
    fn test_static_forward_vs_central_vs_backward() {
        // Quadratic dynamics to compare difference methods
        let dynamics =
            |_t: f64, state: &SVector<f64, 1>, params: &SVector<f64, 1>| -> SVector<f64, 1> {
                SVector::<f64, 1>::new(params[0] * params[0] * state[0])
            };

        let state = SVector::<f64, 1>::new(1.0);
        let params = SVector::<f64, 1>::new(2.0);

        let forward = SNumericalSensitivity::forward(Box::new(dynamics));
        let central = SNumericalSensitivity::central(Box::new(dynamics));
        let backward = SNumericalSensitivity::backward(Box::new(dynamics));

        let sens_forward = forward.compute(0.0, &state, &params);
        let sens_central = central.compute(0.0, &state, &params);
        let sens_backward = backward.compute(0.0, &state, &params);

        // Analytical: ∂f/∂p = 2*p*x = 4
        // Central should be most accurate
        assert_abs_diff_eq!(sens_central[(0, 0)], 4.0, epsilon = 1e-8);
        assert_abs_diff_eq!(sens_forward[(0, 0)], 4.0, epsilon = 1e-4);
        assert_abs_diff_eq!(sens_backward[(0, 0)], 4.0, epsilon = 1e-4);
    }

    #[test]
    fn test_dynamic_backward_vs_central() {
        // Quadratic dynamics to compare backward and central
        let dynamics = |_t: f64, state: &DVector<f64>, params: &DVector<f64>| -> DVector<f64> {
            DVector::from_vec(vec![params[0] * params[0] * state[0]])
        };

        let state = DVector::from_vec(vec![1.0]);
        let params = DVector::from_vec(vec![2.0]);

        let backward = DNumericalSensitivity::backward(Box::new(dynamics));
        let central = DNumericalSensitivity::central(Box::new(dynamics));

        let sens_backward = backward.compute(0.0, &state, &params);
        let sens_central = central.compute(0.0, &state, &params);

        // Analytical: ∂f/∂p = 2*p*x = 4
        // Central should be more accurate
        assert_abs_diff_eq!(sens_central[(0, 0)], 4.0, epsilon = 1e-8);
        assert_abs_diff_eq!(sens_backward[(0, 0)], 4.0, epsilon = 1e-4);
    }
}