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
/*!
 * CCSDS Orbit Parameter Message (OPM) data structures.
 *
 * OPM messages contain a single state vector (position and velocity) with
 * optional Keplerian elements, spacecraft parameters, maneuvers, and
 * covariance data.
 *
 * Reference: CCSDS 502.0-B-3 (Orbit Data Messages), April 2023
 */

use std::path::Path;

use crate::ccsds::common::{
    CCSDSCovariance, CCSDSFormat, CCSDSRefFrame, CCSDSSpacecraftParameters, CCSDSTimeSystem,
    CCSDSUserDefined, ODMHeader,
};
use crate::time::Epoch;
use crate::utils::errors::BraheError;

/// A complete CCSDS Orbit Parameter Message.
#[derive(Debug, Clone)]
pub struct OPM {
    /// Message header
    pub header: ODMHeader,
    /// Metadata
    pub metadata: OPMMetadata,
    /// State vector
    pub state_vector: OPMStateVector,
    /// Optional Keplerian elements
    pub keplerian_elements: Option<OPMKeplerianElements>,
    /// Optional spacecraft parameters
    pub spacecraft_parameters: Option<CCSDSSpacecraftParameters>,
    /// Optional covariance matrix
    pub covariance: Option<CCSDSCovariance>,
    /// Maneuvers
    pub maneuvers: Vec<OPMManeuver>,
    /// Optional user-defined parameters
    pub user_defined: Option<CCSDSUserDefined>,
}

/// OPM metadata.
#[derive(Debug, Clone)]
pub struct OPMMetadata {
    /// Spacecraft name
    pub object_name: String,
    /// International designator
    pub object_id: String,
    /// Center body name
    pub center_name: String,
    /// Reference frame
    pub ref_frame: CCSDSRefFrame,
    /// Optional reference frame epoch
    pub ref_frame_epoch: Option<Epoch>,
    /// Time system
    pub time_system: CCSDSTimeSystem,
    /// Comments
    pub comments: Vec<String>,
}

/// State vector data in an OPM.
///
/// Position and velocity stored in SI units (meters, m/s).
#[derive(Debug, Clone)]
pub struct OPMStateVector {
    /// Epoch of the state vector
    pub epoch: Epoch,
    /// Position [x, y, z]. Units: meters
    pub position: [f64; 3],
    /// Velocity [vx, vy, vz]. Units: m/s
    pub velocity: [f64; 3],
    /// Comments
    pub comments: Vec<String>,
}

/// Keplerian elements in an OPM.
///
/// Semi-major axis stored in meters (SI). Angles in degrees (CCSDS native).
#[derive(Debug, Clone)]
pub struct OPMKeplerianElements {
    /// Semi-major axis. Units: meters (converted from km)
    pub semi_major_axis: f64,
    /// Eccentricity (dimensionless)
    pub eccentricity: f64,
    /// Inclination. Units: degrees
    pub inclination: f64,
    /// Right ascension of ascending node. Units: degrees
    pub ra_of_asc_node: f64,
    /// Argument of pericenter. Units: degrees
    pub arg_of_pericenter: f64,
    /// True anomaly. Units: degrees
    pub true_anomaly: Option<f64>,
    /// Mean anomaly. Units: degrees
    pub mean_anomaly: Option<f64>,
    /// Gravitational parameter. Units: m³/s² (converted from km³/s²)
    pub gm: Option<f64>,
    /// Comments
    pub comments: Vec<String>,
}

impl OPMMetadata {
    /// Create new metadata with required fields.
    pub fn new(
        object_name: String,
        object_id: String,
        center_name: String,
        ref_frame: CCSDSRefFrame,
        time_system: CCSDSTimeSystem,
    ) -> Self {
        Self {
            object_name,
            object_id,
            center_name,
            ref_frame,
            ref_frame_epoch: None,
            time_system,
            comments: Vec::new(),
        }
    }
}

impl OPMStateVector {
    /// Create a new state vector.
    ///
    /// # Arguments
    ///
    /// * `epoch` - Epoch of the state vector
    /// * `position` - Position [x, y, z]. Units: meters
    /// * `velocity` - Velocity [vx, vy, vz]. Units: m/s
    pub fn new(epoch: Epoch, position: [f64; 3], velocity: [f64; 3]) -> Self {
        Self {
            epoch,
            position,
            velocity,
            comments: Vec::new(),
        }
    }
}

