starfield 0.12.1

Astronomical data reduction toolkit with star catalogs, coordinate systems, and star finding algorithms (inspired by skyfield)
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
//! Two-body Keplerian orbit propagation
//!
//! Port of Skyfield's `keplerlib.py` — propagates orbits using two-body
//! (Keplerian) dynamics for comets, asteroids, and other bodies defined
//! by classical orbital elements.
//!
//! # Example
//!
//! ```ignore
//! use starfield::keplerlib::KeplerOrbit;
//! use starfield::time::Timescale;
//! use starfield::constants::GM_SUN;
//!
//! let ts = Timescale::default();
//! let epoch = ts.tt_jd(2451545.0, None);
//!
//! // Halley's Comet (approximate elements)
//! let halley = KeplerOrbit::from_periapsis(
//!     1.058, 0.967, 162.26, 58.42, 111.33,
//!     &epoch, GM_SUN, Some(10), Some("1P/Halley"),
//! );
//!
//! let t = ts.tt_jd(2458000.0, None);
//! let pos = halley.at(&t);
//! ```

#[cfg(all(test, feature = "python-tests"))]
mod python_tests;

use std::f64::consts::PI;

use nalgebra::{Matrix3, Vector3};

#[cfg(test)]
use crate::constants::GM_SUN;
use crate::constants::{AU_KM, DAY_S, DEG2RAD};
use crate::positions::Position;
use crate::time::Time;

/// Convert GM from km³/s² to AU³/day²
const CONVERT_GM: f64 = DAY_S * DAY_S / (AU_KM * AU_KM * AU_KM);

/// A two-body Keplerian orbit
///
/// Stores the state vector (position + velocity) at an epoch and
/// propagates to arbitrary times using universal variable formulation.
#[derive(Debug, Clone)]
pub struct KeplerOrbit {
    /// Position at epoch in AU
    pub position_au: Vector3<f64>,
    /// Velocity at epoch in AU/day
    pub velocity_au_per_day: Vector3<f64>,
    /// Epoch as TT Julian date
    pub epoch_tt: f64,
    /// Gravitational parameter in AU³/day²
    pub mu_au3_d2: f64,
    /// NAIF ID of the central body (10 = Sun)
    pub center: i32,
    /// Name of the orbiting body
    pub target_name: Option<String>,
    /// Optional rotation matrix (e.g., ECLIPJ2000 → equatorial)
    rotation: Option<Matrix3<f64>>,
}

impl KeplerOrbit {
    /// Create a KeplerOrbit from position/velocity state vectors
    ///
    /// Position in AU, velocity in AU/day, GM in AU³/day².
    pub fn new(
        position_au: Vector3<f64>,
        velocity_au_per_day: Vector3<f64>,
        epoch: &Time,
        mu_au3_d2: f64,
        center: Option<i32>,
        target_name: Option<&str>,
    ) -> Self {
        KeplerOrbit {
            position_au,
            velocity_au_per_day,
            epoch_tt: epoch.tt(),
            mu_au3_d2,
            center: center.unwrap_or(10),
            target_name: target_name.map(|s| s.to_string()),
            rotation: None,
        }
    }

    /// Build from periapsis parameters (used for comets)
    ///
    /// # Arguments
    /// * `semilatus_rectum_au` — semi-latus rectum in AU
    /// * `eccentricity` — orbital eccentricity
    /// * `inclination_deg` — inclination in degrees
    /// * `longitude_of_ascending_node_deg` — Ω in degrees
    /// * `argument_of_perihelion_deg` — ω in degrees
    /// * `epoch` — time of periapsis passage
    /// * `gm_km3_s2` — GM in km³/s²
    /// * `center` — NAIF ID of central body
    /// * `target_name` — name of orbiting body
    #[allow(clippy::too_many_arguments)]
    pub fn from_periapsis(
        semilatus_rectum_au: f64,
        eccentricity: f64,
        inclination_deg: f64,
        longitude_of_ascending_node_deg: f64,
        argument_of_perihelion_deg: f64,
        epoch: &Time,
        gm_km3_s2: f64,
        center: Option<i32>,
        target_name: Option<&str>,
    ) -> Self {
        let gm_au3_d2 = gm_km3_s2 * CONVERT_GM;
        let (pos, vel) = ele_to_vec(
            semilatus_rectum_au,
            eccentricity,
            DEG2RAD * inclination_deg,
            DEG2RAD * longitude_of_ascending_node_deg,
            DEG2RAD * argument_of_perihelion_deg,
            0.0, // true anomaly = 0 at periapsis
            gm_au3_d2,
        );
        KeplerOrbit {
            position_au: pos,
            velocity_au_per_day: vel,
            epoch_tt: epoch.tt(),
            mu_au3_d2: gm_au3_d2,
            center: center.unwrap_or(10),
            target_name: target_name.map(|s| s.to_string()),
            rotation: None,
        }
    }

