Skip to main content

map_3d/
lib.rs

1//! map_3d: geographic coordinates conversions
2
3/// Returns the tuple (x,y,z) of coordinates in the ECEF system
4/// with desired reference frame.
5///
6/// ## Inputs:
7/// - lat = latitude [rad]
8/// - lon = longitude [rad]
9/// - alt = altitude [m]
10/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
11///
12/// ## Outputs:
13/// - x = x ECEF coordinate [m]
14/// - y = y ECEF coordinate [m]
15/// - z = z ECEF coordinate [m]
16pub fn geodetic2ecef(lat: f64, lon: f64, alt: f64, r_ellips: Ellipsoid) -> (f64, f64, f64) {
17    let n = get_radius_normal(lat, r_ellips);
18    let (major, minor, _, _) = r_ellips.parameters();
19
20    let x = (n + alt) * lat.cos() * lon.cos();
21    let y = (n + alt) * lat.cos() * lon.sin();
22    let z = (n * (minor / major) * (minor / major) + alt) * lat.sin();
23
24    (x, y, z)
25}
26
27/// Returns the tuple (azimuth,elevation,slant range) of coordinates in the AER system
28///
29/// ## Inputs:
30/// - lat = input latitude [rad]
31/// - lon = input longitude [rad]
32/// - alt = input altitude [m]
33/// - lat0 = reference latitude [rad]
34/// - lon0 = reference longitude [rad]
35/// - alt0 = reference altitude [m]
36/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
37///
38/// ## Outputs:
39/// - az = azimuth angle [rad] of input geodetic location from reference geodetic location
40/// - el = elevation angle [rad] of input geodetic location from reference geodetic location
41/// - slant_range = slant range [m] of input geodetic location from reference geodetic location
42pub fn geodetic2aer(
43    lat: f64,
44    lon: f64,
45    alt: f64,
46    lat0: f64,
47    lon0: f64,
48    alt0: f64,
49    r_ellips: Ellipsoid,
50) -> (f64, f64, f64) {
51    let (e, n, u) = geodetic2enu(lat, lon, alt, lat0, lon0, alt0, r_ellips);
52    let (az, el, slant_range) = enu2aer(e, n, u);
53
54    (az, el, slant_range)
55}
56
57/// Returns the tuple (east,north,up) of coordinates in the ENU system
58///
59///  ## Inputs:
60/// - lat = input latitude [rad]
61/// - lon = input longitude [rad]
62/// - alt = input altitude [m]
63/// - lat0 = reference latitude [rad]
64/// - lon0 = reference longitude [rad]
65/// - alt0 = reference altitude [m]
66/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
67///
68/// ## Outputs:
69/// - e = east coordinate [m] of input geodetic location from reference geodetic location
70/// - n = north coordinate [m] of input geodetic location from reference geodetic location
71/// - u = up coordinate [m] of input geodetic location from reference geodetic location
72pub fn geodetic2enu(
73    lat: f64,
74    lon: f64,
75    alt: f64,
76    lat0: f64,
77    lon0: f64,
78    alt0: f64,
79    r_ellips: Ellipsoid,
80) -> (f64, f64, f64) {
81    let (x1, y1, z1) = geodetic2ecef(lat, lon, alt, r_ellips);
82    let (x2, y2, z2) = geodetic2ecef(lat0, lon0, alt0, r_ellips);
83
84    let (e, n, u) = uvw2enu(x1 - x2, y1 - y2, z1 - z2, lat0, lon0);
85
86    (e, n, u)
87}
88
89/// Returns the tuple (north,east,down) of coordinates in the NED system
90///
91///  ## Inputs:
92/// - lat = input latitude [rad]
93/// - lon = input longitude [rad]
94/// - alt = input altitude [m]
95/// - lat0 = reference latitude [rad]
96/// - lon0 = reference longitude [rad]
97/// - alt0 = reference altitude [m]
98/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
99///
100/// ## Outputs:
101/// - n = north coordinate [m] of input geodetic location from reference geodetic location
102/// - e = east coordinate [m] of input geodetic location from reference geodetic location
103/// - d = down coordinate [m] of input geodetic location from reference geodetic location
104pub fn geodetic2ned(
105    lat: f64,
106    lon: f64,
107    alt: f64,
108    lat0: f64,
109    lon0: f64,
110    alt0: f64,
111    r_ellips: Ellipsoid,
112) -> (f64, f64, f64) {
113    let enu = geodetic2enu(lat, lon, alt, lat0, lon0, alt0, r_ellips);
114    (enu.1, enu.0, -enu.2)
115}
116
117/// Returns the tuple (x,y,z) of coordinates in the ECEF system
118///
119/// ## Inputs:
120/// - az = azimuth angle [rad]
121/// - el = elevation angle [rad]
122/// - slant_range = slant range [m]
123/// - lat0 = reference latitude [rad]
124/// - lon0 = reference longitude [rad]
125/// - alt0 = reference altitude [m]
126/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
127///
128/// ## Outputs:
129/// - x = x ECEF coordinate [m]
130/// - y = y ECEF coordinate [m]
131/// - z = z ECEF coordinate [m]
132pub fn aer2ecef(
133    az: f64,
134    el: f64,
135    slant_range: f64,
136    lat0: f64,
137    lon0: f64,
138    alt0: f64,
139    r_ellips: Ellipsoid,
140) -> (f64, f64, f64) {
141    let (x0, y0, z0) = geodetic2ecef(lat0, lon0, alt0, r_ellips);
142    let (e, n, u) = aer2enu(az, el, slant_range);
143    let (dx, dy, dz) = enu2uvw(e, n, u, lat0, lon0);
144    (x0 + dx, y0 + dy, z0 + dz)
145}
146
147/// Returns the tuple (east,north,up) of coordinates in the ENU system
148///
149/// ## Inputs:
150/// - az = azimuth angle [rad]
151/// - el = elevation angle [rad]
152/// - slant_range = slant range [m]
153///
154/// ## Outputs:
155/// - e = east coordinate [m] of input location from reference geodetic location
156/// - n = north coordinate [m] of input location from reference geodetic location
157/// - u = up coordinate [m] of input location from reference geodetic location
158pub fn aer2enu(az: f64, el: f64, slant_range: f64) -> (f64, f64, f64) {
159    let r = slant_range * el.cos();
160    (r * az.sin(), r * az.cos(), slant_range * el.sin())
161}
162
163/// Returns the tuple (x,y,z) of coordinates in the ECI system
164///
165/// ## Inputs:
166/// - az = azimuth angle [rad]
167/// - el = elevation angle [rad]
168/// - slant_range = slant range [m]
169/// - lat0 = reference latitude [rad]
170/// - lon0 = reference longitude [rad]
171/// - alt0 = reference altitude [m]
172/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
173///
174/// ## Outputs:
175/// - x = x ECI coordinate [m]
176/// - y = y ECI coordinate [m]
177/// - z = z ECI coordinate [m]
178pub fn aer2eci(
179    gst: f64,
180    az: f64,
181    el: f64,
182    slant_range: f64,
183    lat0: f64,
184    lon0: f64,
185    alt0: f64,
186    r_ellips: Ellipsoid,
187) -> (f64, f64, f64) {
188    let (x1, y1, z1) = aer2ecef(az, el, slant_range, lat0, lon0, alt0, r_ellips);
189    ecef2eci(gst, x1, y1, z1)
190}
191
192/// Returns the tuple (north,east,down) of coordinates in the NED system
193///
194/// ## Inputs:
195/// - az = azimuth angle [rad]
196/// - el = elevation angle [rad]
197/// - slant_range = slant range [m]
198///
199/// ## Outputs:
200/// - n = north coordinate [m] of input location from reference geodetic location
201/// - e = east coordinate [m] of input location from reference geodetic location
202/// - d = down coordinate [m] of input location from reference geodetic location
203pub fn aer2ned(az: f64, el: f64, slant_range: f64) -> (f64, f64, f64) {
204    let enu = aer2enu(az, el, slant_range);
205    (enu.1, enu.0, -enu.2)
206}
207
208/// Returns the tuple (latitude,longitude,altitude) of coordinates in the Geodetic system
209///
210/// ## Inputs:
211/// - az = azimuth angle [rad]
212/// - el = elevation angle [rad]
213/// - slant_range = slant range [m]
214/// - lat0 = reference latitude [rad]
215/// - lon0 = reference longitude [rad]
216/// - alt0 = reference altitude [m]
217/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
218///
219/// ## Outputs:
220/// - lat = input latitude [rad]
221/// - lon = input longitude [rad]
222/// - alt = input altitude [m]
223pub fn aer2geodetic(
224    az: f64,
225    el: f64,
226    slant_range: f64,
227    lat0: f64,
228    lon0: f64,
229    alt0: f64,
230    r_ellips: Ellipsoid,
231) -> (f64, f64, f64) {
232    let (x, y, z) = aer2ecef(az, el, slant_range, lat0, lon0, alt0, r_ellips);
233    ecef2geodetic(x, y, z, r_ellips)
234}
235
236/// Returns the tuple (u,v,w) of coordinates in the local vector system
237///
238/// ## Inputs:
239/// - e = east coordinate [m] from reference geodetic location
240/// - n = north coordinate [m] from reference geodetic location
241/// - u = up coordinate [m] from reference geodetic location
242/// - lat0 = reference latitude [rad]
243/// - lon0 = reference longitude [rad]
244///
245/// ## Outputs:
246/// - u = tangent vector component
247/// - v = tangent vector component
248/// - w = tangent vector component
249pub fn enu2uvw(e: f64, n: f64, up: f64, lat0: f64, lon0: f64) -> (f64, f64, f64) {
250    let t = lat0.cos() * up - lat0.sin() * n;
251
252    let u = lon0.cos() * t - lon0.sin() * e;
253    let v = lon0.sin() * t + lon0.cos() * e;
254    let w = lat0.sin() * up + lat0.cos() * n;
255    (u, v, w)
256}
257
258/// Returns the tuple (azimuth,elevation,slant range) of coordinates in the AER system
259///
260/// ## Inputs:
261/// - e = east coordinate [m] from reference geodetic location
262/// - n = north coordinate [m] from reference geodetic location
263/// - u = up coordinate [m] from reference geodetic location
264///
265/// ## Outputs:
266/// - az = azimuth angle [rad] of input location from reference geodetic location
267/// - el = elevation angle [rad] of input location from reference geodetic location
268/// - slant_range = slant range [m] of input location from reference geodetic location
269pub fn enu2aer(e: f64, n: f64, u: f64) -> (f64, f64, f64) {
270    let r = (e * e + n * n).sqrt();
271
272    let slant_range = (r * r + u * u).sqrt();
273    let el = u.atan2(r);
274    let az = e.atan2(n).rem_euclid(2.0 * std::f64::consts::PI);
275
276    (az, el, slant_range)
277}
278
279/// Returns the tuple (x,y,z) of coordinates in the ECEF system
280///
281/// ## Inputs:
282/// - e = east coordinate [m] from reference geodetic location
283/// - n = north coordinate [m] from reference geodetic location
284/// - u = up coordinate [m] from reference geodetic location
285/// - lat0 = reference latitude [rad]
286/// - lon0 = reference longitude [rad]
287/// - alt0 = reference altitude [m]
288/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
289///
290/// ## Outputs:
291/// - x = x ECEF coordinate [m]
292/// - y = y ECEF coordinate [m]
293/// - z = z ECEF coordinate [m]
294pub fn enu2ecef(
295    e: f64,
296    n: f64,
297    u: f64,
298    lat0: f64,
299    lon0: f64,
300    alt0: f64,
301    r_ellips: Ellipsoid,
302) -> (f64, f64, f64) {
303    let (x0, y0, z0) = geodetic2ecef(lat0, lon0, alt0, r_ellips);
304    let (dx, dy, dz) = enu2uvw(e, n, u, lat0, lon0);
305
306    (x0 + dx, y0 + dy, z0 + dz)
307}
308
309/// Returns the tuple (latitude,longitude,altitude) of coordinates in the Geodetic system
310///
311/// ## Inputs:
312/// - e = east coordinate [m]  from reference geodetic location
313/// - n = north coordinate [m]  from reference geodetic location
314/// - u = up coordinate [m]  from reference geodetic location
315/// - lat0 = reference latitude [rad]
316/// - lon0 = reference longitude [rad]
317/// - alt0 = reference altitude [m]
318/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
319/// ## Outputs:
320/// - lat = latitude [rad]
321/// - lon = longitude [rad]
322/// - alt = altitude [m]
323pub fn enu2geodetic(
324    e: f64,
325    n: f64,
326    u: f64,
327    lat0: f64,
328    lon0: f64,
329    alt0: f64,
330    r_ellips: Ellipsoid,
331) -> (f64, f64, f64) {
332    let (x, y, z) = enu2ecef(e, n, u, lat0, lon0, alt0, r_ellips);
333    let (lat, lon, alt) = ecef2geodetic(x, y, z, r_ellips);
334
335    (lat, lon, alt)
336}
337
338/// Returns the tuple (x,y,z) of coordinates in the ECI system
339///
340/// ## Inputs:
341/// - x = x ECEF coordinate [m]
342/// - y = y ECEF coordinate [m]
343/// - z = z ECEF coordinate [m]
344///
345/// ## Outputs:
346/// - x = x ECI coordinate [m]
347/// - y = y ECI coordinate [m]
348/// - z = z ECI coordinate [m]
349pub fn ecef2eci(gst: f64, x: f64, y: f64, z: f64) -> (f64, f64, f64) {
350    let arr = matmul3(transpose3(r3(gst)), [x, y, z]);
351    (arr[0], arr[1], arr[2])
352}
353
354/// Returns the tuple (latitude,longitude,altitude) of coordinates in the Geodetic system
355/// with desired reference frame. Returns NaN for lat, lon, alt if the input ECEF coordinates are at the center of the Earth.
356///
357/// ## Inputs:
358/// - x = x ECEF coordinate [m]
359/// - y = y ECEF coordinate [m]
360/// - z = z ECEF coordinate [m]
361/// - r_ellips = reference ellipsoid, defaults to WGS84
362///
363/// ## Outputs:
364/// - lat = latitude [rad]
365/// - lon = longitude [rad]
366/// - alt = altitude [m]
367pub fn ecef2geodetic(x: f64, y: f64, z: f64, r_ellips: Ellipsoid) -> (f64, f64, f64) {
368    let (major, minor, _, _) = r_ellips.parameters();
369
370    let q = (x * x + y * y).sqrt();
371
372    // Handle poles case where q is zero to avoid division by zero
373    if q <= f64::EPSILON {
374        if z > 0.0 {
375            return (std::f64::consts::FRAC_PI_2, 0.0, z - minor);
376        }
377        if z < 0.0 {
378            return (-std::f64::consts::FRAC_PI_2, 0.0, -z - minor);
379        }
380        return (f64::NAN, f64::NAN, f64::NAN); // Indeterminate case at the center of the Earth
381    }
382
383    let r = (x * x + y * y + z * z).sqrt();
384    let e = (major * major - minor * minor).sqrt();
385    let var = r * r - e * e;
386    let u = (0.5 * var + 0.5 * (var * var + 4.0 * e * e * z * z).sqrt()).sqrt();
387
388    let hu_e = (u * u + e * e).sqrt();
389
390    // FIX: this should improve stability near poles.
391    let mut beta = (hu_e * z).atan2(u * q);
392    let cos_beta = beta.cos();
393    if cos_beta.abs() > 1e-12 {
394        let sin_beta = beta.sin();
395        let eps = ((minor * u - major * hu_e + e * e) * sin_beta)
396            / (major * hu_e / cos_beta - e * e * cos_beta);
397        beta += eps;
398    }
399
400    let lat = (major / minor * beta.tan()).atan();
401    let lon = y.atan2(x);
402
403    let v1 = z - minor * beta.sin();
404    let v2 = q - major * beta.cos();
405
406    let inside = (x * x / major / major) + (y * y / major / major) + (z * z / minor / minor) < 1.0;
407    let alt = if inside {
408        -(v1 * v1 + v2 * v2).sqrt()
409    } else {
410        (v1 * v1 + v2 * v2).sqrt()
411    };
412
413    (lat, lon, alt)
414}
415
416/// Returns the tuple (east,north,up) of coordinates in the ENU system
417///
418/// ## Inputs:
419/// - x = x ECEF coordinate [m]
420/// - y = y ECEF coordinate [m]
421/// - z = z ECEF coordinate [m]
422/// - lat0 = reference latitude [rad]
423/// - lon0 = reference longitude [rad]
424/// - alt0 = reference altitude [m]
425/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
426///
427/// ## Outputs:
428/// - e = east coordinate [m] of input ECEF location from reference geodetic location
429/// - n = north coordinate [m] of input ECEF location from reference geodetic location
430/// - u = up coordinate [m] of input ECEF location from reference geodetic location
431pub fn ecef2enu(
432    x: f64,
433    y: f64,
434    z: f64,
435    lat0: f64,
436    lon0: f64,
437    alt0: f64,
438    r_ellips: Ellipsoid,
439) -> (f64, f64, f64) {
440    let (x0, y0, z0) = geodetic2ecef(lat0, lon0, alt0, r_ellips);
441    let (e, n, u) = uvw2enu(x - x0, y - y0, z - z0, lat0, lon0);
442    (e, n, u)
443}
444
445/// Returns the tuple (north,east,down) of coordinates in the NED system
446///
447/// ## Inputs:
448/// - x = x ECEF coordinate [m]
449/// - y = y ECEF coordinate [m]
450/// - z = z ECEF coordinate [m]
451/// - lat0 = reference latitude [rad]
452/// - lon0 = reference longitude [rad]
453/// - alt0 = reference altitude [m]
454/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
455///
456/// ## Outputs:
457/// - n = north coordinate [m] of input location from reference geodetic location
458/// - e = east coordinate [m] of input location from reference geodetic location
459/// - d = down coordinate [m] of input location from reference geodetic location
460pub fn ecef2ned(
461    x: f64,
462    y: f64,
463    z: f64,
464    lat0: f64,
465    lon0: f64,
466    alt0: f64,
467    r_ellips: Ellipsoid,
468) -> (f64, f64, f64) {
469    let enu = ecef2enu(x, y, z, lat0, lon0, alt0, r_ellips);
470    (enu.1, enu.0, -enu.2)
471}
472
473/// Returns the tuple (east,north,up) of coordinates in the ENU system
474///
475/// ## Inputs:
476/// - u = tangent vector component
477/// - v = tangent vector component
478/// - w = tangent vector component
479/// - lat0 = reference latitude [rad]
480/// - lon0 = reference longitude [rad]
481///
482/// ## Outputs:
483/// - e = east coordinate [m] of input location from reference geodetic location
484/// - n = north coordinate [m] of input location from reference geodetic location
485/// - u = up coordinate [m] of input location from reference geodetic location
486pub fn uvw2enu(u: f64, v: f64, w: f64, lat0: f64, lon0: f64) -> (f64, f64, f64) {
487    let t = lon0.cos() * u + lon0.sin() * v;
488    let e = -lon0.sin() * u + lon0.cos() * v;
489    let n = -lat0.sin() * t + lat0.cos() * w;
490    let up = lat0.cos() * t + lat0.sin() * w;
491    (e, n, up)
492}
493
494/// Returns the tuple (azimuth,elevation,slant range) of coordinates in the AER system
495///
496/// ## Inputs:
497/// - x = x ECEF coordinate [m]
498/// - y = y ECEF coordinate [m]
499/// - z = z ECEF coordinate [m]
500/// - lat0 = reference latitude [rad]
501/// - lon0 = reference longitude [rad]
502/// - alt0 = reference altitude [m]
503/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
504///
505/// ## Outputs:
506/// - az = azimuth angle [rad] of input location from reference geodetic location
507/// - el = elevation angle [rad] of input location from reference geodetic location
508/// - slant_range = slant range [m] of input location from reference geodetic location
509pub fn ecef2aer(
510    x: f64,
511    y: f64,
512    z: f64,
513    lat0: f64,
514    lon0: f64,
515    alt0: f64,
516    r_ellips: Ellipsoid,
517) -> (f64, f64, f64) {
518    let (e, n, u) = ecef2enu(x, y, z, lat0, lon0, alt0, r_ellips);
519    let (az, el, slant_range) = enu2aer(e, n, u);
520
521    (az, el, slant_range)
522}
523
524/// Returns the tuple (azimuth,elevation,slant range) of coordinates in the AER system
525///  
526/// ## Inputs:
527/// - x = x ECI coordinate [m]
528/// - y = y ECI coordinate [m]
529/// - z = z ECI coordinate [m]
530/// - lat = reference latitude [rad]
531/// - lon = reference longitude [rad]
532/// - alt = reference altitude [m]
533/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
534///
535/// ## Outputs:
536/// - az = azimuth angle [rad] of input location from reference geodetic location
537/// - el = elevation angle [rad] of input location from reference geodetic location
538/// - slant_range = slant range [m] of input location from reference geodetic location
539pub fn eci2aer(
540    gst: f64,
541    x: f64,
542    y: f64,
543    z: f64,
544    lat: f64,
545    lon: f64,
546    alt: f64,
547    r_ellips: Ellipsoid,
548) -> (f64, f64, f64) {
549    let (x, y, z) = eci2ecef(gst, x, y, z);
550    let (az, el, slant_range) = ecef2aer(x, y, z, lat, lon, alt, r_ellips);
551    (az, el, slant_range)
552}
553
554/// Returns the tuple (x,y,z) of coordinates in the ECEF system
555///
556/// ## Inputs:
557/// - gst = greenwich sidereal time
558/// - x = x ECI coordinate [m]
559/// - y = y ECI coordinate [m]
560/// - z = z ECI coordinate [m]
561///
562/// ## Outputs:
563/// - x = x ECEF coordinate [m]
564/// - y = y ECEF coordinate [m]
565/// - z = z ECEF coordinate [m]
566pub fn eci2ecef(gst: f64, x: f64, y: f64, z: f64) -> (f64, f64, f64) {
567    let arr = matmul3(r3(gst), [x, y, z]);
568    (arr[0], arr[1], arr[2])
569}
570
571/// Returns the tuple (azimuth,elevation,slant range) of coordinates in the AER system
572///
573/// ## Inputs:
574/// - n = north coordinate [m] of input location from reference geodetic location
575/// - e = east coordinate [m] of input location from reference geodetic location
576/// - d = down coordinate [m] of input location from reference geodetic location
577///
578/// ## Outputs:
579/// - az = azimuth angle [rad] of input location from reference geodetic location
580/// - el = elevation angle [rad] of input location from reference geodetic location
581/// - slant_range = slant range [m] of input location from reference geodetic location
582pub fn ned2aer(n: f64, e: f64, d: f64) -> (f64, f64, f64) {
583    enu2aer(e, n, -d)
584}
585
586/// Returns the tuple (latitude,longitude,altitude) of coordinates in the Geodetic system
587///
588/// ## Inputs:
589/// - n = north coordinate [m] of input location from reference geodetic location
590/// - e = east coordinate [m] of input location from reference geodetic location
591/// - d = down coordinate [m] of input location from reference geodetic location
592/// - lat0 = reference latitude [rad]
593/// - lon0 = reference longitude [rad]
594/// - alt0 = reference altitude [m]
595/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
596///
597/// ## Outputs:
598/// - lat = latitude [rad]
599/// - lon = longitude [rad]
600/// - alt = altitude [m]
601pub fn ned2geodetic(
602    n: f64,
603    e: f64,
604    d: f64,
605    lat0: f64,
606    lon0: f64,
607    alt0: f64,
608    r_ellips: Ellipsoid,
609) -> (f64, f64, f64) {
610    enu2geodetic(e, n, -d, lat0, lon0, alt0, r_ellips)
611}
612
613/// Returns the tuple (x,y,z) of coordinates in the ECEF system
614///
615/// ## Inputs:
616/// - n = north coordinate [m] from reference geodetic location
617/// - e = east coordinate [m] from reference geodetic location
618/// - d = down coordinate [m] from reference geodetic location
619/// - lat0 = reference latitude [rad]
620/// - lon0 = reference longitude [rad]
621/// - alt0 = reference altitude [m]
622/// - r_ellips = reference ellipsoid, defaults ref. frame is WGS84
623///
624/// ## Outputs:
625/// - x = x ECEF coordinate [m]
626/// - y = y ECEF coordinate [m]
627/// - z = z ECEF coordinate [m]
628pub fn ned2ecef(
629    n: f64,
630    e: f64,
631    d: f64,
632    lat0: f64,
633    lon0: f64,
634    alt0: f64,
635    r_ellips: Ellipsoid,
636) -> (f64, f64, f64) {
637    enu2ecef(e, n, -d, lat0, lon0, alt0, r_ellips)
638}
639
640/// Returns the array result of 3-by-3-matrix that multiplies a 3-by-1 column array
641pub fn matmul3(matrix: [f64; 9], col: [f64; 3]) -> [f64; 3] {
642    [
643        matrix[0] * col[0] + matrix[1] * col[1] + matrix[2] * col[2],
644        matrix[3] * col[0] + matrix[4] * col[1] + matrix[5] * col[2],
645        matrix[6] * col[0] + matrix[7] * col[1] + matrix[8] * col[2],
646    ]
647}
648
649/// Returns the array representing a 3-by-3 rotation matrix of the input
650pub fn r3(x: f64) -> [f64; 9] {
651    [x.cos(), x.sin(), 0.0, -x.sin(), x.cos(), 0.0, 0.0, 0.0, 1.0]
652}
653
654/// Returns the array representing the transpose of the input 3-by-3 matrix
655pub fn transpose3(x: [f64; 9]) -> [f64; 9] {
656    [x[0], x[3], x[6], x[1], x[4], x[7], x[2], x[5], x[8]]
657}
658
659#[derive(Debug, Copy, Clone)]
660pub enum Ellipsoid {
661    /// WGS84: GPS Ellipsoid frame  
662    /// semi-major axis: 6378137.0 [m]  
663    /// flattening: 1.0/298.2572235630
664    WGS84,
665    /// WGS72: semi-major axis: 6378135.0 [m]    
666    /// flattening: 1.0/298.26
667    WGS72,
668    /// WGS66: semi-major axis: 6378145.0 [m]    
669    /// flattening: 1.0/298.25
670    WGS66,
671    /// WGS60: semi-major axis: 6378165.0 [m]    
672    /// flattening: 1.0/298.3
673    WGS60,
674    /// PZ90: Glonass Ellipsoid frame  
675    /// semi-major axis: 6378136.0 [m]
676    /// flattening: 1/298.257839303
677    PZ90,
678    /// BDC, also known as CGCS2000,
679    /// is the reference frame used by the
680    /// Beidou constellation.  
681    /// Semi-major axis: 6378137.0 [m]
682    /// flattening: 1/298.257222101
683    BDC,
684    /// GRS80 reference ellipsoid  
685    /// semi-major axis: 6378137.0 [m]  
686    /// flattening: 1.0/298.257222101
687    GRS80,
688    /// Bessel reference ellipsoid   
689    /// semi-major axis: 6377397.155 [m]
690    /// flattening: 1.0/299.1528128
691    Bessel,
692    /// Airy reference ellipsoid   
693    /// semi-major axis: 6377563.396 [m]  
694    /// flattening: 1.0/299.3249646
695    Airy,
696    /// International reference ellipsoid   
697    /// semi-major axis: 6378388.0 [m]  
698    /// flattening: 1.0/297.0
699    International,
700}
701
702impl Default for Ellipsoid {
703    fn default() -> Ellipsoid {
704        Ellipsoid::WGS84
705    }
706}
707
708impl Ellipsoid {
709    /// Returns the tuple representing the Ellipsoid frame.
710    ///
711    /// ## Outputs:
712    /// - tuple.0 = semi-major axis [m]
713    /// - tuple.1 = semi-minor axis [m]
714    /// - tuple.2 = flattening [-]
715    /// - tuple.3 = squared eccentricity [rad^2]
716    pub fn parameters(&self) -> (f64, f64, f64, f64) {
717        let (major, flattening): (f64, f64) = match self {
718            Ellipsoid::WGS84 => (6378137.0, 1.0 / 298.257223563),
719            Ellipsoid::WGS72 => (6378135.0, 1.0 / 298.26),
720            Ellipsoid::WGS66 => (6378145.0, 1.0 / 298.25),
721            Ellipsoid::WGS60 => (6378165.0, 1.0 / 298.3),
722            Ellipsoid::PZ90 => (6378136.0, 1.0 / 298.257839303),
723            Ellipsoid::BDC => (6378137.0, 1.0 / 298.257222101),
724            Ellipsoid::GRS80 => (6378137.0, 1.0 / 298.2572221009),
725            Ellipsoid::Bessel => (6377397.155, 1.0 / 299.1528128),
726            Ellipsoid::Airy => (6377563.396, 1.0 / 299.3249646),
727            Ellipsoid::International => (6378388.0, 1.0 / 297.0),
728        };
729
730        let minor = major * (1.0 - flattening);
731        let ecc_sq = ((major * major) - (minor * minor)) / (major * major);
732        (major, minor, flattening, ecc_sq)
733    }
734}
735
736/// Returns the normal radius based on given latitude
737/// and desired reference frame
738pub fn get_radius_normal(lat: f64, r_ellips: Ellipsoid) -> f64 {
739    let (major, _, _, squared_eccentricity) = r_ellips.parameters();
740    major / ((1.0 - squared_eccentricity * lat.sin() * lat.sin()).sqrt())
741}
742
743/// Returns the radians [rad] value of the decimal degree [deg] input
744#[deprecated(
745    since = "0.1.6",
746    note = "conversion from degrees to radians is natively supported. Use value.to_radians(). See https://doc.rust-lang.org/std/primitive.f64.html#method.to_radians"
747)]
748pub fn deg2rad(x: f64) -> f64 {
749    x.to_radians()
750}
751
752/// Returns the decimal degree [deg] value of the radians [rad] input
753#[deprecated(
754    since = "0.1.6",
755    note = "conversion from radians to degrees is natively supported. Use value.to_degrees(). See https://doc.rust-lang.org/std/primitive.f64.html#method.to_degrees"
756)]
757pub fn rad2deg(x: f64) -> f64 {
758    x.to_degrees()
759}
760
761/// Returns the GST time as f64
762///
763/// ## Input
764/// UTC time defined as: [year,month,day,hour,minute,second]
765///
766/// ## Output
767/// Gst time as f64
768pub fn utc2gst(utc: [i32; 6]) -> f64 {
769    let mut year = utc[0] as f64;
770    let mut month = utc[1] as f64;
771    let day = utc[2] as f64;
772    let h = utc[3] as f64;
773    let m = utc[4] as f64;
774    let s = utc[5] as f64;
775
776    if month < 3.0 {
777        year -= 1.0;
778        month += 12.0;
779    }
780
781    let a = (year / 100.0).trunc();
782
783    let b = 2.0 - a + (a / 4.0).trunc();
784
785    let c = ((s / 60.0 + m) / 60.0 + h) / 24.0;
786
787    let jd = (365.25 * (year + 4716.0)).trunc() + (30.6001 * (month + 1.0)).trunc() + day + b
788        - 1524.5
789        + c;
790
791    let t_ut1 = (jd - 2451545.0) / 36525.0;
792
793    let gmst_sec = 67310.54841 + 3.164400184812866e+09 * t_ut1 + 0.093104 * t_ut1 * t_ut1
794        - 6.2e-6 * t_ut1 * t_ut1 * t_ut1;
795
796    (gmst_sec * 2.0 * std::f64::consts::PI / 86400.0).rem_euclid(2.0 * std::f64::consts::PI)
797}
798
799/// Return the round toward zero value of the input
800#[deprecated(
801    since = "0.1.7",
802    note = "rounding toward zero is natively supported. Use value.trunc(). See https://doc.rust-lang.org/std/primitive.f64.html#method.trunc"
803)]
804pub fn fix(x: f64) -> f64 {
805    let mut out = x;
806    if out < 0.0 {
807        out = x.ceil();
808    } else {
809        out = x.floor();
810    }
811    out
812}
813
814/// Earth radius (m)
815pub const EARTH_RADIUS: f64 = 6371E3_f64;
816
817/// Returns distance (m) between two decimal degrees coordinates::
818/// coord1: (lat,lon), coord2: (lat, lon)
819pub fn distance(coord1: (f64, f64), coord2: (f64, f64)) -> f64 {
820    let dphi = coord2.0.to_radians() - coord1.0.to_radians();
821    let d_lambda = coord2.1.to_radians() - coord1.1.to_radians();
822    let a: f64 = (dphi / 2.0_f64).sin().powf(2.0_f64)
823        + coord1.0.to_radians().cos()
824            * coord2.0.to_radians().cos()
825            * (d_lambda / 2.0_f64).sin().powf(2.0_f64);
826
827    let a = a.clamp(0.0, 1.0); // Clamp a to the range [0, 1] to prevent numerical issues
828    let c = 2.0_f64 * a.sqrt().atan2((1.0 - a).sqrt());
829    EARTH_RADIUS * c
830}
831
832#[cfg(test)]
833mod tests {
834    use super::*;
835
836    #[test]
837    fn test_utc2gst() {
838        let datetime: [i32; 6] = [2020, 5, 12, 18, 2, 10];
839        let gst_ref = 2.469809475597415;
840
841        let gst = utc2gst(datetime);
842
843        assert!((gst - gst_ref).abs() < 1e-8);
844
845        let datetime2: [i32; 6] = [2020, 1, 12, 18, 2, 10];
846        let gst_ref2 = 0.388271658105431;
847
848        let gst2 = utc2gst(datetime2);
849
850        assert!((gst2 - gst_ref2).abs() < 1e-8);
851    }
852
853    #[test]
854    #[allow(deprecated)]
855    fn test_fix() {
856        let x1 = 3.7;
857        let x2 = -4.67;
858
859        assert_eq!(fix(x1), 3.0);
860        assert_eq!(fix(x2), -4.0);
861    }
862
863    #[test]
864    fn test_geodetic2ecef() {
865        let lat = 30.14988205_f64.to_radians();
866        let lon = 91.38733072_f64.to_radians();
867        let alt = 4031.0;
868
869        let (x, y, z) = geodetic2ecef(lat, lon, alt, Ellipsoid::default());
870
871        let xref = -1.337281037300386e+05;
872        let yref = 5.521796910920261e+06;
873        let zref = 3.186776473672415e+06;
874
875        assert!((x - xref).abs() < 1e-3);
876        assert!((y - yref).abs() < 1e-3);
877        assert!((z - zref).abs() < 1e-3);
878    }
879
880    #[test]
881    fn test_geodetic2aer() {
882        let lat0 = 42.0_f64.to_radians();
883        let lon0 = -82.0_f64.to_radians();
884        let alt0 = 200.0;
885
886        let lat = 42.002581974253744_f64.to_radians();
887        let lon = -81.997751960067460_f64.to_radians();
888        let alt = 1.139701799575106e+03;
889
890        let azref = 32.999999999989740_f64.to_radians();
891        let elref = 69.999999999945540_f64.to_radians();
892        let rangeref = 1000.0;
893
894        let (a, e, r) = geodetic2aer(lat, lon, alt, lat0, lon0, alt0, Ellipsoid::default());
895
896        assert!((a - azref).abs() < 1e-3);
897        assert!((e - elref).abs() < 1e-3);
898        assert!((r - rangeref).abs() < 1e-3);
899    }
900
901    #[test]
902    fn test_geodetic2enu() {
903        let lat0 = 42.0_f64.to_radians();
904        let lon0 = -82.0_f64.to_radians();
905        let alt0 = 200.0;
906
907        let lat = 42.002581974253744_f64.to_radians();
908        let lon = -81.997751960067460_f64.to_radians();
909        let alt = 1.139701799575106e+03;
910
911        let eref = 1.862775208168244e+02;
912        let nref = 2.868422278521820e+02;
913        let uref = 9.396926207845534e+02;
914
915        let (e, n, u) = geodetic2enu(lat, lon, alt, lat0, lon0, alt0, Ellipsoid::default());
916
917        assert!((e - eref).abs() < 1e-3);
918        assert!((n - nref).abs() < 1e-3);
919        assert!((u - uref).abs() < 1e-3);
920    }
921
922    #[test]
923    fn test_aer2ecef() {
924        let lat0 = 42.0_f64.to_radians();
925        let lon0 = -82.0_f64.to_radians();
926        let alt0 = 200.0;
927
928        let az = 33.0_f64.to_radians();
929        let el = 70.0_f64.to_radians();
930        let slant_range = 1000.0;
931
932        let (x, y, z) = aer2ecef(az, el, slant_range, lat0, lon0, alt0, Ellipsoid::default());
933
934        let xref = 6.609301927610816e+05;
935        let yref = -4.701424222957011e+06;
936        let zref = 4.246579604632881e+06;
937
938        assert!((x - xref).abs() < 1e-3);
939        assert!((y - yref).abs() < 1e-3);
940        assert!((z - zref).abs() < 1e-3);
941    }
942
943    #[test]
944    fn test_aer2enu() {
945        let az = 33.0_f64.to_radians();
946        let el = 70.0_f64.to_radians();
947        let slant_range = 1000.0;
948
949        let eref = 1.862775208165935e+02;
950        let nref = 2.868422278517140e+02;
951        let uref = 9.396926207859083e+02;
952
953        let (e, n, u) = aer2enu(az, el, slant_range);
954
955        assert!((e - eref).abs() < 1e-3);
956        assert!((n - nref).abs() < 1e-3);
957        assert!((u - uref).abs() < 1e-3);
958    }
959
960    #[test]
961    fn test_aer2eci() {
962        let az = 162.55_f64.to_radians();
963        let el = 55.12_f64.to_radians();
964        let slant_range = 384013940.9;
965        let gst = 4.501012562811752;
966
967        let lat0 = 28.4_f64.to_radians();
968        let lon0 = -80.5_f64.to_radians();
969        let alt0 = 2.7;
970
971        let xref = -3.849714979138141e+08;
972        let yref = -4.836588977863766e+07;
973        let zref = -3.143285462295778e+07;
974
975        let (x, y, z) = aer2eci(
976            gst,
977            az,
978            el,
979            slant_range,
980            lat0,
981            lon0,
982            alt0,
983            Ellipsoid::default(),
984        );
985
986        assert!((x - xref).abs() < 1e-3);
987        assert!((y - yref).abs() < 1e-3);
988        assert!((z - zref).abs() < 1e-3);
989    }
990
991    #[test]
992    fn test_aer2geodetic() {
993        let lat0 = 42.0_f64.to_radians();
994        let lon0 = -82.0_f64.to_radians();
995        let alt0 = 200.0;
996
997        let az = 32.999999999989740_f64.to_radians();
998        let el = 69.999999999945540_f64.to_radians();
999        let slant_range = 1000.0;
1000
1001        let latref = 42.002581974253744_f64.to_radians();
1002        let lonref = -81.997751960067460_f64.to_radians();
1003        let altref = 1.139701799575106e+03;
1004
1005        let (lat, lon, alt) =
1006            aer2geodetic(az, el, slant_range, lat0, lon0, alt0, Ellipsoid::default());
1007
1008        assert!((lat - latref).abs() < 1e-8);
1009        assert!((lon - lonref).abs() < 1e-8);
1010        assert!((alt - altref).abs() < 1e-8);
1011    }
1012
1013    #[test]
1014    fn test_enu2aer() {
1015        let e = 1.862775210000000e+02;
1016        let n = 2.868422200000000e+02;
1017        let u = 9.396926200000000e+02;
1018
1019        let azref = 33.0_f64.to_radians();
1020        let elref = 70.0_f64.to_radians();
1021        let rangeref = 1000.0;
1022
1023        let (az, el, range) = enu2aer(e, n, u);
1024
1025        assert!((az - azref).abs() < 1e-3);
1026        assert!((el - elref).abs() < 1e-3);
1027        assert!((range - rangeref).abs() < 1e-3);
1028    }
1029
1030    #[test]
1031    fn test_enu2ecef() {
1032        let lat0 = 42.0_f64.to_radians();
1033        let lon0 = -82.0_f64.to_radians();
1034        let alt0 = 200.0;
1035        let e = 1.862775210000000e+02;
1036        let n = 2.868422200000000e+02;
1037        let u = 9.396926200000000e+02;
1038
1039        let xref = 6.609301927610815e+05;
1040        let yref = -4.701424222957011e+06;
1041        let zref = 4.246579604632881e+06;
1042
1043        let (x, y, z) = enu2ecef(e, n, u, lat0, lon0, alt0, Ellipsoid::default());
1044
1045        assert!((x - xref).abs() < 1e-3);
1046        assert!((y - yref).abs() < 1e-3);
1047        assert!((z - zref).abs() < 1e-3);
1048    }
1049
1050    #[test]
1051    fn test_enu2geodetic() {
1052        let lat0 = 42.0_f64.to_radians();
1053        let lon0 = -82.0_f64.to_radians();
1054        let alt0 = 200.0;
1055        let e = 0.0;
1056        let n = 0.0;
1057        let u = -1.0;
1058
1059        let latref = 41.999999999999993_f64.to_radians();
1060        let lonref = -82.0_f64.to_radians();
1061        let altref = 1.990000000007368e+02;
1062
1063        let (lat, lon, alt) = enu2geodetic(e, n, u, lat0, lon0, alt0, Ellipsoid::default());
1064
1065        assert!((lat - latref).abs() < 1e-8);
1066        assert!((lon - lonref).abs() < 1e-8);
1067        assert!((alt - altref).abs() < 1e-8);
1068    }
1069
1070    #[test]
1071    fn test_ecef2geodetic() {
1072        let latref = 30.14988205_f64.to_radians();
1073        let lonref = 91.38733072_f64.to_radians();
1074        let altref = 4031.0;
1075
1076        let (x, y, z) = geodetic2ecef(latref, lonref, altref, Ellipsoid::default());
1077        let (lat, lon, alt) = ecef2geodetic(x, y, z, Ellipsoid::default());
1078
1079        assert!((lat - latref).abs() < 1e-8);
1080        assert!((lon - lonref).abs() < 1e-8);
1081        assert!((alt - altref).abs() < 1e-8);
1082
1083        let (x, y, z) = geodetic2ecef(latref, lonref, altref - 5000.0, Ellipsoid::default());
1084        let (lat, lon, alt) = ecef2geodetic(x, y, z, Ellipsoid::default());
1085
1086        assert!((lat - latref).abs() < 1e-8);
1087        assert!((lon - lonref).abs() < 1e-8);
1088        assert!((alt - (altref - 5000.0)).abs() < 1e-8);
1089    }
1090
1091    #[test]
1092    fn test_ecef2enu() {
1093        let lat0 = 42.0_f64.to_radians();
1094        let lon0 = -82.0_f64.to_radians();
1095        let alt0 = 200.0;
1096        let eref = 1.862775210000000e+02;
1097        let nref = 2.868422200000000e+02;
1098        let uref = 9.396926200000000e+02;
1099
1100        let (x, y, z) = enu2ecef(eref, nref, uref, lat0, lon0, alt0, Ellipsoid::default());
1101        let (e, n, u) = ecef2enu(x, y, z, lat0, lon0, alt0, Ellipsoid::default());
1102
1103        assert!((e - eref).abs() < 1e-3);
1104        assert!((n - nref).abs() < 1e-3);
1105        assert!((u - uref).abs() < 1e-3);
1106    }
1107
1108    #[test]
1109    fn test_ecef2aer() {
1110        let lat0 = 42.0_f64.to_radians();
1111        let lon0 = -82.0_f64.to_radians();
1112        let alt0 = 200.0;
1113
1114        let azref = 33.0_f64.to_radians();
1115        let elref = 70.0_f64.to_radians();
1116        let rangeref = 1000.0;
1117
1118        let (x, y, z) = aer2ecef(
1119            azref,
1120            elref,
1121            rangeref,
1122            lat0,
1123            lon0,
1124            alt0,
1125            Ellipsoid::default(),
1126        );
1127        let (az, el, range) = ecef2aer(x, y, z, lat0, lon0, alt0, Ellipsoid::default());
1128
1129        assert!((az - azref).abs() < 1e-3);
1130        assert!((el - elref).abs() < 1e-3);
1131        assert!((range - rangeref).abs() < 1e-3);
1132    }
1133
1134    #[test]
1135    fn test_eci2aer() {
1136        let azref = 162.55_f64.to_radians();
1137        let elref = 55.12_f64.to_radians();
1138        let rangeref = 384013940.9;
1139        let gst = 4.501012562811752;
1140
1141        let lat0 = 28.4_f64.to_radians();
1142        let lon0 = -80.5_f64.to_radians();
1143        let alt0 = 2.7;
1144
1145        let (x, y, z) = aer2eci(
1146            gst,
1147            azref,
1148            elref,
1149            rangeref,
1150            lat0,
1151            lon0,
1152            alt0,
1153            Ellipsoid::default(),
1154        );
1155        let (az, el, range) = eci2aer(gst, x, y, z, lat0, lon0, alt0, Ellipsoid::default());
1156
1157        assert!((az - azref).abs() < 1e-3);
1158        assert!((el - elref).abs() < 1e-3);
1159        assert!((range - rangeref).abs() < 1e-3);
1160    }
1161
1162    #[test]
1163    fn test_ned2geodetic() {
1164        let lat0 = 42.0_f64.to_radians();
1165        let lon0 = -82.0_f64.to_radians();
1166        let alt0 = 200.0;
1167        let e = 0.0;
1168        let n = 0.0;
1169        let d = 1.0;
1170
1171        let latref = 41.999999999999993_f64.to_radians();
1172        let lonref = -82.0_f64.to_radians();
1173        let altref = 1.990000000007368e+02;
1174
1175        let (lat, lon, alt) = ned2geodetic(n, e, d, lat0, lon0, alt0, Ellipsoid::default());
1176
1177        assert!((lat - latref).abs() < 1e-8);
1178        assert!((lon - lonref).abs() < 1e-8);
1179        assert!((alt - altref).abs() < 1e-8);
1180    }
1181
1182    #[test]
1183    fn test_geodetic2ned() {
1184        let lat = 41.999999999999993_f64.to_radians();
1185        let lon = -82.0_f64.to_radians();
1186        let alt = 1.990000000007368e+02;
1187        let lat0 = 42.0_f64.to_radians();
1188        let lon0 = -82.0_f64.to_radians();
1189        let alt0 = 200.0;
1190
1191        let eref = 0.0;
1192        let nref = 0.0;
1193        let dref = 1.0;
1194
1195        let (n, e, d) = geodetic2ned(lat, lon, alt, lat0, lon0, alt0, Ellipsoid::default());
1196
1197        assert!((e - eref).abs() < 1e-3);
1198        assert!((n - nref).abs() < 1e-3);
1199        assert!((d - dref).abs() < 1e-3);
1200    }
1201
1202    #[test]
1203    fn test_aer2ned() {
1204        let az = 33.0_f64.to_radians();
1205        let el = 70.0_f64.to_radians();
1206        let slant_range = 1000.0;
1207
1208        let eref = 1.862775208165935e+02;
1209        let nref = 2.868422278517140e+02;
1210        let dref = -9.396926207859083e+02;
1211
1212        let (n, e, d) = aer2ned(az, el, slant_range);
1213
1214        assert!((e - eref).abs() < 1e-3);
1215        assert!((n - nref).abs() < 1e-3);
1216        assert!((d - dref).abs() < 1e-3);
1217    }
1218
1219    #[test]
1220    fn test_ned2aer() {
1221        let az_ref = 33.0_f64.to_radians();
1222        let el_ref = 70.0_f64.to_radians();
1223        let range_ref = 1000.0;
1224
1225        let e = 1.862775208165935e+02;
1226        let n = 2.868422278517140e+02;
1227        let d = -9.396926207859083e+02;
1228
1229        let (az, el, range) = ned2aer(n, e, d);
1230
1231        assert!((az - az_ref).abs() < 1e-6);
1232        assert!((el - el_ref).abs() < 1e-6);
1233        assert!((range - range_ref).abs() < 1e-3);
1234    }
1235
1236    #[test]
1237    fn test_ned2ecef() {
1238        let lat0 = 42.0_f64.to_radians();
1239        let lon0 = -82.0_f64.to_radians();
1240        let alt0 = 200.0;
1241        let e = 1.862775210000000e+02;
1242        let n = 2.868422200000000e+02;
1243        let d = -9.396926200000000e+02;
1244
1245        let xref = 6.609301927610815e+05;
1246        let yref = -4.701424222957011e+06;
1247        let zref = 4.246579604632881e+06;
1248
1249        let (x, y, z) = ned2ecef(n, e, d, lat0, lon0, alt0, Ellipsoid::default());
1250
1251        assert!((x - xref).abs() < 1e-3);
1252        assert!((y - yref).abs() < 1e-3);
1253        assert!((z - zref).abs() < 1e-3);
1254    }
1255
1256    #[test]
1257    fn test_ellipsoid_references() {
1258        let (a, b, f, e) = Ellipsoid::WGS84.parameters();
1259        assert!((a - 6378137.0).abs() < 1E-6);
1260        assert!((b - 6356752.314245).abs() < 1E-6);
1261        assert!((1.0 / f - 298.257223563).abs() < 1E-6);
1262        assert!((e - 6.6943799E-3_f64).abs() < 1E-6);
1263        let (a, b, f, e) = Ellipsoid::WGS72.parameters();
1264        assert!((b - a * (1.0 - f)).abs() < 1E-6);
1265        assert!((e - (f * (2.0 - f))).abs() < 1E-6);
1266        let (a, b, f, e) = Ellipsoid::WGS66.parameters();
1267        assert!((b - a * (1.0 - f)).abs() < 1E-6);
1268        assert!((e - (f * (2.0 - f))).abs() < 1E-6);
1269        let (a, b, f, e) = Ellipsoid::WGS60.parameters();
1270        assert!((b - a * (1.0 - f)).abs() < 1E-6);
1271        assert!((e - (f * (2.0 - f))).abs() < 1E-6);
1272        let (a, b, f, e) = Ellipsoid::PZ90.parameters();
1273        assert!((a - 6378136.0).abs() < 1E-6);
1274        assert!((b - a * (1.0 - f)).abs() < 1E-6);
1275        assert!((1.0 / f - 298.257839303).abs() < 1E-6);
1276        assert!((e - (f * (2.0 - f))).abs() < 1E-6);
1277        let (a, b, f, e) = Ellipsoid::GRS80.parameters();
1278        assert!((a - 6378137.0).abs() < 1E-6);
1279        assert!((b - a * (1.0 - f)).abs() < 1E-6);
1280        assert!((1.0 / f - 298.257222101).abs() < 1E-6);
1281        assert!((e - (f * (2.0 - f))).abs() < 1E-6);
1282        let (a, b, f, e) = Ellipsoid::BDC.parameters();
1283        assert!((a - 6378137.0).abs() < 1E-6);
1284        assert!((b - a * (1.0 - f)).abs() < 1E-6);
1285        assert!((1.0 / f - 298.257222101).abs() < 1E-6);
1286        assert!((e - (f * (2.0 - f))).abs() < 1E-6);
1287        let (a, b, f, e) = Ellipsoid::Bessel.parameters();
1288        assert!((b - a * (1.0 - f)).abs() < 1E-6);
1289        assert!((e - (f * (2.0 - f))).abs() < 1E-6);
1290        let (a, b, f, e) = Ellipsoid::International.parameters();
1291        assert!((b - a * (1.0 - f)).abs() < 1E-6);
1292        assert!((e - (f * (2.0 - f))).abs() < 1E-6);
1293        let (a, b, f, e) = Ellipsoid::Airy.parameters();
1294        assert!((b - a * (1.0 - f)).abs() < 1E-6);
1295        assert!((e - (f * (2.0 - f))).abs() < 1E-6);
1296    }
1297
1298    #[test]
1299    fn test_ecef2ned() {
1300        let lat0 = 42.0_f64.to_radians();
1301        let lon0 = -82.0_f64.to_radians();
1302        let alt0 = 200.0;
1303        let eref = 1.862775210000000e+02;
1304        let nref = 2.868422200000000e+02;
1305        let dref = -9.396926200000000e+02;
1306
1307        let (x, y, z) = ned2ecef(nref, eref, dref, lat0, lon0, alt0, Ellipsoid::default());
1308        let (n, e, d) = ecef2ned(x, y, z, lat0, lon0, alt0, Ellipsoid::default());
1309
1310        assert!((e - eref).abs() < 1e-3);
1311        assert!((n - nref).abs() < 1e-3);
1312        assert!((d - dref).abs() < 1e-3);
1313    }
1314
1315    #[test]
1316    fn test_distance_calculation() {
1317        let new_york = (40.730610, -73.935242);
1318        let paris = (48.856614, 2.3522219);
1319        let buenos_aires = (-34.603722, -58.381592);
1320        let sydney = (-33.867487, 151.20699);
1321        // TEST 1
1322        let expected_km = 5831.0_f64;
1323        let d_km = distance(new_york, paris) / 1000.0_f64;
1324        assert!((expected_km - d_km).abs() < 1.0);
1325        // TEST2
1326        let expected_km = 8527.0_f64;
1327        let d_km = distance(new_york, buenos_aires) / 1000.0_f64;
1328        assert!((expected_km - d_km).abs() < 1.0);
1329        // TEST3
1330        let expected_km = 15990.0_f64;
1331        let d_km = distance(new_york, sydney) / 1000.0_f64;
1332        assert!((expected_km - d_km).abs() < 10.0);
1333        // TEST4
1334        let expected_km = 11050.0_f64;
1335        let d_km = distance(buenos_aires, paris) / 1000.0_f64;
1336        assert!((expected_km - d_km).abs() < 10.0);
1337    }
1338
1339    #[test]
1340    fn test_ecef2geodetic_poles() {
1341        let (_, minor, _, _) = Ellipsoid::default().parameters();
1342
1343        let (lat, lon, alt) = ecef2geodetic(0.0, 0.0, minor, Ellipsoid::default());
1344        assert!((lat - std::f64::consts::FRAC_PI_2).abs() < 1e-12);
1345        assert!(lon.abs() < 1e-12);
1346        assert!(alt.abs() < 1e-12);
1347
1348        let (lat, lon, alt) = ecef2geodetic(0.0, 0.0, -minor, Ellipsoid::default());
1349        assert!((lat + std::f64::consts::FRAC_PI_2).abs() < 1e-12);
1350        assert!(lon.abs() < 1e-12);
1351        assert!(alt.abs() < 1e-12);
1352    }
1353
1354    #[test]
1355    fn test_ecef2geodetic_near_pole() {
1356        let lat_ref = 89.999999_f64.to_radians();
1357        let lon_ref = 30.0_f64.to_radians();
1358        let alt_ref = 1000.0;
1359
1360        let (x, y, z) = geodetic2ecef(lat_ref, lon_ref, alt_ref, Ellipsoid::default());
1361        let (lat, lon, alt) = ecef2geodetic(x, y, z, Ellipsoid::default());
1362
1363        assert!((lat - lat_ref).abs() < 1e-8);
1364        assert!((lon - lon_ref).abs() < 1e-8);
1365        assert!((alt - alt_ref).abs() < 1e-5);
1366    }
1367
1368    #[test]
1369    fn test_distance_antipodal() {
1370        let d = distance((0.0, 0.0), (0.0, 180.0));
1371        assert!(d.is_finite());
1372    }
1373}