/// A maneuver specification in an OPM.
///
/// Delta-V stored in m/s (SI).
#[derive(Debug, Clone)]
pub struct OPMManeuver {
    /// Epoch of ignition
    pub epoch_ignition: Epoch,
    /// Duration. Units: seconds
    pub duration: f64,
    /// Mass change. Units: kg (negative for mass decrease)
    pub delta_mass: Option<f64>,
    /// Reference frame for the delta-V
    pub ref_frame: CCSDSRefFrame,
    /// Delta-V vector [dv1, dv2, dv3]. Units: m/s
    pub dv: [f64; 3],
    /// Comments
    pub comments: Vec<String>,
}

impl OPMManeuver {
    /// Create a new maneuver.
    ///
    /// # Arguments
    ///
    /// * `epoch_ignition` - Epoch of ignition
    /// * `duration` - Duration in seconds
    /// * `ref_frame` - Reference frame for the delta-V
    /// * `dv` - Delta-V vector [dv1, dv2, dv3]. Units: m/s
    pub fn new(
        epoch_ignition: Epoch,
        duration: f64,
        ref_frame: CCSDSRefFrame,
        dv: [f64; 3],
    ) -> Self {
        Self {
            epoch_ignition,
            duration,
            delta_mass: None,
            ref_frame,
            dv,
            comments: Vec::new(),
        }
    }

    /// Set the delta mass.
    pub fn with_delta_mass(mut self, delta_mass: f64) -> Self {
        self.delta_mass = Some(delta_mass);
        self
    }
}

impl OPM {
    /// Create a new OPM message with required fields.
    ///
    /// # Arguments
    ///
    /// * `originator` - Originator of the message
    /// * `metadata` - OPM metadata
    /// * `state_vector` - State vector data
    pub fn new(originator: String, metadata: OPMMetadata, state_vector: OPMStateVector) -> Self {
        Self {
            header: ODMHeader {
                format_version: 3.0,
                classification: None,
                creation_date: Epoch::now(),
                originator,
                message_id: None,
                comments: Vec::new(),
            },
            metadata,
            state_vector,
            keplerian_elements: None,
            spacecraft_parameters: None,
            covariance: None,
            maneuvers: Vec::new(),
            user_defined: None,
        }
    }

    /// Add a maneuver to the OPM.
    pub fn push_maneuver(&mut self, maneuver: OPMManeuver) {
        self.maneuvers.push(maneuver);
    }

    /// Parse an OPM message from a string, auto-detecting the format.
    #[allow(clippy::should_implement_trait)]
    pub fn from_str(content: &str) -> Result<Self, BraheError> {
        let format = crate::ccsds::common::detect_format(content);
        match format {
            CCSDSFormat::KVN => crate::ccsds::kvn::parse_opm(content),
            CCSDSFormat::XML => crate::ccsds::xml::parse_opm_xml(content),
            CCSDSFormat::JSON => crate::ccsds::json::parse_opm_json(content),
        }
    }

    /// Parse an OPM message from a file, auto-detecting the format.
    pub fn from_file<P: AsRef<Path>>(path: P) -> Result<Self, BraheError> {
        let content = std::fs::read_to_string(path.as_ref())
            .map_err(|e| BraheError::IoError(format!("Failed to read OPM file: {}", e)))?;
        Self::from_str(&content)
    }

    /// Write the OPM message to a string in the specified format.
    pub fn to_string(&self, format: CCSDSFormat) -> Result<String, BraheError> {
        match format {
            CCSDSFormat::KVN => crate::ccsds::kvn::write_opm(self),
            CCSDSFormat::XML => crate::ccsds::xml::write_opm_xml(self),
            CCSDSFormat::JSON => crate::ccsds::json::write_opm_json(
                self,
                crate::ccsds::common::CCSDSJsonKeyCase::Lower,
            ),
        }
    }

    /// Write the OPM message to JSON with explicit key case control.
    pub fn to_json_string(
        &self,
        key_case: crate::ccsds::common::CCSDSJsonKeyCase,
    ) -> Result<String, BraheError> {
        crate::ccsds::json::write_opm_json(self, key_case)
    }

