geographiclib_rs/
geodesic.rs

1#![allow(non_snake_case)]
2#![allow(clippy::excessive_precision)]
3
4use crate::geodesic_capability as caps;
5use crate::geodesic_line;
6use crate::geomath;
7use std::sync;
8
9use std::f64::consts::{FRAC_1_SQRT_2, PI};
10
11pub const WGS84_A: f64 = 6378137.0;
12// Evaluating this as 1000000000.0 / (298257223563f64) reduces the
13// round-off error by about 10%.  However, expressing the flattening as
14// 1/298.257223563 is well ingrained.
15pub const WGS84_F: f64 = 1.0 / ((298257223563f64) / 1000000000.0);
16
17#[derive(Copy, Clone, PartialEq, PartialOrd, Debug)]
18pub struct Geodesic {
19    pub a: f64,
20    pub f: f64,
21    pub _f1: f64,
22    pub _e2: f64,
23    pub _ep2: f64,
24    _n: f64,
25    pub _b: f64,
26    pub _c2: f64,
27    _etol2: f64,
28    _A3x: [f64; GEODESIC_ORDER],
29    _C3x: [f64; _nC3x_],
30    _C4x: [f64; _nC4x_],
31
32    pub GEODESIC_ORDER: usize,
33    _nC3x_: usize,
34    _nC4x_: usize,
35    maxit1_: u64,
36    maxit2_: u64,
37
38    pub tiny_: f64,
39    tol0_: f64,
40    tol1_: f64,
41    _tol2_: f64,
42    tolb_: f64,
43    xthresh_: f64,
44}
45
46static WGS84_GEOD: sync::OnceLock<Geodesic> = sync::OnceLock::new();
47
48impl Geodesic {
49    pub fn wgs84() -> Self {
50        *WGS84_GEOD.get_or_init(|| Geodesic::new(WGS84_A, WGS84_F))
51    }
52
53    pub fn equatorial_radius(&self) -> f64 {
54        self.a
55    }
56
57    pub fn flattening(&self) -> f64 {
58        self.f
59    }
60}
61
62const COEFF_A3: [f64; 18] = [
63    -3.0, 128.0, -2.0, -3.0, 64.0, -1.0, -3.0, -1.0, 16.0, 3.0, -1.0, -2.0, 8.0, 1.0, -1.0, 2.0,
64    1.0, 1.0,
65];
66
67const COEFF_C3: [f64; 45] = [
68    3.0, 128.0, 2.0, 5.0, 128.0, -1.0, 3.0, 3.0, 64.0, -1.0, 0.0, 1.0, 8.0, -1.0, 1.0, 4.0, 5.0,
69    256.0, 1.0, 3.0, 128.0, -3.0, -2.0, 3.0, 64.0, 1.0, -3.0, 2.0, 32.0, 7.0, 512.0, -10.0, 9.0,
70    384.0, 5.0, -9.0, 5.0, 192.0, 7.0, 512.0, -14.0, 7.0, 512.0, 21.0, 2560.0,
71];
72
73const COEFF_C4: [f64; 77] = [
74    97.0, 15015.0, 1088.0, 156.0, 45045.0, -224.0, -4784.0, 1573.0, 45045.0, -10656.0, 14144.0,
75    -4576.0, -858.0, 45045.0, 64.0, 624.0, -4576.0, 6864.0, -3003.0, 15015.0, 100.0, 208.0, 572.0,
76    3432.0, -12012.0, 30030.0, 45045.0, 1.0, 9009.0, -2944.0, 468.0, 135135.0, 5792.0, 1040.0,
77    -1287.0, 135135.0, 5952.0, -11648.0, 9152.0, -2574.0, 135135.0, -64.0, -624.0, 4576.0, -6864.0,
78    3003.0, 135135.0, 8.0, 10725.0, 1856.0, -936.0, 225225.0, -8448.0, 4992.0, -1144.0, 225225.0,
79    -1440.0, 4160.0, -4576.0, 1716.0, 225225.0, -136.0, 63063.0, 1024.0, -208.0, 105105.0, 3584.0,
80    -3328.0, 1144.0, 315315.0, -128.0, 135135.0, -2560.0, 832.0, 405405.0, 128.0, 99099.0,
81];
82
83pub const GEODESIC_ORDER: usize = 6;
84#[allow(non_upper_case_globals)]
85const _nC3x_: usize = 15;
86#[allow(non_upper_case_globals)]
87const _nC4x_: usize = 21;
88
89impl Geodesic {
90    pub fn new(a: f64, f: f64) -> Self {
91        let maxit1_ = 20;
92        let maxit2_ = maxit1_ + geomath::DIGITS + 10;
93        let tiny_ = geomath::get_min_val().sqrt();
94        let tol0_ = geomath::get_epsilon();
95        let tol1_ = 200.0 * tol0_;
96        let _tol2_ = tol0_.sqrt();
97        let tolb_ = tol0_ * _tol2_;
98        let xthresh_ = 1000.0 * _tol2_;
99
100        let _f1 = 1.0 - f;
101        let _e2 = f * (2.0 - f);
102        let _ep2 = _e2 / geomath::sq(_f1);
103        let _n = f / (2.0 - f);
104        let _b = a * _f1;
105        let _c2 = (geomath::sq(a)
106            + geomath::sq(_b)
107                * (if _e2 == 0.0 {
108                    1.0
109                } else {
110                    geomath::eatanhe(1.0, (if f < 0.0 { -1.0 } else { 1.0 }) * _e2.abs().sqrt())
111                        / _e2
112                }))
113            / 2.0;
114        let _etol2 = 0.1 * _tol2_ / (f.abs().max(0.001) * (1.0 - f / 2.0).min(1.0) / 2.0).sqrt();
115
116        let mut _A3x: [f64; GEODESIC_ORDER] = [0.0; GEODESIC_ORDER];
117        let mut _C3x: [f64; _nC3x_] = [0.0; _nC3x_];
118        let mut _C4x: [f64; _nC4x_] = [0.0; _nC4x_];
119
120        // Call a3coeff
121        let mut o: usize = 0;
122        for (k, j) in (0..GEODESIC_ORDER).rev().enumerate() {
123            let m = j.min(GEODESIC_ORDER - j - 1);
124            _A3x[k] = geomath::polyval(m, &COEFF_A3[o..], _n) / COEFF_A3[o + m + 1];
125            o += m + 2;
126        }
127
128        // c3coeff
129        let mut o = 0;
130        let mut k = 0;
131        for l in 1..GEODESIC_ORDER {
132            for j in (l..GEODESIC_ORDER).rev() {
133                let m = j.min(GEODESIC_ORDER - j - 1);
134                _C3x[k] = geomath::polyval(m, &COEFF_C3[o..], _n) / COEFF_C3[o + m + 1];
135                k += 1;
136                o += m + 2;
137            }
138        }
139
140        // c4coeff
141        let mut o = 0;
142        let mut k = 0;
143        for l in 0..GEODESIC_ORDER {
144            for j in (l..GEODESIC_ORDER).rev() {
145                let m = GEODESIC_ORDER - j - 1;
146                _C4x[k] = geomath::polyval(m, &COEFF_C4[o..], _n) / COEFF_C4[o + m + 1];
147                k += 1;
148                o += m + 2;
149            }
150        }
151
152        Geodesic {
153            a,
154            f,
155            _f1,
156            _e2,
157            _ep2,
158            _n,
159            _b,
160            _c2,
161            _etol2,
162            _A3x,
163            _C3x,
164            _C4x,
165
166            GEODESIC_ORDER,
167            _nC3x_,
168            _nC4x_,
169            maxit1_,
170            maxit2_,
171
172            tiny_,
173            tol0_,
174            tol1_,
175            _tol2_,
176            tolb_,
177            xthresh_,
178        }
179    }
180
181    pub fn _A3f(&self, eps: f64) -> f64 {
182        geomath::polyval(self.GEODESIC_ORDER - 1, &self._A3x, eps)
183    }
184
185    pub fn _C3f(&self, eps: f64, c: &mut [f64]) {
186        let mut mult = 1.0;
187        let mut o = 0;
188        // Clippy wants us to turn this into `c.iter_mut().enumerate().take(geodesic_order + 1).skip(1)`
189        // but benching (rust-1.75) shows that it would be slower.
190        #[allow(clippy::needless_range_loop)]
191        for l in 1..GEODESIC_ORDER {
192            let m = GEODESIC_ORDER - l - 1;
193            mult *= eps;
194            c[l] = mult * geomath::polyval(m, &self._C3x[o..], eps);
195            o += m + 1;
196        }
197    }
198
199    pub fn _C4f(&self, eps: f64, c: &mut [f64]) {
200        let mut mult = 1.0;
201        let mut o = 0;
202        // Clippy wants us to turn this into `c.iter_mut().enumerate().take(geodesic_order + 1).skip(1)`
203        // but benching (rust-1.75) shows that it would be slower.
204        #[allow(clippy::needless_range_loop)]
205        for l in 0..GEODESIC_ORDER {
206            let m = GEODESIC_ORDER - l - 1;
207            c[l] = mult * geomath::polyval(m, &self._C4x[o..], eps);
208            o += m + 1;
209            mult *= eps;
210        }
211    }
212
213    #[allow(clippy::too_many_arguments)]
214    pub fn _Lengths(
215        &self,
216        eps: f64,
217        sig12: f64,
218        ssig1: f64,
219        csig1: f64,
220        dn1: f64,
221        ssig2: f64,
222        csig2: f64,
223        dn2: f64,
224        cbet1: f64,
225        cbet2: f64,
226        outmask: u64,
227        C1a: &mut [f64],
228        C2a: &mut [f64],
229    ) -> (f64, f64, f64, f64, f64) {
230        let outmask = outmask & caps::OUT_MASK;
231        let mut s12b = f64::NAN;
232        let mut m12b = f64::NAN;
233        let mut m0 = f64::NAN;
234        let mut M12 = f64::NAN;
235        let mut M21 = f64::NAN;
236
237        let mut A1 = 0.0;
238        let mut A2 = 0.0;
239        let mut m0x = 0.0;
240        let mut J12 = 0.0;
241
242        if outmask & (caps::DISTANCE | caps::REDUCEDLENGTH | caps::GEODESICSCALE) != 0 {
243            A1 = geomath::_A1m1f(eps, self.GEODESIC_ORDER);
244            geomath::_C1f(eps, C1a, self.GEODESIC_ORDER);
245            if outmask & (caps::REDUCEDLENGTH | caps::GEODESICSCALE) != 0 {
246                A2 = geomath::_A2m1f(eps, self.GEODESIC_ORDER);
247                geomath::_C2f(eps, C2a, self.GEODESIC_ORDER);
248                m0x = A1 - A2;
249                A2 += 1.0;
250            }
251            A1 += 1.0;
252        }
253        if outmask & caps::DISTANCE != 0 {
254            let B1 = geomath::sin_cos_series(true, ssig2, csig2, C1a)
255                - geomath::sin_cos_series(true, ssig1, csig1, C1a);
256            s12b = A1 * (sig12 + B1);
257            if outmask & (caps::REDUCEDLENGTH | caps::GEODESICSCALE) != 0 {
258                let B2 = geomath::sin_cos_series(true, ssig2, csig2, C2a)
259                    - geomath::sin_cos_series(true, ssig1, csig1, C2a);
260                J12 = m0x * sig12 + (A1 * B1 - A2 * B2);
261            }
262        } else if outmask & (caps::REDUCEDLENGTH | caps::GEODESICSCALE) != 0 {
263            for l in 1..=self.GEODESIC_ORDER {
264                C2a[l] = A1 * C1a[l] - A2 * C2a[l];
265            }
266            J12 = m0x * sig12
267                + (geomath::sin_cos_series(true, ssig2, csig2, C2a)
268                    - geomath::sin_cos_series(true, ssig1, csig1, C2a));
269        }
270        if outmask & caps::REDUCEDLENGTH != 0 {
271            m0 = m0x;
272            // J12 is wrong
273            m12b = dn2 * (csig1 * ssig2) - dn1 * (ssig1 * csig2) - csig1 * csig2 * J12;
274        }
275        if outmask & caps::GEODESICSCALE != 0 {
276            let csig12 = csig1 * csig2 + ssig1 * ssig2;
277            let t = self._ep2 * (cbet1 - cbet2) * (cbet1 + cbet2) / (dn1 + dn2);
278            M12 = csig12 + (t * ssig2 - csig2 * J12) * ssig1 / dn1;
279            M21 = csig12 - (t * ssig1 - csig1 * J12) * ssig2 / dn2;
280        }
281        (s12b, m12b, m0, M12, M21)
282    }
283
284    #[allow(clippy::too_many_arguments)]
285    pub fn _InverseStart(
286        &self,
287        sbet1: f64,
288        cbet1: f64,
289        dn1: f64,
290        sbet2: f64,
291        cbet2: f64,
292        dn2: f64,
293        lam12: f64,
294        slam12: f64,
295        clam12: f64,
296        C1a: &mut [f64],
297        C2a: &mut [f64],
298    ) -> (f64, f64, f64, f64, f64, f64) {
299        let mut sig12 = -1.0;
300        let mut salp2 = f64::NAN;
301        let mut calp2 = f64::NAN;
302        let mut dnm = f64::NAN;
303
304        let mut somg12: f64;
305        let mut comg12: f64;
306
307        let sbet12 = sbet2 * cbet1 - cbet2 * sbet1;
308        let cbet12 = cbet2 * cbet1 + sbet2 * sbet1;
309
310        let mut sbet12a = sbet2 * cbet1;
311        sbet12a += cbet2 * sbet1;
312
313        let shortline = cbet12 >= 0.0 && sbet12 < 0.5 && cbet2 * lam12 < 0.5;
314        if shortline {
315            let mut sbetm2 = geomath::sq(sbet1 + sbet2);
316            sbetm2 /= sbetm2 + geomath::sq(cbet1 + cbet2);
317            dnm = (1.0 + self._ep2 * sbetm2).sqrt();
318            let omg12 = lam12 / (self._f1 * dnm);
319            somg12 = omg12.sin();
320            comg12 = omg12.cos();
321        } else {
322            somg12 = slam12;
323            comg12 = clam12;
324        }
325
326        let mut salp1 = cbet2 * somg12;
327
328        let mut calp1 = if comg12 >= 0.0 {
329            sbet12 + cbet2 * sbet1 * geomath::sq(somg12) / (1.0 + comg12)
330        } else {
331            sbet12a - cbet2 * sbet1 * geomath::sq(somg12) / (1.0 - comg12)
332        };
333
334        let ssig12 = salp1.hypot(calp1);
335        let csig12 = sbet1 * sbet2 + cbet1 * cbet2 * comg12;
336
337        if shortline && ssig12 < self._etol2 {
338            salp2 = cbet1 * somg12;
339            calp2 = sbet12
340                - cbet1
341                    * sbet2
342                    * (if comg12 >= 0.0 {
343                        geomath::sq(somg12) / (1.0 + comg12)
344                    } else {
345                        1.0 - comg12
346                    });
347            geomath::norm(&mut salp2, &mut calp2);
348            sig12 = ssig12.atan2(csig12);
349        } else if self._n.abs() > 0.1
350            || csig12 >= 0.0
351            || ssig12 >= 6.0 * self._n.abs() * PI * geomath::sq(cbet1)
352        {
353        } else {
354            let x: f64;
355            let y: f64;
356            let betscale: f64;
357            let lamscale: f64;
358            let lam12x = (-slam12).atan2(-clam12);
359            if self.f >= 0.0 {
360                let k2 = geomath::sq(sbet1) * self._ep2;
361                let eps = k2 / (2.0 * (1.0 + (1.0 + k2).sqrt()) + k2);
362                lamscale = self.f * cbet1 * self._A3f(eps) * PI;
363                betscale = lamscale * cbet1;
364                x = lam12x / lamscale;
365                y = sbet12a / betscale;
366            } else {
367                let cbet12a = cbet2 * cbet1 - sbet2 * sbet1;
368                let bet12a = sbet12a.atan2(cbet12a);
369                let (_, m12b, m0, _, _) = self._Lengths(
370                    self._n,
371                    PI + bet12a,
372                    sbet1,
373                    -cbet1,
374                    dn1,
375                    sbet2,
376                    cbet2,
377                    dn2,
378                    cbet1,
379                    cbet2,
380                    caps::REDUCEDLENGTH,
381                    C1a,
382                    C2a,
383                );
384                x = -1.0 + m12b / (cbet1 * cbet2 * m0 * PI);
385                betscale = if x < -0.01 {
386                    sbet12a / x
387                } else {
388                    -self.f * geomath::sq(cbet1) * PI
389                };
390                lamscale = betscale / cbet1;
391                y = lam12x / lamscale;
392            }
393            if y > -self.tol1_ && x > -1.0 - self.xthresh_ {
394                if self.f >= 0.0 {
395                    salp1 = (-x).min(1.0);
396                    calp1 = -(1.0 - geomath::sq(salp1)).sqrt()
397                } else {
398                    calp1 = x.max(if x > -self.tol1_ { 0.0 } else { -1.0 });
399                    salp1 = (1.0 - geomath::sq(calp1)).sqrt();
400                }
401            } else {
402                let k = geomath::astroid(x, y);
403                let omg12a = lamscale
404                    * if self.f >= 0.0 {
405                        -x * k / (1.0 + k)
406                    } else {
407                        -y * (1.0 + k) / k
408                    };
409                somg12 = omg12a.sin();
410                comg12 = -(omg12a.cos());
411                salp1 = cbet2 * somg12;
412                calp1 = sbet12a - cbet2 * sbet1 * geomath::sq(somg12) / (1.0 - comg12);
413            }
414        }
415
416        if salp1 > 0.0 || salp1.is_nan() {
417            geomath::norm(&mut salp1, &mut calp1);
418        } else {
419            salp1 = 1.0;
420            calp1 = 0.0;
421        };
422        (sig12, salp1, calp1, salp2, calp2, dnm)
423    }
424
425    #[allow(clippy::too_many_arguments)]
426    pub fn _Lambda12(
427        &self,
428        sbet1: f64,
429        cbet1: f64,
430        dn1: f64,
431        sbet2: f64,
432        cbet2: f64,
433        dn2: f64,
434        salp1: f64,
435        mut calp1: f64,
436        slam120: f64,
437        clam120: f64,
438        diffp: bool,
439        C1a: &mut [f64],
440        C2a: &mut [f64],
441        C3a: &mut [f64],
442    ) -> (f64, f64, f64, f64, f64, f64, f64, f64, f64, f64, f64) {
443        if sbet1 == 0.0 && calp1 == 0.0 {
444            calp1 = -self.tiny_;
445        }
446        let salp0 = salp1 * cbet1;
447        let calp0 = calp1.hypot(salp1 * sbet1);
448
449        let mut ssig1 = sbet1;
450        let somg1 = salp0 * sbet1;
451        let mut csig1 = calp1 * cbet1;
452        let comg1 = calp1 * cbet1;
453        geomath::norm(&mut ssig1, &mut csig1);
454
455        let salp2 = if cbet2 != cbet1 { salp0 / cbet2 } else { salp1 };
456        let calp2 = if cbet2 != cbet1 || sbet2.abs() != -sbet1 {
457            (geomath::sq(calp1 * cbet1)
458                + if cbet1 < -sbet1 {
459                    (cbet2 - cbet1) * (cbet1 + cbet2)
460                } else {
461                    (sbet1 - sbet2) * (sbet1 + sbet2)
462                })
463            .sqrt()
464                / cbet2
465        } else {
466            calp1.abs()
467        };
468        let mut ssig2 = sbet2;
469        let somg2 = salp0 * sbet2;
470        let mut csig2 = calp2 * cbet2;
471        let comg2 = calp2 * cbet2;
472        geomath::norm(&mut ssig2, &mut csig2);
473
474        let sig12 = ((csig1 * ssig2 - ssig1 * csig2).max(0.0)).atan2(csig1 * csig2 + ssig1 * ssig2);
475        let somg12 = (comg1 * somg2 - somg1 * comg2).max(0.0);
476        let comg12 = comg1 * comg2 + somg1 * somg2;
477        let eta = (somg12 * clam120 - comg12 * slam120).atan2(comg12 * clam120 + somg12 * slam120);
478
479        let k2 = geomath::sq(calp0) * self._ep2;
480        let eps = k2 / (2.0 * (1.0 + (1.0 + k2).sqrt()) + k2);
481        self._C3f(eps, C3a);
482        let B312 = geomath::sin_cos_series(true, ssig2, csig2, C3a)
483            - geomath::sin_cos_series(true, ssig1, csig1, C3a);
484        let domg12 = -self.f * self._A3f(eps) * salp0 * (sig12 + B312);
485        let lam12 = eta + domg12;
486
487        let mut dlam12: f64;
488        if diffp {
489            if calp2 == 0.0 {
490                dlam12 = -2.0 * self._f1 * dn1 / sbet1;
491            } else {
492                let res = self._Lengths(
493                    eps,
494                    sig12,
495                    ssig1,
496                    csig1,
497                    dn1,
498                    ssig2,
499                    csig2,
500                    dn2,
501                    cbet1,
502                    cbet2,
503                    caps::REDUCEDLENGTH,
504                    C1a,
505                    C2a,
506                );
507                dlam12 = res.1;
508                dlam12 *= self._f1 / (calp2 * cbet2);
509            }
510        } else {
511            dlam12 = f64::NAN;
512        }
513        (
514            lam12, salp2, calp2, sig12, ssig1, csig1, ssig2, csig2, eps, domg12, dlam12,
515        )
516    }
517
518    // returns (a12, s12, azi1, azi2, m12, M12, M21, S12)
519    pub fn _gen_inverse_azi(
520        &self,
521        lat1: f64,
522        lon1: f64,
523        lat2: f64,
524        lon2: f64,
525        outmask: u64,
526    ) -> (f64, f64, f64, f64, f64, f64, f64, f64) {
527        let mut azi1 = f64::NAN;
528        let mut azi2 = f64::NAN;
529        let outmask = outmask & caps::OUT_MASK;
530
531        let (a12, s12, salp1, calp1, salp2, calp2, m12, M12, M21, S12) =
532            self._gen_inverse(lat1, lon1, lat2, lon2, outmask);
533        if outmask & caps::AZIMUTH != 0 {
534            azi1 = geomath::atan2d(salp1, calp1);
535            azi2 = geomath::atan2d(salp2, calp2);
536        }
537        (a12, s12, azi1, azi2, m12, M12, M21, S12)
538    }
539
540    // returns (a12, s12, salp1, calp1, salp2, calp2, m12, M12, M21, S12)
541    pub fn _gen_inverse(
542        &self,
543        lat1: f64,
544        lon1: f64,
545        lat2: f64,
546        lon2: f64,
547        outmask: u64,
548    ) -> (f64, f64, f64, f64, f64, f64, f64, f64, f64, f64) {
549        let mut lat1 = lat1;
550        let mut lat2 = lat2;
551        let mut a12 = f64::NAN;
552        let mut s12 = f64::NAN;
553        let mut m12 = f64::NAN;
554        let mut M12 = f64::NAN;
555        let mut M21 = f64::NAN;
556        let mut S12 = f64::NAN;
557        let outmask = outmask & caps::OUT_MASK;
558
559        let (mut lon12, mut lon12s) = geomath::ang_diff(lon1, lon2);
560        let mut lonsign = if lon12 >= 0.0 { 1.0 } else { -1.0 };
561
562        lon12 = lonsign * geomath::ang_round(lon12);
563        lon12s = geomath::ang_round((180.0 - lon12) - lonsign * lon12s);
564        let lam12 = lon12.to_radians();
565        let slam12: f64;
566        let mut clam12: f64;
567        if lon12 > 90.0 {
568            let res = geomath::sincosd(lon12s);
569            slam12 = res.0;
570            clam12 = res.1;
571            clam12 = -clam12;
572        } else {
573            let res = geomath::sincosd(lon12);
574            slam12 = res.0;
575            clam12 = res.1;
576        };
577        lat1 = geomath::ang_round(geomath::lat_fix(lat1));
578        lat2 = geomath::ang_round(geomath::lat_fix(lat2));
579
580        let swapp = if lat1.abs() < lat2.abs() { -1.0 } else { 1.0 };
581        if swapp < 0.0 {
582            lonsign *= -1.0;
583            std::mem::swap(&mut lat2, &mut lat1);
584        }
585        let latsign = if lat1 < 0.0 { 1.0 } else { -1.0 };
586        lat1 *= latsign;
587        lat2 *= latsign;
588
589        let (mut sbet1, mut cbet1) = geomath::sincosd(lat1);
590        sbet1 *= self._f1;
591
592        geomath::norm(&mut sbet1, &mut cbet1);
593        cbet1 = cbet1.max(self.tiny_);
594
595        let (mut sbet2, mut cbet2) = geomath::sincosd(lat2);
596        sbet2 *= self._f1;
597
598        geomath::norm(&mut sbet2, &mut cbet2);
599        cbet2 = cbet2.max(self.tiny_);
600
601        if cbet1 < -sbet1 {
602            if cbet2 == cbet1 {
603                sbet2 = if sbet2 < 0.0 { sbet1 } else { -sbet1 };
604            }
605        } else if sbet2.abs() == -sbet1 {
606            cbet2 = cbet1;
607        }
608
609        let dn1 = (1.0 + self._ep2 * geomath::sq(sbet1)).sqrt();
610        let dn2 = (1.0 + self._ep2 * geomath::sq(sbet2)).sqrt();
611
612        const CARR_SIZE: usize = GEODESIC_ORDER + 1;
613        let mut C1a: [f64; CARR_SIZE] = [0.0; CARR_SIZE];
614        let mut C2a: [f64; CARR_SIZE] = [0.0; CARR_SIZE];
615        let mut C3a: [f64; GEODESIC_ORDER] = [0.0; GEODESIC_ORDER];
616
617        let mut meridian = lat1 == -90.0 || slam12 == 0.0;
618        let mut calp1 = 0.0;
619        let mut salp1 = 0.0;
620        let mut calp2 = 0.0;
621        let mut salp2 = 0.0;
622        let mut ssig1 = 0.0;
623        let mut csig1 = 0.0;
624        let mut ssig2 = 0.0;
625        let mut csig2 = 0.0;
626        let mut sig12: f64;
627        let mut s12x = 0.0;
628        let mut m12x = 0.0;
629
630        if meridian {
631            calp1 = clam12;
632            salp1 = slam12;
633            calp2 = 1.0;
634            salp2 = 0.0;
635
636            ssig1 = sbet1;
637            csig1 = calp1 * cbet1;
638            ssig2 = sbet2;
639            csig2 = calp2 * cbet2;
640
641            sig12 = ((csig1 * ssig2 - ssig1 * csig2).max(0.0)).atan2(csig1 * csig2 + ssig1 * ssig2);
642            let res = self._Lengths(
643                self._n,
644                sig12,
645                ssig1,
646                csig1,
647                dn1,
648                ssig2,
649                csig2,
650                dn2,
651                cbet1,
652                cbet2,
653                outmask | caps::DISTANCE | caps::REDUCEDLENGTH,
654                &mut C1a,
655                &mut C2a,
656            );
657            s12x = res.0;
658            m12x = res.1;
659            M12 = res.3;
660            M21 = res.4;
661
662            if sig12 < 1.0 || m12x >= 0.0 {
663                if sig12 < 3.0 * self.tiny_ {
664                    sig12 = 0.0;
665                    m12x = 0.0;
666                    s12x = 0.0;
667                }
668                m12x *= self._b;
669                s12x *= self._b;
670                a12 = sig12.to_degrees();
671            } else {
672                meridian = false;
673            }
674        }
675
676        let mut somg12 = 2.0;
677        let mut comg12 = 0.0;
678        let mut omg12 = 0.0;
679        let dnm: f64;
680        let mut eps = 0.0;
681        if !meridian && sbet1 == 0.0 && (self.f <= 0.0 || lon12s >= self.f * 180.0) {
682            calp1 = 0.0;
683            calp2 = 0.0;
684            salp1 = 1.0;
685            salp2 = 1.0;
686
687            s12x = self.a * lam12;
688            sig12 = lam12 / self._f1;
689            omg12 = lam12 / self._f1;
690            m12x = self._b * sig12.sin();
691            if outmask & caps::GEODESICSCALE != 0 {
692                M12 = sig12.cos();
693                M21 = sig12.cos();
694            }
695            a12 = lon12 / self._f1;
696        } else if !meridian {
697            let res = self._InverseStart(
698                sbet1, cbet1, dn1, sbet2, cbet2, dn2, lam12, slam12, clam12, &mut C1a, &mut C2a,
699            );
700            sig12 = res.0;
701            salp1 = res.1;
702            calp1 = res.2;
703            salp2 = res.3;
704            calp2 = res.4;
705            dnm = res.5;
706
707            if sig12 >= 0.0 {
708                s12x = sig12 * self._b * dnm;
709                m12x = geomath::sq(dnm) * self._b * (sig12 / dnm).sin();
710                if outmask & caps::GEODESICSCALE != 0 {
711                    M12 = (sig12 / dnm).cos();
712                    M21 = (sig12 / dnm).cos();
713                }
714                a12 = sig12.to_degrees();
715                omg12 = lam12 / (self._f1 * dnm);
716            } else {
717                let mut tripn = false;
718                let mut tripb = false;
719                let mut salp1a = self.tiny_;
720                let mut calp1a = 1.0;
721                let mut salp1b = self.tiny_;
722                let mut calp1b = -1.0;
723                let mut domg12 = 0.0;
724                for numit in 0..self.maxit2_ {
725                    let res = self._Lambda12(
726                        sbet1,
727                        cbet1,
728                        dn1,
729                        sbet2,
730                        cbet2,
731                        dn2,
732                        salp1,
733                        calp1,
734                        slam12,
735                        clam12,
736                        numit < self.maxit1_,
737                        &mut C1a,
738                        &mut C2a,
739                        &mut C3a,
740                    );
741                    let v = res.0;
742                    salp2 = res.1;
743                    calp2 = res.2;
744                    sig12 = res.3;
745                    ssig1 = res.4;
746                    csig1 = res.5;
747                    ssig2 = res.6;
748                    csig2 = res.7;
749                    eps = res.8;
750                    domg12 = res.9;
751                    let dv = res.10;
752
753                    if tripb
754                        || v.abs() < if tripn { 8.0 } else { 1.0 } * self.tol0_
755                        || v.abs().is_nan()
756                    {
757                        break;
758                    };
759                    if v > 0.0 && (numit > self.maxit1_ || calp1 / salp1 > calp1b / salp1b) {
760                        salp1b = salp1;
761                        calp1b = calp1;
762                    } else if v < 0.0 && (numit > self.maxit1_ || calp1 / salp1 < calp1a / salp1a) {
763                        salp1a = salp1;
764                        calp1a = calp1;
765                    }
766                    if numit < self.maxit1_ && dv > 0.0 {
767                        let dalp1 = -v / dv;
768                        let sdalp1 = dalp1.sin();
769                        let cdalp1 = dalp1.cos();
770                        let nsalp1 = salp1 * cdalp1 + calp1 * sdalp1;
771                        if nsalp1 > 0.0 && dalp1.abs() < PI {
772                            calp1 = calp1 * cdalp1 - salp1 * sdalp1;
773                            salp1 = nsalp1;
774                            geomath::norm(&mut salp1, &mut calp1);
775                            tripn = v.abs() <= 16.0 * self.tol0_;
776                            continue;
777                        }
778                    }
779
780                    salp1 = (salp1a + salp1b) / 2.0;
781                    calp1 = (calp1a + calp1b) / 2.0;
782                    geomath::norm(&mut salp1, &mut calp1);
783                    tripn = false;
784                    tripb = (salp1a - salp1).abs() + (calp1a - calp1) < self.tolb_
785                        || (salp1 - salp1b).abs() + (calp1 - calp1b) < self.tolb_;
786                }
787                let lengthmask = outmask
788                    | if outmask & (caps::REDUCEDLENGTH | caps::GEODESICSCALE) != 0 {
789                        caps::DISTANCE
790                    } else {
791                        caps::EMPTY
792                    };
793                let res = self._Lengths(
794                    eps, sig12, ssig1, csig1, dn1, ssig2, csig2, dn2, cbet1, cbet2, lengthmask,
795                    &mut C1a, &mut C2a,
796                );
797                s12x = res.0;
798                m12x = res.1;
799                M12 = res.3;
800                M21 = res.4;
801
802                m12x *= self._b;
803                s12x *= self._b;
804                a12 = sig12.to_degrees();
805                if outmask & caps::AREA != 0 {
806                    let sdomg12 = domg12.sin();
807                    let cdomg12 = domg12.cos();
808                    somg12 = slam12 * cdomg12 - clam12 * sdomg12;
809                    comg12 = clam12 * cdomg12 + slam12 * sdomg12;
810                }
811            }
812        }
813        if outmask & caps::DISTANCE != 0 {
814            s12 = 0.0 + s12x;
815        }
816        if outmask & caps::REDUCEDLENGTH != 0 {
817            m12 = 0.0 + m12x;
818        }
819        if outmask & caps::AREA != 0 {
820            let salp0 = salp1 * cbet1;
821            let calp0 = calp1.hypot(salp1 * sbet1);
822            if calp0 != 0.0 && salp0 != 0.0 {
823                ssig1 = sbet1;
824                csig1 = calp1 * cbet1;
825                ssig2 = sbet2;
826                csig2 = calp2 * cbet2;
827                let k2 = geomath::sq(calp0) * self._ep2;
828                eps = k2 / (2.0 * (1.0 + (1.0 + k2).sqrt()) + k2);
829                let A4 = geomath::sq(self.a) * calp0 * salp0 * self._e2;
830                geomath::norm(&mut ssig1, &mut csig1);
831                geomath::norm(&mut ssig2, &mut csig2);
832                let mut C4a: [f64; GEODESIC_ORDER] = [0.0; GEODESIC_ORDER];
833                self._C4f(eps, &mut C4a);
834                let B41 = geomath::sin_cos_series(false, ssig1, csig1, &C4a);
835                let B42 = geomath::sin_cos_series(false, ssig2, csig2, &C4a);
836                S12 = A4 * (B42 - B41);
837            } else {
838                S12 = 0.0;
839            }
840
841            if !meridian && somg12 > 1.0 {
842                somg12 = omg12.sin();
843                comg12 = omg12.cos();
844            }
845
846            // We're diverging from Karney's implementation here
847            // which uses the hardcoded constant: -0.7071 for FRAC_1_SQRT_2
848            let alp12: f64;
849            if !meridian && comg12 > -FRAC_1_SQRT_2 && sbet2 - sbet1 < 1.75 {
850                let domg12 = 1.0 + comg12;
851                let dbet1 = 1.0 + cbet1;
852                let dbet2 = 1.0 + cbet2;
853                alp12 = 2.0
854                    * (somg12 * (sbet1 * dbet2 + sbet2 * dbet1))
855                        .atan2(domg12 * (sbet1 * sbet2 + dbet1 * dbet2));
856            } else {
857                let mut salp12 = salp2 * calp1 - calp2 * salp1;
858                let mut calp12 = calp2 * calp1 + salp2 * salp1;
859
860                if salp12 == 0.0 && calp12 < 0.0 {
861                    salp12 = self.tiny_ * calp1;
862                    calp12 = -1.0;
863                }
864                alp12 = salp12.atan2(calp12);
865            }
866            S12 += self._c2 * alp12;
867            S12 *= swapp * lonsign * latsign;
868            S12 += 0.0;
869        }
870
871        if swapp < 0.0 {
872            std::mem::swap(&mut salp2, &mut salp1);
873
874            std::mem::swap(&mut calp2, &mut calp1);
875
876            if outmask & caps::GEODESICSCALE != 0 {
877                std::mem::swap(&mut M21, &mut M12);
878            }
879        }
880        salp1 *= swapp * lonsign;
881        calp1 *= swapp * latsign;
882        salp2 *= swapp * lonsign;
883        calp2 *= swapp * latsign;
884        (a12, s12, salp1, calp1, salp2, calp2, m12, M12, M21, S12)
885    }
886
887    ///  returns (a12, lat2, lon2, azi2, s12, m12, M12, M21, S12)
888    pub fn _gen_direct(
889        &self,
890        lat1: f64,
891        lon1: f64,
892        azi1: f64,
893        arcmode: bool,
894        s12_a12: f64,
895        mut outmask: u64,
896    ) -> (f64, f64, f64, f64, f64, f64, f64, f64, f64) {
897        if !arcmode {
898            outmask |= caps::DISTANCE_IN
899        };
900
901        let line =
902            geodesic_line::GeodesicLine::new(self, lat1, lon1, azi1, Some(outmask), None, None);
903        line._gen_position(arcmode, s12_a12, outmask)
904    }
905
906    /// Get the area of the geodesic in square meters
907    pub fn area(&self) -> f64 {
908        self._c2 * 4.0 * std::f64::consts::PI
909    }
910}
911
912/// Place a second point, given the first point, an azimuth, and a distance.
913///
914/// # Arguments
915///   - lat1 - Latitude of 1st point (degrees) [-90.,90.]
916///   - lon1 - Longitude of 1st point (degrees) [-180., 180.]
917///   - azi1 - Azimuth at 1st point (degrees) [-180., 180.]
918///   - s12 - Distance from 1st to 2nd point (meters) Value may be negative
919///
920/// # Returns
921///
922/// There are a variety of outputs associated with this calculation. We save computation by
923/// only calculating the outputs you need. See the following impls which return different subsets of
924/// the following outputs:
925///
926///  - lat2 latitude of point 2 (degrees).
927///  - lon2 longitude of point 2 (degrees).
928///  - azi2 (forward) azimuth at point 2 (degrees).
929///  - m12 reduced length of geodesic (meters).
930///  - M12 geodesic scale of point 2 relative to point 1 (dimensionless).
931///  - M21 geodesic scale of point 1 relative to point 2 (dimensionless).
932///  - S12 area under the geodesic (meters<sup>2</sup>).
933///  - a12 arc length between point 1 and point 2 (degrees).
934///
935///  If either point is at a pole, the azimuth is defined by keeping the
936///  longitude fixed, writing lat = ±(90° − ε), and taking the limit ε → 0+.
937///  An arc length greater that 180° signifies a geodesic which is not a
938///  shortest path. (For a prolate ellipsoid, an additional condition is
939///  necessary for a shortest path: the longitudinal extent must not
940///  exceed of 180°.)
941/// ```rust
942/// // Example, determine the point 10000 km NE of JFK:
943/// use geographiclib_rs::{Geodesic, DirectGeodesic};
944///
945/// let g = Geodesic::wgs84();
946/// let (lat, lon, az) = g.direct(40.64, -73.78, 45.0, 10e6);
947///
948/// use approx::assert_relative_eq;
949/// assert_relative_eq!(lat, 32.621100463725796);
950/// assert_relative_eq!(lon, 49.052487092959836);
951/// assert_relative_eq!(az,  140.4059858768007);
952/// ```
953pub trait DirectGeodesic<T> {
954    fn direct(&self, lat1: f64, lon1: f64, azi1: f64, s12: f64) -> T;
955}
956
957impl DirectGeodesic<(f64, f64)> for Geodesic {
958    /// See the documentation for the DirectGeodesic trait.
959    ///
960    /// # Returns
961    ///  - lat2 latitude of point 2 (degrees).
962    ///  - lon2 longitude of point 2 (degrees).
963    fn direct(&self, lat1: f64, lon1: f64, azi1: f64, s12: f64) -> (f64, f64) {
964        let capabilities = caps::LATITUDE | caps::LONGITUDE;
965        let (_a12, lat2, lon2, _azi2, _s12, _m12, _M12, _M21, _S12) =
966            self._gen_direct(lat1, lon1, azi1, false, s12, capabilities);
967
968        (lat2, lon2)
969    }
970}
971
972impl DirectGeodesic<(f64, f64, f64)> for Geodesic {
973    /// See the documentation for the DirectGeodesic trait.
974    ///
975    /// # Returns
976    ///  - lat2 latitude of point 2 (degrees).
977    ///  - lon2 longitude of point 2 (degrees).
978    ///  - azi2 (forward) azimuth at point 2 (degrees).
979    fn direct(&self, lat1: f64, lon1: f64, azi1: f64, s12: f64) -> (f64, f64, f64) {
980        let capabilities = caps::LATITUDE | caps::LONGITUDE | caps::AZIMUTH;
981        let (_a12, lat2, lon2, azi2, _s12, _m12, _M12, _M21, _S12) =
982            self._gen_direct(lat1, lon1, azi1, false, s12, capabilities);
983
984        (lat2, lon2, azi2)
985    }
986}
987
988impl DirectGeodesic<(f64, f64, f64, f64)> for Geodesic {
989    /// See the documentation for the DirectGeodesic trait.
990    ///
991    /// # Returns
992    ///  - lat2 latitude of point 2 (degrees).
993    ///  - lon2 longitude of point 2 (degrees).
994    ///  - azi2 (forward) azimuth at point 2 (degrees).
995    ///  - m12 reduced length of geodesic (meters).
996    fn direct(&self, lat1: f64, lon1: f64, azi1: f64, s12: f64) -> (f64, f64, f64, f64) {
997        let capabilities = caps::LATITUDE | caps::LONGITUDE | caps::AZIMUTH | caps::REDUCEDLENGTH;
998        let (_a12, lat2, lon2, azi2, _s12, m12, _M12, _M21, _S12) =
999            self._gen_direct(lat1, lon1, azi1, false, s12, capabilities);
1000
1001        (lat2, lon2, azi2, m12)
1002    }
1003}
1004
1005impl DirectGeodesic<(f64, f64, f64, f64, f64)> for Geodesic {
1006    /// See the documentation for the DirectGeodesic trait.
1007    ///
1008    /// # Returns
1009    ///  - lat2 latitude of point 2 (degrees).
1010    ///  - lon2 longitude of point 2 (degrees).
1011    ///  - azi2 (forward) azimuth at point 2 (degrees).
1012    ///  - M12 geodesic scale of point 2 relative to point 1 (dimensionless).
1013    ///  - M21 geodesic scale of point 1 relative to point 2 (dimensionless).
1014    fn direct(&self, lat1: f64, lon1: f64, azi1: f64, s12: f64) -> (f64, f64, f64, f64, f64) {
1015        let capabilities = caps::LATITUDE | caps::LONGITUDE | caps::AZIMUTH | caps::GEODESICSCALE;
1016        let (_a12, lat2, lon2, azi2, _s12, _m12, M12, M21, _S12) =
1017            self._gen_direct(lat1, lon1, azi1, false, s12, capabilities);
1018
1019        (lat2, lon2, azi2, M12, M21)
1020    }
1021}
1022
1023impl DirectGeodesic<(f64, f64, f64, f64, f64, f64)> for Geodesic {
1024    /// See the documentation for the DirectGeodesic trait.
1025    ///
1026    /// # Returns
1027    ///  - lat2 latitude of point 2 (degrees).
1028    ///  - lon2 longitude of point 2 (degrees).
1029    ///  - azi2 (forward) azimuth at point 2 (degrees).
1030    ///  - m12 reduced length of geodesic (meters).
1031    ///  - M12 geodesic scale of point 2 relative to point 1 (dimensionless).
1032    ///  - M21 geodesic scale of point 1 relative to point 2 (dimensionless).
1033    fn direct(&self, lat1: f64, lon1: f64, azi1: f64, s12: f64) -> (f64, f64, f64, f64, f64, f64) {
1034        let capabilities = caps::LATITUDE
1035            | caps::LONGITUDE
1036            | caps::AZIMUTH
1037            | caps::REDUCEDLENGTH
1038            | caps::GEODESICSCALE;
1039        let (_a12, lat2, lon2, azi2, _s12, m12, M12, M21, _S12) =
1040            self._gen_direct(lat1, lon1, azi1, false, s12, capabilities);
1041
1042        (lat2, lon2, azi2, m12, M12, M21)
1043    }
1044}
1045
1046impl DirectGeodesic<(f64, f64, f64, f64, f64, f64, f64, f64)> for Geodesic {
1047    /// See the documentation for the DirectGeodesic trait.
1048    ///
1049    /// # Returns
1050    ///  - lat2 latitude of point 2 (degrees).
1051    ///  - lon2 longitude of point 2 (degrees).
1052    ///  - azi2 (forward) azimuth at point 2 (degrees).
1053    ///  - m12 reduced length of geodesic (meters).
1054    ///  - M12 geodesic scale of point 2 relative to point 1 (dimensionless).
1055    ///  - M21 geodesic scale of point 1 relative to point 2 (dimensionless).
1056    ///  - S12 area under the geodesic (meters<sup>2</sup>).
1057    ///  - a12 arc length between point 1 and point 2 (degrees).
1058    fn direct(
1059        &self,
1060        lat1: f64,
1061        lon1: f64,
1062        azi1: f64,
1063        s12: f64,
1064    ) -> (f64, f64, f64, f64, f64, f64, f64, f64) {
1065        let capabilities = caps::LATITUDE
1066            | caps::LONGITUDE
1067            | caps::AZIMUTH
1068            | caps::REDUCEDLENGTH
1069            | caps::GEODESICSCALE
1070            | caps::AREA;
1071        let (a12, lat2, lon2, azi2, _s12, m12, M12, M21, S12) =
1072            self._gen_direct(lat1, lon1, azi1, false, s12, capabilities);
1073
1074        (lat2, lon2, azi2, m12, M12, M21, S12, a12)
1075    }
1076}
1077
1078/// Measure the distance (and other values) between two points.
1079///
1080/// # Arguments
1081/// - lat1 latitude of point 1 (degrees).
1082/// - lon1 longitude of point 1 (degrees).
1083/// - lat2 latitude of point 2 (degrees).
1084/// - lon2 longitude of point 2 (degrees).
1085///
1086/// # Returns
1087///
1088/// There are a variety of outputs associated with this calculation. We save computation by
1089/// only calculating the outputs you need. See the following impls which return different subsets of
1090/// the following outputs:
1091///
1092/// - s12 distance between point 1 and point 2 (meters).
1093/// - azi1 azimuth at point 1 (degrees).
1094/// - azi2 (forward) azimuth at point 2 (degrees).
1095/// - m12 reduced length of geodesic (meters).
1096/// - M12 geodesic scale of point 2 relative to point 1 (dimensionless).
1097/// - M21 geodesic scale of point 1 relative to point 2 (dimensionless).
1098/// - S12 area under the geodesic (meters<sup>2</sup>).
1099/// - a12 arc length between point 1 and point 2 (degrees).
1100///
1101///  `lat1` and `lat2` should be in the range [&minus;90&deg;, 90&deg;].
1102///  The values of `azi1` and `azi2` returned are in the range
1103///  [&minus;180&deg;, 180&deg;].
1104///
1105/// If either point is at a pole, the azimuth is defined by keeping the
1106/// longitude fixed, writing `lat` = &plusmn;(90&deg; &minus; &epsilon;),
1107/// and taking the limit &epsilon; &rarr; 0+.
1108///
1109/// The solution to the inverse problem is found using Newton's method.  If
1110/// this fails to converge (this is very unlikely in geodetic applications
1111/// but does occur for very eccentric ellipsoids), then the bisection method
1112/// is used to refine the solution.
1113///
1114/// ```rust
1115/// // Example, determine the distance between two points
1116/// use geographiclib_rs::{Geodesic, InverseGeodesic};
1117///
1118/// let g = Geodesic::wgs84();
1119/// let p1 = (34.095925, -118.2884237);
1120/// let p2 = (59.4323439, 24.7341649);
1121/// let s12: f64 = g.inverse(p1.0, p1.1, p2.0, p2.1);
1122///
1123/// use approx::assert_relative_eq;
1124/// assert_relative_eq!(s12, 9094718.72751138);
1125/// ```
1126pub trait InverseGeodesic<T> {
1127    fn inverse(&self, lat1: f64, lon1: f64, lat2: f64, lon2: f64) -> T;
1128}
1129
1130impl InverseGeodesic<f64> for Geodesic {
1131    /// See the documentation for the InverseGeodesic trait.
1132    ///
1133    /// # Returns
1134    /// - s12 distance between point 1 and point 2 (meters).
1135    fn inverse(&self, lat1: f64, lon1: f64, lat2: f64, lon2: f64) -> f64 {
1136        let capabilities = caps::DISTANCE;
1137        let (_a12, s12, _azi1, _azi2, _m12, _M12, _M21, _S12) =
1138            self._gen_inverse_azi(lat1, lon1, lat2, lon2, capabilities);
1139
1140        s12
1141    }
1142}
1143
1144impl InverseGeodesic<(f64, f64)> for Geodesic {
1145    /// See the documentation for the InverseGeodesic trait.
1146    ///
1147    /// # Returns
1148    /// - s12 distance between point 1 and point 2 (meters).
1149    /// - a12 arc length between point 1 and point 2 (degrees).
1150    fn inverse(&self, lat1: f64, lon1: f64, lat2: f64, lon2: f64) -> (f64, f64) {
1151        let capabilities = caps::DISTANCE;
1152        let (a12, s12, _azi1, _azi2, _m12, _M12, _M21, _S12) =
1153            self._gen_inverse_azi(lat1, lon1, lat2, lon2, capabilities);
1154
1155        (s12, a12)
1156    }
1157}
1158
1159impl InverseGeodesic<(f64, f64, f64)> for Geodesic {
1160    /// See the documentation for the InverseGeodesic trait.
1161    ///
1162    /// # Returns
1163    /// - azi1 azimuth at point 1 (degrees).
1164    /// - azi2 (forward) azimuth at point 2 (degrees).
1165    /// - a12 arc length between point 1 and point 2 (degrees).
1166    fn inverse(&self, lat1: f64, lon1: f64, lat2: f64, lon2: f64) -> (f64, f64, f64) {
1167        let capabilities = caps::AZIMUTH;
1168        let (a12, _s12, azi1, azi2, _m12, _M12, _M21, _S12) =
1169            self._gen_inverse_azi(lat1, lon1, lat2, lon2, capabilities);
1170
1171        (azi1, azi2, a12)
1172    }
1173}
1174
1175impl InverseGeodesic<(f64, f64, f64, f64)> for Geodesic {
1176    /// See the documentation for the InverseGeodesic trait.
1177    ///
1178    /// # Returns
1179    /// - s12 distance between point 1 and point 2 (meters).
1180    /// - azi1 azimuth at point 1 (degrees).
1181    /// - azi2 (forward) azimuth at point 2 (degrees).
1182    /// - a12 arc length between point 1 and point 2 (degrees).
1183    fn inverse(&self, lat1: f64, lon1: f64, lat2: f64, lon2: f64) -> (f64, f64, f64, f64) {
1184        let capabilities = caps::DISTANCE | caps::AZIMUTH;
1185        let (a12, s12, azi1, azi2, _m12, _M12, _M21, _S12) =
1186            self._gen_inverse_azi(lat1, lon1, lat2, lon2, capabilities);
1187
1188        (s12, azi1, azi2, a12)
1189    }
1190}
1191
1192impl InverseGeodesic<(f64, f64, f64, f64, f64)> for Geodesic {
1193    /// See the documentation for the InverseGeodesic trait.
1194    ///
1195    /// # Returns
1196    /// - s12 distance between point 1 and point 2 (meters).
1197    /// - azi1 azimuth at point 1 (degrees).
1198    /// - azi2 (forward) azimuth at point 2 (degrees).
1199    /// - m12 reduced length of geodesic (meters).
1200    /// - a12 arc length between point 1 and point 2 (degrees).
1201    fn inverse(&self, lat1: f64, lon1: f64, lat2: f64, lon2: f64) -> (f64, f64, f64, f64, f64) {
1202        let capabilities = caps::DISTANCE | caps::AZIMUTH | caps::REDUCEDLENGTH;
1203        let (a12, s12, azi1, azi2, m12, _M12, _M21, _S12) =
1204            self._gen_inverse_azi(lat1, lon1, lat2, lon2, capabilities);
1205
1206        (s12, azi1, azi2, m12, a12)
1207    }
1208}
1209
1210impl InverseGeodesic<(f64, f64, f64, f64, f64, f64)> for Geodesic {
1211    /// See the documentation for the InverseGeodesic trait.
1212    ///
1213    /// # Returns
1214    /// - s12 distance between point 1 and point 2 (meters).
1215    /// - azi1 azimuth at point 1 (degrees).
1216    /// - azi2 (forward) azimuth at point 2 (degrees).
1217    /// - M12 geodesic scale of point 2 relative to point 1 (dimensionless).
1218    /// - M21 geodesic scale of point 1 relative to point 2 (dimensionless).
1219    /// - a12 arc length between point 1 and point 2 (degrees).
1220    fn inverse(
1221        &self,
1222        lat1: f64,
1223        lon1: f64,
1224        lat2: f64,
1225        lon2: f64,
1226    ) -> (f64, f64, f64, f64, f64, f64) {
1227        let capabilities = caps::DISTANCE | caps::AZIMUTH | caps::GEODESICSCALE;
1228        let (a12, s12, azi1, azi2, _m12, M12, M21, _S12) =
1229            self._gen_inverse_azi(lat1, lon1, lat2, lon2, capabilities);
1230
1231        (s12, azi1, azi2, M12, M21, a12)
1232    }
1233}
1234
1235impl InverseGeodesic<(f64, f64, f64, f64, f64, f64, f64)> for Geodesic {
1236    /// See the documentation for the InverseGeodesic trait.
1237    ///
1238    /// # Returns
1239    /// - s12 distance between point 1 and point 2 (meters).
1240    /// - azi1 azimuth at point 1 (degrees).
1241    /// - azi2 (forward) azimuth at point 2 (degrees).
1242    /// - m12 reduced length of geodesic (meters).
1243    /// - M12 geodesic scale of point 2 relative to point 1 (dimensionless).
1244    /// - M21 geodesic scale of point 1 relative to point 2 (dimensionless).
1245    /// - a12 arc length between point 1 and point 2 (degrees).
1246    fn inverse(
1247        &self,
1248        lat1: f64,
1249        lon1: f64,
1250        lat2: f64,
1251        lon2: f64,
1252    ) -> (f64, f64, f64, f64, f64, f64, f64) {
1253        let capabilities =
1254            caps::DISTANCE | caps::AZIMUTH | caps::REDUCEDLENGTH | caps::GEODESICSCALE;
1255        let (a12, s12, azi1, azi2, m12, M12, M21, _S12) =
1256            self._gen_inverse_azi(lat1, lon1, lat2, lon2, capabilities);
1257
1258        (s12, azi1, azi2, m12, M12, M21, a12)
1259    }
1260}
1261
1262impl InverseGeodesic<(f64, f64, f64, f64, f64, f64, f64, f64)> for Geodesic {
1263    /// See the documentation for the InverseGeodesic trait.
1264    ///
1265    /// # Returns
1266    /// - s12 distance between point 1 and point 2 (meters).
1267    /// - azi1 azimuth at point 1 (degrees).
1268    /// - azi2 (forward) azimuth at point 2 (degrees).
1269    /// - m12 reduced length of geodesic (meters).
1270    /// - M12 geodesic scale of point 2 relative to point 1 (dimensionless).
1271    /// - M21 geodesic scale of point 1 relative to point 2 (dimensionless).
1272    /// - S12 area under the geodesic (meters<sup>2</sup>).
1273    /// - a12 arc length between point 1 and point 2 (degrees).
1274    fn inverse(
1275        &self,
1276        lat1: f64,
1277        lon1: f64,
1278        lat2: f64,
1279        lon2: f64,
1280    ) -> (f64, f64, f64, f64, f64, f64, f64, f64) {
1281        let capabilities =
1282            caps::DISTANCE | caps::AZIMUTH | caps::REDUCEDLENGTH | caps::GEODESICSCALE | caps::AREA;
1283        let (a12, s12, azi1, azi2, m12, M12, M21, S12) =
1284            self._gen_inverse_azi(lat1, lon1, lat2, lon2, capabilities);
1285
1286        (s12, azi1, azi2, m12, M12, M21, S12, a12)
1287    }
1288}
1289
1290#[cfg(test)]
1291mod tests {
1292    use super::*;
1293    use crate::geodesic_line::GeodesicLine;
1294    use approx::assert_relative_eq;
1295    use std::io::BufRead;
1296
1297    #[allow(clippy::type_complexity)]
1298    const TESTCASES: &[(f64, f64, f64, f64, f64, f64, f64, f64, f64, f64, f64, f64)] = &[
1299        (
1300            35.60777,
1301            -139.44815,
1302            111.098748429560326,
1303            -11.17491,
1304            -69.95921,
1305            129.289270889708762,
1306            8935244.5604818305,
1307            80.50729714281974,
1308            6273170.2055303837,
1309            0.16606318447386067,
1310            0.16479116945612937,
1311            12841384694976.432,
1312        ),
1313        (
1314            55.52454,
1315            106.05087,
1316            22.020059880982801,
1317            77.03196,
1318            197.18234,
1319            109.112041110671519,
1320            4105086.1713924406,
1321            36.892740690445894,
1322            3828869.3344387607,
1323            0.80076349608092607,
1324            0.80101006984201008,
1325            61674961290615.615,
1326        ),
1327        (
1328            -21.97856,
1329            142.59065,
1330            -32.44456876433189,
1331            41.84138,
1332            98.56635,
1333            -41.84359951440466,
1334            8394328.894657671,
1335            75.62930491011522,
1336            6161154.5773110616,
1337            0.24816339233950381,
1338            0.24930251203627892,
1339            -6637997720646.717,
1340        ),
1341        (
1342            -66.99028,
1343            112.2363,
1344            173.73491240878403,
1345            -12.70631,
1346            285.90344,
1347            2.512956620913668,
1348            11150344.2312080241,
1349            100.278634181155759,
1350            6289939.5670446687,
1351            -0.17199490274700385,
1352            -0.17722569526345708,
1353            -121287239862139.744,
1354        ),
1355        (
1356            -17.42761,
1357            173.34268,
1358            -159.033557661192928,
1359            -15.84784,
1360            5.93557,
1361            -20.787484651536988,
1362            16076603.1631180673,
1363            144.640108810286253,
1364            3732902.1583877189,
1365            -0.81273638700070476,
1366            -0.81299800519154474,
1367            97825992354058.708,
1368        ),
1369        (
1370            32.84994,
1371            48.28919,
1372            150.492927788121982,
1373            -56.28556,
1374            202.29132,
1375            48.113449399816759,
1376            16727068.9438164461,
1377            150.565799985466607,
1378            3147838.1910180939,
1379            -0.87334918086923126,
1380            -0.86505036767110637,
1381            -72445258525585.010,
1382        ),
1383        (
1384            6.96833,
1385            52.74123,
1386            92.581585386317712,
1387            -7.39675,
1388            206.17291,
1389            90.721692165923907,
1390            17102477.2496958388,
1391            154.147366239113561,
1392            2772035.6169917581,
1393            -0.89991282520302447,
1394            -0.89986892177110739,
1395            -1311796973197.995,
1396        ),
1397        (
1398            -50.56724,
1399            -16.30485,
1400            -105.439679907590164,
1401            -33.56571,
1402            -94.97412,
1403            -47.348547835650331,
1404            6455670.5118668696,
1405            58.083719495371259,
1406            5409150.7979815838,
1407            0.53053508035997263,
1408            0.52988722644436602,
1409            41071447902810.047,
1410        ),
1411        (
1412            -58.93002,
1413            -8.90775,
1414            140.965397902500679,
1415            -8.91104,
1416            133.13503,
1417            19.255429433416599,
1418            11756066.0219864627,
1419            105.755691241406877,
1420            6151101.2270708536,
1421            -0.26548622269867183,
1422            -0.27068483874510741,
1423            -86143460552774.735,
1424        ),
1425        (
1426            -68.82867,
1427            -74.28391,
1428            93.774347763114881,
1429            -50.63005,
1430            -8.36685,
1431            34.65564085411343,
1432            3956936.926063544,
1433            35.572254987389284,
1434            3708890.9544062657,
1435            0.81443963736383502,
1436            0.81420859815358342,
1437            -41845309450093.787,
1438        ),
1439        (
1440            -10.62672,
1441            -32.0898,
1442            -86.426713286747751,
1443            5.883,
1444            -134.31681,
1445            -80.473780971034875,
1446            11470869.3864563009,
1447            103.387395634504061,
1448            6184411.6622659713,
1449            -0.23138683500430237,
1450            -0.23155097622286792,
1451            4198803992123.548,
1452        ),
1453        (
1454            -21.76221,
1455            166.90563,
1456            29.319421206936428,
1457            48.72884,
1458            213.97627,
1459            43.508671946410168,
1460            9098627.3986554915,
1461            81.963476716121964,
1462            6299240.9166992283,
1463            0.13965943368590333,
1464            0.14152969707656796,
1465            10024709850277.476,
1466        ),
1467        (
1468            -19.79938,
1469            -174.47484,
1470            71.167275780171533,
1471            -11.99349,
1472            -154.35109,
1473            65.589099775199228,
1474            2319004.8601169389,
1475            20.896611684802389,
1476            2267960.8703918325,
1477            0.93427001867125849,
1478            0.93424887135032789,
1479            -3935477535005.785,
1480        ),
1481        (
1482            -11.95887,
1483            -116.94513,
1484            92.712619830452549,
1485            4.57352,
1486            7.16501,
1487            78.64960934409585,
1488            13834722.5801401374,
1489            124.688684161089762,
1490            5228093.177931598,
1491            -0.56879356755666463,
1492            -0.56918731952397221,
1493            -9919582785894.853,
1494        ),
1495        (
1496            -87.85331,
1497            85.66836,
1498            -65.120313040242748,
1499            66.48646,
1500            16.09921,
1501            -4.888658719272296,
1502            17286615.3147144645,
1503            155.58592449699137,
1504            2635887.4729110181,
1505            -0.90697975771398578,
1506            -0.91095608883042767,
1507            42667211366919.534,
1508        ),
1509        (
1510            1.74708,
1511            128.32011,
1512            -101.584843631173858,
1513            -11.16617,
1514            11.87109,
1515            -86.325793296437476,
1516            12942901.1241347408,
1517            116.650512484301857,
1518            5682744.8413270572,
1519            -0.44857868222697644,
1520            -0.44824490340007729,
1521            10763055294345.653,
1522        ),
1523        (
1524            -25.72959,
1525            -144.90758,
1526            -153.647468693117198,
1527            -57.70581,
1528            -269.17879,
1529            -48.343983158876487,
1530            9413446.7452453107,
1531            84.664533838404295,
1532            6356176.6898881281,
1533            0.09492245755254703,
1534            0.09737058264766572,
1535            74515122850712.444,
1536        ),
1537        (
1538            -41.22777,
1539            122.32875,
1540            14.285113402275739,
1541            -7.57291,
1542            130.37946,
1543            10.805303085187369,
1544            3812686.035106021,
1545            34.34330804743883,
1546            3588703.8812128856,
1547            0.82605222593217889,
1548            0.82572158200920196,
1549            -2456961531057.857,
1550        ),
1551        (
1552            11.01307,
1553            138.25278,
1554            79.43682622782374,
1555            6.62726,
1556            247.05981,
1557            103.708090215522657,
1558            11911190.819018408,
1559            107.341669954114577,
1560            6070904.722786735,
1561            -0.29767608923657404,
1562            -0.29785143390252321,
1563            17121631423099.696,
1564        ),
1565        (
1566            -29.47124,
1567            95.14681,
1568            -163.779130441688382,
1569            -27.46601,
1570            -69.15955,
1571            -15.909335945554969,
1572            13487015.8381145492,
1573            121.294026715742277,
1574            5481428.9945736388,
1575            -0.51527225545373252,
1576            -0.51556587964721788,
1577            104679964020340.318,
1578        ),
1579    ];
1580
1581    #[test]
1582    fn test_inverse_and_direct() -> Result<(), String> {
1583        // See python/test_geodesic.py
1584        let geod = Geodesic::wgs84();
1585        let (_a12, s12, _azi1, _azi2, _m12, _M12, _M21, _S12) =
1586            geod._gen_inverse_azi(0.0, 0.0, 1.0, 1.0, caps::STANDARD);
1587        assert_eq!(s12, 156899.56829134026);
1588
1589        // Test inverse
1590        for (lat1, lon1, azi1, lat2, lon2, azi2, s12, a12, m12, M12, M21, S12) in TESTCASES.iter() {
1591            let (
1592                computed_a12,
1593                computed_s12,
1594                computed_azi1,
1595                computed_azi2,
1596                computed_m12,
1597                computed_M12,
1598                computed_M21,
1599                computed_S12,
1600            ) = geod._gen_inverse_azi(*lat1, *lon1, *lat2, *lon2, caps::ALL | caps::LONG_UNROLL);
1601            assert_relative_eq!(computed_azi1, azi1, epsilon = 1e-13f64);
1602            assert_relative_eq!(computed_azi2, azi2, epsilon = 1e-13f64);
1603            assert_relative_eq!(computed_s12, s12, epsilon = 1e-8f64);
1604            assert_relative_eq!(computed_a12, a12, epsilon = 1e-13f64);
1605            assert_relative_eq!(computed_m12, m12, epsilon = 1e-8f64);
1606            assert_relative_eq!(computed_M12, M12, epsilon = 1e-15f64);
1607            assert_relative_eq!(computed_M21, M21, epsilon = 1e-15f64);
1608            assert_relative_eq!(computed_S12, S12, epsilon = 0.1f64);
1609        }
1610
1611        // Test direct
1612        for (lat1, lon1, azi1, lat2, lon2, azi2, s12, a12, m12, M12, M21, S12) in TESTCASES.iter() {
1613            let (
1614                computed_a12,
1615                computed_lat2,
1616                computed_lon2,
1617                computed_azi2,
1618                _computed_s12,
1619                computed_m12,
1620                computed_M12,
1621                computed_M21,
1622                computed_S12,
1623            ) = geod._gen_direct(
1624                *lat1,
1625                *lon1,
1626                *azi1,
1627                false,
1628                *s12,
1629                caps::ALL | caps::LONG_UNROLL,
1630            );
1631            assert_relative_eq!(computed_lat2, lat2, epsilon = 1e-13f64);
1632            assert_relative_eq!(computed_lon2, lon2, epsilon = 1e-13f64);
1633            assert_relative_eq!(computed_azi2, azi2, epsilon = 1e-13f64);
1634            assert_relative_eq!(computed_a12, a12, epsilon = 1e-13f64);
1635            assert_relative_eq!(computed_m12, m12, epsilon = 1e-8f64);
1636            assert_relative_eq!(computed_M12, M12, epsilon = 1e-15f64);
1637            assert_relative_eq!(computed_M21, M21, epsilon = 1e-15f64);
1638            assert_relative_eq!(computed_S12, S12, epsilon = 0.1f64);
1639        }
1640        Ok(())
1641    }
1642
1643    #[test]
1644    fn test_arcdirect() {
1645        // Corresponds with ArcDirectCheck from Java, or test_arcdirect from Python
1646        let geod = Geodesic::wgs84();
1647        for (lat1, lon1, azi1, lat2, lon2, azi2, s12, a12, m12, M12, M21, S12) in TESTCASES.iter() {
1648            let (
1649                _computed_a12,
1650                computed_lat2,
1651                computed_lon2,
1652                computed_azi2,
1653                computed_s12,
1654                computed_m12,
1655                computed_M12,
1656                computed_M21,
1657                computed_S12,
1658            ) = geod._gen_direct(
1659                *lat1,
1660                *lon1,
1661                *azi1,
1662                true,
1663                *a12,
1664                caps::ALL | caps::LONG_UNROLL,
1665            );
1666            assert_relative_eq!(computed_lat2, lat2, epsilon = 1e-13);
1667            assert_relative_eq!(computed_lon2, lon2, epsilon = 1e-13);
1668            assert_relative_eq!(computed_azi2, azi2, epsilon = 1e-13);
1669            assert_relative_eq!(computed_s12, s12, epsilon = 1e-8);
1670            assert_relative_eq!(computed_m12, m12, epsilon = 1e-8);
1671            assert_relative_eq!(computed_M12, M12, epsilon = 1e-15);
1672            assert_relative_eq!(computed_M21, M21, epsilon = 1e-15);
1673            assert_relative_eq!(computed_S12, S12, epsilon = 0.1);
1674        }
1675    }
1676
1677    #[test]
1678    fn test_geninverse() {
1679        let geod = Geodesic::wgs84();
1680        let res = geod._gen_inverse(0.0, 0.0, 1.0, 1.0, caps::STANDARD);
1681        assert_eq!(res.0, 1.4141938478710363);
1682        assert_eq!(res.1, 156899.56829134026);
1683        assert_eq!(res.2, 0.7094236375834774);
1684        assert_eq!(res.3, 0.7047823085448635);
1685        assert_eq!(res.4, 0.7095309793242709);
1686        assert_eq!(res.5, 0.7046742434480923);
1687        assert!(res.6.is_nan());
1688        assert!(res.7.is_nan());
1689        assert!(res.8.is_nan());
1690        assert!(res.9.is_nan());
1691    }
1692
1693    #[test]
1694    fn test_inverse_start() {
1695        let geod = Geodesic::wgs84();
1696        let res = geod._InverseStart(
1697            -0.017393909556108908,
1698            0.9998487145115275,
1699            1.0000010195104125,
1700            -0.0,
1701            1.0,
1702            1.0,
1703            0.017453292519943295,
1704            0.01745240643728351,
1705            0.9998476951563913,
1706            &mut [0.0, 1.0, 2.0, 3.0, 4.0, 5.0, 6.0],
1707            &mut [0.0, 1.0, 2.0, 3.0, 4.0, 5.0, 6.0],
1708        );
1709        assert_eq!(res.0, -1.0);
1710        assert_relative_eq!(res.1, 0.7095310092765433, epsilon = 1e-13);
1711        assert_relative_eq!(res.2, 0.7046742132893822, epsilon = 1e-13);
1712        assert!(res.3.is_nan());
1713        assert!(res.4.is_nan());
1714        assert_eq!(res.5, 1.0000002548969817);
1715
1716        let res = geod._InverseStart(
1717            -0.017393909556108908,
1718            0.9998487145115275,
1719            1.0000010195104125,
1720            -0.0,
1721            1.0,
1722            1.0,
1723            0.017453292519943295,
1724            0.01745240643728351,
1725            0.9998476951563913,
1726            &mut [0.0, 1.0, 2.0, 3.0, 4.0, 5.0, 6.0],
1727            &mut [0.0, 1.0, 2.0, 3.0, 4.0, 5.0, 6.0],
1728        );
1729        assert_eq!(res.0, -1.0);
1730        assert_relative_eq!(res.1, 0.7095310092765433, epsilon = 1e-13);
1731        assert_relative_eq!(res.2, 0.7046742132893822, epsilon = 1e-13);
1732        assert!(res.3.is_nan());
1733        assert!(res.4.is_nan());
1734        assert_eq!(res.5, 1.0000002548969817);
1735    }
1736
1737    #[test]
1738    fn test_lambda12() {
1739        let geod = Geodesic::wgs84();
1740        let res1 = geod._Lambda12(
1741            -0.017393909556108908,
1742            0.9998487145115275,
1743            1.0000010195104125,
1744            -0.0,
1745            1.0,
1746            1.0,
1747            0.7095310092765433,
1748            0.7046742132893822,
1749            0.01745240643728351,
1750            0.9998476951563913,
1751            true,
1752            &mut [0.0, 1.0, 2.0, 3.0, 4.0, 5.0, 6.0],
1753            &mut [0.0, 1.0, 2.0, 3.0, 4.0, 5.0, 6.0],
1754            &mut [0.0, 1.0, 2.0, 3.0, 4.0, 5.0],
1755        );
1756        assert_eq!(res1.0, 1.4834408705897495e-09);
1757        assert_eq!(res1.1, 0.7094236675312185);
1758        assert_eq!(res1.2, 0.7047822783999007);
1759        assert_eq!(res1.3, 0.024682339962725352);
1760        assert_eq!(res1.4, -0.024679833885152578);
1761        assert_eq!(res1.5, 0.9996954065111039);
1762        assert_eq!(res1.6, -0.0);
1763        assert_eq!(res1.7, 1.0);
1764        assert_relative_eq!(res1.8, 0.0008355095326524276, epsilon = 1e-13);
1765        assert_eq!(res1.9, -5.8708496511415445e-05);
1766        assert_eq!(res1.10, 0.034900275148485);
1767
1768        let res2 = geod._Lambda12(
1769            -0.017393909556108908,
1770            0.9998487145115275,
1771            1.0000010195104125,
1772            -0.0,
1773            1.0,
1774            1.0,
1775            0.7095309793242709,
1776            0.7046742434480923,
1777            0.01745240643728351,
1778            0.9998476951563913,
1779            true,
1780            &mut [
1781                0.0,
1782                -0.00041775465696698233,
1783                -4.362974596862037e-08,
1784                -1.2151022357848552e-11,
1785                -4.7588881620421004e-15,
1786                -2.226614930167366e-18,
1787                -1.1627237498131586e-21,
1788            ],
1789            &mut [
1790                0.0,
1791                -0.0008355098973052918,
1792                -1.7444619952659748e-07,
1793                -7.286557795511902e-11,
1794                -3.80472772706481e-14,
1795                -2.2251271876594078e-17,
1796                1.2789961247944744e-20,
1797            ],
1798            &mut [
1799                0.0,
1800                0.00020861391868413911,
1801                4.3547247296823945e-08,
1802                1.515432276542012e-11,
1803                6.645637323698485e-15,
1804                3.3399223952510497e-18,
1805            ],
1806        );
1807        assert_eq!(res2.0, 6.046459990680098e-17);
1808        assert_eq!(res2.1, 0.7094236375834774);
1809        assert_eq!(res2.2, 0.7047823085448635);
1810        assert_eq!(res2.3, 0.024682338906797385);
1811        assert_eq!(res2.4, -0.02467983282954624);
1812        assert_eq!(res2.5, 0.9996954065371639);
1813        assert_eq!(res2.6, -0.0);
1814        assert_eq!(res2.7, 1.0);
1815        assert_relative_eq!(res2.8, 0.0008355096040059597, epsilon = 1e-18);
1816        assert_eq!(res2.9, -5.870849152149326e-05);
1817        assert_eq!(res2.10, 0.03490027216297455);
1818    }
1819
1820    #[test]
1821    fn test_lengths() {
1822        // Results taken from the python implementation
1823        let geod = Geodesic::wgs84();
1824        let mut c1a = vec![0.0, 1.0, 2.0, 3.0, 4.0, 5.0, 6.0];
1825        let mut c2a = vec![0.0, 1.0, 2.0, 3.0, 4.0, 5.0, 6.0];
1826        let res1 = geod._Lengths(
1827            0.0008355095326524276,
1828            0.024682339962725352,
1829            -0.024679833885152578,
1830            0.9996954065111039,
1831            1.0000010195104125,
1832            -0.0,
1833            1.0,
1834            1.0,
1835            0.9998487145115275,
1836            1.0,
1837            4101,
1838            &mut c1a,
1839            &mut c2a,
1840        );
1841        assert!(res1.0.is_nan());
1842        assert_eq!(res1.1, 0.024679842274314294);
1843        assert_eq!(res1.2, 0.0016717180169067588);
1844        assert!(res1.3.is_nan());
1845        assert!(res1.4.is_nan());
1846
1847        let res2 = geod._Lengths(
1848            0.0008355096040059597,
1849            0.024682338906797385,
1850            -0.02467983282954624,
1851            0.9996954065371639,
1852            1.0000010195104125,
1853            -0.0,
1854            1.0,
1855            1.0,
1856            0.9998487145115275,
1857            1.0,
1858            4101,
1859            &mut [
1860                0.0,
1861                -0.00041775465696698233,
1862                -4.362974596862037e-08,
1863                -1.2151022357848552e-11,
1864                -4.7588881620421004e-15,
1865                -2.226614930167366e-18,
1866                -1.1627237498131586e-21,
1867            ],
1868            &mut [
1869                0.0,
1870                -0.0008355098973052918,
1871                -1.7444619952659748e-07,
1872                -7.286557795511902e-11,
1873                -3.80472772706481e-14,
1874                -2.2251271876594078e-17,
1875                1.2789961247944744e-20,
1876            ],
1877        );
1878        assert!(res2.0.is_nan());
1879        assert_eq!(res2.1, 0.02467984121870759);
1880        assert_eq!(res2.2, 0.0016717181597332804);
1881        assert!(res2.3.is_nan());
1882        assert!(res2.4.is_nan());
1883
1884        let res3 = geod._Lengths(
1885            0.0008355096040059597,
1886            0.024682338906797385,
1887            -0.02467983282954624,
1888            0.9996954065371639,
1889            1.0000010195104125,
1890            -0.0,
1891            1.0,
1892            1.0,
1893            0.9998487145115275,
1894            1.0,
1895            1920,
1896            &mut [
1897                0.0,
1898                -0.00041775469264372037,
1899                -4.362975342068502e-08,
1900                -1.215102547098435e-11,
1901                -4.758889787701359e-15,
1902                -2.2266158809456692e-18,
1903                -1.1627243456014359e-21,
1904            ],
1905            &mut [
1906                0.0,
1907                -0.0008355099686589174,
1908                -1.744462293162189e-07,
1909                -7.286559662008413e-11,
1910                -3.804729026574989e-14,
1911                -2.2251281376754273e-17,
1912                1.2789967801615795e-20,
1913            ],
1914        );
1915        assert_eq!(res3.0, 0.024682347295447677);
1916        assert!(res3.1.is_nan());
1917        assert!(res3.2.is_nan());
1918        assert!(res3.3.is_nan());
1919        assert!(res3.4.is_nan());
1920
1921        let res = geod._Lengths(
1922            0.0007122620325664751,
1923            1.405117407023628,
1924            -0.8928657853278468,
1925            0.45032287238256896,
1926            1.0011366173804046,
1927            0.2969032234925426,
1928            0.9549075745221299,
1929            1.0001257451360057,
1930            0.8139459053827204,
1931            0.9811634781422108,
1932            1920,
1933            &mut [
1934                0.0,
1935                -0.0003561309485314716,
1936                -3.170731714689771e-08,
1937                -7.527972480734327e-12,
1938                -2.5133854116682488e-15,
1939                -1.0025061462383107e-18,
1940                -4.462794158625518e-22,
1941            ],
1942            &mut [
1943                0.0,
1944                -0.0007122622584701569,
1945                -1.2678416507678478e-07,
1946                -4.514641118748122e-11,
1947                -2.0096353119518367e-14,
1948                -1.0019350865558619e-17,
1949                4.90907357448807e-21,
1950            ],
1951        );
1952        assert_eq!(res.0, 1.4056304412645388);
1953        assert!(res.1.is_nan());
1954        assert!(res.2.is_nan());
1955        assert!(res.3.is_nan());
1956        assert!(res.4.is_nan());
1957    }
1958
1959    #[test]
1960    fn test_goed__C4f() {
1961        let geod = Geodesic::wgs84();
1962        let mut c = vec![1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0];
1963        geod._C4f(0.12, &mut c);
1964        assert_eq!(
1965            c,
1966            [
1967                0.6420952961066771,
1968                0.0023680700061156517,
1969                9.96704067834604e-05,
1970                5.778187189466089e-06,
1971                3.9979026199316593e-07,
1972                3.2140078103714466e-08,
1973                7.0
1974            ]
1975        );
1976    }
1977
1978    #[test]
1979    fn test_goed__C3f() {
1980        let geod = Geodesic::wgs84();
1981        let mut c = vec![1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0];
1982        geod._C3f(0.12, &mut c);
1983
1984        assert_eq!(
1985            c,
1986            [
1987                1.0,
1988                0.031839442894193756,
1989                0.0009839921354137713,
1990                5.0055242248766214e-05,
1991                3.1656788204092044e-06,
1992                2.0412e-07,
1993                7.0
1994            ]
1995        );
1996    }
1997
1998    #[test]
1999    fn test_goed__A3f() {
2000        let geod = Geodesic::wgs84();
2001        assert_eq!(geod._A3f(0.12), 0.9363788874000158);
2002    }
2003
2004    #[test]
2005    fn test_geod_init() {
2006        // Check that after the init the variables are correctly set.
2007        // Actual values are taken from the python implementation
2008        let geod = Geodesic::wgs84();
2009        assert_eq!(geod.a, 6378137.0, "geod.a wrong");
2010        assert_eq!(geod.f, 0.0033528106647474805, "geod.f wrong");
2011        assert_eq!(geod._f1, 0.9966471893352525, "geod._f1 wrong");
2012        assert_eq!(geod._e2, 0.0066943799901413165, "geod._e2 wrong");
2013        assert_eq!(geod._ep2, 0.006739496742276434, "geod._ep2 wrong");
2014        assert_eq!(geod._n, 0.0016792203863837047, "geod._n wrong");
2015        assert_eq!(geod._b, 6356752.314245179, "geod._b wrong");
2016        assert_eq!(geod._c2, 40589732499314.76, "geod._c2 wrong");
2017        assert_eq!(geod._etol2, 3.6424611488788524e-08, "geod._etol2 wrong");
2018        assert_eq!(
2019            geod._A3x,
2020            [
2021                -0.0234375,
2022                -0.046927475637074494,
2023                -0.06281503005876607,
2024                -0.2502088451303832,
2025                -0.49916038980680816,
2026                1.0
2027            ],
2028            "geod._A3x wrong"
2029        );
2030
2031        assert_eq!(
2032            geod._C3x,
2033            [
2034                0.0234375,
2035                0.03908873781853724,
2036                0.04695366939653196,
2037                0.12499964752736174,
2038                0.24958019490340408,
2039                0.01953125,
2040                0.02345061890926862,
2041                0.046822392185686165,
2042                0.062342661206936094,
2043                0.013671875,
2044                0.023393770302437927,
2045                0.025963026642854565,
2046                0.013671875,
2047                0.01362595881755982,
2048                0.008203125
2049            ],
2050            "geod._C3x wrong"
2051        );
2052        assert_eq!(
2053            geod._C4x,
2054            [
2055                0.00646020646020646,
2056                0.0035037627212872787,
2057                0.034742279454780166,
2058                -0.01921732223244865,
2059                -0.19923321555984239,
2060                0.6662190894642603,
2061                0.000111000111000111,
2062                0.003426620602971002,
2063                -0.009510765372597735,
2064                -0.01893413691235592,
2065                0.0221370239510936,
2066                0.0007459207459207459,
2067                -0.004142006291321442,
2068                -0.00504225176309005,
2069                0.007584982177746079,
2070                -0.0021565735851450138,
2071                -0.001962613370670692,
2072                0.0036104265913438913,
2073                -0.0009472009472009472,
2074                0.0020416649913317735,
2075                0.0012916376552740189
2076            ],
2077            "geod._C4x wrong"
2078        );
2079    }
2080
2081    // The test_std_geodesic_* tests below are based on Karney's GeodSolve unit
2082    // tests, found in many geographiclib variants.
2083    // The versions below are mostly adapted from their Java counterparts,
2084    // which use a testing structure more similar to Rust than do the C++ versions.
2085    // Note that the Java tests often incorporate more than one of the C++ tests,
2086    // and take their name from the lowest-numbered test in the set.
2087    // These tests use that convention as well.
2088
2089    #[test]
2090    fn test_std_geodesic_geodsolve0() {
2091        let geod = Geodesic::wgs84();
2092        let (s12, azi1, azi2, _a12) = geod.inverse(40.6, -73.8, 49.01666667, 2.55);
2093        assert_relative_eq!(azi1, 53.47022, epsilon = 0.5e-5);
2094        assert_relative_eq!(azi2, 111.59367, epsilon = 0.5e-5);
2095        assert_relative_eq!(s12, 5853226.0, epsilon = 0.5);
2096    }
2097
2098    #[test]
2099    fn test_std_geodesic_geodsolve1() {
2100        let geod = Geodesic::wgs84();
2101        let (lat2, lon2, azi2) = geod.direct(40.63972222, -73.77888889, 53.5, 5850e3);
2102        assert_relative_eq!(lat2, 49.01467, epsilon = 0.5e-5);
2103        assert_relative_eq!(lon2, 2.56106, epsilon = 0.5e-5);
2104        assert_relative_eq!(azi2, 111.62947, epsilon = 0.5e-5);
2105    }
2106
2107    #[test]
2108    fn test_std_geodesic_geodsolve2() {
2109        // Check fix for antipodal prolate bug found 2010-09-04
2110        let geod = Geodesic::new(6.4e6, -1f64 / 150.0);
2111        let (s12, azi1, azi2, _a12) = geod.inverse(0.07476, 0.0, -0.07476, 180.0);
2112        assert_relative_eq!(azi1, 90.00078, epsilon = 0.5e-5);
2113        assert_relative_eq!(azi2, 90.00078, epsilon = 0.5e-5);
2114        assert_relative_eq!(s12, 20106193.0, epsilon = 0.5);
2115        let (s12, azi1, azi2, _a12) = geod.inverse(0.1, 0.0, -0.1, 180.0);
2116        assert_relative_eq!(azi1, 90.00105, epsilon = 0.5e-5);
2117        assert_relative_eq!(azi2, 90.00105, epsilon = 0.5e-5);
2118        assert_relative_eq!(s12, 20106193.0, epsilon = 0.5);
2119    }
2120
2121    #[test]
2122    fn test_std_geodesic_geodsolve4() {
2123        // Check fix for short line bug found 2010-05-21
2124        let geod = Geodesic::wgs84();
2125        let s12: f64 = geod.inverse(36.493349428792, 0.0, 36.49334942879201, 0.0000008);
2126        assert_relative_eq!(s12, 0.072, epsilon = 0.5e-3);
2127    }
2128
2129    #[test]
2130    fn test_std_geodesic_geodsolve5() {
2131        // Check fix for point2=pole bug found 2010-05-03
2132        let geod = Geodesic::wgs84();
2133        let (lat2, lon2, azi2) = geod.direct(0.01777745589997, 30.0, 0.0, 10e6);
2134        assert_relative_eq!(lat2, 90.0, epsilon = 0.5e-5);
2135        if lon2 < 0.0 {
2136            assert_relative_eq!(lon2, -150.0, epsilon = 0.5e-5);
2137            assert_relative_eq!(azi2.abs(), 180.0, epsilon = 0.5e-5);
2138        } else {
2139            assert_relative_eq!(lon2, 30.0, epsilon = 0.5e-5);
2140            assert_relative_eq!(azi2, 0.0, epsilon = 0.5e-5);
2141        }
2142    }
2143
2144    #[test]
2145    fn test_std_geodesic_geodsolve6() {
2146        // Check fix for volatile sbet12a bug found 2011-06-25 (gcc 4.4.4
2147        // x86 -O3).  Found again on 2012-03-27 with tdm-mingw32 (g++ 4.6.1).
2148        let geod = Geodesic::wgs84();
2149        let s12: f64 = geod.inverse(
2150            88.202499451857,
2151            0.0,
2152            -88.202499451857,
2153            179.981022032992859592,
2154        );
2155        assert_relative_eq!(s12, 20003898.214, epsilon = 0.5e-3);
2156        let s12: f64 = geod.inverse(
2157            89.333123580033,
2158            0.0,
2159            -89.333123580032997687,
2160            179.99295812360148422,
2161        );
2162        assert_relative_eq!(s12, 20003926.881, epsilon = 0.5e-3);
2163    }
2164
2165    #[test]
2166    fn test_std_geodesic_geodsolve9() {
2167        // Check fix for volatile x bug found 2011-06-25 (gcc 4.4.4 x86 -O3)
2168        let geod = Geodesic::wgs84();
2169        let s12: f64 = geod.inverse(
2170            56.320923501171,
2171            0.0,
2172            -56.320923501171,
2173            179.664747671772880215,
2174        );
2175        assert_relative_eq!(s12, 19993558.287, epsilon = 0.5e-3);
2176    }
2177
2178    #[test]
2179    fn test_std_geodesic_geodsolve10() {
2180        // Check fix for adjust tol1_ bug found 2011-06-25 (Visual Studio
2181        // 10 rel + debug)
2182        let geod = Geodesic::wgs84();
2183        let s12: f64 = geod.inverse(
2184            52.784459512564,
2185            0.0,
2186            -52.784459512563990912,
2187            179.634407464943777557,
2188        );
2189        assert_relative_eq!(s12, 19991596.095, epsilon = 0.5e-3);
2190    }
2191
2192    #[test]
2193    fn test_std_geodesic_geodsolve11() {
2194        // Check fix for bet2 = -bet1 bug found 2011-06-25 (Visual Studio
2195        // 10 rel + debug)
2196        let geod = Geodesic::wgs84();
2197        let s12: f64 = geod.inverse(
2198            48.522876735459,
2199            0.0,
2200            -48.52287673545898293,
2201            179.599720456223079643,
2202        );
2203        assert_relative_eq!(s12, 19989144.774, epsilon = 0.5e-3);
2204    }
2205
2206    #[test]
2207    fn test_std_geodesic_geodsolve12() {
2208        // Check fix for inverse geodesics on extreme prolate/oblate
2209        // ellipsoids Reported 2012-08-29 Stefan Guenther
2210        // <stefan.gunther@embl.de>; fixed 2012-10-07
2211        let geod = Geodesic::new(89.8, -1.83);
2212        let (s12, azi1, azi2, _a12) = geod.inverse(0.0, 0.0, -10.0, 160.0);
2213        assert_relative_eq!(azi1, 120.27, epsilon = 1e-2);
2214        assert_relative_eq!(azi2, 105.15, epsilon = 1e-2);
2215        assert_relative_eq!(s12, 266.7, epsilon = 1e-1);
2216    }
2217
2218    #[test]
2219    fn test_std_geodesic_geodsolve14() {
2220        // Check fix for inverse ignoring lon12 = nan
2221        let geod = Geodesic::wgs84();
2222        let (s12, azi1, azi2, _a12) = geod.inverse(0.0, 0.0, 1.0, f64::NAN);
2223        assert!(azi1.is_nan());
2224        assert!(azi2.is_nan());
2225        assert!(s12.is_nan());
2226    }
2227
2228    #[test]
2229    fn test_std_geodesic_geodsolve15() {
2230        // Initial implementation of Math::eatanhe was wrong for e^2 < 0.  This
2231        // checks that this is fixed.
2232        let geod = Geodesic::new(6.4e6, -1f64 / 150.0);
2233        let (_lat2, _lon2, _azi2, _m12, _M12, _M21, S12, _a12) = geod.direct(1.0, 2.0, 3.0, 4.0);
2234        assert_relative_eq!(S12, 23700.0, epsilon = 0.5);
2235    }
2236
2237    #[test]
2238    fn test_std_geodesic_geodsolve17() {
2239        // Check fix for LONG_UNROLL bug found on 2015-05-07
2240        let geod = Geodesic::new(6.4e6, -1f64 / 150.0);
2241        let (_a12, lat2, lon2, azi2, _s12, _m12, _M12, _M21, _S12) = geod._gen_direct(
2242            40.0,
2243            -75.0,
2244            -10.0,
2245            false,
2246            2e7,
2247            caps::STANDARD | caps::LONG_UNROLL,
2248        );
2249        assert_relative_eq!(lat2, -39.0, epsilon = 1.0);
2250        assert_relative_eq!(lon2, -254.0, epsilon = 1.0);
2251        assert_relative_eq!(azi2, -170.0, epsilon = 1.0);
2252
2253        let line = GeodesicLine::new(&geod, 40.0, -75.0, -10.0, None, None, None);
2254        let (_a12, lat2, lon2, azi2, _s12, _m12, _M12, _M21, _S12) =
2255            line._gen_position(false, 2e7, caps::STANDARD | caps::LONG_UNROLL);
2256        assert_relative_eq!(lat2, -39.0, epsilon = 1.0);
2257        assert_relative_eq!(lon2, -254.0, epsilon = 1.0);
2258        assert_relative_eq!(azi2, -170.0, epsilon = 1.0);
2259
2260        let (lat2, lon2, azi2) = geod.direct(40.0, -75.0, -10.0, 2e7);
2261        assert_relative_eq!(lat2, -39.0, epsilon = 1.0);
2262        assert_relative_eq!(lon2, 105.0, epsilon = 1.0);
2263        assert_relative_eq!(azi2, -170.0, epsilon = 1.0);
2264
2265        let (_a12, lat2, lon2, azi2, _s12, _m12, _M12, _M21, _S12) =
2266            line._gen_position(false, 2e7, caps::STANDARD);
2267        assert_relative_eq!(lat2, -39.0, epsilon = 1.0);
2268        assert_relative_eq!(lon2, 105.0, epsilon = 1.0);
2269        assert_relative_eq!(azi2, -170.0, epsilon = 1.0);
2270    }
2271
2272    #[test]
2273    fn test_std_geodesic_geodsolve26() {
2274        // Check 0/0 problem with area calculation on sphere 2015-09-08
2275        let geod = Geodesic::new(6.4e6, 0.0);
2276        let (_a12, _s12, _salp1, _calp1, _salp2, _calp2, _m12, _M12, _M21, S12) =
2277            geod._gen_inverse(1.0, 2.0, 3.0, 4.0, caps::AREA);
2278        assert_relative_eq!(S12, 49911046115.0, epsilon = 0.5);
2279    }
2280
2281    #[test]
2282    fn test_std_geodesic_geodsolve28() {
2283        // Check for bad placement of assignment of r.a12 with |f| > 0.01 (bug in
2284        // Java implementation fixed on 2015-05-19).
2285        let geod = Geodesic::new(6.4e6, 0.1);
2286        let (a12, _lat2, _lon2, _azi2, _s12, _m12, _M12, _M21, _S12) =
2287            geod._gen_direct(1.0, 2.0, 10.0, false, 5e6, caps::STANDARD);
2288        assert_relative_eq!(a12, 48.55570690, epsilon = 0.5e-8);
2289    }
2290
2291    #[test]
2292    fn test_std_geodesic_geodsolve29() {
2293        // Check longitude unrolling with inverse calculation 2015-09-16
2294        let geod = Geodesic::wgs84();
2295        let (_a12, s12, _salp1, _calp1, _salp2, _calp2, _m12, _M12, _M21, _S12) =
2296            geod._gen_inverse(0.0, 539.0, 0.0, 181.0, caps::STANDARD);
2297        // Note: This is also supposed to check adjusted longitudes, but geographiclib-rs
2298        //       doesn't seem to support that as of 2021/01/18.
2299        // assert_relative_eq!(lon1, 179, epsilon = 1e-10);
2300        // assert_relative_eq!(lon2, -179, epsilon = 1e-10);
2301        assert_relative_eq!(s12, 222639.0, epsilon = 0.5);
2302        let (_a12, s12, _salp1, _calp1, _salp2, _calp2, _m12, _M12, _M21, _S12) =
2303            geod._gen_inverse(0.0, 539.0, 0.0, 181.0, caps::STANDARD | caps::LONG_UNROLL);
2304        // assert_relative_eq!(lon1, 539, epsilon = 1e-10);
2305        // assert_relative_eq!(lon2, 541, epsilon = 1e-10);
2306        assert_relative_eq!(s12, 222639.0, epsilon = 0.5);
2307    }
2308
2309    #[test]
2310    fn test_std_geodesic_geodsolve33() {
2311        // Check max(-0.0,+0.0) issues 2015-08-22 (triggered by bugs in Octave --
2312        // sind(-0.0) = +0.0 -- and in some version of Visual Studio --
2313        // fmod(-0.0, 360.0) = +0.0.
2314        let geod = Geodesic::wgs84();
2315        let (s12, azi1, azi2, _a12) = geod.inverse(0.0, 0.0, 0.0, 179.0);
2316        assert_relative_eq!(azi1, 90.0, epsilon = 0.5e-5);
2317        assert_relative_eq!(azi2, 90.0, epsilon = 0.5e-5);
2318        assert_relative_eq!(s12, 19926189.0, epsilon = 0.5);
2319        let (s12, azi1, azi2, _a12) = geod.inverse(0.0, 0.0, 0.0, 179.5);
2320        assert_relative_eq!(azi1, 55.96650, epsilon = 0.5e-5);
2321        assert_relative_eq!(azi2, 124.03350, epsilon = 0.5e-5);
2322        assert_relative_eq!(s12, 19980862.0, epsilon = 0.5);
2323        let (s12, azi1, azi2, _a12) = geod.inverse(0.0, 0.0, 0.0, 180.0);
2324        assert_relative_eq!(azi1, 0.0, epsilon = 0.5e-5);
2325        assert_relative_eq!(azi2.abs(), 180.0, epsilon = 0.5e-5);
2326        assert_relative_eq!(s12, 20003931.0, epsilon = 0.5);
2327        let (s12, azi1, azi2, _a12) = geod.inverse(0.0, 0.0, 1.0, 180.0);
2328        assert_relative_eq!(azi1, 0.0, epsilon = 0.5e-5);
2329        assert_relative_eq!(azi2.abs(), 180.0, epsilon = 0.5e-5);
2330        assert_relative_eq!(s12, 19893357.0, epsilon = 0.5);
2331
2332        let geod = Geodesic::new(6.4e6, 0.0);
2333        let (s12, azi1, azi2, _a12) = geod.inverse(0.0, 0.0, 0.0, 179.0);
2334        assert_relative_eq!(azi1, 90.0, epsilon = 0.5e-5);
2335        assert_relative_eq!(azi2, 90.0, epsilon = 0.5e-5);
2336        assert_relative_eq!(s12, 19994492.0, epsilon = 0.5);
2337        let (s12, azi1, azi2, _a12) = geod.inverse(0.0, 0.0, 0.0, 180.0);
2338        assert_relative_eq!(azi1, 0.0, epsilon = 0.5e-5);
2339        assert_relative_eq!(azi2.abs(), 180.0, epsilon = 0.5e-5);
2340        assert_relative_eq!(s12, 20106193.0, epsilon = 0.5);
2341        let (s12, azi1, azi2, _a12) = geod.inverse(0.0, 0.0, 1.0, 180.0);
2342        assert_relative_eq!(azi1, 0.0, epsilon = 0.5e-5);
2343        assert_relative_eq!(azi2.abs(), 180.0, epsilon = 0.5e-5);
2344        assert_relative_eq!(s12, 19994492.0, epsilon = 0.5);
2345
2346        let geod = Geodesic::new(6.4e6, -1.0 / 300.0);
2347        let (s12, azi1, azi2, _a12) = geod.inverse(0.0, 0.0, 0.0, 179.0);
2348        assert_relative_eq!(azi1, 90.0, epsilon = 0.5e-5);
2349        assert_relative_eq!(azi2, 90.0, epsilon = 0.5e-5);
2350        assert_relative_eq!(s12, 19994492.0, epsilon = 0.5);
2351        let (s12, azi1, azi2, _a12) = geod.inverse(0.0, 0.0, 0.0, 180.0);
2352        assert_relative_eq!(azi1, 90.0, epsilon = 0.5e-5);
2353        assert_relative_eq!(azi2, 90.0, epsilon = 0.5e-5);
2354        assert_relative_eq!(s12, 20106193.0, epsilon = 0.5);
2355        let (s12, azi1, azi2, _a12) = geod.inverse(0.0, 0.0, 0.5, 180.0);
2356        assert_relative_eq!(azi1, 33.02493, epsilon = 0.5e-5);
2357        assert_relative_eq!(azi2, 146.97364, epsilon = 0.5e-5);
2358        assert_relative_eq!(s12, 20082617.0, epsilon = 0.5);
2359        let (s12, azi1, azi2, _a12) = geod.inverse(0.0, 0.0, 1.0, 180.0);
2360        assert_relative_eq!(azi1, 0.0, epsilon = 0.5e-5);
2361        assert_relative_eq!(azi2.abs(), 180.0, epsilon = 0.5e-5);
2362        assert_relative_eq!(s12, 20027270.0, epsilon = 0.5);
2363    }
2364
2365    #[test]
2366    fn test_std_geodesic_geodsolve55() {
2367        // Check fix for nan + point on equator or pole not returning all nans in
2368        // Geodesic::Inverse, found 2015-09-23.
2369        let geod = Geodesic::wgs84();
2370        let (s12, azi1, azi2, _a12) = geod.inverse(f64::NAN, 0.0, 0.0, 90.0);
2371        assert!(azi1.is_nan());
2372        assert!(azi2.is_nan());
2373        assert!(s12.is_nan());
2374        let (s12, azi1, azi2, _a12) = geod.inverse(f64::NAN, 0.0, 90.0, 3.0);
2375        assert!(azi1.is_nan());
2376        assert!(azi2.is_nan());
2377        assert!(s12.is_nan());
2378    }
2379
2380    #[test]
2381    fn test_std_geodesic_geodsolve59() {
2382        // Check for points close with longitudes close to 180 deg apart.
2383        let geod = Geodesic::wgs84();
2384        let (s12, azi1, azi2, _a12) = geod.inverse(5.0, 0.00000000000001, 10.0, 180.0);
2385        assert_relative_eq!(azi1, 0.000000000000035, epsilon = 1.5e-14);
2386        assert_relative_eq!(azi2, 179.99999999999996, epsilon = 1.5e-14);
2387        assert_relative_eq!(s12, 18345191.174332713, epsilon = 5e-9);
2388    }
2389
2390    #[test]
2391    fn test_std_geodesic_geodsolve61() {
2392        // Make sure small negative azimuths are west-going
2393        let geod = Geodesic::wgs84();
2394        let (_a12, lat2, lon2, azi2, _s12, _m12, _M12, _M21, _S12) = geod._gen_direct(
2395            45.0,
2396            0.0,
2397            -0.000000000000000003,
2398            false,
2399            1e7,
2400            caps::STANDARD | caps::LONG_UNROLL,
2401        );
2402        assert_relative_eq!(lat2, 45.30632, epsilon = 0.5e-5);
2403        assert_relative_eq!(lon2, -180.0, epsilon = 0.5e-5);
2404        assert_relative_eq!(azi2.abs(), 180.0, epsilon = 0.5e-5);
2405        // geographiclib-rs does not appear to support Geodesic.inverse_line or
2406        // or GeodesicLine.position as of 2021/01/18.
2407        // let line = geod.inverse_line(45, 0, 80, -0.000000000000000003);
2408        // let res = line.position(1e7, caps::STANDARD | caps::LONG_UNROLL);
2409        // assert_relative_eq!(lat2, 45.30632, epsilon = 0.5e-5);
2410        // assert_relative_eq!(lon2, -180, epsilon = 0.5e-5);
2411        // assert_relative_eq!(azi2.abs(), 180, epsilon = 0.5e-5);
2412    }
2413
2414    // #[test]
2415    // fn test_std_geodesic_geodsolve65() {
2416    //     // Check for bug in east-going check in GeodesicLine (needed to check for
2417    //     // sign of 0) and sign error in area calculation due to a bogus override
2418    //     // of the code for alp12.  Found/fixed on 2015-12-19.
2419    //     // These tests rely on Geodesic.inverse_line, which is not supported by
2420    //     // geographiclib-rs as of 2021/01/18.
2421    // }
2422
2423    // #[test]
2424    // fn test_std_geodesic_geodsolve69() {
2425    //     // Check for InverseLine if line is slightly west of S and that s13 is
2426    //     // correctly set.
2427    //     // These tests rely on Geodesic.inverse_line, which is not supported by
2428    //     // geographiclib-rs as of 2021/01/18.
2429    // }
2430
2431    // #[test]
2432    // fn test_std_geodesic_geodsolve71() {
2433    //     // Check that DirectLine sets s13.
2434    //     // These tests rely on Geodesic.direct_line, which is not supported by
2435    //     // geographiclib-rs as of 2021/01/18.
2436    // }
2437
2438    #[test]
2439    fn test_std_geodesic_geodsolve73() {
2440        // Check for backwards from the pole bug reported by Anon on 2016-02-13.
2441        // This only affected the Java implementation.  It was introduced in Java
2442        // version 1.44 and fixed in 1.46-SNAPSHOT on 2016-01-17.
2443        // Also the + sign on azi2 is a check on the normalizing of azimuths
2444        // (converting -0.0 to +0.0).
2445        let geod = Geodesic::wgs84();
2446        let (lat2, lon2, azi2) = geod.direct(90.0, 10.0, 180.0, -1e6);
2447        assert_relative_eq!(lat2, 81.04623, epsilon = 0.5e-5);
2448        assert_relative_eq!(lon2, -170.0, epsilon = 0.5e-5);
2449        assert_relative_eq!(azi2, 0.0, epsilon = 0.5e-5);
2450        assert!(azi2.is_sign_positive());
2451    }
2452
2453    #[test]
2454    fn test_std_geodesic_geodsolve74() {
2455        // Check fix for inaccurate areas, bug introduced in v1.46, fixed
2456        // 2015-10-16.
2457        let geod = Geodesic::wgs84();
2458        let (a12, s12, azi1, azi2, m12, M12, M21, S12) =
2459            geod._gen_inverse_azi(54.1589, 15.3872, 54.1591, 15.3877, caps::ALL);
2460        assert_relative_eq!(azi1, 55.723110355, epsilon = 5e-9);
2461        assert_relative_eq!(azi2, 55.723515675, epsilon = 5e-9);
2462        assert_relative_eq!(s12, 39.527686385, epsilon = 5e-9);
2463        assert_relative_eq!(a12, 0.000355495, epsilon = 5e-9);
2464        assert_relative_eq!(m12, 39.527686385, epsilon = 5e-9);
2465        assert_relative_eq!(M12, 0.999999995, epsilon = 5e-9);
2466        assert_relative_eq!(M21, 0.999999995, epsilon = 5e-9);
2467        assert_relative_eq!(S12, 286698586.30197, epsilon = 5e-4);
2468    }
2469
2470    #[test]
2471    fn test_std_geodesic_geodsolve76() {
2472        // The distance from Wellington and Salamanca (a classic failure of
2473        // Vincenty)
2474        let geod = Geodesic::wgs84();
2475        let (s12, azi1, azi2, _a12) = geod.inverse(
2476            -(41.0 + 19.0 / 60.0),
2477            174.0 + 49.0 / 60.0,
2478            40.0 + 58.0 / 60.0,
2479            -(5.0 + 30.0 / 60.0),
2480        );
2481        assert_relative_eq!(azi1, 160.39137649664, epsilon = 0.5e-11);
2482        assert_relative_eq!(azi2, 19.50042925176, epsilon = 0.5e-11);
2483        assert_relative_eq!(s12, 19960543.857179, epsilon = 0.5e-6);
2484    }
2485
2486    #[test]
2487    fn test_std_geodesic_geodsolve78() {
2488        // An example where the NGS calculator fails to converge
2489        let geod = Geodesic::wgs84();
2490        let (s12, azi1, azi2, _a12) = geod.inverse(27.2, 0.0, -27.1, 179.5);
2491        assert_relative_eq!(azi1, 45.82468716758, epsilon = 0.5e-11);
2492        assert_relative_eq!(azi2, 134.22776532670, epsilon = 0.5e-11);
2493        assert_relative_eq!(s12, 19974354.765767, epsilon = 0.5e-6);
2494    }
2495
2496    #[test]
2497    fn test_std_geodesic_geodsolve80() {
2498        // Some tests to add code coverage: computing scale in special cases + zero
2499        // length geodesic (includes GeodSolve80 - GeodSolve83).
2500        let geod = Geodesic::wgs84();
2501        let (_a12, _s12, _salp1, _calp1, _salp2, _calp2, _m12, M12, M21, _S12) =
2502            geod._gen_inverse(0.0, 0.0, 0.0, 90.0, caps::GEODESICSCALE);
2503        assert_relative_eq!(M12, -0.00528427534, epsilon = 0.5e-10);
2504        assert_relative_eq!(M21, -0.00528427534, epsilon = 0.5e-10);
2505
2506        let (_a12, _s12, _salp1, _calp1, _salp2, _calp2, _m12, M12, M21, _S12) =
2507            geod._gen_inverse(0.0, 0.0, 1e-6, 1e-6, caps::GEODESICSCALE);
2508        assert_relative_eq!(M12, 1.0, epsilon = 0.5e-10);
2509        assert_relative_eq!(M21, 1.0, epsilon = 0.5e-10);
2510
2511        let (a12, s12, azi1, azi2, m12, M12, M21, S12) =
2512            geod._gen_inverse_azi(20.001, 0.0, 20.001, 0.0, caps::ALL);
2513        assert_relative_eq!(a12, 0.0, epsilon = 1e-13);
2514        assert_relative_eq!(s12, 0.0, epsilon = 1e-8);
2515        assert_relative_eq!(azi1, 180.0, epsilon = 1e-13);
2516        assert_relative_eq!(azi2, 180.0, epsilon = 1e-13);
2517        assert_relative_eq!(m12, 0.0, epsilon = 1e-8);
2518        assert_relative_eq!(M12, 1.0, epsilon = 1e-15);
2519        assert_relative_eq!(M21, 1.0, epsilon = 1e-15);
2520        assert_relative_eq!(S12, 0.0, epsilon = 1e-10);
2521
2522        let (a12, s12, azi1, azi2, m12, M12, M21, S12) =
2523            geod._gen_inverse_azi(90.0, 0.0, 90.0, 180.0, caps::ALL);
2524        assert_relative_eq!(a12, 0.0, epsilon = 1e-13);
2525        assert_relative_eq!(s12, 0.0, epsilon = 1e-8);
2526        assert_relative_eq!(azi1, 0.0, epsilon = 1e-13);
2527        assert_relative_eq!(azi2, 180.0, epsilon = 1e-13);
2528        assert_relative_eq!(m12, 0.0, epsilon = 1e-8);
2529        assert_relative_eq!(M12, 1.0, epsilon = 1e-15);
2530        assert_relative_eq!(M21, 1.0, epsilon = 1e-15);
2531        assert_relative_eq!(S12, 127516405431022.0, epsilon = 0.5);
2532
2533        // An incapable line which can't take distance as input
2534        let line = GeodesicLine::new(&geod, 1.0, 2.0, 90.0, Some(caps::LATITUDE), None, None);
2535        let (a12, _lat2, _lon2, _azi2, _s12, _m12, _M12, _M21, _S12) =
2536            line._gen_position(false, 1000.0, caps::CAP_NONE);
2537        assert!(a12.is_nan());
2538    }
2539
2540    #[test]
2541    fn test_std_geodesic_geodsolve84() {
2542        // Tests for python implementation to check fix for range errors with
2543        // {fmod,sin,cos}(inf) (includes GeodSolve84 - GeodSolve91).
2544        let geod = Geodesic::wgs84();
2545        let (lat2, lon2, azi2) = geod.direct(0.0, 0.0, 90.0, f64::INFINITY);
2546        assert!(lat2.is_nan());
2547        assert!(lon2.is_nan());
2548        assert!(azi2.is_nan());
2549        let (lat2, lon2, azi2) = geod.direct(0.0, 0.0, 90.0, f64::NAN);
2550        assert!(lat2.is_nan());
2551        assert!(lon2.is_nan());
2552        assert!(azi2.is_nan());
2553        let (lat2, lon2, azi2) = geod.direct(0.0, 0.0, f64::INFINITY, 1000.0);
2554        assert!(lat2.is_nan());
2555        assert!(lon2.is_nan());
2556        assert!(azi2.is_nan());
2557        let (lat2, lon2, azi2) = geod.direct(0.0, 0.0, f64::NAN, 1000.0);
2558        assert!(lat2.is_nan());
2559        assert!(lon2.is_nan());
2560        assert!(azi2.is_nan());
2561        let (lat2, lon2, azi2) = geod.direct(0.0, f64::INFINITY, 90.0, 1000.0);
2562        assert_eq!(lat2, 0.0);
2563        assert!(lon2.is_nan());
2564        assert_eq!(azi2, 90.0);
2565        let (lat2, lon2, azi2) = geod.direct(0.0, f64::NAN, 90.0, 1000.0);
2566        assert_eq!(lat2, 0.0);
2567        assert!(lon2.is_nan());
2568        assert_eq!(azi2, 90.0);
2569        let (lat2, lon2, azi2) = geod.direct(f64::INFINITY, 0.0, 90.0, 1000.0);
2570        assert!(lat2.is_nan());
2571        assert!(lon2.is_nan());
2572        assert!(azi2.is_nan());
2573        let (lat2, lon2, azi2) = geod.direct(f64::NAN, 0.0, 90.0, 1000.0);
2574        assert!(lat2.is_nan());
2575        assert!(lon2.is_nan());
2576        assert!(azi2.is_nan());
2577    }
2578
2579    // *_geodtest_* tests are based on Karney's GeodTest*.dat test datasets.
2580    // A description of these files' content can be found at:
2581    //     https://geographiclib.sourceforge.io/html/geodesic.html#testgeod
2582    // Here are some key excerpts...
2583    //    This consists of a set of geodesics for the WGS84 ellipsoid.
2584    //     Each line of the test set gives 10 space delimited numbers
2585    //         latitude at point 1, lat1 (degrees, exact)
2586    //         longitude at point 1, lon1 (degrees, always 0)
2587    //         azimuth at point 1, azi1 (clockwise from north in degrees, exact)
2588    //         latitude at point 2, lat2 (degrees, accurate to 10−18 deg)
2589    //         longitude at point 2, lon2 (degrees, accurate to 10−18 deg)
2590    //         azimuth at point 2, azi2 (degrees, accurate to 10−18 deg)
2591    //         geodesic distance from point 1 to point 2, s12 (meters, exact)
2592    //         arc distance on the auxiliary sphere, a12 (degrees, accurate to 10−18 deg)
2593    //         reduced length of the geodesic, m12 (meters, accurate to 0.1 pm)
2594    //         the area under the geodesic, S12 (m2, accurate to 1 mm2)
2595
2596    static FULL_TEST_PATH: &str = "test_fixtures/test_data_unzipped/GeodTest.dat";
2597    static SHORT_TEST_PATH: &str = "test_fixtures/test_data_unzipped/GeodTest-short.dat";
2598    static BUILTIN_TEST_PATH: &str = "test_fixtures/GeodTest-100.dat";
2599    fn test_input_path() -> &'static str {
2600        if cfg!(feature = "test_full") {
2601            FULL_TEST_PATH
2602        } else if cfg!(feature = "test_short") {
2603            SHORT_TEST_PATH
2604        } else {
2605            BUILTIN_TEST_PATH
2606        }
2607    }
2608
2609    fn geodtest_basic<T>(path: &str, f: T)
2610    where
2611        T: Fn(usize, &(f64, f64, f64, f64, f64, f64, f64, f64, f64, f64)),
2612    {
2613        let dir_base = std::env::current_dir().expect("Failed to determine current directory");
2614        let path_base = dir_base.as_path();
2615        let pathbuf = std::path::Path::new(path_base).join(path);
2616        let path = pathbuf.as_path();
2617        let file = match std::fs::File::open(path) {
2618            Ok(val) => val,
2619            Err(_error) => {
2620                let path_str = path
2621                    .to_str()
2622                    .expect("Failed to convert GeodTest path to string during error reporting");
2623                panic!("Failed to open test input file. Run `script/download-test-data.sh` to download test input to: {}\nFor details see https://geographiclib.sourceforge.io/html/geodesic.html#testgeod", path_str)
2624            }
2625        };
2626        let reader = std::io::BufReader::new(file);
2627        reader.lines().enumerate().for_each(|(i, line)| {
2628            let line_safe = line.expect("Failed to read line");
2629            let items: Vec<f64> = line_safe
2630                .split(' ')
2631                .enumerate()
2632                .map(|(j, item)| match item.parse::<f64>() {
2633                    Ok(parsed) => parsed,
2634                    Err(_error) => {
2635                        panic!("Error parsing item {} on line {}: {}", j + 1, i + 1, item)
2636                    }
2637                })
2638                .collect();
2639            assert_eq!(items.len(), 10);
2640            let tuple = (
2641                items[0], items[1], items[2], items[3], items[4], items[5], items[6], items[7],
2642                items[8], items[9],
2643            );
2644            f(i + 1, &tuple); // report 1-based line number rather than 0-based
2645        });
2646    }
2647
2648    #[test]
2649    fn test_geodtest_geodesic_direct12() {
2650        let g = std::sync::Arc::new(std::sync::Mutex::new(Geodesic::wgs84()));
2651
2652        geodtest_basic(
2653            test_input_path(),
2654            |_line_num, &(lat1, lon1, azi1, lat2, lon2, azi2, s12, a12, m12, S12)| {
2655                let g = g.lock().unwrap();
2656                let (lat2_out, lon2_out, azi2_out, m12_out, _M12_out, _M21_out, S12_out, a12_out) =
2657                    g.direct(lat1, lon1, azi1, s12);
2658                assert_relative_eq!(lat2, lat2_out, epsilon = 1e-13);
2659                assert_relative_eq!(lon2, lon2_out, epsilon = 2e-8);
2660                assert_relative_eq!(azi2, azi2_out, epsilon = 2e-8);
2661                assert_relative_eq!(m12, m12_out, epsilon = 9e-9);
2662                assert_relative_eq!(S12, S12_out, epsilon = 2e4); // Note: unreasonable tolerance
2663                assert_relative_eq!(a12, a12_out, epsilon = 9e-14);
2664            },
2665        );
2666    }
2667
2668    #[test]
2669    fn test_geodtest_geodesic_direct21() {
2670        let g = std::sync::Arc::new(std::sync::Mutex::new(Geodesic::wgs84()));
2671
2672        geodtest_basic(
2673            test_input_path(),
2674            |_line_num, &(lat1, lon1, azi1, lat2, lon2, azi2, s12, a12, m12, S12)| {
2675                let g = g.lock().unwrap();
2676                // Reverse some values for 2->1 instead of 1->2
2677                let (lat1, lon1, azi1, lat2, lon2, azi2, s12, a12, m12, S12) =
2678                    (lat2, lon2, azi2, lat1, lon1, azi1, -s12, -a12, -m12, -S12);
2679                let (lat2_out, lon2_out, azi2_out, m12_out, _M12_out, _M21_out, S12_out, a12_out) =
2680                    g.direct(lat1, lon1, azi1, s12);
2681                assert_relative_eq!(lat2, lat2_out, epsilon = 8e-14);
2682                assert_relative_eq!(lon2, lon2_out, epsilon = 4e-6);
2683                assert_relative_eq!(azi2, azi2_out, epsilon = 4e-6);
2684                assert_relative_eq!(m12, m12_out, epsilon = 1e-8);
2685                assert_relative_eq!(S12, S12_out, epsilon = 3e6); // Note: unreasonable tolerance
2686                assert_relative_eq!(a12, a12_out, epsilon = 9e-14);
2687            },
2688        );
2689    }
2690
2691    #[test]
2692    fn test_geodtest_geodesic_inverse12() {
2693        let g = std::sync::Arc::new(std::sync::Mutex::new(Geodesic::wgs84()));
2694
2695        geodtest_basic(
2696            test_input_path(),
2697            |line_num, &(lat1, lon1, azi1, lat2, lon2, azi2, s12, a12, m12, S12)| {
2698                let g = g.lock().unwrap();
2699                let (s12_out, azi1_out, azi2_out, m12_out, _M12_out, _M21_out, S12_out, a12_out) =
2700                    g.inverse(lat1, lon1, lat2, lon2);
2701                assert_relative_eq!(s12, s12_out, epsilon = 8e-9);
2702                assert_relative_eq!(azi1, azi1_out, epsilon = 2e-2);
2703                assert_relative_eq!(azi2, azi2_out, epsilon = 2e-2);
2704                assert_relative_eq!(m12, m12_out, epsilon = 5e-5);
2705                // Our area calculation differs significantly (~1e7) from the value in GeodTest.dat for
2706                // line 400001, BUT our value also perfectly matches the value returned by GeographicLib
2707                // (C++) 1.51. Here's the problem line, for reference:
2708                // 4.199535552987 0 90 -4.199535552987 179.398106343454992238 90 19970505.608097404994 180 0
2709                if line_num != 400001 {
2710                    assert_relative_eq!(S12, S12_out, epsilon = 3e10); // Note: unreasonable tolerance
2711                }
2712                assert_relative_eq!(a12, a12_out, epsilon = 2e-10);
2713            },
2714        );
2715    }
2716
2717    #[test]
2718    fn test_turnaround() {
2719        let g = Geodesic::wgs84();
2720
2721        let start = (0.0, 0.0);
2722        let destination = (0.0, 1.0);
2723
2724        let (distance, azi1, _, _) = g.inverse(start.0, start.1, destination.0, destination.1);
2725
2726        // Confirm that we've gone due-east
2727        assert_eq!(azi1, 90.0);
2728
2729        // Turn around by adding 180 degrees to the azimuth
2730        let turn_around = azi1 + 180.0;
2731
2732        // Confirm that turn around is due west
2733        assert_eq!(turn_around, 270.0);
2734
2735        // Test that we can turn around and get back to the starting point.
2736        let (lat, lon) = g.direct(destination.0, destination.1, turn_around, distance);
2737        assert_relative_eq!(lat, start.0, epsilon = 1.0e-3);
2738        assert_relative_eq!(lon, start.1, epsilon = 1.0e-3);
2739    }
2740}