Skip to main content

lox_core/
coords.rs

1// SPDX-FileCopyrightText: 2025 Helge Eichhorn <git@helgeeichhorn.de>
2//
3// SPDX-License-Identifier: MPL-2.0
4
5//! Coordinate types for representing positions, velocities, and trajectories.
6
7use core::f64::consts::{FRAC_PI_2, PI, TAU};
8use std::{
9    ops::{Add, AddAssign, Div, DivAssign, Mul, MulAssign, Neg, Sub, SubAssign},
10    sync::Arc,
11};
12
13use glam::{DMat3, DVec3};
14use lox_test_utils::ApproxEq;
15use thiserror::Error;
16
17use crate::{
18    math::{
19        roots::{FindRoot, RootFinderError, Secant},
20        series::{InterpolationType, Series},
21    },
22    time::deltas::TimeDelta,
23    units::{Angle, Distance, Velocity},
24};
25
26/// Azimuth-elevation pair for representing direction in a topocentric frame.
27#[derive(Copy, Clone, Debug, Default, PartialEq)]
28#[cfg_attr(feature = "serde", derive(serde::Serialize, serde::Deserialize))]
29pub struct AzEl(Angle, Angle);
30
31impl AzEl {
32    /// Returns a new [`AzElBuilder`].
33    pub fn builder() -> AzElBuilder {
34        AzElBuilder::default()
35    }
36
37    /// Returns the azimuth angle.
38    pub fn azimuth(&self) -> Angle {
39        self.0
40    }
41
42    /// Returns the elevation angle.
43    pub fn elevation(&self) -> Angle {
44        self.1
45    }
46}
47
48/// Error returned when constructing an [`AzEl`] with invalid angles.
49#[derive(Copy, Clone, Debug, Error, PartialEq)]
50pub enum AzElError {
51    /// The azimuth angle is outside the valid range [0°, 360°].
52    #[error("azimuth must be between 0 deg and 360 deg but was {0}")]
53    InvalidAzimuth(Angle),
54    /// The elevation angle is outside the valid range [0°, 360°].
55    #[error("elevation must be between 0 deg and 360 deg but was {0}")]
56    InvalidElevation(Angle),
57}
58
59/// Builder for constructing validated [`AzEl`] instances.
60#[derive(Copy, Clone, Debug)]
61pub struct AzElBuilder {
62    azimuth: Result<Angle, AzElError>,
63    elevation: Result<Angle, AzElError>,
64}
65
66impl Default for AzElBuilder {
67    fn default() -> Self {
68        Self::new()
69    }
70}
71
72impl AzElBuilder {
73    /// Creates a new builder with default (zero) angles.
74    pub fn new() -> Self {
75        Self {
76            azimuth: Ok(Angle::default()),
77            elevation: Ok(Angle::default()),
78        }
79    }
80
81    /// Sets the azimuth angle. Must be between 0° and 360°.
82    pub fn azimuth(&mut self, azimuth: Angle) -> &mut Self {
83        self.azimuth = match azimuth.to_radians() {
84            lon if lon < 0.0 => Err(AzElError::InvalidAzimuth(azimuth)),
85            lon if lon > TAU => Err(AzElError::InvalidAzimuth(azimuth)),
86            _ => Ok(azimuth),
87        };
88        self
89    }
90
91    /// Sets the elevation angle. Must be between 0° and 360°.
92    pub fn elevation(&mut self, elevation: Angle) -> &mut Self {
93        self.elevation = match elevation.to_radians() {
94            lat if lat < 0.0 => Err(AzElError::InvalidElevation(elevation)),
95            lat if lat > TAU => Err(AzElError::InvalidElevation(elevation)),
96            _ => Ok(elevation),
97        };
98        self
99    }
100
101    /// Builds the [`AzEl`], returning an error if any angle is invalid.
102    pub fn build(&self) -> Result<AzEl, AzElError> {
103        Ok(AzEl(self.azimuth?, self.elevation?))
104    }
105}
106
107/// Geodetic coordinates: longitude, latitude, and altitude.
108#[derive(Copy, Clone, Debug, Default, PartialEq)]
109#[cfg_attr(feature = "serde", derive(serde::Serialize, serde::Deserialize))]
110pub struct LonLatAlt(Angle, Angle, Distance);
111
112impl LonLatAlt {
113    /// Creates a `LonLatAlt` from longitude and latitude in degrees and altitude in meters.
114    pub fn from_degrees(lon_deg: f64, lat_deg: f64, alt_m: f64) -> Result<Self, LonLatAltError> {
115        LonLatAltBuilder::new()
116            .longitude(Angle::degrees(lon_deg))
117            .latitude(Angle::degrees(lat_deg))
118            .altitude(Distance::meters(alt_m))
119            .build()
120    }
121
122    /// Returns a new [`LonLatAltBuilder`].
123    pub fn builder() -> LonLatAltBuilder {
124        LonLatAltBuilder::default()
125    }
126
127    /// Returns the longitude.
128    pub fn lon(&self) -> Angle {
129        self.0
130    }
131
132    /// Returns the latitude.
133    pub fn lat(&self) -> Angle {
134        self.1
135    }
136
137    /// Returns the altitude.
138    pub fn alt(&self) -> Distance {
139        self.2
140    }
141
142    /// Converts geodetic coordinates (LLA) to body-fixed Cartesian position (meters).
143    pub fn to_body_fixed(&self, equatorial_radius: Distance, flattening: f64) -> DVec3 {
144        let alt = self.alt().to_meters();
145        let (lon_sin, lon_cos) = self.lon().sin_cos();
146        let (lat_sin, lat_cos) = self.lat().sin_cos();
147        let r_eq = equatorial_radius.to_meters();
148        let e = (2.0 * flattening - flattening.powi(2)).sqrt();
149        let c = r_eq / (1.0 - e.powi(2) * lat_sin.powi(2)).sqrt();
150        let s = c * (1.0 - e.powi(2));
151        let r_delta = (c + alt) * lat_cos;
152        let r_kappa = (s + alt) * lat_sin;
153        DVec3::new(r_delta * lon_cos, r_delta * lon_sin, r_kappa)
154    }
155
156    /// Converts a body-fixed Cartesian position (meters) to geodetic coordinates (LLA).
157    ///
158    /// Returns [`FromBodyFixedError::ZeroPosition`] if the position vector has
159    /// zero length. Polar positions (where the equatorial projection is zero)
160    /// are handled as a special case without root-finding.
161    pub fn from_body_fixed(
162        pos: DVec3,
163        equatorial_radius: Distance,
164        flattening: f64,
165    ) -> Result<Self, FromBodyFixedError> {
166        let r_eq = equatorial_radius.to_meters();
167        let rm = pos.length();
168
169        if rm < 1e-10 {
170            return Err(FromBodyFixedError::ZeroPosition);
171        }
172
173        let r_delta = (pos.x.powi(2) + pos.y.powi(2)).sqrt();
174
175        // Polar special case: r_delta ≈ 0 means we're on or near a pole.
176        // The iterative solver divides by r_delta so we handle this directly.
177        if r_delta < 1e-10 {
178            let lat = if pos.z >= 0.0 { FRAC_PI_2 } else { -FRAC_PI_2 };
179            let e = (2.0 * flattening - flattening.powi(2)).sqrt();
180            let r_polar = r_eq * (1.0 - e.powi(2)).sqrt();
181            let alt = pos.z.abs() - r_polar;
182            return Ok(LonLatAlt(
183                Angle::radians(0.0),
184                Angle::radians(lat),
185                Distance::meters(alt),
186            ));
187        }
188
189        let mut lon = pos.y.atan2(pos.x);
190
191        if lon.abs() >= PI {
192            if lon < 0.0 {
193                lon += TAU;
194            } else {
195                lon -= TAU;
196            }
197        }
198
199        let delta = (pos.z / rm).asin();
200
201        let root_finder = Secant::default();
202
203        let f = flattening;
204        let lat = root_finder.find(
205            |lat: f64| {
206                let e = (2.0 * f - f.powi(2)).sqrt();
207                let c = r_eq / (1.0 - e.powi(2) * lat.sin().powi(2)).sqrt();
208                Ok(pos.z + c * e.powi(2) * lat.sin()).map(|v| v / r_delta - lat.tan())
209            },
210            delta,
211        )?;
212
213        let e = (2.0 * f - f.powi(2)).sqrt();
214        let c = r_eq / (1.0 - e.powi(2) * lat.sin().powi(2)).sqrt();
215        let alt = r_delta / lat.cos() - c;
216
217        Ok(LonLatAlt(
218            Angle::radians(lon),
219            Angle::radians(lat),
220            Distance::meters(alt),
221        ))
222    }
223
224    /// Returns the rotation matrix from body-fixed to topocentric (SEZ) frame.
225    pub fn rotation_to_topocentric(&self) -> DMat3 {
226        let rot1 = DMat3::from_rotation_z(self.lon().to_radians()).transpose();
227        let rot2 = DMat3::from_rotation_y(FRAC_PI_2 - self.lat().to_radians()).transpose();
228        rot2 * rot1
229    }
230}
231
232/// Error returned when constructing a [`LonLatAlt`] with invalid values.
233#[derive(Copy, Clone, Debug, Error, PartialEq)]
234pub enum LonLatAltError {
235    /// The longitude is outside the valid range [-180°, 180°].
236    #[error("longitude must be between -180 deg and 180 deg but was {0}")]
237    InvalidLongitude(Angle),
238    /// The latitude is outside the valid range [-90°, 90°].
239    #[error("latitude must between -90 deg and 90 deg but was {0}")]
240    InvalidLatitude(Angle),
241    /// The altitude is not a finite value.
242    #[error("invalid altitude {0}")]
243    InvalidAltitude(Distance),
244}
245
246/// Error returned when converting from body-fixed coordinates to geodetic.
247#[derive(Debug, Error)]
248pub enum FromBodyFixedError {
249    /// The position vector has zero length.
250    #[error("position vector has zero length")]
251    ZeroPosition,
252    /// The root finder failed to converge.
253    #[error(transparent)]
254    RootFinder(#[from] RootFinderError),
255}
256
257/// Builder for constructing validated [`LonLatAlt`] instances.
258#[derive(Copy, Clone, Debug)]
259pub struct LonLatAltBuilder {
260    longitude: Result<Angle, LonLatAltError>,
261    latitude: Result<Angle, LonLatAltError>,
262    altitude: Result<Distance, LonLatAltError>,
263}
264
265impl Default for LonLatAltBuilder {
266    fn default() -> Self {
267        Self::new()
268    }
269}
270
271impl LonLatAltBuilder {
272    /// Creates a new builder with default (zero) values.
273    pub fn new() -> Self {
274        Self {
275            longitude: Ok(Angle::default()),
276            latitude: Ok(Angle::default()),
277            altitude: Ok(Distance::default()),
278        }
279    }
280
281    /// Sets the longitude. Must be between -180° and 180°.
282    pub fn longitude(&mut self, longitude: Angle) -> &mut Self {
283        self.longitude = match longitude.to_degrees() {
284            lon if lon < -180.0 => Err(LonLatAltError::InvalidLongitude(longitude)),
285            lon if lon > 180.0 => Err(LonLatAltError::InvalidLongitude(longitude)),
286            _ => Ok(longitude),
287        };
288        self
289    }
290
291    /// Sets the latitude. Must be between -90° and 90°.
292    pub fn latitude(&mut self, latitude: Angle) -> &mut Self {
293        self.latitude = match latitude.to_degrees() {
294            lat if lat < -90.0 => Err(LonLatAltError::InvalidLatitude(latitude)),
295            lat if lat > 90.0 => Err(LonLatAltError::InvalidLatitude(latitude)),
296            _ => Ok(latitude),
297        };
298        self
299    }
300
301    /// Sets the altitude. Must be a finite value.
302    pub fn altitude(&mut self, altitude: Distance) -> &mut Self {
303        self.altitude = if !altitude.to_meters().is_finite() {
304            Err(LonLatAltError::InvalidAltitude(altitude))
305        } else {
306            Ok(altitude)
307        };
308        self
309    }
310
311    /// Builds the [`LonLatAlt`], returning an error if any value is invalid.
312    pub fn build(&self) -> Result<LonLatAlt, LonLatAltError> {
313        Ok(LonLatAlt(self.longitude?, self.latitude?, self.altitude?))
314    }
315}
316
317/// Cartesian state vector with position and velocity in meters and m/s.
318#[derive(Copy, Clone, Debug, Default, PartialEq, ApproxEq)]
319#[cfg_attr(feature = "serde", derive(serde::Serialize, serde::Deserialize))]
320pub struct Cartesian {
321    pos: DVec3,
322    vel: DVec3,
323}
324
325impl Cartesian {
326    /// Creates a new Cartesian state from individual position and velocity components.
327    pub const fn new(
328        x: Distance,
329        y: Distance,
330        z: Distance,
331        vx: Velocity,
332        vy: Velocity,
333        vz: Velocity,
334    ) -> Self {
335        Self {
336            pos: DVec3::new(x.to_meters(), y.to_meters(), z.to_meters()),
337            vel: DVec3::new(
338                vx.to_meters_per_second(),
339                vy.to_meters_per_second(),
340                vz.to_meters_per_second(),
341            ),
342        }
343    }
344
345    /// Creates a new Cartesian state from position and velocity vectors in meters and m/s.
346    #[inline]
347    pub const fn from_vecs(pos: DVec3, vel: DVec3) -> Self {
348        Self { pos, vel }
349    }
350
351    /// Creates a new Cartesian state from a `[x, y, z, vx, vy, vz]` array in meters and m/s.
352    pub const fn from_array([x, y, z, vx, vy, vz]: [f64; 6]) -> Self {
353        Self {
354            pos: DVec3::new(x, y, z),
355            vel: DVec3::new(vx, vy, vz),
356        }
357    }
358
359    /// Returns a new [`CartesianBuilder`].
360    pub const fn builder() -> CartesianBuilder {
361        CartesianBuilder::new()
362    }
363
364    /// Returns the position vector in meters.
365    #[inline]
366    pub fn position(&self) -> DVec3 {
367        self.pos
368    }
369
370    /// Sets the position vector in meters.
371    pub fn set_position(&mut self, position: DVec3) {
372        self.pos = position
373    }
374
375    /// Returns the velocity vector in m/s.
376    #[inline]
377    pub fn velocity(&self) -> DVec3 {
378        self.vel
379    }
380
381    /// Sets the velocity vector in m/s.
382    pub fn set_velocity(&mut self, velocity: DVec3) {
383        self.vel = velocity
384    }
385
386    /// Returns the x position component.
387    pub fn x(&self) -> Distance {
388        Distance::meters(self.pos.x)
389    }
390
391    /// Sets the `N`-th component (0–5 for x, y, z, vx, vy, vz) to `value`.
392    pub fn set<const N: usize>(&mut self, value: f64) {
393        const { assert!(N < 6, "index out of bounds") }
394
395        match N {
396            0 => {
397                self.pos.x = value;
398            }
399            1 => {
400                self.pos.y = value;
401            }
402            2 => {
403                self.pos.z = value;
404            }
405            3 => {
406                self.vel.x = value;
407            }
408            4 => {
409                self.vel.y = value;
410            }
411            5 => {
412                self.vel.z = value;
413            }
414            _ => unreachable!(),
415        }
416    }
417
418    /// Returns the y position component.
419    pub fn y(&self) -> Distance {
420        Distance::meters(self.pos.y)
421    }
422
423    /// Returns the z position component.
424    pub fn z(&self) -> Distance {
425        Distance::meters(self.pos.z)
426    }
427
428    /// Returns the x velocity component.
429    pub fn vx(&self) -> Velocity {
430        Velocity::meters_per_second(self.vel.x)
431    }
432
433    /// Returns the y velocity component.
434    pub fn vy(&self) -> Velocity {
435        Velocity::meters_per_second(self.vel.y)
436    }
437
438    /// Returns the z velocity component.
439    pub fn vz(&self) -> Velocity {
440        Velocity::meters_per_second(self.vel.z)
441    }
442}
443
444/// Builder for constructing [`Cartesian`] states from unitful components.
445#[derive(Debug, Default, Clone, Copy)]
446pub struct CartesianBuilder {
447    pos: Option<DVec3>,
448    vel: Option<DVec3>,
449}
450
451impl CartesianBuilder {
452    /// Creates a new builder with no position or velocity set.
453    pub const fn new() -> Self {
454        Self {
455            pos: None,
456            vel: None,
457        }
458    }
459
460    /// Sets the position components.
461    pub const fn position(&mut self, x: Distance, y: Distance, z: Distance) -> &mut Self {
462        self.pos = Some(DVec3::new(x.to_meters(), y.to_meters(), z.to_meters()));
463        self
464    }
465
466    /// Sets the velocity components.
467    pub const fn velocity(&mut self, vx: Velocity, vy: Velocity, vz: Velocity) -> &mut Self {
468        self.vel = Some(DVec3::new(
469            vx.to_meters_per_second(),
470            vy.to_meters_per_second(),
471            vz.to_meters_per_second(),
472        ));
473        self
474    }
475
476    /// Builds the [`Cartesian`] state. Unset components default to zero.
477    pub const fn build(&self) -> Cartesian {
478        Cartesian {
479            pos: match self.pos {
480                Some(pos) => pos,
481                None => DVec3::ZERO,
482            },
483            vel: match self.vel {
484                Some(vel) => vel,
485                None => DVec3::ZERO,
486            },
487        }
488    }
489}
490
491impl Add for Cartesian {
492    type Output = Self;
493
494    fn add(self, rhs: Self) -> Self::Output {
495        Self::from_vecs(self.pos + rhs.pos, self.vel + rhs.vel)
496    }
497}
498
499impl AddAssign for Cartesian {
500    fn add_assign(&mut self, rhs: Self) {
501        self.pos += rhs.pos;
502        self.vel += rhs.vel;
503    }
504}
505
506impl Sub for Cartesian {
507    type Output = Self;
508
509    fn sub(self, rhs: Self) -> Self::Output {
510        Self::from_vecs(self.pos - rhs.pos, self.vel - rhs.vel)
511    }
512}
513
514impl SubAssign for Cartesian {
515    fn sub_assign(&mut self, rhs: Self) {
516        self.pos -= rhs.pos;
517        self.vel -= rhs.vel;
518    }
519}
520
521impl Mul<f64> for Cartesian {
522    type Output = Self;
523
524    fn mul(self, rhs: f64) -> Self::Output {
525        Self {
526            pos: self.pos * rhs,
527            vel: self.vel * rhs,
528        }
529    }
530}
531
532impl MulAssign<f64> for Cartesian {
533    fn mul_assign(&mut self, rhs: f64) {
534        self.pos *= rhs;
535        self.vel *= rhs;
536    }
537}
538
539impl Div<f64> for Cartesian {
540    type Output = Self;
541
542    fn div(self, rhs: f64) -> Self::Output {
543        Self {
544            pos: self.pos / rhs,
545            vel: self.vel / rhs,
546        }
547    }
548}
549
550impl DivAssign<f64> for Cartesian {
551    fn div_assign(&mut self, rhs: f64) {
552        self.pos /= rhs;
553        self.vel /= rhs;
554    }
555}
556
557impl Neg for Cartesian {
558    type Output = Self;
559
560    fn neg(self) -> Self::Output {
561        Self::from_vecs(-self.pos, -self.vel)
562    }
563}
564
565/// Generic interpolated time series data with `N` components.
566#[derive(Debug, Clone)]
567pub struct TrajectoryData<const N: usize> {
568    epoch: TimeDelta,
569    time_steps: Arc<[f64]>,
570    data: [Arc<[f64]>; N],
571    series: [Series; N],
572}
573
574impl<const N: usize> TrajectoryData<N> {
575    /// Creates trajectory data from fixed-size arrays of time steps and component data.
576    pub fn from_arrays<const M: usize>(
577        epoch: TimeDelta,
578        time_steps: [TimeDelta; M],
579        data: &[[f64; M]; N],
580    ) -> Self {
581        let time_steps: Arc<[f64]> = Arc::from_iter(
582            time_steps
583                .into_iter()
584                .map(|t| (t - epoch).to_seconds().to_f64()),
585        );
586        let data: [Arc<[f64]>; N] = data.map(Arc::from);
587        let series = data.clone().map(|d| {
588            Series::new(
589                time_steps.clone(),
590                d.clone(),
591                InterpolationType::CubicSpline,
592            )
593        });
594        Self {
595            epoch,
596            time_steps,
597            data,
598            series,
599        }
600    }
601
602    /// Returns the time steps relative to the epoch in seconds.
603    pub fn time_steps(&self) -> Arc<[f64]> {
604        self.time_steps.clone()
605    }
606
607    /// Interpolates the `M`-th component at time `t` (seconds since epoch).
608    #[inline]
609    pub fn interpolate<const M: usize>(&self, t: f64) -> f64 {
610        const { assert!(M < N, "index is out-of-bounds") }
611
612        self.series[M].interpolate(t)
613    }
614
615    /// Interpolates all `N` components at time `t` (seconds since epoch).
616    #[inline]
617    pub fn interpolate_all(&self, t: f64) -> [f64; N] {
618        let idx = self.series[0].find_index(t);
619        self.series
620            .each_ref()
621            .map(|s| s.interpolate_at_index(t, idx))
622    }
623}
624
625/// A Cartesian state paired with a timestamp.
626#[cfg_attr(feature = "serde", derive(serde::Serialize, serde::Deserialize))]
627pub struct TimeStampedCartesian {
628    /// The timestamp of the state.
629    pub time: TimeDelta,
630    /// The Cartesian state vector.
631    pub state: Cartesian,
632}
633
634/// A 6-component trajectory (x, y, z, vx, vy, vz) with cubic spline interpolation.
635pub type CartesianTrajectory = TrajectoryData<6>;
636
637impl CartesianTrajectory {
638    /// Creates a trajectory from an iterator of timestamped Cartesian states.
639    pub fn from_states(states: impl IntoIterator<Item = TimeStampedCartesian>) -> Self {
640        let mut iter = states.into_iter().peekable();
641        let epoch = iter.peek().expect("should have at least two states").time;
642        let _ = iter.peek().expect("should have at least two states");
643        let (n, _) = iter.size_hint();
644
645        let mut time_steps: Vec<f64> = Vec::with_capacity(n);
646        let mut x: Vec<f64> = Vec::with_capacity(n);
647        let mut y: Vec<f64> = Vec::with_capacity(n);
648        let mut z: Vec<f64> = Vec::with_capacity(n);
649        let mut vx: Vec<f64> = Vec::with_capacity(n);
650        let mut vy: Vec<f64> = Vec::with_capacity(n);
651        let mut vz: Vec<f64> = Vec::with_capacity(n);
652
653        iter.for_each(|TimeStampedCartesian { time, state }| {
654            time_steps.push((time - epoch).to_seconds().to_f64());
655            x.push(state.x().as_f64());
656            y.push(state.y().as_f64());
657            z.push(state.z().as_f64());
658            vx.push(state.vx().as_f64());
659            vy.push(state.vy().as_f64());
660            vz.push(state.vz().as_f64());
661        });
662
663        let time_steps: Arc<[f64]> = Arc::from(time_steps);
664
665        let x: Arc<[f64]> = Arc::from(x);
666        let y: Arc<[f64]> = Arc::from(y);
667        let z: Arc<[f64]> = Arc::from(z);
668        let vx: Arc<[f64]> = Arc::from(vx);
669        let vy: Arc<[f64]> = Arc::from(vy);
670        let vz: Arc<[f64]> = Arc::from(vz);
671
672        let data = [
673            x.clone(),
674            y.clone(),
675            z.clone(),
676            vx.clone(),
677            vy.clone(),
678            vz.clone(),
679        ];
680
681        let series = data.clone().map(|d| {
682            Series::new(
683                time_steps.clone(),
684                d.clone(),
685                InterpolationType::CubicSpline,
686            )
687        });
688
689        Self {
690            epoch,
691            time_steps,
692            data,
693            series,
694        }
695    }
696
697    /// Returns the x position data points.
698    pub fn x(&self) -> Arc<[f64]> {
699        self.data[0].clone()
700    }
701
702    /// Returns the y position data points.
703    pub fn y(&self) -> Arc<[f64]> {
704        self.data[1].clone()
705    }
706
707    /// Returns the z position data points.
708    pub fn z(&self) -> Arc<[f64]> {
709        self.data[2].clone()
710    }
711
712    /// Returns the x velocity data points.
713    pub fn vx(&self) -> Arc<[f64]> {
714        self.data[3].clone()
715    }
716
717    /// Returns the y velocity data points.
718    pub fn vy(&self) -> Arc<[f64]> {
719        self.data[4].clone()
720    }
721
722    /// Returns the z velocity data points.
723    pub fn vz(&self) -> Arc<[f64]> {
724        self.data[5].clone()
725    }
726
727    /// Interpolates the x position at time `t` (seconds since epoch).
728    #[inline]
729    pub fn interpolate_x(&self, t: f64) -> f64 {
730        self.interpolate::<0>(t)
731    }
732
733    /// Interpolates the y position at time `t` (seconds since epoch).
734    #[inline]
735    pub fn interpolate_y(&self, t: f64) -> f64 {
736        self.interpolate::<1>(t)
737    }
738
739    /// Interpolates the z position at time `t` (seconds since epoch).
740    #[inline]
741    pub fn interpolate_z(&self, t: f64) -> f64 {
742        self.interpolate::<2>(t)
743    }
744
745    /// Interpolates the x velocity at time `t` (seconds since epoch).
746    #[inline]
747    pub fn interpolate_vx(&self, t: f64) -> f64 {
748        self.interpolate::<3>(t)
749    }
750
751    /// Interpolates the y velocity at time `t` (seconds since epoch).
752    #[inline]
753    pub fn interpolate_vy(&self, t: f64) -> f64 {
754        self.interpolate::<4>(t)
755    }
756
757    /// Interpolates the z velocity at time `t` (seconds since epoch).
758    #[inline]
759    pub fn interpolate_vz(&self, t: f64) -> f64 {
760        self.interpolate::<5>(t)
761    }
762
763    /// Interpolates the position vector at time `t` (seconds since epoch).
764    #[inline]
765    pub fn position(&self, t: f64) -> DVec3 {
766        let idx = self.series[0].find_index(t);
767        DVec3::new(
768            self.series[0].interpolate_at_index(t, idx),
769            self.series[1].interpolate_at_index(t, idx),
770            self.series[2].interpolate_at_index(t, idx),
771        )
772    }
773
774    /// Interpolates the velocity vector at time `t` (seconds since epoch).
775    #[inline]
776    pub fn velocity(&self, t: f64) -> DVec3 {
777        let idx = self.series[3].find_index(t);
778        DVec3::new(
779            self.series[3].interpolate_at_index(t, idx),
780            self.series[4].interpolate_at_index(t, idx),
781            self.series[5].interpolate_at_index(t, idx),
782        )
783    }
784
785    /// Interpolates the full Cartesian state at time `t` (seconds since epoch).
786    #[inline]
787    pub fn at(&self, t: f64) -> Cartesian {
788        let vals = self.interpolate_all(t);
789        Cartesian::from_array(vals)
790    }
791}
792
793/// Iterator over the discrete states in a [`CartesianTrajectory`].
794pub struct CartesianTrajectoryIterator {
795    data: CartesianTrajectory,
796    curr: usize,
797}
798
799impl CartesianTrajectoryIterator {
800    fn new(data: CartesianTrajectory) -> Self {
801        Self { data, curr: 0 }
802    }
803
804    fn len(&self) -> usize {
805        self.data.time_steps.len()
806    }
807
808    fn get_item(&self, idx: usize) -> Option<TimeStampedCartesian> {
809        let n = self.len();
810        if idx >= n {
811            return None;
812        }
813
814        let time = self.data.time_steps[idx];
815        let state = Cartesian::from_array(self.data.data.clone().map(|d| d[idx]));
816        Some(TimeStampedCartesian {
817            time: self.data.epoch + TimeDelta::from_seconds_f64(time),
818            state,
819        })
820    }
821}
822
823impl Iterator for CartesianTrajectoryIterator {
824    type Item = TimeStampedCartesian;
825
826    fn next(&mut self) -> Option<Self::Item> {
827        let item = self.get_item(self.curr);
828        self.curr += 1;
829        item
830    }
831}
832
833impl IntoIterator for CartesianTrajectory {
834    type Item = TimeStampedCartesian;
835
836    type IntoIter = CartesianTrajectoryIterator;
837
838    fn into_iter(self) -> Self::IntoIter {
839        Self::IntoIter::new(self)
840    }
841}
842
843impl FromIterator<TimeStampedCartesian> for CartesianTrajectory {
844    fn from_iter<T: IntoIterator<Item = TimeStampedCartesian>>(iter: T) -> Self {
845        TrajectoryData::from_states(iter)
846    }
847}
848
849#[cfg(test)]
850mod tests {
851    use rstest::rstest;
852
853    use crate::units::{AngleUnits, DistanceUnits, VelocityUnits};
854
855    use super::*;
856
857    #[test]
858    fn test_azel_api() {
859        let azel = AzEl::builder()
860            .azimuth(45.0.deg())
861            .elevation(45.0.deg())
862            .build()
863            .unwrap();
864        assert_eq!(azel.azimuth(), 45.0.deg());
865        assert_eq!(azel.elevation(), 45.0.deg());
866    }
867
868    #[rstest]
869    #[case(0.0.deg(), 0.0.deg(), Ok(AzEl(0.0.deg(), 0.0.deg())))]
870    #[case(-1.0.deg(), 0.0.deg(), Err(AzElError::InvalidAzimuth(-1.0.deg())))]
871    #[case(361.0.deg(), 0.0.deg(), Err(AzElError::InvalidAzimuth(361.0.deg())))]
872    #[case(0.0.deg(), -1.0.deg(), Err(AzElError::InvalidElevation(-1.0.deg())))]
873    #[case(0.0.deg(), 361.0.deg(), Err(AzElError::InvalidElevation(361.0.deg())))]
874    fn test_azel(#[case] az: Angle, #[case] el: Angle, #[case] exp: Result<AzEl, AzElError>) {
875        let act = AzEl::builder().azimuth(az).elevation(el).build();
876        assert_eq!(act, exp)
877    }
878
879    #[test]
880    fn test_lla_api() {
881        let lla = LonLatAlt::builder()
882            .longitude(45.0.deg())
883            .latitude(45.0.deg())
884            .altitude(100.0.m())
885            .build()
886            .unwrap();
887        assert_eq!(lla.lon(), 45.0.deg());
888        assert_eq!(lla.lat(), 45.0.deg());
889        assert_eq!(lla.alt(), 100.0.m());
890    }
891
892    #[rstest]
893    #[case(0.0.deg(), 0.0.deg(), 0.0.m(), Ok(LonLatAlt(0.0.deg(), 0.0.deg(), 0.0.m())))]
894    #[case(-181.0.deg(), 0.0.deg(), 0.0.m(), Err(LonLatAltError::InvalidLongitude(-181.0.deg())))]
895    #[case(181.0.deg(), 0.0.deg(), 0.0.m(), Err(LonLatAltError::InvalidLongitude(181.0.deg())))]
896    #[case(0.0.deg(), -91.0.deg(), 0.0.m(), Err(LonLatAltError::InvalidLatitude(-91.0.deg())))]
897    #[case(0.0.deg(), 91.0.deg(), 0.0.m(), Err(LonLatAltError::InvalidLatitude(91.0.deg())))]
898    #[case(0.0.deg(), 0.0.deg(), f64::INFINITY.m(), Err(LonLatAltError::InvalidAltitude(f64::INFINITY.m())))]
899    fn test_lla(
900        #[case] lon: Angle,
901        #[case] lat: Angle,
902        #[case] alt: Distance,
903        #[case] exp: Result<LonLatAlt, LonLatAltError>,
904    ) {
905        let act = LonLatAlt::builder()
906            .longitude(lon)
907            .latitude(lat)
908            .altitude(alt)
909            .build();
910        assert_eq!(act, exp)
911    }
912
913    #[test]
914    fn test_cartesian() {
915        let c = Cartesian::builder()
916            .position(1000.0.km(), 1000.0.km(), 1000.0.km())
917            .velocity(1.0.kps(), 1.0.kps(), 1.0.kps())
918            .build();
919        assert_eq!(c.position(), DVec3::new(1e6, 1e6, 1e6));
920        assert_eq!(c.velocity(), DVec3::new(1e3, 1e3, 1e3));
921        assert_eq!(c.x(), 1e6.m());
922        assert_eq!(c.y(), 1e6.m());
923        assert_eq!(c.z(), 1e6.m());
924        assert_eq!(c.vx(), 1e3.mps());
925        assert_eq!(c.vy(), 1e3.mps());
926        assert_eq!(c.vz(), 1e3.mps());
927    }
928
929    // Earth constants matching lox-bodies generated values
930    const EARTH_R_EQ: f64 = 6378136.6; // meters (6378.1366 km)
931    const EARTH_F: f64 = (6378.1366 - 6356.7519) / 6378.1366;
932
933    #[test]
934    fn test_lla_to_body_fixed() {
935        let coords = LonLatAlt::from_degrees(-4.3676, 40.4527, 0.0).unwrap();
936        let r_eq = Distance::meters(EARTH_R_EQ);
937        let result = coords.to_body_fixed(r_eq, EARTH_F);
938        let expected = DVec3::new(4846130.017870638, -370132.8551351891, 4116364.272747229);
939        assert!((result - expected).length() < 1e-3);
940    }
941
942    #[test]
943    fn test_lla_from_body_fixed_roundtrip() {
944        let coords = LonLatAlt::from_degrees(-4.3676, 40.4527, 100.0).unwrap();
945        let r_eq = Distance::meters(EARTH_R_EQ);
946        let body_fixed = coords.to_body_fixed(r_eq, EARTH_F);
947        let roundtrip = LonLatAlt::from_body_fixed(body_fixed, r_eq, EARTH_F).unwrap();
948        assert!((roundtrip.lon().to_degrees() - coords.lon().to_degrees()).abs() < 1e-6);
949        assert!((roundtrip.lat().to_degrees() - coords.lat().to_degrees()).abs() < 1e-6);
950        assert!((roundtrip.alt().to_meters() - coords.alt().to_meters()).abs() < 1e-3);
951    }
952
953    #[test]
954    fn test_lla_rotation_to_topocentric() {
955        let coords = LonLatAlt::from_degrees(-4.3676, 40.4527, 0.0).unwrap();
956        let act = coords.rotation_to_topocentric();
957        let exp = DMat3::from_cols(
958            DVec3::new(0.6469358921661584, 0.07615519584215287, 0.7587320591443464),
959            DVec3::new(
960                -0.049411020334552434,
961                0.9970959763965771,
962                -0.05794967578213965,
963            ),
964            DVec3::new(-0.7609418522440956, 0.0, 0.6488200809957448),
965        );
966        // Check element-wise within tolerance
967        for i in 0..3 {
968            assert!((act.col(i) - exp.col(i)).length() < 1e-10);
969        }
970    }
971
972    #[test]
973    fn test_from_body_fixed_north_pole() {
974        let r_eq = Distance::meters(EARTH_R_EQ);
975        let e = (2.0 * EARTH_F - EARTH_F.powi(2)).sqrt();
976        let r_polar = EARTH_R_EQ * (1.0 - e.powi(2)).sqrt();
977        let pos = DVec3::new(0.0, 0.0, r_polar);
978        let result = LonLatAlt::from_body_fixed(pos, r_eq, EARTH_F).unwrap();
979        assert!((result.lat().to_degrees() - 90.0).abs() < 1e-10);
980        assert!(result.alt().to_meters().abs() < 1e-3);
981    }
982
983    #[test]
984    fn test_from_body_fixed_south_pole() {
985        let r_eq = Distance::meters(EARTH_R_EQ);
986        let e = (2.0 * EARTH_F - EARTH_F.powi(2)).sqrt();
987        let r_polar = EARTH_R_EQ * (1.0 - e.powi(2)).sqrt();
988        let pos = DVec3::new(0.0, 0.0, -r_polar - 1000.0);
989        let result = LonLatAlt::from_body_fixed(pos, r_eq, EARTH_F).unwrap();
990        assert!((result.lat().to_degrees() + 90.0).abs() < 1e-10);
991        assert!((result.alt().to_meters() - 1000.0).abs() < 1e-3);
992    }
993
994    #[test]
995    fn test_from_body_fixed_zero_position() {
996        let r_eq = Distance::meters(EARTH_R_EQ);
997        let pos = DVec3::ZERO;
998        let result = LonLatAlt::from_body_fixed(pos, r_eq, EARTH_F);
999        assert!(matches!(result, Err(FromBodyFixedError::ZeroPosition)));
1000    }
1001}