    /// Write the OPM message to a file in the specified format.
    pub fn to_file<P: AsRef<Path>>(&self, path: P, format: CCSDSFormat) -> Result<(), BraheError> {
        let content = self.to_string(format)?;
        std::fs::write(path.as_ref(), content)
            .map_err(|e| BraheError::IoError(format!("Failed to write OPM file: {}", e)))
    }
}

#[cfg(test)]
#[cfg_attr(coverage_nightly, coverage(off))]
mod tests {
    use super::*;
    use crate::ccsds::common::CCSDSJsonKeyCase;
    use crate::time::TimeSystem;

    #[test]
    fn test_opm_builder() {
        let metadata = OPMMetadata::new(
            "SAT1".to_string(),
            "2024-001A".to_string(),
            "EARTH".to_string(),
            CCSDSRefFrame::GCRF,
            CCSDSTimeSystem::UTC,
        );
        let sv = OPMStateVector::new(Epoch::now(), [7000e3, 0.0, 0.0], [0.0, 7500.0, 0.0]);
        let mut opm = OPM::new("TEST_ORG".to_string(), metadata, sv);
        assert_eq!(opm.header.originator, "TEST_ORG");
        assert_eq!(opm.maneuvers.len(), 0);

        let m = OPMManeuver::new(Epoch::now(), 120.0, CCSDSRefFrame::RTN, [10.0, 0.0, 0.0])
            .with_delta_mass(-15.0);
        opm.push_maneuver(m);
        assert_eq!(opm.maneuvers.len(), 1);
        assert_eq!(opm.maneuvers[0].delta_mass, Some(-15.0));
    }

    #[test]
    fn test_opm_json_round_trip_via_dispatch() {
        let opm = OPM::from_file("test_assets/ccsds/opm/OPMExample1.txt").unwrap();
        let json_str = opm.to_string(CCSDSFormat::JSON).unwrap();
        assert!(json_str.contains("object_name") || json_str.contains("OBJECT_NAME"));
        let opm2 = OPM::from_str(&json_str).unwrap();
        assert_eq!(opm2.metadata.object_name, opm.metadata.object_name);
        assert_eq!(opm2.metadata.object_id, opm.metadata.object_id);
        assert!((opm2.state_vector.position[0] - opm.state_vector.position[0]).abs() < 1.0);
        assert!((opm2.state_vector.velocity[0] - opm.state_vector.velocity[0]).abs() < 0.001);
    }

    #[test]
    fn test_opm_from_file_nonexistent() {
        let result = OPM::from_file("nonexistent_file.txt");
        assert!(result.is_err());
    }

    #[test]
    fn test_opm_metadata_new() {
        let meta = OPMMetadata::new(
            "ISS".to_string(),
            "1998-067A".to_string(),
            "EARTH".to_string(),
            CCSDSRefFrame::ITRF2000,
            CCSDSTimeSystem::UTC,
        );
        assert_eq!(meta.object_name, "ISS");
        assert_eq!(meta.object_id, "1998-067A");
        assert_eq!(meta.center_name, "EARTH");
        assert!(matches!(meta.ref_frame, CCSDSRefFrame::ITRF2000));
        assert!(matches!(meta.time_system, CCSDSTimeSystem::UTC));
        assert!(meta.ref_frame_epoch.is_none());
        assert!(meta.comments.is_empty());
    }

    #[test]
    fn test_opm_state_vector_new() {
        let epoch = Epoch::from_datetime(2024, 3, 1, 0, 0, 0.0, 0.0, TimeSystem::UTC);
        let pos = [6503.514e3, 1239.647e3, -717.490e3];
        let vel = [-873.160, 8740.420, -4191.076];
        let sv = OPMStateVector::new(epoch, pos, vel);

        assert!((sv.position[0] - 6503.514e3).abs() < 1e-6);
        assert!((sv.position[1] - 1239.647e3).abs() < 1e-6);
        assert!((sv.position[2] - (-717.490e3)).abs() < 1e-6);
        assert!((sv.velocity[0] - (-873.160)).abs() < 1e-6);
        assert!((sv.velocity[1] - 8740.420).abs() < 1e-6);
        assert!((sv.velocity[2] - (-4191.076)).abs() < 1e-6);
        assert!(sv.comments.is_empty());
    }

