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
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
/*!
 * Batch Least Squares (BLS) estimator for orbit determination.
 *
 * The BLS estimator processes all observations simultaneously through an
 * iterative Gauss-Newton algorithm. At each iteration it propagates the
 * current state estimate through all observation epochs, accumulates
 * measurement residuals and their Jacobians mapped back to the reference
 * epoch via the State Transition Matrix (STM), and solves for a state
 * correction that minimizes the weighted sum of squared residuals.
 *
 * # Algorithm (Gauss-Newton iteration)
 *
 * 1. Reinitialize the propagator at (t₀, x_k)
 * 2. Propagate sequentially through each observation epoch t_i, collecting:
 *    - Residual: δy_i = y_i - h(x_k(t_i))
 *    - Mapped Jacobian: H_i = H̃_i * Φ(t_i, t₀) (solve-for partition)
 * 3. Solve for state correction δx via either:
 *    - **Normal Equations**: Λ = P̄₀⁻¹ + Σ Hᵢᵀ Rᵢ⁻¹ Hᵢ, N = Σ Hᵢᵀ Rᵢ⁻¹ δyᵢ
 *    - **Stacked Observation Matrix**: Build full H and δy, solve HᵀWH δx = HᵀW δy
 * 4. Apply correction: x_{k+1} = x_k + δx
 * 5. Check convergence (state correction norm and/or cost function change)
 *
 * Two solver formulations are available via [`BLSSolverMethod`]:
 * - [`NormalEquations`]: Memory-efficient O(n²), standard Tapley/Schutz/Born
 * - [`StackedObservationMatrix`]: Better numerical conditioning, O(m*n) memory
 *
 * Optional **consider parameters** allow partitioning the state into solve-for
 * and consider subsets, with the consider covariance contribution computed via
 * cross-correlation terms.
 *
 * [`BLSSolverMethod`]: crate::estimation::BLSSolverMethod
 * [`NormalEquations`]: crate::estimation::BLSSolverMethod::NormalEquations
 * [`StackedObservationMatrix`]: crate::estimation::BLSSolverMethod::StackedObservationMatrix
 */

use nalgebra::{DMatrix, DVector};

use crate::integrators::traits::{DControlInput, DStateDynamics};
use crate::propagators::force_model_config::ForceModelConfig;
use crate::propagators::{DNumericalOrbitPropagator, NumericalPropagationConfig};
use crate::time::Epoch;
use crate::utils::errors::BraheError;

use super::config::{BLSConfig, BLSSolverMethod};
use super::dynamics_source::DynamicsSource;
use super::traits::MeasurementModel;
use super::types::{BLSIterationRecord, BLSObservationResidual, Observation};

/// Batch Least Squares estimator for orbit determination.
///
/// Processes all observations simultaneously through an iterative
/// Gauss-Newton algorithm. The state correction at each iteration
/// minimizes the weighted sum of squared residuals subject to the
/// a priori constraint.
///
/// # Examples
///
/// ```no_run
/// use brahe::estimation::*;
/// use brahe::propagators::*;
/// use brahe::time::{Epoch, TimeSystem};
/// use nalgebra::{DMatrix, DVector};
///
/// let epoch = Epoch::from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0, TimeSystem::UTC);
/// let state = DVector::from_vec(vec![6878e3, 0.0, 0.0, 0.0, 7612.0, 0.0]);
/// let p0 = DMatrix::from_diagonal(&DVector::from_vec(vec![1e6, 1e6, 1e6, 1e2, 1e2, 1e2]));
///
/// // let mut bls = BatchLeastSquares::new(
/// //     epoch, state, p0,
/// //     NumericalPropagationConfig::default(),
/// //     ForceModelConfig::two_body_gravity(),
/// //     None, None, None,
/// //     models, BLSConfig::default(),
/// // )?;
/// // bls.solve(&observations)?;
/// ```
pub struct BatchLeastSquares {
    /// Dynamics source (propagator).
    dynamics: DynamicsSource,

    /// Measurement models (supports multiple types)
    measurement_models: Vec<Box<dyn MeasurementModel>>,

    /// BLS configuration
    config: BLSConfig,

    /// A priori epoch (reference epoch for the batch)
    apriori_epoch: Epoch,

    /// A priori state estimate
    apriori_state: DVector<f64>,

    /// A priori covariance matrix
    apriori_covariance: DMatrix<f64>,

    /// Current state estimate (updated each iteration)
    current_state: DVector<f64>,

    /// Current covariance (formal, solve-for partition)
    current_covariance: DMatrix<f64>,

    /// Current epoch (same as apriori_epoch for BLS)
    current_epoch: Epoch,

    /// Whether the solver has converged
    converged: bool,

    /// Number of iterations completed
    iterations_completed: usize,

    /// Final cost function value J
    final_cost: f64,

    /// Per-iteration diagnostic records
    iteration_records: Vec<BLSIterationRecord>,

    /// Per-observation residuals for each iteration
    observation_residuals: Vec<Vec<BLSObservationResidual>>,

    /// Cached consider parameter information (Λ_ss_inv, Λ_sc) for
    /// computing the consider covariance contribution.
    consider_info: Option<(DMatrix<f64>, DMatrix<f64>)>,
}

impl std::fmt::Debug for BatchLeastSquares {
    fn fmt(&self, f: &mut std::fmt::Formatter<'_>) -> std::fmt::Result {
        f.debug_struct("BatchLeastSquares")
            .field("apriori_epoch", &self.apriori_epoch)
            .field("state_dim", &self.current_state.len())
            .field("converged", &self.converged)
            .field("iterations_completed", &self.iterations_completed)
            .field("final_cost", &self.final_cost)
            .field("num_models", &self.measurement_models.len())
            .finish()
    }
}

/// Internal result type for solver formulations.
struct SolverResult {
    /// State correction vector δx
    state_correction: DVector<f64>,
    /// Formal covariance (solve-for partition)
    formal_covariance: DMatrix<f64>,
    /// Consider info: (Λ_ss_inv, Λ_sc) if consider params configured
    consider_info: Option<(DMatrix<f64>, DMatrix<f64>)>,
}