    /// Build from mean anomaly (used for asteroids / MPC data)
    ///
    /// # Arguments
    /// * `semilatus_rectum_au` — semi-latus rectum in AU
    /// * `eccentricity` — orbital eccentricity
    /// * `inclination_deg` — inclination in degrees
    /// * `longitude_of_ascending_node_deg` — Ω in degrees
    /// * `argument_of_perihelion_deg` — ω in degrees
    /// * `mean_anomaly_deg` — mean anomaly M in degrees
    /// * `epoch` — epoch of elements
    /// * `gm_km3_s2` — GM in km³/s²
    /// * `center` — NAIF ID of central body
    /// * `target_name` — name of orbiting body
    #[allow(clippy::too_many_arguments)]
    pub fn from_mean_anomaly(
        semilatus_rectum_au: f64,
        eccentricity: f64,
        inclination_deg: f64,
        longitude_of_ascending_node_deg: f64,
        argument_of_perihelion_deg: f64,
        mean_anomaly_deg: f64,
        epoch: &Time,
        gm_km3_s2: f64,
        center: Option<i32>,
        target_name: Option<&str>,
    ) -> Self {
        let m = DEG2RAD * mean_anomaly_deg;
        let gm_au3_d2 = gm_km3_s2 * CONVERT_GM;

        let v = if eccentricity < 1.0 {
            let ea = eccentric_anomaly(eccentricity, m);
            true_anomaly_closed(eccentricity, ea)
        } else if eccentricity > 1.0 {
            let ea = eccentric_anomaly(eccentricity, m);
            true_anomaly_hyperbolic(eccentricity, ea)
        } else {
            true_anomaly_parabolic(semilatus_rectum_au, gm_au3_d2, m)
        };

        let (pos, vel) = ele_to_vec(
            semilatus_rectum_au,
            eccentricity,
            DEG2RAD * inclination_deg,
            DEG2RAD * longitude_of_ascending_node_deg,
            DEG2RAD * argument_of_perihelion_deg,
            v,
            gm_au3_d2,
        );
        KeplerOrbit {
            position_au: pos,
            velocity_au_per_day: vel,
            epoch_tt: epoch.tt(),
            mu_au3_d2: gm_au3_d2,
            center: center.unwrap_or(10),
            target_name: target_name.map(|s| s.to_string()),
            rotation: None,
        }
    }

    /// Set the rotation matrix to convert from ecliptic (ECLIPJ2000) to equatorial
    ///
    /// The ECLIPJ2000 frame matrix rotates equatorial → ecliptic.
    /// Its transpose rotates ecliptic → equatorial, which is what we need
    /// for orbits defined in the ecliptic plane (comets, asteroids from MPC data).
    pub fn set_ecliptic_rotation(&mut self) {
        // ECLIPJ2000 rotation matrix (equatorial → ecliptic) from SPICE/Skyfield.
        // This is Rx(obliquity) where obliquity = 23.4392911 degrees.
        #[rustfmt::skip]
        let eclip = Matrix3::new(
            1.0,  0.0,                    0.0,
            0.0,  0.917_482_062_069_181_8,  0.397_777_155_931_913_7,
            0.0, -0.397_777_155_931_913_7,  0.917_482_062_069_181_8,
        );
        // Transpose: ecliptic → equatorial
        self.rotation = Some(eclip.transpose());
    }

    /// Propagate the orbit to the given time
    ///
    /// Returns a Barycentric `Position` in the ICRF.
    pub fn at(&self, time: &Time) -> Position {
        let (mut pos, mut vel) = propagate(
            &self.position_au,
            &self.velocity_au_per_day,
            self.epoch_tt,
            time.tt(),
            self.mu_au3_d2,
        );

        if let Some(rot) = &self.rotation {
            pos = rot * pos;
            vel = rot * vel;
        }

        Position::barycentric(pos, vel, self.center)
    }
}