    #[test]
    fn test_opm_new() {
        let meta = OPMMetadata::new(
            "SAT1".to_string(),
            "2024-001A".to_string(),
            "EARTH".to_string(),
            CCSDSRefFrame::GCRF,
            CCSDSTimeSystem::UTC,
        );
        let epoch = Epoch::from_datetime(2024, 3, 1, 0, 0, 0.0, 0.0, TimeSystem::UTC);
        let sv = OPMStateVector::new(epoch, [7000e3, 0.0, 0.0], [0.0, 7500.0, 0.0]);
        let opm = OPM::new("TEST_ORG".to_string(), meta, sv);

        assert_eq!(opm.header.originator, "TEST_ORG");
        assert!((opm.header.format_version - 3.0).abs() < 1e-15);
        assert!(opm.header.classification.is_none());
        assert!(opm.header.message_id.is_none());
        assert_eq!(opm.metadata.object_name, "SAT1");
        assert_eq!(opm.metadata.object_id, "2024-001A");
        assert!(opm.keplerian_elements.is_none());
        assert!(opm.spacecraft_parameters.is_none());
        assert!(opm.covariance.is_none());
        assert!(opm.maneuvers.is_empty());
        assert!(opm.user_defined.is_none());
    }

    #[test]
    fn test_opm_kvn_parse_example1() {
        let opm = OPM::from_file("test_assets/ccsds/opm/OPMExample1.txt").unwrap();
        assert_eq!(opm.metadata.object_name, "GODZILLA 5");
        assert_eq!(opm.metadata.object_id, "1998-999A");
        assert_eq!(opm.metadata.center_name, "EARTH");
        assert!(matches!(opm.metadata.ref_frame, CCSDSRefFrame::ITRF2000));
        assert!(matches!(opm.metadata.time_system, CCSDSTimeSystem::UTC));
        // Position in meters (converted from km in file)
        assert!((opm.state_vector.position[0] - 6503.514e3).abs() < 1.0);
        assert!((opm.state_vector.position[1] - 1239.647e3).abs() < 1.0);
        assert!((opm.state_vector.position[2] - (-717.490e3)).abs() < 1.0);
        // Velocity in m/s (converted from km/s in file)
        assert!((opm.state_vector.velocity[0] - (-873.160)).abs() < 0.01);
        assert!((opm.state_vector.velocity[1] - 8740.420).abs() < 0.01);
        assert!((opm.state_vector.velocity[2] - (-4191.076)).abs() < 0.01);
    }

    #[test]
    fn test_opm_with_keplerian_elements() {
        let opm = OPM::from_file("test_assets/ccsds/opm/OPMExample5.txt").unwrap();
        assert!(
            opm.keplerian_elements.is_some(),
            "OPMExample5 should have Keplerian elements"
        );
        let ke = opm.keplerian_elements.as_ref().unwrap();
        assert!((ke.eccentricity - 0.020842611).abs() < 1e-9);
        assert!((ke.inclination - 0.117746).abs() < 1e-6);
        assert!((ke.ra_of_asc_node - 17.604721).abs() < 1e-6);
        assert!((ke.arg_of_pericenter - 218.242943).abs() < 1e-6);
        assert!(ke.true_anomaly.is_some());
        assert!((ke.true_anomaly.unwrap() - 41.922339).abs() < 1e-6);
    }

    #[test]
    fn test_opm_with_maneuvers() {
        // OPMExample5 has 3 maneuvers
        let opm = OPM::from_file("test_assets/ccsds/opm/OPMExample5.txt").unwrap();
        assert!(
            opm.maneuvers.len() >= 2,
            "OPMExample5 should have maneuvers"
        );
        // Verify first maneuver has expected delta_mass
        assert!(opm.maneuvers[0].delta_mass.is_some());
        assert!((opm.maneuvers[0].delta_mass.unwrap() - (-18.418)).abs() < 0.001);
        assert!((opm.maneuvers[0].duration - 132.60).abs() < 0.01);
    }