impl BatchLeastSquares {
    /// Create a Batch Least Squares estimator with orbit propagator dynamics.
    ///
    /// Builds a numerical orbit propagator internally with STM enabled for
    /// the state transition matrix computation required by BLS.
    ///
    /// # Arguments
    ///
    /// * `epoch` - Initial (reference) epoch
    /// * `state` - Initial state vector (meters, m/s)
    /// * `initial_covariance` - A priori covariance matrix
    /// * `propagation_config` - Numerical propagation configuration
    /// * `force_config` - Force model configuration
    /// * `params` - Optional parameter vector for force models
    /// * `additional_dynamics` - Optional additional dynamics function
    /// * `control_input` - Optional control input function
    /// * `measurement_models` - List of measurement models
    /// * `config` - BLS configuration
    #[allow(clippy::too_many_arguments)]
    pub fn new(
        epoch: Epoch,
        state: DVector<f64>,
        initial_covariance: DMatrix<f64>,
        propagation_config: NumericalPropagationConfig,
        force_config: ForceModelConfig,
        params: Option<DVector<f64>>,
        additional_dynamics: Option<DStateDynamics>,
        control_input: DControlInput,
        measurement_models: Vec<Box<dyn MeasurementModel>>,
        config: BLSConfig,
    ) -> Result<Self, BraheError> {
        // Force STM enabled for BLS Jacobian mapping
        let mut prop_config = propagation_config;
        prop_config.variational.enable_stm = true;

        let prop = DNumericalOrbitPropagator::new(
            epoch,
            state,
            prop_config,
            force_config,
            params,
            additional_dynamics,
            control_input,
            Some(initial_covariance),
        )
        .map_err(|e| BraheError::Error(format!("Failed to create propagator: {}", e)))?;

        Self::from_propagator(
            DynamicsSource::OrbitPropagator(prop),
            measurement_models,
            config,
        )
    }

    /// Create a Batch Least Squares estimator from an existing propagator.
    ///
    /// Use this when you have a pre-built propagator with custom dynamics
    /// or need full control over propagator configuration.
    ///
    /// # Arguments
    ///
    /// * `propagator` - Pre-built dynamics source (must have STM enabled and covariance set)
    /// * `measurement_models` - List of measurement models
    /// * `config` - BLS configuration
    ///
    /// # Errors
    ///
    /// Returns error if:
    /// - STM propagation is not enabled
    /// - No measurement models are provided
    /// - Covariance dimensions do not match state dimensions
    /// - No convergence threshold is configured
    /// - Consider parameter dimensions are inconsistent
    pub fn from_propagator(
        propagator: DynamicsSource,
        measurement_models: Vec<Box<dyn MeasurementModel>>,
        config: BLSConfig,
    ) -> Result<Self, BraheError> {
        // Validate STM enabled
        if !propagator.has_stm() {
            return Err(BraheError::Error(
                "BatchLeastSquares requires STM propagation to be enabled. \
                 Set propagation_config.variational.enable_stm = true or provide \
                 an initial_covariance to the propagator constructor."
                    .to_string(),
            ));
        }

        // Validate at least one measurement model
        if measurement_models.is_empty() {
            return Err(BraheError::Error(
                "At least one measurement model is required".to_string(),
            ));
        }

        // Validate covariance is set and dimensions match
        let covariance = propagator
            .current_covariance()
            .ok_or_else(|| {
                BraheError::Error(
                    "BatchLeastSquares requires an initial covariance on the propagator. \
                     Provide initial_covariance when constructing the propagator."
                        .to_string(),
                )
            })?
            .clone();

        let state_dim = propagator.state_dim();
        if covariance.nrows() != state_dim || covariance.ncols() != state_dim {
            return Err(BraheError::Error(format!(
                "Covariance dimensions ({}x{}) do not match state dimension ({})",
                covariance.nrows(),
                covariance.ncols(),
                state_dim
            )));
        }

        // Validate at least one convergence threshold is set
        if config.state_correction_threshold.is_none()
            && config.cost_convergence_threshold.is_none()
        {
            return Err(BraheError::Error(
                "At least one convergence threshold must be set: \
                 state_correction_threshold and/or cost_convergence_threshold"
                    .to_string(),
            ));
        }

        // Validate consider parameter dimensions if configured
        if let Some(ref consider) = config.consider_params {
            if consider.n_solve >= state_dim {
                return Err(BraheError::Error(format!(
                    "n_solve ({}) must be less than state_dim ({}) when consider \
                     parameters are configured",
                    consider.n_solve, state_dim
                )));
            }
            let n_consider = state_dim - consider.n_solve;
            if consider.consider_covariance.nrows() != n_consider
                || consider.consider_covariance.ncols() != n_consider
            {
                return Err(BraheError::Error(format!(
                    "Consider covariance dimensions ({}x{}) do not match expected \
                     consider parameter count ({})",
                    consider.consider_covariance.nrows(),
                    consider.consider_covariance.ncols(),
                    n_consider
                )));
            }
        }

        let epoch = propagator.current_epoch();
        let state = propagator.current_state();

        Ok(Self {
            dynamics: propagator,
            measurement_models,
            config,
            apriori_epoch: epoch,
            apriori_state: state.clone(),
            apriori_covariance: covariance.clone(),
            current_state: state,
            current_covariance: covariance,
            current_epoch: epoch,
            converged: false,
            iterations_completed: 0,
            final_cost: f64::INFINITY,
            iteration_records: Vec::new(),
            observation_residuals: Vec::new(),
            consider_info: None,
        })
    }

    // =========================================================================
    // Accessor methods
    // =========================================================================

    /// Get current state estimate at the reference epoch.
    pub fn current_state(&self) -> &DVector<f64> {
        &self.current_state
    }

    /// Get current formal covariance (solve-for partition).
    pub fn current_covariance(&self) -> &DMatrix<f64> {
        &self.current_covariance
    }

    /// Get current epoch (reference epoch for the batch).
    pub fn current_epoch(&self) -> Epoch {
        self.current_epoch
    }

    /// Whether the solver has converged.
    pub fn converged(&self) -> bool {
        self.converged
    }

    /// Number of Gauss-Newton iterations completed.
    pub fn iterations_completed(&self) -> usize {
        self.iterations_completed
    }