impl std::fmt::Display for KeplerOrbit {
    fn fmt(&self, f: &mut std::fmt::Formatter<'_>) -> std::fmt::Result {
        if let Some(name) = &self.target_name {
            write!(f, "KeplerOrbit {}{}", self.center, name)
        } else {
            write!(
                f,
                "KeplerOrbit {} → e={:.3} a={:.2} AU",
                self.center,
                0.0, // placeholder
                0.0,
            )
        }
    }
}

/// Build a comet orbit from MPC-style parameters
///
/// This is the Rust equivalent of Skyfield's `comet_orbit()`.
#[allow(clippy::too_many_arguments)]
pub fn comet_orbit(
    perihelion_distance_au: f64,
    eccentricity: f64,
    inclination_deg: f64,
    longitude_of_ascending_node_deg: f64,
    argument_of_perihelion_deg: f64,
    t_perihelion: &Time,
    gm_km3_s2: f64,
    target_name: Option<&str>,
) -> KeplerOrbit {
    let p = if eccentricity == 1.0 {
        perihelion_distance_au * 2.0
    } else {
        let a = perihelion_distance_au / (1.0 - eccentricity);
        a * (1.0 - eccentricity * eccentricity)
    };

    let mut orbit = KeplerOrbit::from_periapsis(
        p,
        eccentricity,
        inclination_deg,
        longitude_of_ascending_node_deg,
        argument_of_perihelion_deg,
        t_perihelion,
        gm_km3_s2,
        Some(10),
        target_name,
    );
    orbit.set_ecliptic_rotation();
    orbit
}

/// Build a minor planet orbit from MPC-style parameters
///
/// This is the Rust equivalent of Skyfield's `mpcorb_orbit()`.
#[allow(clippy::too_many_arguments)]
pub fn mpcorb_orbit(
    semimajor_axis_au: f64,
    eccentricity: f64,
    inclination_deg: f64,
    longitude_of_ascending_node_deg: f64,
    argument_of_perihelion_deg: f64,
    mean_anomaly_deg: f64,
    epoch: &Time,
    gm_km3_s2: f64,
    target_name: Option<&str>,
) -> KeplerOrbit {
    let p = semimajor_axis_au * (1.0 - eccentricity * eccentricity);

    let mut orbit = KeplerOrbit::from_mean_anomaly(
        p,
        eccentricity,
        inclination_deg,
        longitude_of_ascending_node_deg,
        argument_of_perihelion_deg,
        mean_anomaly_deg,
        epoch,
        gm_km3_s2,
        Some(10),
        target_name,
    );
    orbit.set_ecliptic_rotation();
    orbit
}

// ---------------------------------------------------------------------------
// Core orbital mechanics functions
// ---------------------------------------------------------------------------

/// Normalize angle to [-π, π]
fn normpi(m: f64) -> f64 {
    let mut x = m % (2.0 * PI);
    if x > PI {
        x -= 2.0 * PI;
    }
    if x < -PI {
        x += 2.0 * PI;
    }
    x
}

/// Solve Kepler's equation for eccentric anomaly
///
/// Iterative solver following arXiv:2108.03215.
/// Works for both elliptic (e < 1) and hyperbolic (e > 1) orbits.
pub(crate) fn eccentric_anomaly(e: f64, m: f64) -> f64 {
    let m = normpi(m);
    let sign_m = m.signum();
    let m = m * sign_m;

    let ebar = 0.25 * PI / e - 1.0;
    let mut ea = 0.5 * PI * ebar * (ebar.signum() * (1.0 + m / (e * ebar * ebar)).sqrt() - 1.0);

    for _ in 0..10 {
        let f1 = 1.0 - e * ea.cos();
        let f2 = e * ea.sin();
        let f = ea - f2 - m;
        let d_ea = f * f1 / (f1 * f1 - 0.5 * f * f2);
        ea -= d_ea;
        if d_ea.abs() < 1e-14 {
            return ea * sign_m;
        }
    }

    ea * sign_m
}

/// True anomaly from eccentric anomaly for closed (elliptic) orbits
fn true_anomaly_closed(e: f64, ea: f64) -> f64 {
    2.0 * (((1.0 + e) / (1.0 - e)).sqrt() * (ea / 2.0).tan()).atan()
}