    #[test]
    fn test_opm_kvn_round_trip() {
        let opm = OPM::from_file("test_assets/ccsds/opm/OPMExample1.txt").unwrap();
        let kvn_str = opm.to_string(CCSDSFormat::KVN).unwrap();
        let opm2 = OPM::from_str(&kvn_str).unwrap();
        assert_eq!(opm2.metadata.object_name, opm.metadata.object_name);
        assert_eq!(opm2.metadata.object_id, opm.metadata.object_id);
        // Position round-trip: m → km → m
        assert!((opm2.state_vector.position[0] - opm.state_vector.position[0]).abs() < 1.0);
        assert!((opm2.state_vector.velocity[0] - opm.state_vector.velocity[0]).abs() < 0.001);
        // Spacecraft parameters
        assert!(opm2.spacecraft_parameters.is_some());
        let sp2 = opm2.spacecraft_parameters.as_ref().unwrap();
        assert!((sp2.mass.unwrap() - 3000.0).abs() < 0.01);
    }

    #[test]
    fn test_opm_kvn_round_trip_with_keplerian_and_maneuvers() {
        let opm = OPM::from_file("test_assets/ccsds/opm/OPMExample5.txt").unwrap();
        let kvn_str = opm.to_string(CCSDSFormat::KVN).unwrap();
        let opm2 = OPM::from_str(&kvn_str).unwrap();
        // Keplerian elements
        assert!(opm2.keplerian_elements.is_some());
        let ke1 = opm.keplerian_elements.as_ref().unwrap();
        let ke2 = opm2.keplerian_elements.as_ref().unwrap();
        assert!((ke2.eccentricity - ke1.eccentricity).abs() < 1e-9);
        assert!((ke2.semi_major_axis - ke1.semi_major_axis).abs() < 1.0);
        // Maneuvers
        assert_eq!(opm2.maneuvers.len(), opm.maneuvers.len());
        assert!((opm2.maneuvers[0].duration - opm.maneuvers[0].duration).abs() < 0.01);
        assert!((opm2.maneuvers[0].dv[0] - opm.maneuvers[0].dv[0]).abs() < 0.01);
    }

    #[test]
    fn test_opm_xml_round_trip() {
        let opm = OPM::from_file("test_assets/ccsds/opm/OPMExample3.xml").unwrap();
        let xml_str = opm.to_string(CCSDSFormat::XML).unwrap();
        let opm2 = OPM::from_str(&xml_str).unwrap();
        assert_eq!(opm2.metadata.object_name, opm.metadata.object_name);
        assert_eq!(opm2.metadata.object_id, opm.metadata.object_id);
        assert!((opm2.state_vector.position[0] - opm.state_vector.position[0]).abs() < 1.0);
        assert!((opm2.state_vector.velocity[0] - opm.state_vector.velocity[0]).abs() < 0.001);
        // Covariance
        assert!(opm2.covariance.is_some());
        let cov1 = opm.covariance.as_ref().unwrap();
        let cov2 = opm2.covariance.as_ref().unwrap();
        assert!((cov2.matrix[(0, 0)] - cov1.matrix[(0, 0)]).abs() < 1.0);
    }

    #[test]
    fn test_opm_xml_parse_example3() {
        let opm = OPM::from_file("test_assets/ccsds/opm/OPMExample3.xml").unwrap();
        assert_eq!(opm.metadata.object_name, "OSPREY 5");
        assert_eq!(opm.metadata.object_id, "1998-999A");
        assert_eq!(opm.metadata.center_name, "EARTH");
        assert!(matches!(opm.metadata.ref_frame, CCSDSRefFrame::TOD));
        assert!(opm.metadata.ref_frame_epoch.is_some());
        // Position: 6503.514 km → m
        assert!((opm.state_vector.position[0] - 6503514.0).abs() < 1.0);
        assert!((opm.state_vector.velocity[0] - (-873.16)).abs() < 0.01);
        // Spacecraft parameters
        assert!(opm.spacecraft_parameters.is_some());
        let sp = opm.spacecraft_parameters.as_ref().unwrap();
        assert!((sp.mass.unwrap() - 3000.0).abs() < 0.01);
        // Covariance
        assert!(opm.covariance.is_some());
        let cov = opm.covariance.as_ref().unwrap();
        assert_eq!(cov.cov_ref_frame.as_ref().unwrap(), &CCSDSRefFrame::ITRF97);
        // CX_X = 0.316 km² = 316000 m²
        assert!((cov.matrix[(0, 0)] - 316000.0).abs() < 1.0);
    }