    /// Final cost function value J.
    pub fn final_cost(&self) -> f64 {
        self.final_cost
    }

    /// Per-iteration diagnostic records.
    ///
    /// Only populated when `config.store_iteration_records` is true.
    pub fn iteration_records(&self) -> &[BLSIterationRecord] {
        &self.iteration_records
    }

    /// Per-observation residuals for each iteration.
    ///
    /// Only populated when `config.store_observation_residuals` is true.
    /// Outer index is iteration, inner index is observation.
    pub fn observation_residuals(&self) -> &[Vec<BLSObservationResidual>] {
        &self.observation_residuals
    }

    /// Formal covariance (same as `current_covariance()`).
    ///
    /// This is the inverse of the information matrix, representing the
    /// uncertainty from the solve-for parameters only.
    pub fn formal_covariance(&self) -> &DMatrix<f64> {
        &self.current_covariance
    }

    /// Consider covariance contribution P_consider.
    ///
    /// Computes P_consider = Λ_ss⁻¹ Λ_sc P̄_c Λ_cs Λ_ss⁻¹ from cached
    /// consider information. Returns None if consider parameters are not
    /// configured or solve() has not been called.
    pub fn consider_covariance(&self) -> Option<DMatrix<f64>> {
        let (lambda_ss_inv, lambda_sc) = self.consider_info.as_ref()?;
        let consider_cov = &self.config.consider_params.as_ref()?.consider_covariance;

        // P_consider = Λ_ss⁻¹ * Λ_sc * P̄_c * Λ_csᵀ * Λ_ss⁻¹
        // Note: Λ_cs = Λ_scᵀ
        let p_consider =
            lambda_ss_inv * lambda_sc * consider_cov * lambda_sc.transpose() * lambda_ss_inv;
        Some(p_consider)
    }

    /// Total covariance: formal + consider contribution.
    ///
    /// P_total = P_formal + P_consider
    ///
    /// Returns formal covariance if consider parameters are not configured.
    pub fn total_covariance(&self) -> DMatrix<f64> {
        match self.consider_covariance() {
            Some(p_consider) => &self.current_covariance + p_consider,
            None => self.current_covariance.clone(),
        }
    }

    // =========================================================================
    // solve() — main entry point
    // =========================================================================