/// True anomaly from eccentric anomaly for hyperbolic orbits
fn true_anomaly_hyperbolic(e: f64, ea: f64) -> f64 {
    2.0 * (((e + 1.0) / (e - 1.0)).sqrt() * (ea / 2.0).tanh()).atan()
}

/// True anomaly for parabolic orbits
fn true_anomaly_parabolic(p: f64, gm: f64, m: f64) -> f64 {
    let delta_t = (2.0 * p.powi(3) / gm).sqrt() * m;
    let periapsis_distance = p / 2.0;
    let a = 1.5 * (gm / (2.0 * periapsis_distance.powi(3))).sqrt() * delta_t;
    let b = (a + (a * a + 1.0)).cbrt();
    2.0 * (b - 1.0 / b).atan()
}

/// Convert orbital elements to position and velocity state vectors
///
/// Based on CCAR equations:
/// <https://web.archive.org/web/*/http://ccar.colorado.edu/asen5070/handouts/kep2cart_2002.doc>
///
/// # Arguments
/// * `p` — semi-latus rectum (AU)
/// * `e` — eccentricity
/// * `i` — inclination (radians)
/// * `om` — longitude of ascending node Ω (radians)
/// * `w` — argument of periapsis ω (radians)
/// * `v` — true anomaly ν (radians)
/// * `mu` — gravitational parameter (AU³/day²)
pub(crate) fn ele_to_vec(
    p: f64,
    e: f64,
    i: f64,
    om: f64,
    w: f64,
    v: f64,
    mu: f64,
) -> (Vector3<f64>, Vector3<f64>) {
    let r = p / (1.0 + e * v.cos());
    let h = (p * mu).sqrt();
    let u = v + w;

    let (sin_om, cos_om) = om.sin_cos();
    let (sin_u, cos_u) = u.sin_cos();
    let cos_i = i.cos();
    let sin_i = i.sin();

    let x = r * (cos_om * cos_u - sin_om * sin_u * cos_i);
    let y = r * (sin_om * cos_u + cos_om * sin_u * cos_i);
    let z = r * (sin_i * sin_u);

    let he_rp = h * e / (r * p) * v.sin();
    let h_r = h / r;

    let x_dot = x * he_rp - h_r * (cos_om * sin_u + sin_om * cos_u * cos_i);
    let y_dot = y * he_rp - h_r * (sin_om * sin_u - cos_om * cos_u * cos_i);
    let z_dot = z * he_rp + h_r * sin_i * cos_u;

    (Vector3::new(x, y, z), Vector3::new(x_dot, y_dot, z_dot))
}

/// Stumpff functions c0, c1, c2, c3
///
/// Based on SPICE toolkit's `stmp03.f`.
fn stumpff(x: f64) -> (f64, f64, f64, f64) {
    let z = x.abs().sqrt();

    if x < -1.0 {
        // Hyperbolic regime
        let c0 = z.cosh();
        let c1 = z.sinh() / z;
        let c2 = (1.0 - c0) / x;
        let c3 = (1.0 - c1) / x;
        (c0, c1, c2, c3)
    } else if x > 1.0 {
        // Elliptic regime
        let c0 = z.cos();
        let c1 = z.sin() / z;
        let c2 = (1.0 - c0) / x;
        let c3 = (1.0 - c1) / x;
        (c0, c1, c2, c3)
    } else {
        // Near-zero: use series expansion for c2 and c3
        // c2 = 1/2! - x/4! + x²/6! - x³/8! + ...
        // c3 = 1/3! - x/5! + x²/7! - x³/9! + ...
        let mut c2 = 0.0;
        let mut c3 = 0.0;
        let mut xn = 1.0;
        // Even factorials for c2: 2!, 4!, 6!, 8!, ...
        let c2_facts: [f64; 8] = [
            2.0,
            24.0,
            720.0,
            40320.0,
            3628800.0,
            479001600.0,
            87178291200.0,
            20922789888000.0,
        ];
        // Odd factorials for c3: 3!, 5!, 7!, 9!, ...
        let c3_facts: [f64; 8] = [
            6.0,
            120.0,
            5040.0,
            362880.0,
            39916800.0,
            6227020800.0,
            1307674368000.0,
            355687428096000.0,
        ];
        let mut sign = 1.0;
        for k in 0..8 {
            c2 += sign * xn / c2_facts[k];
            c3 += sign * xn / c3_facts[k];
            xn *= x;
            sign = -sign;
        }
        let c1 = 1.0 - x * c3;
        let c0 = 1.0 - x * c2;
        (c0, c1, c2, c3)
    }
}