    #[test]
    fn test_opm_to_file_kvn() {
        let opm = OPM::from_file("test_assets/ccsds/opm/OPMExample1.txt").unwrap();
        let dir = std::env::temp_dir();
        let path = dir.join("brahe_test_opm.txt");
        opm.to_file(&path, CCSDSFormat::KVN).unwrap();
        let opm2 = OPM::from_file(&path).unwrap();
        assert_eq!(opm2.metadata.object_name, opm.metadata.object_name);
        std::fs::remove_file(&path).ok();
    }

    #[test]
    fn test_opm_maneuver_new() {
        let epoch = Epoch::from_datetime(2024, 6, 1, 9, 0, 0.0, 0.0, TimeSystem::UTC);
        let m = OPMManeuver::new(epoch, 132.6, CCSDSRefFrame::EME2000, [10.0, -5.0, 2.0]);
        assert!((m.duration - 132.6).abs() < 1e-15);
        assert!(matches!(m.ref_frame, CCSDSRefFrame::EME2000));
        assert!((m.dv[0] - 10.0).abs() < 1e-15);
        assert!((m.dv[1] - (-5.0)).abs() < 1e-15);
        assert!((m.dv[2] - 2.0).abs() < 1e-15);
        assert!(m.delta_mass.is_none());
        assert!(m.comments.is_empty());
    }

    #[test]
    fn test_opm_maneuver_with_delta_mass() {
        let epoch = Epoch::from_datetime(2024, 6, 1, 9, 0, 0.0, 0.0, TimeSystem::UTC);
        let m = OPMManeuver::new(epoch, 100.0, CCSDSRefFrame::RTN, [1.0, 0.0, 0.0])
            .with_delta_mass(-10.5);
        assert_eq!(m.delta_mass, Some(-10.5));
    }

    #[test]
    fn test_opm_push_maneuver() {
        let meta = OPMMetadata::new(
            "SAT".to_string(),
            "2024-001A".to_string(),
            "EARTH".to_string(),
            CCSDSRefFrame::GCRF,
            CCSDSTimeSystem::UTC,
        );
        let sv = OPMStateVector::new(Epoch::now(), [7000e3, 0.0, 0.0], [0.0, 7500.0, 0.0]);
        let mut opm = OPM::new("ORG".to_string(), meta, sv);
        assert_eq!(opm.maneuvers.len(), 0);

        let epoch = Epoch::from_datetime(2024, 6, 1, 9, 0, 0.0, 0.0, TimeSystem::UTC);
        opm.push_maneuver(OPMManeuver::new(
            epoch,
            60.0,
            CCSDSRefFrame::RTN,
            [1.0, 0.0, 0.0],
        ));
        opm.push_maneuver(OPMManeuver::new(
            epoch,
            30.0,
            CCSDSRefFrame::RTN,
            [0.0, 1.0, 0.0],
        ));
        assert_eq!(opm.maneuvers.len(), 2);
        assert!((opm.maneuvers[0].duration - 60.0).abs() < 1e-15);
        assert!((opm.maneuvers[1].duration - 30.0).abs() < 1e-15);
    }

    #[test]
    fn test_opm_kvn_parse_spacecraft_params() {
        // OPMExample1 has spacecraft parameters
        let opm = OPM::from_file("test_assets/ccsds/opm/OPMExample1.txt").unwrap();
        assert!(opm.spacecraft_parameters.is_some());
        let sp = opm.spacecraft_parameters.as_ref().unwrap();
        assert!(sp.mass.is_some());
        assert!((sp.mass.unwrap() - 3000.0).abs() < 0.01);
    }

    #[test]
    fn test_opm_to_json_string_upper_key_case() {
        let opm = OPM::from_file("test_assets/ccsds/opm/OPMExample1.txt").unwrap();
        let json_str = opm.to_json_string(CCSDSJsonKeyCase::Upper).unwrap();
        assert!(json_str.contains("OBJECT_NAME"));
        assert!(json_str.contains("OBJECT_ID"));
    }