    /// Solve the batch least squares problem.
    ///
    /// Iteratively processes all observations to find the state that minimizes
    /// the weighted sum of squared residuals subject to the a priori constraint.
    ///
    /// # Arguments
    ///
    /// * `observations` - Slice of observations to process
    ///
    /// # Errors
    ///
    /// Returns error if:
    /// - No observations are provided
    /// - Any observation's model_index is out of bounds
    /// - A numerical error occurs (singular matrix, propagation failure)
    pub fn solve(&mut self, observations: &[Observation]) -> Result<(), BraheError> {
        // Validate observations non-empty
        if observations.is_empty() {
            return Err(BraheError::Error(
                "No observations provided to BLS solver".to_string(),
            ));
        }

        // Validate model indices
        for (i, obs) in observations.iter().enumerate() {
            if obs.model_index >= self.measurement_models.len() {
                return Err(BraheError::Error(format!(
                    "Observation {} has model_index {} but only {} models are configured",
                    i,
                    obs.model_index,
                    self.measurement_models.len()
                )));
            }
        }

        // Sort observations by epoch
        let mut sorted_obs: Vec<&Observation> = observations.iter().collect();
        sorted_obs.sort_by(|a, b| {
            let dt: f64 = a.epoch - b.epoch;
            dt.partial_cmp(&0.0)
                .unwrap_or(std::cmp::Ordering::Equal)
                .then(std::cmp::Ordering::Equal)
        });

        // Determine solve-for dimension
        let state_dim = self.current_state.len();
        let n_solve = self
            .config
            .consider_params
            .as_ref()
            .map(|c| c.n_solve)
            .unwrap_or(state_dim);

        // Compute a priori information matrix (solve-for partition)
        let p0_solve = self.apriori_covariance.view((0, 0), (n_solve, n_solve));
        let p0_solve_inv = p0_solve.clone_owned().try_inverse().ok_or_else(|| {
            BraheError::NumericalError(
                "A priori covariance (solve-for partition) is singular".to_string(),
            )
        })?;

        let mut prev_cost = f64::INFINITY;

        for iteration in 0..self.config.max_iterations {
            // Reinitialize propagator at reference epoch with current state
            self.dynamics.reinitialize(
                self.apriori_epoch,
                self.current_state.clone(),
                Some(self.current_covariance.clone()),
            );

            // Collect residuals and mapped Jacobians at each observation epoch
            let mut residuals: Vec<DVector<f64>> = Vec::with_capacity(sorted_obs.len());
            let mut h_mapped_solve: Vec<DMatrix<f64>> = Vec::with_capacity(sorted_obs.len());
            let mut h_mapped_full: Vec<DMatrix<f64>> = Vec::with_capacity(sorted_obs.len());
            let mut r_matrices: Vec<DMatrix<f64>> = Vec::with_capacity(sorted_obs.len());
            let mut obs_epochs: Vec<Epoch> = Vec::with_capacity(sorted_obs.len());
            let mut obs_model_names: Vec<String> = Vec::with_capacity(sorted_obs.len());

            for obs in &sorted_obs {
                // Propagate to observation epoch
                self.dynamics.propagate_to(obs.epoch);

                let state_at_obs = self.dynamics.current_state();
                let model = &self.measurement_models[obs.model_index];

                // Compute predicted measurement
                let z_predicted = model.predict(&obs.epoch, &state_at_obs, None)?;

                // Compute residual: δy = y - h(x_k(t_i))
                let residual = &obs.measurement - &z_predicted;

                // Get local Jacobian H̃ (m x n_state)
                let h_local = model.jacobian(&obs.epoch, &state_at_obs, None)?;

                // Get STM Φ(t_i, t_0) from propagator
                let stm = self
                    .dynamics
                    .stm()
                    .ok_or_else(|| {
                        BraheError::NumericalError(
                            "STM not available from propagator during BLS solve".to_string(),
                        )
                    })?
                    .clone();

                // Mapped Jacobian: H_i = H̃_i * Φ(t_i, t_0)
                let h_full = &h_local * &stm;

                // Extract solve-for columns
                let h_solve = h_full.columns(0, n_solve).into_owned();

                residuals.push(residual);
                h_mapped_solve.push(h_solve);
                if self.config.consider_params.is_some() {
                    h_mapped_full.push(h_full);
                }
                r_matrices.push(model.noise_covariance());
                obs_epochs.push(obs.epoch);
                obs_model_names.push(model.name().to_string());
            }

            // Dispatch to formulation-specific solver
            let solver_result = match self.config.solver_method {
                BLSSolverMethod::NormalEquations => self.solve_normal_equations(
                    &residuals,
                    &h_mapped_solve,
                    &h_mapped_full,
                    &r_matrices,
                    &p0_solve_inv,
                    n_solve,
                )?,
                BLSSolverMethod::StackedObservationMatrix => self.solve_stacked_matrix(
                    &residuals,
                    &h_mapped_solve,
                    &h_mapped_full,
                    &r_matrices,
                    &p0_solve_inv,
                    n_solve,
                )?,
            };

            // Compute cost function J = Σ δy_iᵀ R_i⁻¹ δy_i + (x̄₀ - x_k)ᵀ P̄₀⁻¹ (x̄₀ - x_k)
            let mut cost = 0.0;
            for (i, residual) in residuals.iter().enumerate() {
                let r_inv = r_matrices[i].clone().try_inverse().ok_or_else(|| {
                    BraheError::NumericalError(format!(
                        "Measurement noise covariance R is singular for observation {}",
                        i
                    ))
                })?;
                cost += (residual.transpose() * &r_inv * residual)[(0, 0)];
            }
            let apriori_diff = self.apriori_state.rows(0, n_solve).into_owned()
                - self.current_state.rows(0, n_solve).into_owned();
            cost += (apriori_diff.transpose() * &p0_solve_inv * &apriori_diff)[(0, 0)];

            // Compute pre-fit RMS
            let total_meas: usize = residuals.iter().map(|r| r.len()).sum();
            let prefit_rms = if total_meas > 0 {
                let sum_sq: f64 = residuals.iter().map(|r| r.dot(r)).sum();
                (sum_sq / total_meas as f64).sqrt()
            } else {
                0.0
            };

            // Apply state correction to solve-for parameters
            let dx = &solver_result.state_correction;
            let dx_norm = dx.norm();
            for i in 0..n_solve {
                self.current_state[i] += dx[i];
            }
            self.current_covariance = solver_result.formal_covariance;
            self.consider_info = solver_result.consider_info;

            // Compute post-fit residuals by re-propagating with updated state
            self.dynamics.reinitialize(
                self.apriori_epoch,
                self.current_state.clone(),
                Some(self.current_covariance.clone()),
            );

            let mut postfit_residuals: Vec<DVector<f64>> = Vec::with_capacity(sorted_obs.len());
            for obs in &sorted_obs {
                self.dynamics.propagate_to(obs.epoch);
                let state_at_obs = self.dynamics.current_state();
                let model = &self.measurement_models[obs.model_index];
                let z_predicted = model.predict(&obs.epoch, &state_at_obs, None)?;
                postfit_residuals.push(&obs.measurement - &z_predicted);
            }

            let postfit_rms = if total_meas > 0 {
                let sum_sq: f64 = postfit_residuals.iter().map(|r| r.dot(r)).sum();
                (sum_sq / total_meas as f64).sqrt()
            } else {
                0.0
            };

            // Store iteration record
            if self.config.store_iteration_records {
                self.iteration_records.push(BLSIterationRecord {
                    iteration,
                    epoch: self.apriori_epoch,
                    state: self.current_state.clone(),
                    covariance: self.current_covariance.clone(),
                    state_correction: dx.clone(),
                    state_correction_norm: dx_norm,
                    cost,
                    rms_prefit_residual: prefit_rms,
                    rms_postfit_residual: postfit_rms,
                });
            }

            // Store per-observation residuals
            if self.config.store_observation_residuals {
                let obs_residuals: Vec<BLSObservationResidual> = (0..sorted_obs.len())
                    .map(|i| BLSObservationResidual {
                        epoch: obs_epochs[i],
                        model_name: obs_model_names[i].clone(),
                        prefit_residual: residuals[i].clone(),
                        postfit_residual: postfit_residuals[i].clone(),
                    })
                    .collect();
                self.observation_residuals.push(obs_residuals);
            }

            self.iterations_completed = iteration + 1;
            self.final_cost = cost;

            // Check convergence
            let state_converged = self
                .config
                .state_correction_threshold
                .map(|thresh| dx_norm < thresh)
                .unwrap_or(false);

            let cost_converged = self
                .config
                .cost_convergence_threshold
                .map(|thresh| {
                    if prev_cost.is_infinite() {
                        false
                    } else {
                        let rel_change = (prev_cost - cost).abs() / prev_cost.abs().max(1e-30);
                        rel_change < thresh
                    }
                })
                .unwrap_or(false);

            if state_converged || cost_converged {
                self.converged = true;
                break;
            }

            prev_cost = cost;
        }

        self.current_epoch = self.apriori_epoch;
        Ok(())
    }

    // =========================================================================
    // Normal Equations solver
    // =========================================================================