/// Propagate an orbit from t0 to t1 using universal variable formulation
///
/// Based on SPICE toolkit's `prop2b.f` via Skyfield's `propagate()`.
///
/// # Arguments
/// * `position` — position at epoch (AU)
/// * `velocity` — velocity at epoch (AU/day)
/// * `t0` — epoch TT Julian date
/// * `t1` — target TT Julian date
/// * `gm` — gravitational parameter (AU³/day²)
///
/// # Returns
/// (position, velocity) at time t1
pub(crate) fn propagate(
    position: &Vector3<f64>,
    velocity: &Vector3<f64>,
    t0: f64,
    t1: f64,
    gm: f64,
) -> (Vector3<f64>, Vector3<f64>) {
    let r0 = position.norm();
    let rv = position.dot(velocity);

    let hvec = position.cross(velocity);
    let h2 = hvec.dot(&hvec);

    let eqvec = velocity.cross(&hvec) / gm - position / r0;
    let e = eqvec.norm();
    let q = h2 / (gm * (1.0 + e));

    let f = 1.0 - e;
    let b = (q / gm).sqrt();

    let br0 = b * r0;
    let b2rv = b * b * rv;
    let bq = b * q;
    let qovr0 = q / r0;

    let maxc = br0
        .abs()
        .max(b2rv.abs())
        .max(bq.abs())
        .max((qovr0 / bq).abs());

    // Compute upper bound for x
    let bound = if f < 0.0 {
        // Hyperbolic
        let fixed = (f64::MAX / 2.0).ln() - maxc.ln();
        let rootf = (-f).sqrt();
        let logf = (-f).ln();
        (fixed / rootf).min((fixed + 1.5 * logf) / rootf)
    } else {
        // Elliptic or parabolic
        let logbound = (1.5_f64.ln() + (f64::MAX).ln() - maxc.ln()) / 3.0;
        logbound.exp()
    };

    let dt = t1 - t0;

    // Kepler function: x*(br0*c1 + x*(b2rv*c2 + x*bq*c3))
    let kepler = |x: f64| -> f64 {
        let (_, c1, c2, c3) = stumpff(f * x * x);
        x * (br0 * c1 + x * (b2rv * c2 + x * bq * c3))
    };

    // Initial guess
    let mut x = (dt / bq).clamp(-bound, bound);
    let mut kfun = kepler(x);

    // Bracket the root
    let mut lower;
    let mut upper;

    if dt < 0.0 {
        upper = 0.0;
        lower = x;
        while kfun > dt {
            upper = lower;
            lower *= 2.0;
            let old_x = x;
            x = lower.clamp(-bound, bound);
            if x == old_x {
                break;
            }
            kfun = kepler(x);
        }
    } else if dt > 0.0 {
        lower = 0.0;
        upper = x;
        while kfun < dt {
            lower = upper;
            upper *= 2.0;
            let old_x = x;
            x = upper.clamp(-bound, bound);
            if x == old_x {
                break;
            }
            kfun = kepler(x);
        }
    } else {
        // dt == 0: no propagation needed
        return (*position, *velocity);
    }

    // Bisection refinement
    x = if lower <= upper {
        (upper + lower) / 2.0
    } else {
        upper
    };

    let mut lcount = 0;
    let mut mostc = 1000;

    while lower < x && x < upper && lcount < mostc {
        kfun = kepler(x);

        if kfun > dt {
            upper = x;
        } else if kfun < dt {
            lower = x;
        } else {
            upper = x;
            lower = x;
        }

        if mostc > 64 && upper != 0.0 && lower != 0.0 {
            mostc = 64;
            lcount = 0;
        }

        x = if lower > upper {
            upper
        } else {
            (upper + lower) / 2.0
        };

        lcount += 1;
    }

    // Compute final state vectors using Lagrange coefficients
    let (c0, c1, c2, c3) = stumpff(f * x * x);
    let br = br0 * c0 + x * (b2rv * c1 + x * bq * c2);

    let pc = 1.0 - qovr0 * x * x * c2;
    let vc = dt - bq * x.powi(3) * c3;
    let pcdot = -qovr0 / br * x * c1;
    let vcdot = 1.0 - bq / br * x * x * c2;

    let pos = pc * position + vc * velocity;
    let vel = pcdot * position + vcdot * velocity;

    (pos, vel)
}

