unit_sphere/
lib.rs

1// Copyright (c) 2024-2025 Ken Barker
2
3// Permission is hereby granted, free of charge, to any person obtaining a copy
4// of this software and associated documentation files (the "Software"),
5// to deal in the Software without restriction, including without limitation the
6// rights to use, copy, modify, merge, publish, distribute, sublicense, and/or
7// sell copies of the Software, and to permit persons to whom the Software is
8// furnished to do so, subject to the following conditions:
9
10// The above copyright notice and this permission notice shall be included in
11// all copies or substantial portions of the Software.
12
13// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
14// IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
15// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
16// AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
17// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
18// OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
19// THE SOFTWARE.
20
21//! # unit-sphere
22//!
23//! [![crates.io](https://img.shields.io/crates/v/unit-sphere.svg)](https://crates.io/crates/unit-sphere)
24//! [![docs.io](https://docs.rs/unit-sphere/badge.svg)](https://docs.rs/unit-sphere/)
25//! [![License](https://img.shields.io/badge/License-MIT-blue)](https://opensource.org/license/mit/)
26//! [![Rust](https://github.com/kenba/unit-sphere-rs/actions/workflows/rust.yml/badge.svg)](https://github.com/kenba/unit-sphere-rs/actions)
27//! [![codecov](https://codecov.io/gh/kenba/unit-sphere-rs/graph/badge.svg?token=G1H1XINERW)](https://codecov.io/gh/kenba/unit-sphere-rs)
28//!
29//! A library for performing geometric calculations on the surface of a sphere.
30//!
31//! The library uses a combination of spherical trigonometry and vector geometry
32//! to perform [great-circle navigation](https://en.wikipedia.org/wiki/Great-circle_navigation)
33//! on the surface of a unit sphere, see *Figure 1*.
34//!
35//! ![great circle arc](https://via-technology.aero/img/navigation/sphere/great_circle_arc.svg)
36//!
37//! *Figure 1 A Great Circle Arc*
38//!
39//! A [great circle](https://en.wikipedia.org/wiki/Great_circle) is the
40//! shortest path between positions on the surface of a sphere.
41//! It is the spherical equivalent of a straight line in planar geometry.
42//!
43//! ## Spherical trigonometry
44//!
45//! A great circle path between positions may be found using
46//! [spherical trigonometry](https://en.wikipedia.org/wiki/Spherical_trigonometry).
47//!
48//! The [course](https://en.wikipedia.org/wiki/Great-circle_navigation#Course)
49//! (initial azimuth) of a great circle can be calculated from the
50//! latitudes and longitudes of the start and end points.
51//! While great circle distance can also be calculated from the latitudes and
52//! longitudes of the start and end points using the
53//! [haversine formula](https://en.wikipedia.org/wiki/Haversine_formula).
54//! The resulting distance in `Radians` can be converted to the required units by
55//! multiplying the distance by the Earth radius measured in the required units.
56//!
57//! ## Vector geometry
58//!
59//! Points on the surface of a sphere and great circle poles may be represented
60//! by 3D [vectors](https://www.movable-type.co.uk/scripts/latlong-vectors.html).
61//! Many calculations are simpler using vectors than spherical trigonometry.
62//!
63//! ![Spherical Vector Coordinates](https://via-technology.aero/img/navigation/sphere/ecef_coordinates.svg)
64//!
65//! *Figure 2 Spherical Vector Coordinates*
66//!
67//! For example, the across track distance of a point from a great circle can
68//! be calculated from the [dot product](https://en.wikipedia.org/wiki/Dot_product)
69//! of the point and the great circle pole vectors.
70//! While intersection points of great circles can simply be calculated from
71//! the [cross product](https://en.wikipedia.org/wiki/Cross_product) of their
72//! pole vectors.
73//!
74//! ## Design
75//!
76//! The `great_circle` module performs spherical trigonometric calculations
77//! and the `vector` module performs vector geometry calculations.
78//! See: [spherical vector geometry](https://via-technology.aero/navigation/spherical-vector-geometry/).
79//!
80//! The software uses types: `Angle`, `Degrees` and `Radians` from the
81//! [angle-sc](https://crates.io/crates/angle-sc) crate.
82//!
83//! The library is declared [no_std](https://docs.rust-embedded.org/book/intro/no-std.html)
84//! so it can be used in embedded applications.
85//!
86//! ## Example
87//!
88//! The following example calculates the intersection between two Great Circle `Arc`s
89//! it is taken from Charles Karney's original solution to
90//! [Intersection between two geodesic lines](https://sourceforge.net/p/geographiclib/discussion/1026621/thread/21aaff9f/#fe0a).
91//!
92//! ```rust
93//! use unit_sphere::{Arc, Degrees, LatLong, calculate_intersection_point};
94//! use angle_sc::is_within_tolerance;
95//!
96//! let istanbul = LatLong::new(Degrees(42.0), Degrees(29.0));
97//! let washington = LatLong::new(Degrees(39.0), Degrees(-77.0));
98//! let reyjavik = LatLong::new(Degrees(64.0), Degrees(-22.0));
99//! let accra = LatLong::new(Degrees(6.0), Degrees(0.0));
100//!
101//! let arc1 = Arc::try_from((&istanbul, &washington)).unwrap();
102//! let arc2 = Arc::try_from((&reyjavik, &accra)).unwrap();
103//!
104//! let intersection_point = calculate_intersection_point(&arc1, &arc2).unwrap();
105//! let lat_long = LatLong::from(&intersection_point);
106//! // Geodesic intersection latitude is 54.7170296089477
107//! assert!(is_within_tolerance(54.72, lat_long.lat().0, 0.05));
108//! // Geodesic intersection longitude is -14.56385574430775
109//! assert!(is_within_tolerance(-14.56, lat_long.lon().0, 0.02));
110//! ```
111
112#![cfg_attr(not(test), no_std)]
113
114extern crate angle_sc;
115extern crate nalgebra as na;
116
117pub mod great_circle;
118pub mod vector;
119
120use angle_sc::trig;
121pub use angle_sc::{Angle, Degrees, Radians, Validate};
122
123/// Test whether a latitude in degrees is a valid latitude.
124///
125/// I.e. whether it lies in the range: -90.0 <= degrees <= 90.0
126#[must_use]
127pub fn is_valid_latitude(degrees: f64) -> bool {
128    (-90.0..=90.0).contains(&degrees)
129}
130
131/// Test whether a longitude in degrees is a valid longitude.
132///
133/// I.e. whether it lies in the range: -180.0 <= degrees <= 180.0
134#[must_use]
135pub fn is_valid_longitude(degrees: f64) -> bool {
136    (-180.0..=180.0).contains(&degrees)
137}
138
139/// A position as a latitude and longitude pair of `Degrees`.
140#[derive(Clone, Copy, Debug, PartialEq)]
141pub struct LatLong {
142    lat: Degrees,
143    lon: Degrees,
144}
145
146impl Validate for LatLong {
147    /// Test whether a `LatLong` is valid.
148    ///
149    /// I.e. whether the latitude lies in the range: -90.0 <= lat <= 90.0
150    /// and the longitude lies in the range: -180.0 <= lon <= 180.0
151    #[must_use]
152    fn is_valid(&self) -> bool {
153        is_valid_latitude(self.lat.0) && is_valid_longitude(self.lon.0)
154    }
155}
156
157impl LatLong {
158    #[must_use]
159    pub const fn new(lat: Degrees, lon: Degrees) -> Self {
160        Self { lat, lon }
161    }
162
163    #[must_use]
164    pub const fn lat(&self) -> Degrees {
165        self.lat
166    }
167
168    #[must_use]
169    pub const fn lon(&self) -> Degrees {
170        self.lon
171    }
172}
173
174impl TryFrom<(f64, f64)> for LatLong {
175    type Error = &'static str;
176
177    /// Attempt to convert a pair of f64 values in latitude, longitude order.
178    ///
179    /// return a valid `LatLong`.
180    fn try_from(lat_long: (f64, f64)) -> Result<Self, Self::Error> {
181        if !is_valid_latitude(lat_long.0) {
182            Err("invalid latitude")
183        } else if !is_valid_longitude(lat_long.1) {
184            Err("invalid longitude")
185        } else {
186            Ok(Self::new(Degrees(lat_long.0), Degrees(lat_long.1)))
187        }
188    }
189}
190
191/// Calculate the azimuth and distance along the great circle of point b from
192/// point a.
193/// * `a`, `b` - the start and end positions
194///
195/// returns the great-circle azimuth relative to North and distance of point b
196/// from point a.
197#[must_use]
198pub fn calculate_azimuth_and_distance(a: &LatLong, b: &LatLong) -> (Angle, Radians) {
199    let a_lat = Angle::from(a.lat);
200    let b_lat = Angle::from(b.lat);
201    let delta_long = Angle::from((b.lon, a.lon));
202    (
203        great_circle::calculate_gc_azimuth(a_lat, b_lat, delta_long),
204        great_circle::calculate_gc_distance(a_lat, b_lat, delta_long),
205    )
206}
207
208/// Calculate the distance along the great circle of point b from point a.
209///
210/// See: [Haversine formula](https://en.wikipedia.org/wiki/Haversine_formula).
211/// This function is less accurate than `calculate_azimuth_and_distance`.
212/// * `a`, `b` - the start and end positions
213///
214/// returns the great-circle distance of point b from point a in `Radians`.
215#[must_use]
216pub fn haversine_distance(a: &LatLong, b: &LatLong) -> Radians {
217    let a_lat = Angle::from(a.lat);
218    let b_lat = Angle::from(b.lat);
219    let delta_lat = Angle::from((b.lat, a.lat));
220    let delta_long = Angle::from(b.lon - a.lon);
221    great_circle::calculate_haversine_distance(a_lat, b_lat, delta_long, delta_lat)
222}
223
224/// A `Vector3d` is a [nalgebra](https://crates.io/crates/nalgebra) `Vector3<f64>`.
225#[allow(clippy::module_name_repetitions)]
226pub type Vector3d = na::Vector3<f64>;
227
228impl From<&LatLong> for Vector3d {
229    /// Convert a `LatLong` to a point on the unit sphere.
230    ///
231    /// @pre |lat| <= 90.0 degrees.
232    /// * `lat` - the latitude.
233    /// * `lon` - the longitude.
234    ///
235    /// returns a `Vector3d` of the point on the unit sphere.
236    #[must_use]
237    fn from(a: &LatLong) -> Self {
238        vector::to_point(Angle::from(a.lat), Angle::from(a.lon))
239    }
240}
241
242/// Calculate the latitude of a Point.
243#[must_use]
244pub fn latitude(a: &Vector3d) -> Angle {
245    let sin_a = trig::UnitNegRange(a.z);
246    Angle::new(sin_a, trig::swap_sin_cos(sin_a))
247}
248
249/// Calculate the longitude of a Point.
250#[must_use]
251pub fn longitude(a: &Vector3d) -> Angle {
252    Angle::from_y_x(a.y, a.x)
253}
254
255impl From<&Vector3d> for LatLong {
256    /// Convert a point to a `LatLong`
257    #[must_use]
258    fn from(value: &Vector3d) -> Self {
259        Self::new(
260            Degrees::from(latitude(value)),
261            Degrees::from(longitude(value)),
262        )
263    }
264}
265
266/// An `Arc` of a Great Circle on a unit sphere.
267#[derive(Clone, Copy, Debug, PartialEq)]
268pub struct Arc {
269    /// The start point of the `Arc`.
270    a: Vector3d,
271    /// The right hand pole of the Great Circle of the `Arc`.
272    pole: Vector3d,
273    /// The length of the `Arc`.
274    length: Radians,
275    /// The half width of the `Arc`.
276    half_width: Radians,
277}
278
279impl Validate for Arc {
280    /// Test whether an `Arc` is valid.
281    ///
282    /// I.e. both a and pole are on the unit sphere and are orthogonal and
283    /// both length and `half_width` are not negative.
284    fn is_valid(&self) -> bool {
285        vector::is_unit(&self.a)
286            && vector::is_unit(&self.pole)
287            && vector::are_orthogonal(&self.a, &self.pole)
288            && !self.length.0.is_sign_negative()
289            && !self.half_width.0.is_sign_negative()
290    }
291}
292
293impl Arc {
294    /// Construct an `Arc`
295    ///
296    /// * `a` - the start point of the `Arc`.
297    /// * `pole` - the right hand pole of the Great Circle of the `Arc`.
298    /// * `length` - the length of the `Arc`.
299    /// * `half_width` - the half width of the `Arc`.
300    #[must_use]
301    pub const fn new(a: Vector3d, pole: Vector3d, length: Radians, half_width: Radians) -> Self {
302        Self {
303            a,
304            pole,
305            length,
306            half_width,
307        }
308    }
309
310    /// Construct an `Arc`
311    ///
312    /// * `a` - the start position
313    /// * `azimuth` - the azimuth at a.
314    /// * `length` - the length of the `Arc`.
315    #[must_use]
316    pub fn from_lat_lon_azi_length(a: &LatLong, azimuth: Angle, length: Radians) -> Self {
317        Self::new(
318            Vector3d::from(a),
319            vector::calculate_pole(Angle::from(a.lat()), Angle::from(a.lon()), azimuth),
320            length,
321            Radians(0.0),
322        )
323    }
324
325    /// Construct an `Arc` from the start and end positions.
326    ///
327    /// Note: if the points are the same or antipodal, the pole will be invalid.
328    /// * `a`, `b` - the start and end positions
329    #[must_use]
330    pub fn between_positions(a: &LatLong, b: &LatLong) -> Self {
331        let (azimuth, length) = calculate_azimuth_and_distance(a, b);
332        let a_lat = Angle::from(a.lat());
333        // if a is at the North or South pole
334        if a_lat.cos().0 < great_circle::MIN_VALUE {
335            // use b's longitude
336            Self::from_lat_lon_azi_length(&LatLong::new(a.lat(), b.lon()), azimuth, length)
337        } else {
338            Self::from_lat_lon_azi_length(a, azimuth, length)
339        }
340    }
341
342    /// Set the `half_width` of an `Arc`.
343    ///
344    /// * `half_width` - the half width of the `Arc`.
345    #[must_use]
346    pub const fn set_half_width(&mut self, half_width: Radians) -> &mut Self {
347        self.half_width = half_width;
348        self
349    }
350
351    /// The start point of the `Arc`.
352    #[must_use]
353    pub const fn a(&self) -> Vector3d {
354        self.a
355    }
356
357    /// The right hand pole of the Great Circle at the start point of the `Arc`.
358    #[must_use]
359    pub const fn pole(&self) -> Vector3d {
360        self.pole
361    }
362
363    /// The length of the `Arc`.
364    #[must_use]
365    pub const fn length(&self) -> Radians {
366        self.length
367    }
368
369    /// The half width of the `Arc`.
370    #[must_use]
371    pub const fn half_width(&self) -> Radians {
372        self.half_width
373    }
374
375    /// The azimuth at the start point.
376    #[must_use]
377    pub fn azimuth(&self) -> Angle {
378        vector::calculate_azimuth(&self.a, &self.pole)
379    }
380
381    /// The direction vector of the `Arc` at the start point.
382    #[must_use]
383    pub fn direction(&self) -> Vector3d {
384        vector::direction(&self.a, &self.pole)
385    }
386
387    /// A position vector at distance along the `Arc`.
388    #[must_use]
389    pub fn position(&self, distance: Radians) -> Vector3d {
390        vector::position(&self.a, &self.direction(), Angle::from(distance))
391    }
392
393    /// The end point of the `Arc`.
394    #[must_use]
395    pub fn b(&self) -> Vector3d {
396        self.position(self.length)
397    }
398
399    /// The mid point of the `Arc`.
400    #[must_use]
401    pub fn mid_point(&self) -> Vector3d {
402        self.position(Radians(0.5 * self.length.0))
403    }
404
405    /// The position of a perpendicular point at distance from the `Arc`.
406    ///
407    /// * `point` a point on the `Arc`'s great circle.
408    /// * `distance` the perpendicular distance from the `Arc`'s great circle.
409    ///
410    /// returns the point at perpendicular distance from point.
411    #[must_use]
412    pub fn perp_position(&self, point: &Vector3d, distance: Radians) -> Vector3d {
413        vector::position(point, &self.pole, Angle::from(distance))
414    }
415
416    /// The position of a point at angle from the `Arc` start, at `Arc` length.
417    ///
418    /// * `angle` the angle from the `Arc` start.
419    ///
420    /// returns the point at angle from the `Arc` start, at `Arc` length.
421    #[must_use]
422    pub fn angle_position(&self, angle: Angle) -> Vector3d {
423        vector::rotate_position(&self.a, &self.pole, angle, Angle::from(self.length))
424    }
425
426    /// The `Arc` at the end of an `Arc`, just the point if `half_width` is zero.
427    ///
428    /// @param `at_b` if true the `Arc` at b, else the `Arc` at a.
429    ///
430    /// @return the end `Arc` at a or b.
431    #[must_use]
432    pub fn end_arc(&self, at_b: bool) -> Self {
433        let p = if at_b { self.b() } else { self.a };
434        let pole = vector::direction(&p, &self.pole);
435        if self.half_width.0 < great_circle::MIN_VALUE {
436            Self::new(p, pole, Radians(0.0), Radians(0.0))
437        } else {
438            let a = self.perp_position(&p, self.half_width);
439            Self::new(a, pole, self.half_width + self.half_width, Radians(0.0))
440        }
441    }
442
443    /// Calculate great-circle along and across track distances of point from
444    /// the `Arc`.
445    ///
446    /// * `point` - the point.
447    ///
448    /// returns the along and across track distances of the point in Radians.
449    #[must_use]
450    pub fn calculate_atd_and_xtd(&self, point: &Vector3d) -> (Radians, Radians) {
451        vector::calculate_atd_and_xtd(&self.a, &self.pole(), point)
452    }
453}
454
455impl TryFrom<(&LatLong, &LatLong)> for Arc {
456    type Error = &'static str;
457
458    /// Construct an `Arc` from a pair of positions.
459    ///
460    /// * `params` - the start and end positions
461    fn try_from(params: (&LatLong, &LatLong)) -> Result<Self, Self::Error> {
462        // Convert positions to vectors
463        let a = Vector3d::from(params.0);
464        let b = Vector3d::from(params.1);
465        // Calculate the great circle pole
466        vector::normalise(&a.cross(&b), vector::MIN_SQ_NORM).map_or_else(
467            || {
468                if vector::sq_distance(&a, &b) < 1.0 {
469                    Err("Positions are too close")
470                } else {
471                    Err("Positions are antipodal")
472                }
473            },
474            |pole| {
475                Ok(Self::new(
476                    a,
477                    pole,
478                    great_circle::e2gc_distance(vector::distance(&a, &b)),
479                    Radians(0.0),
480                ))
481            },
482        )
483    }
484}
485
486/// Calculate the great-circle distances along a pair of `Arc`s to their
487/// closest intersection point or their coincident arc distances if the
488/// `Arc`s are on coincident Great Circles.
489///
490/// * `arc1`, `arc2` the `Arc`s.
491///
492/// returns the distances along the first `Arc` and second `Arc` to the intersection
493/// point or to their coincident arc distances if the `Arc`s do not intersect.
494#[must_use]
495pub fn calculate_intersection_distances(arc1: &Arc, arc2: &Arc) -> (Radians, Radians) {
496    vector::intersection::calculate_intersection_point_distances(
497        &arc1.a,
498        &arc1.pole,
499        arc1.length(),
500        &arc2.a,
501        &arc2.pole,
502        arc2.length(),
503        &(0.5 * (arc1.mid_point() + arc2.mid_point())),
504    )
505}
506
507/// Calculate whether a pair of `Arc`s intersect and (if so) where.
508///
509/// * `arc1`, `arc2` the `Arc`s.
510///
511/// returns the distance along the first `Arc` to the second `Arc` or None if they
512/// don't intersect.
513///
514/// # Examples
515/// ```
516/// use unit_sphere::{Arc, Degrees, LatLong, calculate_intersection_point};
517/// use angle_sc::is_within_tolerance;
518///
519/// let istanbul = LatLong::new(Degrees(42.0), Degrees(29.0));
520/// let washington = LatLong::new(Degrees(39.0), Degrees(-77.0));
521/// let reyjavik = LatLong::new(Degrees(64.0), Degrees(-22.0));
522/// let accra = LatLong::new(Degrees(6.0), Degrees(0.0));
523///
524/// let arc1 = Arc::try_from((&istanbul, &washington)).unwrap();
525/// let arc2 = Arc::try_from((&reyjavik, &accra)).unwrap();
526///
527/// // Calculate the intersection point position
528/// let intersection_point = calculate_intersection_point(&arc1, &arc2).unwrap();
529/// let lat_long = LatLong::from(&intersection_point);
530///
531/// // The expected latitude and longitude are from:
532/// // <https://sourceforge.net/p/geographiclib/discussion/1026621/thread/21aaff9f/#fe0a>
533///
534/// // Geodesic intersection latitude is 54.7170296089477
535/// assert!(is_within_tolerance(54.72, lat_long.lat().0, 0.05));
536/// // Geodesic intersection longitude is -14.56385574430775
537/// assert!(is_within_tolerance(-14.56, lat_long.lon().0, 0.02));
538/// ```
539#[must_use]
540pub fn calculate_intersection_point(arc1: &Arc, arc2: &Arc) -> Option<Vector3d> {
541    let (distance1, distance2) = calculate_intersection_distances(arc1, arc2);
542
543    // Determine whether both distances are within both arcs
544    if vector::intersection::is_within(distance1.0, arc1.length().0)
545        && vector::intersection::is_within(distance2.0, arc2.length().0)
546    {
547        Some(arc1.position(distance1))
548    } else {
549        None
550    }
551}
552
553#[cfg(test)]
554mod tests {
555    use super::*;
556    use angle_sc::{is_within_tolerance, Degrees};
557
558    #[test]
559    fn test_is_valid_latitude() {
560        // value < -90
561        assert!(!is_valid_latitude(-90.0001));
562        // value = -90
563        assert!(is_valid_latitude(-90.0));
564        // value = 90
565        assert!(is_valid_latitude(90.0));
566        // value > 90
567        assert!(!is_valid_latitude(90.0001));
568    }
569
570    #[test]
571    fn test_is_valid_longitude() {
572        // value < -180
573        assert!(!is_valid_longitude(-180.0001));
574        // value = -180
575        assert!(is_valid_longitude(-180.0));
576        // value = 180
577        assert!(is_valid_longitude(180.0));
578        // value > 180
579        assert!(!is_valid_longitude(180.0001));
580    }
581
582    #[test]
583    fn test_latlong_traits() {
584        let a = LatLong::try_from((0.0, 90.0)).unwrap();
585
586        assert!(a.is_valid());
587
588        let a_clone = a.clone();
589        assert!(a_clone == a);
590
591        assert_eq!(Degrees(0.0), a.lat());
592        assert_eq!(Degrees(90.0), a.lon());
593
594        print!("LatLong: {:?}", a);
595
596        let invalid_lat = LatLong::try_from((91.0, 0.0));
597        assert_eq!(Err("invalid latitude"), invalid_lat);
598
599        let invalid_lon = LatLong::try_from((0.0, 181.0));
600        assert_eq!(Err("invalid longitude"), invalid_lon);
601    }
602
603    #[test]
604    fn test_vector3d_traits() {
605        let a = LatLong::try_from((0.0, 90.0)).unwrap();
606        let point = Vector3d::from(&a);
607
608        assert_eq!(0.0, point.x);
609        assert_eq!(1.0, point.y);
610        assert_eq!(0.0, point.z);
611
612        assert_eq!(Degrees(0.0), Degrees::from(latitude(&point)));
613        assert_eq!(Degrees(90.0), Degrees::from(longitude(&point)));
614
615        let result = LatLong::from(&point);
616        assert_eq!(a, result);
617    }
618
619    #[test]
620    fn test_great_circle_90n_0n_0e() {
621        let a = LatLong::new(Degrees(90.0), Degrees(0.0));
622        let b = LatLong::new(Degrees(0.0), Degrees(0.0));
623        let (azimuth, dist) = calculate_azimuth_and_distance(&a, &b);
624
625        assert!(is_within_tolerance(
626            core::f64::consts::FRAC_PI_2,
627            dist.0,
628            f64::EPSILON
629        ));
630        assert_eq!(180.0, Degrees::from(azimuth).0);
631
632        let dist = haversine_distance(&a, &b);
633        assert!(is_within_tolerance(
634            core::f64::consts::FRAC_PI_2,
635            dist.0,
636            f64::EPSILON
637        ));
638    }
639
640    #[test]
641    fn test_great_circle_90s_0n_50e() {
642        let a = LatLong::new(Degrees(-90.0), Degrees(0.0));
643        let b = LatLong::new(Degrees(0.0), Degrees(50.0));
644        let (azimuth, dist) = calculate_azimuth_and_distance(&a, &b);
645
646        assert!(is_within_tolerance(
647            core::f64::consts::FRAC_PI_2,
648            dist.0,
649            f64::EPSILON
650        ));
651        assert_eq!(0.0, Degrees::from(azimuth).0);
652
653        let dist = haversine_distance(&a, &b);
654        assert!(is_within_tolerance(
655            core::f64::consts::FRAC_PI_2,
656            dist.0,
657            f64::EPSILON
658        ));
659    }
660
661    #[test]
662    fn test_great_circle_0n_60e_0n_60w() {
663        let a = LatLong::new(Degrees(0.0), Degrees(60.0));
664        let b = LatLong::new(Degrees(0.0), Degrees(-60.0));
665        let (azimuth, dist) = calculate_azimuth_and_distance(&a, &b);
666
667        assert!(is_within_tolerance(
668            2.0 * core::f64::consts::FRAC_PI_3,
669            dist.0,
670            2.0 * f64::EPSILON
671        ));
672        assert_eq!(-90.0, Degrees::from(azimuth).0);
673
674        let dist = haversine_distance(&a, &b);
675        assert!(is_within_tolerance(
676            2.0 * core::f64::consts::FRAC_PI_3,
677            dist.0,
678            2.0 * f64::EPSILON
679        ));
680    }
681
682    #[test]
683    fn test_arc() {
684        // Greenwich equator
685        let g_eq = LatLong::new(Degrees(0.0), Degrees(0.0));
686
687        // 90 degrees East on the equator
688        let e_eq = LatLong::new(Degrees(0.0), Degrees(90.0));
689
690        let mut arc = Arc::between_positions(&g_eq, &e_eq);
691        let arc = arc.set_half_width(Radians(0.01));
692        assert!(arc.is_valid());
693        assert_eq!(Radians(0.01), arc.half_width());
694
695        assert_eq!(Vector3d::from(&g_eq), arc.a());
696        assert_eq!(Vector3d::new(0.0, 0.0, 1.0), arc.pole());
697        assert!(is_within_tolerance(
698            core::f64::consts::FRAC_PI_2,
699            arc.length().0,
700            f64::EPSILON
701        ));
702        assert_eq!(Angle::from(Degrees(90.0)), arc.azimuth());
703        let b = Vector3d::from(&e_eq);
704        assert!(is_within_tolerance(
705            0.0,
706            vector::distance(&b, &arc.b()),
707            f64::EPSILON
708        ));
709
710        let mid_point = arc.mid_point();
711        assert_eq!(0.0, mid_point.z);
712        assert!(is_within_tolerance(
713            45.0,
714            Degrees::from(longitude(&mid_point)).0,
715            32.0 * f64::EPSILON
716        ));
717
718        let start_arc = arc.end_arc(false);
719        assert_eq!(0.02, start_arc.length().0);
720
721        let start_arc_a = start_arc.a();
722        assert_eq!(start_arc_a, arc.perp_position(&arc.a(), Radians(0.01)));
723
724        let angle_90 = Angle::from(Degrees(90.0));
725        let pole_0 = Vector3d::new(0.0, 0.0, 1.0);
726        assert!(vector::distance(&pole_0, &arc.angle_position(angle_90)) <= f64::EPSILON);
727
728        let end_arc = arc.end_arc(true);
729        assert_eq!(0.02, end_arc.length().0);
730
731        let end_arc_a = end_arc.a();
732        assert_eq!(end_arc_a, arc.perp_position(&arc.b(), Radians(0.01)));
733    }
734
735    #[test]
736    fn test_north_and_south_poles() {
737        let north_pole = LatLong::new(Degrees(90.0), Degrees(0.0));
738        let south_pole = LatLong::new(Degrees(-90.0), Degrees(0.0));
739
740        let (azimuth, distance) = calculate_azimuth_and_distance(&south_pole, &north_pole);
741        assert_eq!(0.0, Degrees::from(azimuth).0);
742        assert_eq!(core::f64::consts::PI, distance.0);
743
744        let (azimuth, distance) = calculate_azimuth_and_distance(&north_pole, &south_pole);
745        assert_eq!(180.0, Degrees::from(azimuth).0);
746        assert_eq!(core::f64::consts::PI, distance.0);
747
748        // 90 degrees East on the equator
749        let e_eq = LatLong::new(Degrees(0.0), Degrees(50.0));
750
751        let arc = Arc::between_positions(&north_pole, &e_eq);
752        assert!(is_within_tolerance(
753            e_eq.lat().0,
754            LatLong::from(&arc.b()).lat().abs().0,
755            1e-13
756        ));
757        assert!(is_within_tolerance(
758            e_eq.lon().0,
759            LatLong::from(&arc.b()).lon().0,
760            50.0 * f64::EPSILON
761        ));
762
763        let arc = Arc::between_positions(&south_pole, &e_eq);
764        assert!(is_within_tolerance(
765            e_eq.lat().0,
766            LatLong::from(&arc.b()).lat().abs().0,
767            1e-13
768        ));
769        assert!(is_within_tolerance(
770            e_eq.lon().0,
771            LatLong::from(&arc.b()).lon().0,
772            50.0 * f64::EPSILON
773        ));
774
775        let w_eq = LatLong::new(Degrees(0.0), Degrees(-140.0));
776
777        let arc = Arc::between_positions(&north_pole, &w_eq);
778        assert!(is_within_tolerance(
779            w_eq.lat().0,
780            LatLong::from(&arc.b()).lat().abs().0,
781            1e-13
782        ));
783        assert!(is_within_tolerance(
784            w_eq.lon().0,
785            LatLong::from(&arc.b()).lon().0,
786            256.0 * f64::EPSILON
787        ));
788
789        let arc = Arc::between_positions(&south_pole, &w_eq);
790        assert!(is_within_tolerance(
791            w_eq.lat().0,
792            LatLong::from(&arc.b()).lat().abs().0,
793            1e-13
794        ));
795        assert!(is_within_tolerance(
796            w_eq.lon().0,
797            LatLong::from(&arc.b()).lon().0,
798            256.0 * f64::EPSILON
799        ));
800
801        let invalid_arc = Arc::try_from((&north_pole, &north_pole));
802        assert_eq!(Err("Positions are too close"), invalid_arc);
803
804        let arc = Arc::between_positions(&north_pole, &north_pole);
805        assert_eq!(north_pole, LatLong::from(&arc.b()));
806
807        let invalid_arc = Arc::try_from((&north_pole, &south_pole));
808        assert_eq!(Err("Positions are antipodal"), invalid_arc);
809
810        let arc = Arc::between_positions(&north_pole, &south_pole);
811        assert_eq!(south_pole, LatLong::from(&arc.b()));
812
813        let arc = Arc::between_positions(&south_pole, &north_pole);
814        assert_eq!(north_pole, LatLong::from(&arc.b()));
815
816        let arc = Arc::between_positions(&south_pole, &south_pole);
817        assert_eq!(south_pole, LatLong::from(&arc.b()));
818    }
819
820    #[test]
821    fn test_arc_atd_and_xtd() {
822        // Greenwich equator
823        let g_eq = LatLong::new(Degrees(0.0), Degrees(0.0));
824
825        // 90 degrees East on the equator
826        let e_eq = LatLong::new(Degrees(0.0), Degrees(90.0));
827
828        let arc = Arc::try_from((&g_eq, &e_eq)).unwrap();
829        assert!(arc.is_valid());
830
831        let start_arc = arc.end_arc(false);
832        assert_eq!(0.0, start_arc.length().0);
833
834        let start_arc_a = start_arc.a();
835        assert_eq!(arc.a(), start_arc_a);
836
837        let longitude = Degrees(1.0);
838
839        // Test across track distance
840        // Accuracy drops off outside of this range
841        for lat in -83..84 {
842            let latitude = Degrees(lat as f64);
843            let latlong = LatLong::new(latitude, longitude);
844            let point = Vector3d::from(&latlong);
845
846            let expected = (lat as f64).to_radians();
847            let (atd, xtd) = arc.calculate_atd_and_xtd(&point);
848            assert!(is_within_tolerance(1_f64.to_radians(), atd.0, f64::EPSILON));
849            assert!(is_within_tolerance(expected, xtd.0, 2.0 * f64::EPSILON));
850        }
851    }
852
853    #[test]
854    fn test_arc_intersection_point() {
855        // Karney's example:
856        // Istanbul, Washington, Reyjavik and Accra
857        // from: <https://sourceforge.net/p/geographiclib/discussion/1026621/thread/21aaff9f/#fe0a>
858        let istanbul = LatLong::new(Degrees(42.0), Degrees(29.0));
859        let washington = LatLong::new(Degrees(39.0), Degrees(-77.0));
860        let reyjavik = LatLong::new(Degrees(64.0), Degrees(-22.0));
861        let accra = LatLong::new(Degrees(6.0), Degrees(0.0));
862
863        let arc1 = Arc::try_from((&istanbul, &washington)).unwrap();
864        let arc2 = Arc::try_from((&reyjavik, &accra)).unwrap();
865
866        let intersection_point = calculate_intersection_point(&arc1, &arc2).unwrap();
867        let lat_long = LatLong::from(&intersection_point);
868        // Geodesic intersection latitude is 54.7170296089477
869        assert!(is_within_tolerance(54.72, lat_long.lat().0, 0.05));
870        // Geodesic intersection longitude is -14.56385574430775
871        assert!(is_within_tolerance(-14.56, lat_long.lon().0, 0.02));
872
873        // Switch arcs
874        let intersection_point = calculate_intersection_point(&arc2, &arc1).unwrap();
875        let lat_long = LatLong::from(&intersection_point);
876        // Geodesic intersection latitude is 54.7170296089477
877        assert!(is_within_tolerance(54.72, lat_long.lat().0, 0.05));
878        // Geodesic intersection longitude is -14.56385574430775
879        assert!(is_within_tolerance(-14.56, lat_long.lon().0, 0.02));
880    }
881
882    #[test]
883    fn test_arc_intersection_same_great_circles() {
884        let south_pole_1 = LatLong::new(Degrees(-88.0), Degrees(-180.0));
885        let south_pole_2 = LatLong::new(Degrees(-87.0), Degrees(0.0));
886
887        let arc1 = Arc::try_from((&south_pole_1, &south_pole_2)).unwrap();
888
889        let intersection_lengths = calculate_intersection_distances(&arc1, &arc1);
890        assert_eq!(Radians(0.0), intersection_lengths.0);
891        assert_eq!(Radians(0.0), intersection_lengths.1);
892
893        let intersection_point = calculate_intersection_point(&arc1, &arc1).unwrap();
894        assert_eq!(0.0, vector::sq_distance(&arc1.a(), &intersection_point));
895
896        let south_pole_3 = LatLong::new(Degrees(-85.0), Degrees(0.0));
897        let south_pole_4 = LatLong::new(Degrees(-86.0), Degrees(0.0));
898        let arc2 = Arc::try_from((&south_pole_3, &south_pole_4)).unwrap();
899        let intersection_point = calculate_intersection_point(&arc1, &arc2);
900        assert!(intersection_point.is_none());
901    }
902}