    /// Solve using the Normal Equations formulation.
    ///
    /// Accumulates:
    /// - Λ = P̄₀⁻¹ + Σ Hᵢᵀ Rᵢ⁻¹ Hᵢ (information matrix)
    /// - N = P̄₀⁻¹(x̄₀ - xₖ) + Σ Hᵢᵀ Rᵢ⁻¹ δyᵢ (normal vector)
    ///
    /// Solves δx = Λ⁻¹N via Cholesky decomposition.
    fn solve_normal_equations(
        &self,
        residuals: &[DVector<f64>],
        h_solve: &[DMatrix<f64>],
        h_full: &[DMatrix<f64>],
        r_matrices: &[DMatrix<f64>],
        p0_solve_inv: &DMatrix<f64>,
        n_solve: usize,
    ) -> Result<SolverResult, BraheError> {
        let state_dim = self.current_state.len();

        // Initialize information matrix and normal vector with a priori
        let mut lambda = p0_solve_inv.clone();
        let apriori_diff = self.apriori_state.rows(0, n_solve).into_owned()
            - self.current_state.rows(0, n_solve).into_owned();
        let mut n_vec = p0_solve_inv * &apriori_diff;

        // Also accumulate cross-term if consider params are configured
        let has_consider = self.config.consider_params.is_some();
        let n_consider = state_dim - n_solve;
        let mut lambda_sc = if has_consider {
            DMatrix::zeros(n_solve, n_consider)
        } else {
            DMatrix::zeros(0, 0)
        };

        // Accumulate over observations
        for i in 0..residuals.len() {
            let r_inv = r_matrices[i].clone().try_inverse().ok_or_else(|| {
                BraheError::NumericalError(format!(
                    "Measurement noise covariance R is singular for observation {}",
                    i
                ))
            })?;

            let h_t = h_solve[i].transpose();
            let h_t_r_inv = &h_t * &r_inv;

            lambda += &h_t_r_inv * &h_solve[i];
            n_vec += &h_t_r_inv * &residuals[i];

            // Cross-term for consider parameters
            if has_consider {
                let h_c = h_full[i].columns(n_solve, n_consider);
                lambda_sc += &h_t_r_inv * h_c;
            }
        }

        // Solve via Cholesky decomposition
        let chol = lambda.clone().cholesky().ok_or_else(|| {
            BraheError::NumericalError(
                "Cholesky decomposition of information matrix Λ failed. \
                 The information matrix may not be positive definite."
                    .to_string(),
            )
        })?;

        let dx = chol.solve(&n_vec);
        let lambda_inv = chol.inverse();

        let consider_info = if has_consider {
            Some((lambda_inv.clone(), lambda_sc))
        } else {
            None
        };

        Ok(SolverResult {
            state_correction: dx,
            formal_covariance: lambda_inv,
            consider_info,
        })
    }

    // =========================================================================
    // Stacked Observation Matrix solver
    // =========================================================================

    /// Solve using the Stacked Observation Matrix formulation.
    ///
    /// Builds the augmented system:
    /// - Stack all Hᵢ and δyᵢ into H (m_total x n_solve) and δy (m_total)
    /// - Augment with a priori: append √(P̄₀⁻¹) rows and √(P̄₀⁻¹)(x̄₀ - xₖ)
    /// - Apply weights √(Rᵢ⁻¹) to observation rows
    /// - Solve via Cholesky on HᵀWH
    fn solve_stacked_matrix(
        &self,
        residuals: &[DVector<f64>],
        h_solve: &[DMatrix<f64>],
        h_full: &[DMatrix<f64>],
        r_matrices: &[DMatrix<f64>],
        p0_solve_inv: &DMatrix<f64>,
        n_solve: usize,
    ) -> Result<SolverResult, BraheError> {
        let state_dim = self.current_state.len();

        // Total measurement dimension
        let m_total: usize = residuals.iter().map(|r| r.len()).sum();

        // Build weighted H and δy matrices
        // Total rows = m_total (observations) + n_solve (a priori)
        let total_rows = m_total + n_solve;
        let mut h_stacked = DMatrix::zeros(total_rows, n_solve);
        let mut dy_stacked = DVector::zeros(total_rows);

        let mut row = 0;
        for i in 0..residuals.len() {
            let m_i = residuals[i].len();

            // Compute √(R⁻¹) via diagonal (valid for diagonal R)
            let r_inv = r_matrices[i].clone().try_inverse().ok_or_else(|| {
                BraheError::NumericalError(format!(
                    "Measurement noise covariance R is singular for observation {}",
                    i
                ))
            })?;

            // For diagonal R, √(R⁻¹) is element-wise sqrt of diagonal
            let sqrt_r_inv =
                DMatrix::from_fn(
                    m_i,
                    m_i,
                    |r, c| {
                        if r == c { r_inv[(r, c)].sqrt() } else { 0.0 }
                    },
                );

            // Weighted H row block: √(R⁻¹) * H_i
            let weighted_h = &sqrt_r_inv * &h_solve[i];
            h_stacked
                .view_mut((row, 0), (m_i, n_solve))
                .copy_from(&weighted_h);

            // Weighted residual: √(R⁻¹) * δy_i
            let weighted_dy = &sqrt_r_inv * &residuals[i];
            dy_stacked.rows_mut(row, m_i).copy_from(&weighted_dy);

            row += m_i;
        }

        // Augment with a priori: √(P̄₀⁻¹)
        let p0_chol = p0_solve_inv.clone().cholesky().ok_or_else(|| {
            BraheError::NumericalError(
                "Cholesky decomposition of a priori information matrix failed".to_string(),
            )
        })?;
        // Use L^T (upper factor) so that (L^T)^T * L^T = L * L^T = P0_inv
        // when accumulated via H^T H in the stacked system
        let sqrt_p0_inv = p0_chol.l().transpose();

        h_stacked
            .view_mut((row, 0), (n_solve, n_solve))
            .copy_from(&sqrt_p0_inv);

        let apriori_diff = self.apriori_state.rows(0, n_solve).into_owned()
            - self.current_state.rows(0, n_solve).into_owned();
        let weighted_apriori = &sqrt_p0_inv * &apriori_diff;
        dy_stacked
            .rows_mut(row, n_solve)
            .copy_from(&weighted_apriori);

        // Solve via HᵀH δx = Hᵀ δy (the weighting is already in H and δy)
        let hth = h_stacked.transpose() * &h_stacked;
        let htdy = h_stacked.transpose() * &dy_stacked;

        let chol = hth.clone().cholesky().ok_or_else(|| {
            BraheError::NumericalError(
                "Cholesky decomposition of HᵀWH failed in stacked matrix formulation".to_string(),
            )
        })?;

        let dx = chol.solve(&htdy);
        let lambda_inv = chol.inverse();

        // Compute cross-term for consider parameters if configured
        let has_consider = self.config.consider_params.is_some();
        let n_consider = state_dim - n_solve;
        let consider_info = if has_consider {
            let mut lambda_sc = DMatrix::zeros(n_solve, n_consider);
            for i in 0..residuals.len() {
                let r_inv = r_matrices[i].clone().try_inverse().ok_or_else(|| {
                    BraheError::NumericalError(
                        "Measurement noise covariance R is singular (consider cross-term)"
                            .to_string(),
                    )
                })?;
                let h_s_t = h_solve[i].transpose();
                let h_c = h_full[i].columns(n_solve, n_consider);
                lambda_sc += &h_s_t * &r_inv * h_c;
            }
            Some((lambda_inv.clone(), lambda_sc))
        } else {
            None
        };

        Ok(SolverResult {
            state_correction: dx,
            formal_covariance: lambda_inv,
            consider_info,
        })
    }
}