#[cfg(test)]
mod tests {
    use super::*;
    use crate::time::Timescale;
    use approx::assert_relative_eq;

    fn ts() -> Timescale {
        Timescale::default()
    }

    #[test]
    fn test_stumpff_zero() {
        let (c0, c1, c2, c3) = stumpff(0.0);
        assert_relative_eq!(c0, 1.0, epsilon = 1e-14);
        assert_relative_eq!(c1, 1.0, epsilon = 1e-14);
        assert_relative_eq!(c2, 0.5, epsilon = 1e-14);
        assert_relative_eq!(c3, 1.0 / 6.0, epsilon = 1e-14);
    }

    #[test]
    fn test_stumpff_positive() {
        // x = PI² → c0 = cos(PI) = -1, c1 = sin(PI)/PI ≈ 0
        let x = PI * PI;
        let (c0, c1, _c2, _c3) = stumpff(x);
        assert_relative_eq!(c0, -1.0, epsilon = 1e-10);
        assert_relative_eq!(c1, 0.0, epsilon = 1e-10);
    }

    #[test]
    fn test_stumpff_negative() {
        // x = -1 → hyperbolic: c0 = cosh(1), c1 = sinh(1)
        let (c0, c1, _c2, _c3) = stumpff(-1.0);
        assert_relative_eq!(c0, 1.0_f64.cosh(), epsilon = 1e-10);
        assert_relative_eq!(c1, 1.0_f64.sinh(), epsilon = 1e-10);
    }

    #[test]
    fn test_stumpff_continuity() {
        // Verify series matches trig at boundaries
        let (c0a, c1a, c2a, c3a) = stumpff(0.99);
        let (c0b, c1b, c2b, c3b) = stumpff(1.01);
        assert!((c0a - c0b).abs() < 0.01);
        assert!((c1a - c1b).abs() < 0.01);
        assert!((c2a - c2b).abs() < 0.01);
        assert!((c3a - c3b).abs() < 0.01);
    }

    #[test]
    fn test_eccentric_anomaly_circular() {
        // For e ≈ 0, E ≈ M
        let e = 0.001;
        let m = 1.0;
        let ea = eccentric_anomaly(e, m);
        assert_relative_eq!(ea, m, epsilon = 0.01);
    }

    #[test]
    fn test_eccentric_anomaly_moderate() {
        // e = 0.5, M = 1.0 → known solution
        let e = 0.5;
        let m = 1.0;
        let ea = eccentric_anomaly(e, m);
        // Verify: E - e*sin(E) = M
        let residual = ea - e * ea.sin() - m;
        assert!(residual.abs() < 1e-13, "Kepler residual = {residual}");
    }

    #[test]
    fn test_eccentric_anomaly_high_e() {
        // e = 0.99, M = 0.5
        let e = 0.99;
        let m = 0.5;
        let ea = eccentric_anomaly(e, m);
        let residual = ea - e * ea.sin() - m;
        assert!(residual.abs() < 1e-12, "Kepler residual = {residual}");
    }

    #[test]
    fn test_ele_to_vec_circular() {
        // Circular orbit: e=0, p=a=1 AU, i=0, face-on
        let gm = GM_SUN * CONVERT_GM;
        let (pos, vel) = ele_to_vec(1.0, 0.0, 0.0, 0.0, 0.0, 0.0, gm);
        // At periapsis of circular orbit: r = p = 1 AU, v = 0
        assert_relative_eq!(pos.norm(), 1.0, epsilon = 1e-10);
        // Velocity should be √(GM/a) ≈ √(GM) for a=1
        let expected_v = gm.sqrt();
        assert_relative_eq!(vel.norm(), expected_v, epsilon = 1e-6);
    }

    #[test]
    fn test_propagate_stationary() {
        // dt = 0 should return original state
        let pos = Vector3::new(1.0, 0.0, 0.0);
        let vel = Vector3::new(0.0, 0.01, 0.0);
        let gm = GM_SUN * CONVERT_GM;
        let (p2, v2) = propagate(&pos, &vel, 0.0, 0.0, gm);
        assert_relative_eq!(p2, pos, epsilon = 1e-14);
        assert_relative_eq!(v2, vel, epsilon = 1e-14);
    }