    /// Helper: assert all OPM fields match between original and round-tripped.
    fn assert_opm_fields_match(opm1: &OPM, opm2: &OPM) {
        // Header
        assert_eq!(opm1.header.format_version, opm2.header.format_version);
        assert_eq!(opm1.header.originator, opm2.header.originator);
        assert_eq!(opm1.header.classification, opm2.header.classification);
        assert_eq!(opm1.header.message_id, opm2.header.message_id);

        // Metadata
        assert_eq!(opm1.metadata.object_name, opm2.metadata.object_name);
        assert_eq!(opm1.metadata.object_id, opm2.metadata.object_id);
        assert_eq!(opm1.metadata.center_name, opm2.metadata.center_name);
        assert_eq!(opm1.metadata.ref_frame, opm2.metadata.ref_frame);
        assert_eq!(opm1.metadata.time_system, opm2.metadata.time_system);

        // State vector
        for i in 0..3 {
            assert!(
                (opm1.state_vector.position[i] - opm2.state_vector.position[i]).abs() < 1.0,
                "position[{}] mismatch: {} vs {}",
                i,
                opm1.state_vector.position[i],
                opm2.state_vector.position[i]
            );
            assert!(
                (opm1.state_vector.velocity[i] - opm2.state_vector.velocity[i]).abs() < 0.001,
                "velocity[{}] mismatch: {} vs {}",
                i,
                opm1.state_vector.velocity[i],
                opm2.state_vector.velocity[i]
            );
        }

        // Keplerian elements
        assert_eq!(
            opm1.keplerian_elements.is_some(),
            opm2.keplerian_elements.is_some()
        );
        if let (Some(ke1), Some(ke2)) = (&opm1.keplerian_elements, &opm2.keplerian_elements) {
            assert!((ke1.semi_major_axis - ke2.semi_major_axis).abs() < 1.0);
            assert!((ke1.eccentricity - ke2.eccentricity).abs() < 1e-9);
            assert!((ke1.inclination - ke2.inclination).abs() < 1e-6);
            assert!((ke1.ra_of_asc_node - ke2.ra_of_asc_node).abs() < 1e-6);
            assert!((ke1.arg_of_pericenter - ke2.arg_of_pericenter).abs() < 1e-6);
            assert_eq!(ke1.true_anomaly.is_some(), ke2.true_anomaly.is_some());
            if let (Some(ta1), Some(ta2)) = (ke1.true_anomaly, ke2.true_anomaly) {
                assert!((ta1 - ta2).abs() < 1e-6);
            }
            assert_eq!(ke1.mean_anomaly.is_some(), ke2.mean_anomaly.is_some());
            assert_eq!(ke1.gm.is_some(), ke2.gm.is_some());
            if let (Some(gm1), Some(gm2)) = (ke1.gm, ke2.gm) {
                assert!((gm1 - gm2).abs() < 1e3);
            }
        }

        // Spacecraft parameters
        assert_eq!(
            opm1.spacecraft_parameters.is_some(),
            opm2.spacecraft_parameters.is_some()
        );
        if let (Some(sp1), Some(sp2)) = (&opm1.spacecraft_parameters, &opm2.spacecraft_parameters) {
            assert_eq!(sp1.mass.is_some(), sp2.mass.is_some());
            if let (Some(m1), Some(m2)) = (sp1.mass, sp2.mass) {
                assert!((m1 - m2).abs() < 0.01);
            }
            assert_eq!(sp1.solar_rad_area.is_some(), sp2.solar_rad_area.is_some());
            if let (Some(a1), Some(a2)) = (sp1.solar_rad_area, sp2.solar_rad_area) {
                assert!((a1 - a2).abs() < 0.01);
            }
            assert_eq!(sp1.solar_rad_coeff.is_some(), sp2.solar_rad_coeff.is_some());
            if let (Some(c1), Some(c2)) = (sp1.solar_rad_coeff, sp2.solar_rad_coeff) {
                assert!((c1 - c2).abs() < 0.01);
            }
            assert_eq!(sp1.drag_area.is_some(), sp2.drag_area.is_some());
            if let (Some(a1), Some(a2)) = (sp1.drag_area, sp2.drag_area) {
                assert!((a1 - a2).abs() < 0.01);
            }
            assert_eq!(sp1.drag_coeff.is_some(), sp2.drag_coeff.is_some());
            if let (Some(c1), Some(c2)) = (sp1.drag_coeff, sp2.drag_coeff) {
                assert!((c1 - c2).abs() < 0.01);
            }
        }

        // Covariance
        assert_eq!(opm1.covariance.is_some(), opm2.covariance.is_some());
        if let (Some(cov1), Some(cov2)) = (&opm1.covariance, &opm2.covariance) {
            assert_eq!(cov1.cov_ref_frame, cov2.cov_ref_frame);
            for i in 0..6 {
                for j in 0..6 {
                    let rel = if cov1.matrix[(i, j)].abs() > 1e-20 {
                        ((cov1.matrix[(i, j)] - cov2.matrix[(i, j)]) / cov1.matrix[(i, j)]).abs()
                    } else {
                        (cov1.matrix[(i, j)] - cov2.matrix[(i, j)]).abs()
                    };
                    assert!(
                        rel < 1e-4,
                        "cov({},{}) mismatch: {} vs {}",
                        i,
                        j,
                        cov1.matrix[(i, j)],
                        cov2.matrix[(i, j)]
                    );
                }
            }
        }

        // Maneuvers
        assert_eq!(opm1.maneuvers.len(), opm2.maneuvers.len());
        for (m1, m2) in opm1.maneuvers.iter().zip(opm2.maneuvers.iter()) {
            assert!((m1.duration - m2.duration).abs() < 0.01);
            assert_eq!(m1.ref_frame, m2.ref_frame);
            for i in 0..3 {
                assert!((m1.dv[i] - m2.dv[i]).abs() < 0.01);
            }
            assert_eq!(m1.delta_mass.is_some(), m2.delta_mass.is_some());
            if let (Some(dm1), Some(dm2)) = (m1.delta_mass, m2.delta_mass) {
                assert!((dm1 - dm2).abs() < 0.01);
            }
        }

        // User-defined parameters
        assert_eq!(opm1.user_defined.is_some(), opm2.user_defined.is_some());
        if let (Some(ud1), Some(ud2)) = (&opm1.user_defined, &opm2.user_defined) {
            assert_eq!(ud1.parameters.len(), ud2.parameters.len());
            for (k, v) in &ud1.parameters {
                assert_eq!(
                    ud2.parameters.get(k),
                    Some(v),
                    "user_defined key {} mismatch",
                    k
                );
            }
        }
    }