#[cfg(test)]
mod tests {
    use super::*;
    use crate::constants::physical::GM_EARTH;
    use crate::eop::{EOPExtrapolation, FileEOPProvider, set_global_eop_provider};
    use crate::estimation::{BLSConfig, BLSSolverMethod, InertialPositionMeasurementModel};
    use crate::propagators::DNumericalOrbitPropagator;
    use crate::propagators::NumericalPropagationConfig;
    use crate::propagators::force_model_config::ForceModelConfig;
    use crate::propagators::traits::DStatePropagator;
    use crate::time::TimeSystem;
    use serial_test::serial;

    fn setup_global_test_eop() {
        let eop = FileEOPProvider::from_default_standard(true, EOPExtrapolation::Hold).unwrap();
        set_global_eop_provider(eop);
    }

    // =========================================================================
    // Test helpers — two-body point-mass scenario (exact dynamics)
    // =========================================================================

    /// Two-body LEO circular orbit. Point-mass gravity only.
    fn two_body_leo() -> (Epoch, DVector<f64>) {
        let epoch = Epoch::from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0, TimeSystem::UTC);
        let r = 6878.0e3;
        let v = (GM_EARTH / r).sqrt();
        (epoch, DVector::from_vec(vec![r, 0.0, 0.0, 0.0, v, 0.0]))
    }

    /// Generate noise-free inertial position observations from a truth trajectory.
    fn generate_position_observations(
        epoch: Epoch,
        true_state: &DVector<f64>,
        num_obs: usize,
        interval_s: f64,
    ) -> Vec<Observation> {
        let mut prop = DNumericalOrbitPropagator::new(
            epoch,
            true_state.clone(),
            NumericalPropagationConfig::default(),
            ForceModelConfig::two_body_gravity(),
            None,
            None,
            None,
            None,
        )
        .unwrap();

        (1..=num_obs)
            .map(|i| {
                let t = epoch + (i as f64) * interval_s;
                prop.propagate_to(t);
                Observation::new(t, prop.current_state().rows(0, 3).into_owned(), 0)
            })
            .collect()
    }

    fn default_p0() -> DMatrix<f64> {
        DMatrix::from_diagonal(&DVector::from_vec(vec![1e6, 1e6, 1e6, 1e2, 1e2, 1e2]))
    }

    /// Create a BLS estimator with two-body point-mass gravity.
    fn create_two_body_bls(
        epoch: Epoch,
        state: DVector<f64>,
        p0: DMatrix<f64>,
        models: Vec<Box<dyn MeasurementModel>>,
        config: BLSConfig,
    ) -> BatchLeastSquares {
        BatchLeastSquares::new(
            epoch,
            state,
            p0,
            NumericalPropagationConfig::default(),
            ForceModelConfig::two_body_gravity(),
            None,
            None,
            None,
            models,
            config,
        )
        .unwrap()
    }

    // =========================================================================
    // Construction tests
    // =========================================================================

    #[test]
    #[serial]
    fn test_bls_construction_valid() {
        setup_global_test_eop();
        let (epoch, state) = two_body_leo();
        let bls = create_two_body_bls(
            epoch,
            state.clone(),
            default_p0(),
            vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
            BLSConfig::default(),
        );
        assert_eq!(bls.current_state().len(), 6);
        assert_eq!(bls.current_covariance().nrows(), 6);
        assert!(!bls.converged());
        assert_eq!(bls.iterations_completed(), 0);
    }

    #[test]
    #[serial]
    fn test_bls_construction_no_stm_errors() {
        setup_global_test_eop();
        let (epoch, state) = two_body_leo();

        let prop = DNumericalOrbitPropagator::new(
            epoch,
            state,
            NumericalPropagationConfig::default(), // no STM
            ForceModelConfig::two_body_gravity(),
            None,
            None,
            None,
            None,
        )
        .unwrap();

        let result = BatchLeastSquares::from_propagator(
            DynamicsSource::OrbitPropagator(prop),
            vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
            BLSConfig::default(),
        );
        assert!(result.is_err());
        assert!(result.unwrap_err().to_string().contains("STM"));
    }

    #[test]
    #[serial]
    fn test_bls_construction_no_models_errors() {
        setup_global_test_eop();
        let (epoch, state) = two_body_leo();

        let mut prop_config = NumericalPropagationConfig::default();
        prop_config.variational.enable_stm = true;

        let prop = DNumericalOrbitPropagator::new(
            epoch,
            state,
            prop_config,
            ForceModelConfig::two_body_gravity(),
            None,
            None,
            None,
            Some(default_p0()),
        )
        .unwrap();

        let result = BatchLeastSquares::from_propagator(
            DynamicsSource::OrbitPropagator(prop),
            vec![], // no models
            BLSConfig::default(),
        );
        assert!(result.is_err());
        assert!(
            result
                .unwrap_err()
                .to_string()
                .contains("measurement model")
        );
    }

    #[test]
    #[serial]
    fn test_bls_construction_no_convergence_threshold_errors() {
        setup_global_test_eop();
        let (epoch, state) = two_body_leo();

        let config = BLSConfig {
            state_correction_threshold: None,
            cost_convergence_threshold: None,
            ..BLSConfig::default()
        };

        let result = BatchLeastSquares::new(
            epoch,
            state,
            default_p0(),
            NumericalPropagationConfig::default(),
            ForceModelConfig::two_body_gravity(),
            None,
            None,
            None,
            vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
            config,
        );
        assert!(result.is_err());
        assert!(
            result
                .unwrap_err()
                .to_string()
                .contains("convergence threshold")
        );
    }

    #[test]
    #[serial]
    fn test_bls_construction_bad_covariance_dims_errors() {
        setup_global_test_eop();
        let (epoch, state) = two_body_leo();

        // The propagator validates covariance dimensions internally, so this
        // will fail at the propagator level. We verify the error propagates up.
        let p_wrong = DMatrix::identity(4, 4);
        let result = std::panic::catch_unwind(|| {
            BatchLeastSquares::new(
                epoch,
                state,
                p_wrong,
                NumericalPropagationConfig::default(),
                ForceModelConfig::two_body_gravity(),
                None,
                None,
                None,
                vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
                BLSConfig::default(),
            )
        });
        // The propagator panics on dimension mismatch, so we expect a panic or error
        assert!(
            result.is_err() || result.unwrap().is_err(),
            "Should fail with mismatched covariance dimensions"
        );
    }

    // =========================================================================
    // Solve tests — Normal Equations
    // =========================================================================

    #[test]
    #[serial]
    fn test_bls_normal_equations_converges() {
        setup_global_test_eop();

        let (epoch, true_state) = two_body_leo();
        let mut perturbed = true_state.clone();
        perturbed[0] += 1000.0; // 1km offset in x
        perturbed[1] += 500.0; // 500m offset in y

        let obs = generate_position_observations(epoch, &true_state, 20, 30.0);

        let config = BLSConfig {
            solver_method: BLSSolverMethod::NormalEquations,
            max_iterations: 10,
            state_correction_threshold: Some(1e-6),
            cost_convergence_threshold: None,
            consider_params: None,
            store_iteration_records: true,
            store_observation_residuals: true,
        };

        let mut bls = create_two_body_bls(
            epoch,
            perturbed.clone(),
            default_p0(),
            vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
            config,
        );

        bls.solve(&obs).unwrap();

        // Verify convergence
        assert!(bls.converged(), "BLS should converge");

        // Verify position error reduced significantly
        let final_state = bls.current_state();
        let pos_error = (final_state.rows(0, 3) - true_state.rows(0, 3)).norm();
        let initial_error = (perturbed.rows(0, 3) - true_state.rows(0, 3)).norm();

        assert!(
            pos_error < initial_error * 0.001,
            "Position error should decrease >1000x: initial={:.1}m, final={:.6}m",
            initial_error,
            pos_error
        );
    }

    #[test]
    #[serial]
    fn test_bls_empty_observations_errors() {
        setup_global_test_eop();

        let (epoch, state) = two_body_leo();
        let mut bls = create_two_body_bls(
            epoch,
            state,
            default_p0(),
            vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
            BLSConfig::default(),
        );

        let result = bls.solve(&[]);
        assert!(result.is_err());
        assert!(result.unwrap_err().to_string().contains("No observations"));
    }

    #[test]
    #[serial]
    fn test_bls_cost_decreases_monotonically() {
        setup_global_test_eop();

        let (epoch, true_state) = two_body_leo();
        let mut perturbed = true_state.clone();
        perturbed[0] += 500.0;

        let obs = generate_position_observations(epoch, &true_state, 15, 30.0);

        let config = BLSConfig {
            solver_method: BLSSolverMethod::NormalEquations,
            max_iterations: 10,
            state_correction_threshold: Some(1e-12), // tight threshold to get multiple iterations
            cost_convergence_threshold: None,
            consider_params: None,
            store_iteration_records: true,
            store_observation_residuals: false,
        };

        let mut bls = create_two_body_bls(
            epoch,
            perturbed,
            default_p0(),
            vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
            config,
        );

        bls.solve(&obs).unwrap();

        let records = bls.iteration_records();
        assert!(
            records.len() >= 2,
            "Need at least 2 iterations to check monotonicity, got {}",
            records.len()
        );

        for i in 1..records.len() {
            assert!(
                records[i].cost <= records[i - 1].cost + 1e-6,
                "Cost should not increase: iteration {} cost={:.6}, iteration {} cost={:.6}",
                i - 1,
                records[i - 1].cost,
                i,
                records[i].cost
            );
        }
    }

    #[test]
    #[serial]
    fn test_bls_iteration_records_stored() {
        setup_global_test_eop();

        let (epoch, true_state) = two_body_leo();
        let mut perturbed = true_state.clone();
        perturbed[0] += 200.0;

        let obs = generate_position_observations(epoch, &true_state, 10, 30.0);

        let config = BLSConfig {
            solver_method: BLSSolverMethod::NormalEquations,
            max_iterations: 5,
            state_correction_threshold: Some(1e-12), // tight to force iterations
            cost_convergence_threshold: None,
            consider_params: None,
            store_iteration_records: true,
            store_observation_residuals: false,
        };

        let mut bls = create_two_body_bls(
            epoch,
            perturbed,
            default_p0(),
            vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
            config,
        );

        bls.solve(&obs).unwrap();

        let records = bls.iteration_records();
        assert!(
            !records.is_empty(),
            "Should have at least one iteration record"
        );
        assert_eq!(records.len(), bls.iterations_completed());

        // Verify iteration numbers
        for (i, record) in records.iter().enumerate() {
            assert_eq!(record.iteration, i);
            assert_eq!(record.state.len(), 6);
            assert_eq!(record.covariance.nrows(), 6);
            assert_eq!(record.state_correction.len(), 6);
        }
    }

    #[test]
    #[serial]
    fn test_bls_observation_residuals_stored() {
        setup_global_test_eop();

        let (epoch, true_state) = two_body_leo();
        let mut perturbed = true_state.clone();
        perturbed[0] += 200.0;

        let num_obs = 10;
        let obs = generate_position_observations(epoch, &true_state, num_obs, 30.0);

        let config = BLSConfig {
            solver_method: BLSSolverMethod::NormalEquations,
            max_iterations: 3,
            state_correction_threshold: Some(1e-12),
            cost_convergence_threshold: None,
            consider_params: None,
            store_iteration_records: true,
            store_observation_residuals: true,
        };

        let mut bls = create_two_body_bls(
            epoch,
            perturbed,
            default_p0(),
            vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
            config,
        );

        bls.solve(&obs).unwrap();

        let obs_residuals = bls.observation_residuals();
        assert_eq!(obs_residuals.len(), bls.iterations_completed());

        for iteration_residuals in obs_residuals {
            assert_eq!(iteration_residuals.len(), num_obs);
            for residual in iteration_residuals {
                assert_eq!(residual.prefit_residual.len(), 3);
                assert_eq!(residual.postfit_residual.len(), 3);
                assert_eq!(residual.model_name, "InertialPosition");
            }
        }
    }

    #[test]
    #[serial]
    fn test_bls_max_iterations_respected() {
        setup_global_test_eop();

        let (epoch, true_state) = two_body_leo();
        let mut perturbed = true_state.clone();
        perturbed[0] += 1000.0;

        let obs = generate_position_observations(epoch, &true_state, 10, 30.0);

        let config = BLSConfig {
            solver_method: BLSSolverMethod::NormalEquations,
            max_iterations: 2,
            state_correction_threshold: Some(1e-30), // impossible to reach
            cost_convergence_threshold: None,
            consider_params: None,
            store_iteration_records: true,
            store_observation_residuals: false,
        };

        let mut bls = create_two_body_bls(
            epoch,
            perturbed,
            default_p0(),
            vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
            config,
        );

        bls.solve(&obs).unwrap();

        assert_eq!(
            bls.iterations_completed(),
            2,
            "Should stop at max_iterations=2"
        );
        assert!(
            !bls.converged(),
            "Should not converge with impossible threshold"
        );
    }

    // =========================================================================
    // Solve tests — Stacked Observation Matrix
    // =========================================================================

    #[test]
    #[serial]
    fn test_bls_stacked_matrix_converges() {
        setup_global_test_eop();

        let (epoch, true_state) = two_body_leo();
        let mut perturbed = true_state.clone();
        perturbed[0] += 1000.0;
        perturbed[1] += 500.0;

        let obs = generate_position_observations(epoch, &true_state, 20, 30.0);

        let config = BLSConfig {
            solver_method: BLSSolverMethod::StackedObservationMatrix,
            max_iterations: 10,
            state_correction_threshold: Some(1e-6),
            cost_convergence_threshold: None,
            consider_params: None,
            store_iteration_records: true,
            store_observation_residuals: false,
        };

        let mut bls = create_two_body_bls(
            epoch,
            perturbed.clone(),
            default_p0(),
            vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
            config,
        );

        bls.solve(&obs).unwrap();

        assert!(bls.converged(), "Stacked matrix BLS should converge");

        let final_state = bls.current_state();
        let pos_error = (final_state.rows(0, 3) - true_state.rows(0, 3)).norm();
        let initial_error = (perturbed.rows(0, 3) - true_state.rows(0, 3)).norm();

        assert!(
            pos_error < initial_error * 0.001,
            "Position error should decrease >1000x: initial={:.1}m, final={:.6}m",
            initial_error,
            pos_error
        );
    }

    #[test]
    #[serial]
    fn test_bls_formulation_equivalence() {
        setup_global_test_eop();

        let (epoch, true_state) = two_body_leo();
        let mut perturbed = true_state.clone();
        perturbed[0] += 500.0;
        perturbed[1] += 200.0;

        let obs = generate_position_observations(epoch, &true_state, 15, 30.0);

        // Normal Equations
        let config_ne = BLSConfig {
            solver_method: BLSSolverMethod::NormalEquations,
            max_iterations: 10,
            state_correction_threshold: Some(1e-8),
            cost_convergence_threshold: None,
            consider_params: None,
            store_iteration_records: true,
            store_observation_residuals: false,
        };

        let mut bls_ne = create_two_body_bls(
            epoch,
            perturbed.clone(),
            default_p0(),
            vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
            config_ne,
        );
        bls_ne.solve(&obs).unwrap();

        // Stacked Matrix
        let config_sm = BLSConfig {
            solver_method: BLSSolverMethod::StackedObservationMatrix,
            max_iterations: 10,
            state_correction_threshold: Some(1e-8),
            cost_convergence_threshold: None,
            consider_params: None,
            store_iteration_records: true,
            store_observation_residuals: false,
        };

        let mut bls_sm = create_two_body_bls(
            epoch,
            perturbed,
            default_p0(),
            vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
            config_sm,
        );
        bls_sm.solve(&obs).unwrap();

        // Both should converge
        assert!(bls_ne.converged());
        assert!(bls_sm.converged());

        // States should match within tolerance
        let state_diff = (bls_ne.current_state() - bls_sm.current_state()).norm();
        assert!(
            state_diff < 1e-3,
            "Formulations should produce same state: diff={:.6}m",
            state_diff
        );

        // Costs should match within tolerance
        let cost_diff = (bls_ne.final_cost() - bls_sm.final_cost()).abs();
        assert!(
            cost_diff < 1e-3,
            "Formulations should produce same cost: NE={:.6}, SM={:.6}",
            bls_ne.final_cost(),
            bls_sm.final_cost()
        );
    }

    #[test]
    #[serial]
    fn test_bls_cost_convergence_criterion() {
        setup_global_test_eop();

        let (epoch, true_state) = two_body_leo();
        let mut perturbed = true_state.clone();
        perturbed[0] += 500.0;

        let obs = generate_position_observations(epoch, &true_state, 15, 30.0);

        // Use only cost convergence (no state correction threshold)
        let config = BLSConfig {
            solver_method: BLSSolverMethod::NormalEquations,
            max_iterations: 20,
            state_correction_threshold: None,
            cost_convergence_threshold: Some(1e-6),
            consider_params: None,
            store_iteration_records: true,
            store_observation_residuals: false,
        };

        let mut bls = create_two_body_bls(
            epoch,
            perturbed,
            default_p0(),
            vec![Box::new(InertialPositionMeasurementModel::new(10.0))],
            config,
        );

        bls.solve(&obs).unwrap();

        assert!(
            bls.converged(),
            "Should converge via cost criterion. Iterations: {}",
            bls.iterations_completed()
        );

        // Cost should be very small (noise-free measurements with exact dynamics)
        assert!(
            bls.final_cost() < 1.0,
            "Final cost should be small: {:.6}",
            bls.final_cost()
        );
    }
}