    #[test]
    fn test_propagate_half_period() {
        // Circular orbit: after half period, should be on opposite side
        let gm = GM_SUN * CONVERT_GM;
        let a = 1.0; // 1 AU
        let v_circ = (gm / a).sqrt();
        let pos = Vector3::new(a, 0.0, 0.0);
        let vel = Vector3::new(0.0, v_circ, 0.0);

        let period = 2.0 * PI * (a.powi(3) / gm).sqrt();
        let (p2, _v2) = propagate(&pos, &vel, 0.0, period / 2.0, gm);

        // Should be at roughly (-1, 0, 0)
        assert_relative_eq!(p2.x, -a, epsilon = 1e-6);
        assert!(p2.y.abs() < 1e-6);
    }

    #[test]
    fn test_propagate_full_period() {
        // After one full period, should return to starting position
        let gm = GM_SUN * CONVERT_GM;
        let a = 1.0;
        let v_circ = (gm / a).sqrt();
        let pos = Vector3::new(a, 0.0, 0.0);
        let vel = Vector3::new(0.0, v_circ, 0.0);

        let period = 2.0 * PI * (a.powi(3) / gm).sqrt();
        let (p2, v2) = propagate(&pos, &vel, 0.0, period, gm);

        assert_relative_eq!(p2.x, pos.x, epsilon = 1e-5);
        assert_relative_eq!(p2.y, pos.y, epsilon = 1e-5);
        assert_relative_eq!(v2.x, vel.x, epsilon = 1e-5);
        assert_relative_eq!(v2.y, vel.y, epsilon = 1e-5);
    }

    #[test]
    fn test_comet_orbit_basic() {
        let ts = ts();
        let epoch = ts.tt_jd(2451545.0, None);

        // Simple comet: q=1 AU, e=0.5
        let orbit = comet_orbit(1.0, 0.5, 0.0, 0.0, 0.0, &epoch, GM_SUN, Some("TestComet"));

        // At periapsis (epoch), distance from Sun should be ~1 AU
        let pos = orbit.at(&epoch);
        let dist = pos.position.norm();
        assert!(
            dist > 0.9 && dist < 1.1,
            "Periapsis distance should be ~1 AU, got {dist}"
        );
    }

    #[test]
    fn test_mpcorb_orbit_basic() {
        let ts = ts();
        let epoch = ts.tt_jd(2451545.0, None);

        // Earth-like orbit: a=1, e=0.01
        let orbit = mpcorb_orbit(
            1.0,
            0.01,
            0.0,
            0.0,
            0.0,
            0.0,
            &epoch,
            GM_SUN,
            Some("TestAsteroid"),
        );

        let pos = orbit.at(&epoch);
        let dist = pos.position.norm();
        assert!(
            dist > 0.9 && dist < 1.1,
            "Distance should be ~1 AU, got {dist}"
        );
    }

    #[test]
    fn test_normpi() {
        assert_relative_eq!(normpi(0.0), 0.0, epsilon = 1e-14);
        assert_relative_eq!(normpi(PI), PI, epsilon = 1e-14);
        assert_relative_eq!(normpi(-PI), -PI, epsilon = 1e-14);
        assert_relative_eq!(normpi(3.0 * PI), PI, epsilon = 1e-10);
        assert_relative_eq!(normpi(-3.0 * PI), -PI, epsilon = 1e-10);
    }

    #[test]
    fn test_true_anomaly_closed_at_zero() {
        let v = true_anomaly_closed(0.5, 0.0);
        assert_relative_eq!(v, 0.0, epsilon = 1e-14);
    }

    #[test]
    fn test_true_anomaly_closed_at_pi() {
        let v = true_anomaly_closed(0.5, PI);
        assert_relative_eq!(v, PI, epsilon = 1e-10);
    }

    #[test]
    fn test_display() {
        let ts = ts();
        let epoch = ts.tt_jd(2451545.0, None);
        let orbit = comet_orbit(1.0, 0.5, 0.0, 0.0, 0.0, &epoch, GM_SUN, Some("1P/Halley"));
        let s = format!("{orbit}");
        assert!(s.contains("Halley"));
    }
}