    #[test]
    fn test_opm_kvn_full_round_trip() {
        // OPMExample4 has: state, Keplerian, spacecraft params, covariance, user_defined
        let opm1 = OPM::from_file("test_assets/ccsds/opm/OPMExample4.txt").unwrap();
        let kvn = opm1.to_string(CCSDSFormat::KVN).unwrap();
        let opm2 = OPM::from_str(&kvn).unwrap();
        assert_opm_fields_match(&opm1, &opm2);
    }

    #[test]
    fn test_opm_xml_full_round_trip() {
        // Start from KVN (richer fields), write to XML, re-parse
        let opm1 = OPM::from_file("test_assets/ccsds/opm/OPMExample4.txt").unwrap();
        let xml = opm1.to_string(CCSDSFormat::XML).unwrap();
        let opm2 = OPM::from_str(&xml).unwrap();
        assert_opm_fields_match(&opm1, &opm2);
    }

    #[test]
    fn test_opm_json_full_round_trip() {
        let opm1 = OPM::from_file("test_assets/ccsds/opm/OPMExample4.txt").unwrap();
        let json = opm1.to_string(CCSDSFormat::JSON).unwrap();
        let opm2 = OPM::from_str(&json).unwrap();
        assert_opm_fields_match(&opm1, &opm2);
    }

    #[test]
    fn test_opm_kvn_full_round_trip_with_maneuvers() {
        // OPMExample5 has Keplerian elements and maneuvers
        let opm1 = OPM::from_file("test_assets/ccsds/opm/OPMExample5.txt").unwrap();
        let kvn = opm1.to_string(CCSDSFormat::KVN).unwrap();
        let opm2 = OPM::from_str(&kvn).unwrap();
        assert_opm_fields_match(&opm1, &opm2);
    }
}