Coverage for pygeodesy/formy.py: 99%
362 statements
« prev ^ index » next coverage.py v7.2.2, created at 2023-04-21 13:14 -0400
« prev ^ index » next coverage.py v7.2.2, created at 2023-04-21 13:14 -0400
2# -*- coding: utf-8 -*-
4u'''Formulary of basic geodesy functions and approximations.
5'''
6# make sure int/int division yields float quotient, see .basics
7from __future__ import division as _; del _ # PYCHOK semicolon
9# from pygeodesy.basics import isscalar # from .fsums
10from pygeodesy.constants import EPS, EPS0, EPS1, PI, PI2, PI3, PI_2, R_M, \
11 _umod_PI2, float0_, isnon0, remainder, \
12 _0_0, _0_125, _0_25, _0_5, _1_0, _2_0, \
13 _4_0, _32_0, _90_0, _180_0, _360_0
14from pygeodesy.datums import Datum, Ellipsoid, _ellipsoidal_datum, \
15 _mean_radius, _spherical_datum, _WGS84
16# from pygeodesy.ellipsoids import Ellipsoid # from .datums
17from pygeodesy.errors import _AssertionError, IntersectionError, LimitError, \
18 limiterrors, _ValueError, _xError
19from pygeodesy.fmath import euclid, hypot, hypot2, sqrt0
20from pygeodesy.fsums import fsum_, isscalar, unstr
21from pygeodesy.interns import NN, _distant_, _too_
22from pygeodesy.lazily import _ALL_LAZY, _ALL_MODS as _MODS
23from pygeodesy.named import _NamedTuple, _xnamed
24from pygeodesy.namedTuples import Bearing2Tuple, Distance4Tuple, \
25 Intersection3Tuple, LatLon2Tuple, \
26 PhiLam2Tuple, Vector3Tuple
27# from pygeodesy.streprs import unstr # from .fsums
28from pygeodesy.units import Bearing, Degrees_, Distance, Distance_, Height, \
29 Lam_, Lat, Lon, Meter_, Phi_, Radians, Radians_, \
30 Radius, Radius_, Scalar, _100km
31from pygeodesy.utily import acos1, atan2b, atan2d, degrees2m, m2degrees, tan_2, \
32 sincos2, sincos2_, sincos2d_, unroll180, unrollPI
34from math import atan, atan2, cos, degrees, fabs, radians, sin, sqrt # pow
36__all__ = _ALL_LAZY.formy
37__version__ = '23.04.14'
39_ratio_ = 'ratio'
40_xline_ = 'xline'
43def _anti2(a, b, n_2, n, n2):
44 '''(INTERNAL) Helper for C{antipode} and C{antipode_}.
45 '''
46 r = remainder(a, n) if fabs(a) > n_2 else a
47 if r == a:
48 r = -r
49 b += n
50 if fabs(b) > n:
51 b = remainder(b, n2)
52 return float0_(r, b)
55def antipode(lat, lon, name=NN):
56 '''Return the antipode, the point diametrically opposite
57 to a given point in C{degrees}.
59 @arg lat: Latitude (C{degrees}).
60 @arg lon: Longitude (C{degrees}).
61 @kwarg name: Optional name (C{str}).
63 @return: A L{LatLon2Tuple}C{(lat, lon)}.
65 @see: Functions L{antipode_} and L{normal} and U{Geosphere
66 <https://CRAN.R-Project.org/web/packages/geosphere/geosphere.pdf>}.
67 '''
68 return LatLon2Tuple(*_anti2(lat, lon, _90_0, _180_0, _360_0), name=name)
71def antipode_(phi, lam, name=NN):
72 '''Return the antipode, the point diametrically opposite
73 to a given point in C{radians}.
75 @arg phi: Latitude (C{radians}).
76 @arg lam: Longitude (C{radians}).
77 @kwarg name: Optional name (C{str}).
79 @return: A L{PhiLam2Tuple}C{(phi, lam)}.
81 @see: Functions L{antipode} and L{normal_} and U{Geosphere
82 <https://CRAN.R-Project.org/web/packages/geosphere/geosphere.pdf>}.
83 '''
84 return PhiLam2Tuple(*_anti2(phi, lam, PI_2, PI, PI2), name=name)
87def _area_or_(excess_, lat1, lat2, radius, d_lon, unused):
88 '''(INTERNAL) Helper for area and spherical excess.
89 '''
90 r = excess_(Phi_(lat2=lat2),
91 Phi_(lat1=lat1), radians(d_lon))
92 if radius:
93 r *= _mean_radius(radius, lat1, lat2)**2
94 return r
97def bearing(lat1, lon1, lat2, lon2, **options):
98 '''Compute the initial or final bearing (forward or reverse
99 azimuth) between a (spherical) start and end point.
101 @arg lat1: Start latitude (C{degrees}).
102 @arg lon1: Start longitude (C{degrees}).
103 @arg lat2: End latitude (C{degrees}).
104 @arg lon2: End longitude (C{degrees}).
105 @kwarg options: Optional keyword arguments for function
106 L{pygeodesy.bearing_}.
108 @return: Initial or final bearing (compass C{degrees360}) or
109 zero if start and end point coincide.
110 '''
111 r = bearing_(Phi_(lat1=lat1), Lam_(lon1=lon1),
112 Phi_(lat2=lat2), Lam_(lon2=lon2), **options)
113 return degrees(r)
116def bearing_(phi1, lam1, phi2, lam2, final=False, wrap=False):
117 '''Compute the initial or final bearing (forward or reverse azimuth)
118 between a (spherical) start and end point.
120 @arg phi1: Start latitude (C{radians}).
121 @arg lam1: Start longitude (C{radians}).
122 @arg phi2: End latitude (C{radians}).
123 @arg lam2: End longitude (C{radians}).
124 @kwarg final: Return final bearing if C{True}, initial otherwise (C{bool}).
125 @kwarg wrap: Wrap and L{pygeodesy.unrollPI} longitudes (C{bool}).
127 @return: Initial or final bearing (compass C{radiansPI2}) or zero if start
128 and end point coincide.
130 @see: U{Bearing<https://www.Movable-Type.co.UK/scripts/latlong.html>}, U{Course
131 between two points<https://www.EdWilliams.org/avform147.htm#Crs>} and
132 U{Bearing Between Two Points<https://web.Archive.org/web/20020630205931/
133 https://MathForum.org/library/drmath/view/55417.html>}.
134 '''
135 if final: # swap plus PI
136 phi1, lam1, phi2, lam2 = phi2, lam2, phi1, lam1
137 r = PI3
138 else:
139 r = PI2
141 db, _ = unrollPI(lam1, lam2, wrap=wrap)
142 sa1, ca1, sa2, ca2, sdb, cdb = sincos2_(phi1, phi2, db)
144 x = ca1 * sa2 - sa1 * ca2 * cdb
145 y = sdb * ca2
146 return _umod_PI2(atan2(y, x) + r) # .utily.wrapPI2
149def _bearingTo2(p1, p2, wrap=False): # for points.ispolar, sphericalTrigonometry.areaOf
150 '''(INTERNAL) Compute initial and final bearing.
151 '''
152 try: # for LatLon_ and ellipsoidal LatLon
153 return p1.bearingTo2(p2, wrap=wrap)
154 except AttributeError:
155 pass
156 # XXX spherical version, OK for ellipsoidal ispolar?
157 a1, b1 = p1.philam
158 a2, b2 = p2.philam
159 return Bearing2Tuple(degrees(bearing_(a1, b1, a2, b2, final=False, wrap=wrap)),
160 degrees(bearing_(a1, b1, a2, b2, final=True, wrap=wrap)),
161 name=_bearingTo2.__name__)
164def compassAngle(lat1, lon1, lat2, lon2, adjust=True, wrap=False):
165 '''Return the angle from North for the direction vector
166 M{(lon2 - lon1, lat2 - lat1)} between two points.
168 Suitable only for short, not near-polar vectors up to a few hundred
169 Km or Miles. Use function L{pygeodesy.bearing} for longer vectors.
171 @arg lat1: From latitude (C{degrees}).
172 @arg lon1: From longitude (C{degrees}).
173 @arg lat2: To latitude (C{degrees}).
174 @arg lon2: To longitude (C{degrees}).
175 @kwarg adjust: Adjust the longitudinal delta by the cosine of the
176 mean latitude (C{bool}).
177 @kwarg wrap: Wrap and L{pygeodesy.unroll180} longitudes (C{bool}).
179 @return: Compass angle from North (C{degrees360}).
181 @note: Courtesy of Martin Schultz.
183 @see: U{Local, flat earth approximation
184 <https://www.EdWilliams.org/avform.htm#flat>}.
185 '''
186 d_lon, _ = unroll180(lon1, lon2, wrap=wrap)
187 if adjust: # scale delta lon
188 d_lon *= _scale_deg(lat1, lat2)
189 return atan2b(d_lon, lat2 - lat1)
192def cosineAndoyerLambert(lat1, lon1, lat2, lon2, datum=_WGS84, wrap=False):
193 '''Compute the distance between two (ellipsoidal) points using the
194 U{Andoyer-Lambert correction<https://navlib.net/wp-content/uploads/
195 2013/10/admiralty-manual-of-navigation-vol-1-1964-english501c.pdf>} of the
196 U{Law of Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>}
197 fromula.
199 @arg lat1: Start latitude (C{degrees}).
200 @arg lon1: Start longitude (C{degrees}).
201 @arg lat2: End latitude (C{degrees}).
202 @arg lon2: End longitude (C{degrees}).
203 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid},
204 L{Ellipsoid2} or L{a_f2Tuple}).
205 @kwarg wrap: Wrap and L{pygeodesy.unroll180} longitudes (C{bool}).
207 @return: Distance (C{meter}, same units as the B{C{datum}}'s
208 ellipsoid axes or C{radians} if B{C{datum}} is C{None}).
210 @raise TypeError: Invalid B{C{datum}}.
212 @see: Functions L{cosineAndoyerLambert_}, L{cosineForsytheAndoyerLambert},
213 L{cosineLaw}, L{equirectangular}, L{euclidean}, L{flatLocal}/L{hubeny},
214 L{flatPolar}, L{haversine}, L{thomas} and L{vincentys} and method
215 L{Ellipsoid.distance2}.
216 '''
217 return _distanceToE(cosineAndoyerLambert_, lat1, lat2, datum,
218 *unroll180(lon1, lon2, wrap=wrap))
221def cosineAndoyerLambert_(phi2, phi1, lam21, datum=_WGS84):
222 '''Compute the I{angular} distance between two (ellipsoidal) points using the
223 U{Andoyer-Lambert correction<https://navlib.net/wp-content/uploads/2013/10/
224 admiralty-manual-of-navigation-vol-1-1964-english501c.pdf>} of the U{Law
225 of Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>}
226 fromula.
228 @arg phi2: End latitude (C{radians}).
229 @arg phi1: Start latitude (C{radians}).
230 @arg lam21: Longitudinal delta, M{end-start} (C{radians}).
231 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid},
232 L{Ellipsoid2} or L{a_f2Tuple}).
234 @return: Angular distance (C{radians}).
236 @raise TypeError: Invalid B{C{datum}}.
238 @see: Functions L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert_},
239 L{cosineLaw_}, L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_},
240 L{flatPolar_}, L{haversine_}, L{thomas_} and L{vincentys_} and U{Geodesy-PHP
241 <https://GitHub.com/jtejido/geodesy-php/blob/master/src/Geodesy/
242 Distance/AndoyerLambert.php>}.
243 '''
244 s2, c2, s1, c1, r, c21 = _sincosa6(phi2, phi1, lam21)
245 if isnon0(c1) and isnon0(c2):
246 E = _ellipsoidal(datum, cosineAndoyerLambert_)
247 if E.f: # ellipsoidal
248 r2 = atan2(E.b_a * s2, c2)
249 r1 = atan2(E.b_a * s1, c1)
250 s2, c2, s1, c1 = sincos2_(r2, r1)
251 r = acos1(s1 * s2 + c1 * c2 * c21)
252 if r:
253 sr, _, sr_2, cr_2 = sincos2_(r, r * _0_5)
254 if isnon0(sr_2) and isnon0(cr_2):
255 s = (sr + r) * ((s1 - s2) / sr_2)**2
256 c = (sr - r) * ((s1 + s2) / cr_2)**2
257 r += (c - s) * E.f * _0_125
258 return r
261def cosineForsytheAndoyerLambert(lat1, lon1, lat2, lon2, datum=_WGS84, wrap=False):
262 '''Compute the distance between two (ellipsoidal) points using the
263 U{Forsythe-Andoyer-Lambert correction<https://www2.UNB.Ca/gge/Pubs/TR77.pdf>} of
264 the U{Law of Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>}
265 formula.
267 @arg lat1: Start latitude (C{degrees}).
268 @arg lon1: Start longitude (C{degrees}).
269 @arg lat2: End latitude (C{degrees}).
270 @arg lon2: End longitude (C{degrees}).
271 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid},
272 L{Ellipsoid2} or L{a_f2Tuple}).
273 @kwarg wrap: Wrap and L{pygeodesy.unroll180} longitudes (C{bool}).
275 @return: Distance (C{meter}, same units as the B{C{datum}}'s
276 ellipsoid axes or C{radians} if B{C{datum}} is C{None}).
278 @raise TypeError: Invalid B{C{datum}}.
280 @see: Functions L{cosineForsytheAndoyerLambert_}, L{cosineAndoyerLambert},
281 L{cosineLaw}, L{equirectangular}, L{euclidean}, L{flatLocal}/L{hubeny},
282 L{flatPolar}, L{haversine}, L{thomas} and L{vincentys} and method
283 L{Ellipsoid.distance2}.
284 '''
285 return _distanceToE(cosineForsytheAndoyerLambert_, lat1, lat2, datum,
286 *unroll180(lon1, lon2, wrap=wrap))
289def cosineForsytheAndoyerLambert_(phi2, phi1, lam21, datum=_WGS84):
290 '''Compute the I{angular} distance between two (ellipsoidal) points using the
291 U{Forsythe-Andoyer-Lambert correction<https://www2.UNB.Ca/gge/Pubs/TR77.pdf>} of
292 the U{Law of Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>}
293 formula.
295 @arg phi2: End latitude (C{radians}).
296 @arg phi1: Start latitude (C{radians}).
297 @arg lam21: Longitudinal delta, M{end-start} (C{radians}).
298 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid},
299 L{Ellipsoid2} or L{a_f2Tuple}).
301 @return: Angular distance (C{radians}).
303 @raise TypeError: Invalid B{C{datum}}.
305 @see: Functions L{cosineForsytheAndoyerLambert}, L{cosineAndoyerLambert_},
306 L{cosineLaw_}, L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_},
307 L{flatPolar_}, L{haversine_}, L{thomas_} and L{vincentys_} and U{Geodesy-PHP
308 <https://GitHub.com/jtejido/geodesy-php/blob/master/src/Geodesy/
309 Distance/ForsytheCorrection.php>}.
310 '''
311 s2, c2, s1, c1, r, _ = _sincosa6(phi2, phi1, lam21)
312 if r and isnon0(c1) and isnon0(c2):
313 E = _ellipsoidal(datum, cosineForsytheAndoyerLambert_)
314 if E.f: # ellipsoidal
315 sr, cr, s2r, _ = sincos2_(r, r * _2_0)
316 if isnon0(sr) and fabs(cr) < EPS1:
317 s = (s1 + s2)**2 / (1 + cr)
318 t = (s1 - s2)**2 / (1 - cr)
319 x = s + t
320 y = s - t
322 s = 8 * r**2 / sr
323 a = 64 * r + _2_0 * s * cr # 16 * r**2 / tan(r)
324 d = 48 * sr + s # 8 * r**2 / tan(r)
325 b = -2 * d
326 e = 30 * s2r
327 c = fsum_(30 * r, e * _0_5, s * cr) # 8 * r**2 / tan(r)
329 t = fsum_( a * x, b * y, -c * x**2, d * x * y, e * y**2)
330 r += fsum_(-r * x, 3 * y * sr, t * E.f / _32_0) * E.f * _0_25
331 return r
334def cosineLaw(lat1, lon1, lat2, lon2, radius=R_M, wrap=False):
335 '''Compute the distance between two points using the
336 U{spherical Law of Cosines
337 <https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>}
338 formula.
340 @arg lat1: Start latitude (C{degrees}).
341 @arg lon1: Start longitude (C{degrees}).
342 @arg lat2: End latitude (C{degrees}).
343 @arg lon2: End longitude (C{degrees}).
344 @kwarg radius: Mean earth radius, ellipsoid or datum (C{meter},
345 L{Ellipsoid}, L{Ellipsoid2}, L{Datum} or L{a_f2Tuple}).
346 @kwarg wrap: Wrap and L{pygeodesy.unroll180} longitudes (C{bool}).
348 @return: Distance (C{meter}, same units as B{C{radius}} or the
349 ellipsoid or datum axes).
351 @raise TypeError: Invalid B{C{radius}}.
353 @see: Functions L{cosineLaw_}, L{cosineAndoyerLambert},
354 L{cosineForsytheAndoyerLambert}, L{equirectangular}, L{euclidean},
355 L{flatLocal}/L{hubeny}, L{flatPolar}, L{haversine}, L{thomas} and
356 L{vincentys} and method L{Ellipsoid.distance2}.
358 @note: See note at function L{vincentys_}.
359 '''
360 return _distanceToS(cosineLaw_, lat1, lat2, radius,
361 *unroll180(lon1, lon2, wrap=wrap))
364def cosineLaw_(phi2, phi1, lam21):
365 '''Compute the I{angular} distance between two points using the
366 U{spherical Law of Cosines
367 <https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>}
368 formula.
370 @arg phi2: End latitude (C{radians}).
371 @arg phi1: Start latitude (C{radians}).
372 @arg lam21: Longitudinal delta, M{end-start} (C{radians}).
374 @return: Angular distance (C{radians}).
376 @see: Functions L{cosineLaw}, L{cosineAndoyerLambert_},
377 L{cosineForsytheAndoyerLambert_}, L{equirectangular_},
378 L{euclidean_}, L{flatLocal_}/L{hubeny_}, L{flatPolar_},
379 L{haversine_}, L{thomas_} and L{vincentys_}.
381 @note: See note at function L{vincentys_}.
382 '''
383 return _sincosa6(phi2, phi1, lam21)[4]
386def _distanceToE(func_, lat1, lat2, earth, d_lon, unused):
387 '''(INTERNAL) Helper for ellipsoidal distances.
388 '''
389 E = _ellipsoidal(earth, func_)
390 r = func_(Phi_(lat2=lat2),
391 Phi_(lat1=lat1), radians(d_lon), datum=E)
392 return r * E.a
395def _distanceToS(func_, lat1, lat2, earth, d_lon, unused, **adjust):
396 '''(INTERNAL) Helper for spherical distances.
397 '''
398 r = func_(Phi_(lat2=lat2),
399 Phi_(lat1=lat1), radians(d_lon), **adjust)
400 return r * _mean_radius(earth, lat1, lat2)
403def _ellipsoidal(earth, where):
404 '''(INTERNAL) Helper for distances.
405 '''
406 return earth if isinstance(earth, Ellipsoid) else (
407 earth if isinstance(earth, Datum) else
408 _ellipsoidal_datum(earth, name=where.__name__)).ellipsoid # PYCHOK indent
411def equirectangular(lat1, lon1, lat2, lon2, radius=R_M, **options):
412 '''Compute the distance between two points using
413 the U{Equirectangular Approximation / Projection
414 <https://www.Movable-Type.co.UK/scripts/latlong.html#equirectangular>}.
416 @arg lat1: Start latitude (C{degrees}).
417 @arg lon1: Start longitude (C{degrees}).
418 @arg lat2: End latitude (C{degrees}).
419 @arg lon2: End longitude (C{degrees}).
420 @kwarg radius: Mean earth radius, ellipsoid or datum
421 (C{meter}, L{Ellipsoid}, L{Ellipsoid2},
422 L{Datum} or L{a_f2Tuple}).
423 @kwarg options: Optional keyword arguments for function
424 L{equirectangular_}.
426 @return: Distance (C{meter}, same units as B{C{radius}} or
427 the ellipsoid or datum axes).
429 @raise TypeError: Invalid B{C{radius}}.
431 @see: Function L{equirectangular_} for more details, the
432 available B{C{options}}, errors, restrictions and other,
433 approximate or accurate distance functions.
434 '''
435 d = sqrt(equirectangular_(Lat(lat1=lat1), Lon(lon1=lon1),
436 Lat(lat2=lat2), Lon(lon2=lon2),
437 **options).distance2) # PYCHOK 4 vs 2-3
438 return degrees2m(d, radius=_mean_radius(radius, lat1, lat2))
441def equirectangular_(lat1, lon1, lat2, lon2,
442 adjust=True, limit=45, wrap=False):
443 '''Compute the distance between two points using
444 the U{Equirectangular Approximation / Projection
445 <https://www.Movable-Type.co.UK/scripts/latlong.html#equirectangular>}.
447 This approximation is valid for short distance of several
448 hundred Km or Miles, see the B{C{limit}} keyword argument and
449 the L{LimitError}.
451 @arg lat1: Start latitude (C{degrees}).
452 @arg lon1: Start longitude (C{degrees}).
453 @arg lat2: End latitude (C{degrees}).
454 @arg lon2: End longitude (C{degrees}).
455 @kwarg adjust: Adjust the wrapped, unrolled longitudinal delta
456 by the cosine of the mean latitude (C{bool}).
457 @kwarg limit: Optional limit for lat- and longitudinal deltas
458 (C{degrees}) or C{None} or C{0} for unlimited.
459 @kwarg wrap: Wrap and L{pygeodesy.unroll180} longitudes (C{bool}).
461 @return: A L{Distance4Tuple}C{(distance2, delta_lat, delta_lon,
462 unroll_lon2)}.
464 @raise LimitError: If the lat- and/or longitudinal delta exceeds the
465 B{C{-limit..+limit}} range and L{pygeodesy.limiterrors}
466 set to C{True}.
468 @see: U{Local, flat earth approximation
469 <https://www.EdWilliams.org/avform.htm#flat>}, functions
470 L{equirectangular}, L{cosineAndoyerLambert},
471 L{cosineForsytheAndoyerLambert}, L{cosineLaw}, L{euclidean},
472 L{flatLocal}/L{hubeny}, L{flatPolar}, L{haversine}, L{thomas}
473 and L{vincentys} and methods L{Ellipsoid.distance2},
474 C{LatLon.distanceTo*} and C{LatLon.equirectangularTo}.
475 '''
476 d_lat = lat2 - lat1
477 d_lon, ulon2 = unroll180(lon1, lon2, wrap=wrap)
479 if limit and limit > 0 and limiterrors() and (fabs(d_lat) > limit or
480 fabs(d_lon) > limit):
481 t = unstr(equirectangular_, lat1, lon1, lat2, lon2, limit=limit)
482 raise LimitError('delta exceeds limit', txt=t)
484 if adjust: # scale delta lon
485 d_lon *= _scale_deg(lat1, lat2)
487 d2 = hypot2(d_lat, d_lon) # degrees squared!
488 return Distance4Tuple(d2, d_lat, d_lon, ulon2 - lon2)
491def euclidean(lat1, lon1, lat2, lon2, radius=R_M, adjust=True, wrap=False):
492 '''Approximate the C{Euclidean} distance between two (spherical) points.
494 @arg lat1: Start latitude (C{degrees}).
495 @arg lon1: Start longitude (C{degrees}).
496 @arg lat2: End latitude (C{degrees}).
497 @arg lon2: End longitude (C{degrees}).
498 @kwarg radius: Mean earth radius, ellipsoid or datum (C{meter},
499 L{Ellipsoid}, L{Ellipsoid2}, L{Datum} or L{a_f2Tuple}).
500 @kwarg adjust: Adjust the longitudinal delta by the cosine of the
501 mean latitude (C{bool}).
502 @kwarg wrap: Wrap and L{pygeodesy.unroll180} longitudes (C{bool}).
504 @return: Distance (C{meter}, same units as B{C{radius}} or the
505 ellipsoid or datum axes).
507 @raise TypeError: Invalid B{C{radius}}.
509 @see: U{Distance between two (spherical) points
510 <https://www.EdWilliams.org/avform.htm#Dist>}, functions L{euclid},
511 L{euclidean_}, L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert},
512 L{cosineLaw}, L{equirectangular}, L{flatLocal}/L{hubeny}, L{flatPolar},
513 L{haversine}, L{thomas} and L{vincentys} and methods L{Ellipsoid.distance2},
514 C{LatLon.distanceTo*} and C{LatLon.equirectangularTo}.
515 '''
516 return _distanceToS(euclidean_, lat1, lat2, radius,
517 *unroll180(lon1, lon2, wrap=wrap),
518 adjust=adjust)
521def euclidean_(phi2, phi1, lam21, adjust=True):
522 '''Approximate the I{angular} C{Euclidean} distance between two
523 (spherical) points.
525 @arg phi2: End latitude (C{radians}).
526 @arg phi1: Start latitude (C{radians}).
527 @arg lam21: Longitudinal delta, M{end-start} (C{radians}).
528 @kwarg adjust: Adjust the longitudinal delta by the cosine
529 of the mean latitude (C{bool}).
531 @return: Angular distance (C{radians}).
533 @see: Functions L{euclid}, L{euclidean}, L{cosineAndoyerLambert_},
534 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, L{equirectangular_},
535 L{flatLocal_}/L{hubeny_}, L{flatPolar_}, L{haversine_}, L{thomas_}
536 and L{vincentys_}.
537 '''
538 if adjust:
539 lam21 *= _scale_rad(phi2, phi1)
540 return euclid(phi2 - phi1, lam21)
543def excessAbc_(A, b, c):
544 '''Compute the I{spherical excess} C{E} of a (spherical) triangle
545 from two sides and the included angle.
547 @arg A: An interior triangle angle (C{radians}).
548 @arg b: Frist adjacent triangle side (C{radians}).
549 @arg c: Second adjacent triangle side (C{radians}).
551 @return: Spherical excess (C{radians}).
553 @raise UnitError: Invalid B{C{A}}, B{C{b}} or B{C{c}}.
555 @see: Functions L{excessGirard_}, L{excessLHuilier_} and U{Spherical
556 trigonometry<https://WikiPedia.org/wiki/Spherical_trigonometry>}.
557 '''
558 sA, cA, sb, cb, sc, cc = sincos2_(Radians_(A=A), Radians_(b=b) * _0_5,
559 Radians_(c=c) * _0_5)
560 return atan2(sA * sb * sc, cb * cc + cA * sb * sc) * _2_0
563def excessGirard_(A, B, C):
564 '''Compute the I{spherical excess} C{E} of a (spherical) triangle using
565 U{Girard's<https://MathWorld.Wolfram.com/GirardsSphericalExcessFormula.html>}
566 formula.
568 @arg A: First interior triangle angle (C{radians}).
569 @arg B: Second interior triangle angle (C{radians}).
570 @arg C: Third interior triangle angle (C{radians}).
572 @return: Spherical excess (C{radians}).
574 @raise UnitError: Invalid B{C{A}}, B{C{B}} or B{C{C}}.
576 @see: Function L{excessLHuilier_} and U{Spherical trigonometry
577 <https://WikiPedia.org/wiki/Spherical_trigonometry>}.
578 '''
579 return Radians(Girard=fsum_(Radians_(A=A),
580 Radians_(B=B),
581 Radians_(C=C), -PI))
584def excessLHuilier_(a, b, c):
585 '''Compute the I{spherical excess} C{E} of a (spherical) triangle using
586 U{L'Huilier's<https://MathWorld.Wolfram.com/LHuiliersTheorem.html>}
587 Theorem.
589 @arg a: First triangle side (C{radians}).
590 @arg b: Second triangle side (C{radians}).
591 @arg c: Third triangle side (C{radians}).
593 @return: Spherical excess (C{radians}).
595 @raise UnitError: Invalid B{C{a}}, B{C{b}} or B{C{c}}.
597 @see: Function L{excessGirard_} and U{Spherical trigonometry
598 <https://WikiPedia.org/wiki/Spherical_trigonometry>}.
599 '''
600 a = Radians_(a=a)
601 b = Radians_(b=b)
602 c = Radians_(c=c)
604 s = fsum_(a, b, c) * _0_5
605 r = tan_2(s) * tan_2(s - a) * tan_2(s - b) * tan_2(s - c)
606 r = atan(sqrt(r)) if r > 0 else _0_0
607 return Radians(LHuilier=r * _4_0)
610def excessKarney(lat1, lon1, lat2, lon2, radius=R_M, wrap=False):
611 '''Compute the surface area of a (spherical) quadrilateral bounded by a
612 segment of a great circle, two meridians and the equator using U{Karney's
613 <https://MathOverflow.net/questions/97711/the-area-of-spherical-polygons>}
614 method.
616 @arg lat1: Start latitude (C{degrees}).
617 @arg lon1: Start longitude (C{degrees}).
618 @arg lat2: End latitude (C{degrees}).
619 @arg lon2: End longitude (C{degrees}).
620 @kwarg radius: Mean earth radius, ellipsoid or datum (C{meter}, L{Ellipsoid},
621 L{Ellipsoid2}, L{Datum} or L{a_f2Tuple}) or C{None}.
622 @kwarg wrap: Wrap and L{pygeodesy.unroll180} longitudes (C{bool}).
624 @return: Surface area, I{signed} (I{square} C{meter} or the same units as
625 B{C{radius}} I{squared}) or the I{spherical excess} (C{radians})
626 if C{B{radius}=0} or C{None}.
628 @raise TypeError: Invalid B{C{radius}}.
630 @raise UnitError: Invalid B{C{lat2}} or B{C{lat1}}.
632 @raise ValueError: Semi-circular longitudinal delta.
634 @see: Functions L{excessKarney_} and L{excessQuad}.
635 '''
636 return _area_or_(excessKarney_, lat1, lat2, radius,
637 *unroll180(lon1, lon2, wrap=wrap))
640def excessKarney_(phi2, phi1, lam21):
641 '''Compute the I{spherical excess} C{E} of a (spherical) quadrilateral bounded
642 by a segment of a great circle, two meridians and the equator using U{Karney's
643 <https://MathOverflow.net/questions/97711/the-area-of-spherical-polygons>}
644 method.
646 @arg phi2: End latitude (C{radians}).
647 @arg phi1: Start latitude (C{radians}).
648 @arg lam21: Longitudinal delta, M{end-start} (C{radians}).
650 @return: Spherical excess, I{signed} (C{radians}).
652 @raise ValueError: Semi-circular longitudinal delta B{C{lam21}}.
654 @see: Function L{excessKarney} and U{Area of a spherical polygon
655 <https://MathOverflow.net/questions/97711/the-area-of-spherical-polygons>}.
656 '''
657 # from: Veness <https://www.Movable-Type.co.UK/scripts/latlong.html> Area
658 # method due to Karney: for each edge of the polygon,
659 #
660 # tan(Δλ / 2) · (tan(φ1 / 2) + tan(φ2 / 2))
661 # tan(E / 2) = -----------------------------------------
662 # 1 + tan(φ1 / 2) · tan(φ2 / 2)
663 #
664 # where E is the spherical excess of the trapezium obtained by extending
665 # the edge to the equator-circle vector for each edge (see also ***).
666 t2 = tan_2(phi2)
667 t1 = tan_2(phi1)
668 t = tan_2(lam21, lam21=None)
669 return Radians(Karney=atan2(t * (t1 + t2),
670 _1_0 + (t1 * t2)) * _2_0)
673# ***) Original post no longer available, following is a copy of the main part
674# <http://OSGeo-org.1560.x6.Nabble.com/Area-of-a-spherical-polygon-td3841625.html>
675#
676# The area of a polygon on a (unit) sphere is given by the spherical excess
677#
678# A = 2 * pi - sum(exterior angles)
679#
680# However this is badly conditioned if the polygon is small. In this case, use
681#
682# A = sum(S12{i, i+1}) over the edges of the polygon
683#
684# where S12 is the area of the quadrilateral bounded by an edge of the polygon,
685# two meridians and the equator, i.e. with vertices (phi1, lambda1), (phi2,
686# lambda2), (0, lambda1) and (0, lambda2). S12 is given by
687#
688# tan(S12 / 2) = tan(lambda21 / 2) * (tan(phi1 / 2) + tan(phi2 / 2)) /
689# (tan(phi1 / 2) * tan(phi2 / 2) + 1)
690#
691# = tan(lambda21 / 2) * tanh((Lambertian(phi1) +
692# Lambertian(phi2)) / 2)
693#
694# where lambda21 = lambda2 - lambda1 and lamb(x) is the Lambertian (or
695# inverse Gudermannian) function
696#
697# Lambertian(x) = asinh(tan(x)) = atanh(sin(x)) = 2 * atanh(tan(x / 2))
698#
699# Notes: The formula for S12 is exact, except that...
700# - it is indeterminate if an edge is a semi-circle
701# - the formula for A applies only if the polygon does not include a pole
702# (if it does, then add +/- 2 * pi to the result)
703# - in the limit of small phi and lambda, S12 reduces to the trapezoidal
704# formula, S12 = (lambda2 - lambda1) * (phi1 + phi2) / 2
705# - I derived this result from the equation for the area of a spherical
706# triangle in terms of two edges and the included angle given by, e.g.
707# U{Todhunter, I. - Spherical Trigonometry (1871), Sec. 103, Eq. (2)
708# <http://Books.Google.com/books?id=3uBHAAAAIAAJ&pg=PA71>}
709# - I would be interested to know if this formula for S12 is already known
710# - Charles Karney
713def excessQuad(lat1, lon1, lat2, lon2, radius=R_M, wrap=False):
714 '''Compute the surface area of a (spherical) quadrilateral bounded by a segment
715 of a great circle, two meridians and the equator.
717 @arg lat1: Start latitude (C{degrees}).
718 @arg lon1: Start longitude (C{degrees}).
719 @arg lat2: End latitude (C{degrees}).
720 @arg lon2: End longitude (C{degrees}).
721 @kwarg radius: Mean earth radius, ellipsoid or datum (C{meter},
722 L{Ellipsoid}, L{Ellipsoid2}, L{Datum} or L{a_f2Tuple}) or C{None}.
723 @kwarg wrap: Wrap and L{pygeodesy.unroll180} longitudes (C{bool}).
725 @return: Surface area, I{signed} (I{square} C{meter} or the same units as
726 B{C{radius}} I{squared}) or the I{spherical excess} (C{radians})
727 if C{B{radius}=0} or C{None}.
729 @raise TypeError: Invalid B{C{radius}}.
731 @raise UnitError: Invalid B{C{lat2}} or B{C{lat1}}.
733 @see: Function L{excessQuad_} and L{excessKarney}.
734 '''
735 return _area_or_(excessQuad_, lat1, lat2, radius,
736 *unroll180(lon1, lon2, wrap=wrap))
739def excessQuad_(phi2, phi1, lam21):
740 '''Compute the I{spherical excess} C{E} of a (spherical) quadrilateral bounded
741 by a segment of a great circle, two meridians and the equator.
743 @arg phi2: End latitude (C{radians}).
744 @arg phi1: Start latitude (C{radians}).
745 @arg lam21: Longitudinal delta, M{end-start} (C{radians}).
747 @return: Spherical excess, I{signed} (C{radians}).
749 @see: Function L{excessQuad}, U{Spherical trigonometry
750 <https://WikiPedia.org/wiki/Spherical_trigonometry>}.
751 '''
752 s = sin((phi2 + phi1) * _0_5)
753 c = cos((phi2 - phi1) * _0_5)
754 return Radians(Quad=atan2(tan_2(lam21) * s, c) * _2_0)
757def flatLocal(lat1, lon1, lat2, lon2, datum=_WGS84, wrap=False):
758 '''Compute the distance between two (ellipsoidal) points using
759 the U{ellipsoidal Earth to plane projection<https://WikiPedia.org/
760 wiki/Geographical_distance#Ellipsoidal_Earth_projected_to_a_plane>}
761 aka U{Hubeny<https://www.OVG.AT/de/vgi/files/pdf/3781/>} formula.
763 @arg lat1: Start latitude (C{degrees}).
764 @arg lon1: Start longitude (C{degrees}).
765 @arg lat2: End latitude (C{degrees}).
766 @arg lon2: End longitude (C{degrees}).
767 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid},
768 L{Ellipsoid2} or L{a_f2Tuple}).
769 @kwarg wrap: Wrap and L{pygeodesy.unroll180} longitudes (C{bool}).
771 @return: Distance (C{meter}, same units as the B{C{datum}}'s
772 ellipsoid axes).
774 @raise TypeError: Invalid B{C{datum}}.
776 @note: The meridional and prime_vertical radii of curvature
777 are taken and scaled at the mean of both latitude.
779 @see: Functions L{flatLocal_} or L{hubeny_}, L{cosineLaw}, L{flatPolar},
780 L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert},
781 L{equirectangular}, L{euclidean}, L{haversine}, L{thomas},
782 L{vincentys}, method L{Ellipsoid.distance2} and U{local, flat
783 earth approximation<https://www.EdWilliams.org/avform.htm#flat>}.
784 '''
785 d, _ = unroll180(lon1, lon2, wrap=wrap)
786 return flatLocal_(Phi_(lat2=lat2),
787 Phi_(lat1=lat1), radians(d), datum=datum)
789hubeny = flatLocal # PYCHOK for Karl Hubeny
792def flatLocal_(phi2, phi1, lam21, datum=_WGS84):
793 '''Compute the I{angular} distance between two (ellipsoidal) points using
794 the U{ellipsoidal Earth to plane projection<https://WikiPedia.org/
795 wiki/Geographical_distance#Ellipsoidal_Earth_projected_to_a_plane>}
796 aka U{Hubeny<https://www.OVG.AT/de/vgi/files/pdf/3781/>} formula.
798 @arg phi2: End latitude (C{radians}).
799 @arg phi1: Start latitude (C{radians}).
800 @arg lam21: Longitudinal delta, M{end-start} (C{radians}).
801 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid},
802 L{Ellipsoid2} or L{a_f2Tuple}).
804 @return: Angular distance (C{radians}).
806 @raise TypeError: Invalid B{C{datum}}.
808 @note: The meridional and prime_vertical radii of curvature
809 are taken and scaled I{at the mean of both latitude}.
811 @see: Functions L{flatLocal} or L{hubeny}, L{cosineAndoyerLambert_},
812 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, L{flatPolar_},
813 L{equirectangular_}, L{euclidean_}, L{haversine_}, L{thomas_}
814 and L{vincentys_} and U{local, flat earth approximation
815 <https://www.EdWilliams.org/avform.htm#flat>}.
816 '''
817 E = _ellipsoidal(datum, flatLocal_)
818 m, n = E.roc2_((phi2 + phi1) * _0_5, scaled=True)
819 return hypot(m * (phi2 - phi1), n * lam21)
821hubeny_ = flatLocal_ # PYCHOK for Karl Hubeny
824def flatPolar(lat1, lon1, lat2, lon2, radius=R_M, wrap=False):
825 '''Compute the distance between two (spherical) points using
826 the U{polar coordinate flat-Earth <https://WikiPedia.org/wiki/
827 Geographical_distance#Polar_coordinate_flat-Earth_formula>}
828 formula.
830 @arg lat1: Start latitude (C{degrees}).
831 @arg lon1: Start longitude (C{degrees}).
832 @arg lat2: End latitude (C{degrees}).
833 @arg lon2: End longitude (C{degrees}).
834 @kwarg radius: Mean earth radius, ellipsoid or datum (C{meter},
835 L{Ellipsoid}, L{Ellipsoid2}, L{Datum} or L{a_f2Tuple}).
836 @kwarg wrap: Wrap and L{pygeodesy.unroll180} longitudes (C{bool}).
838 @return: Distance (C{meter}, same units as B{C{radius}} or the
839 ellipsoid or datum axes).
841 @raise TypeError: Invalid B{C{radius}}.
843 @see: Functions L{flatPolar_}, L{cosineAndoyerLambert},
844 L{cosineForsytheAndoyerLambert},L{cosineLaw},
845 L{flatLocal}/L{hubeny}, L{equirectangular},
846 L{euclidean}, L{haversine}, L{thomas} and
847 L{vincentys}.
848 '''
849 return _distanceToS(flatPolar_, lat1, lat2, radius,
850 *unroll180(lon1, lon2, wrap=wrap))
853def flatPolar_(phi2, phi1, lam21):
854 '''Compute the I{angular} distance between two (spherical) points
855 using the U{polar coordinate flat-Earth<https://WikiPedia.org/wiki/
856 Geographical_distance#Polar_coordinate_flat-Earth_formula>}
857 formula.
859 @arg phi2: End latitude (C{radians}).
860 @arg phi1: Start latitude (C{radians}).
861 @arg lam21: Longitudinal delta, M{end-start} (C{radians}).
863 @return: Angular distance (C{radians}).
865 @see: Functions L{flatPolar}, L{cosineAndoyerLambert_},
866 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_},
867 L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_},
868 L{haversine_}, L{thomas_} and L{vincentys_}.
869 '''
870 a = fabs(PI_2 - phi1) # co-latitude
871 b = fabs(PI_2 - phi2) # co-latitude
872 if a < b:
873 a, b = b, a
874 if a < EPS0:
875 a = _0_0
876 elif b > 0:
877 b = b / a # /= chokes PyChecker
878 c = b * cos(lam21) * _2_0
879 c = fsum_(_1_0, b**2, -fabs(c))
880 a *= sqrt0(c)
881 return a
884def hartzell(pov, los=None, earth=_WGS84, name=NN, **LatLon_and_kwds):
885 '''Compute the intersection of the earth's surface and a Line-Of-Sight
886 from a Point-Of-View in space.
888 @arg pov: Point-Of-View outside the earth (C{Cartesian}, L{Ecef9Tuple}
889 or L{Vector3d}).
890 @kwarg los: Line-Of-Sight, I{direction} to earth (L{Vector3d}) or
891 C{None} to point to the earth' center.
892 @kwarg earth: The earth model (L{Datum}, L{Ellipsoid}, L{Ellipsoid2},
893 L{a_f2Tuple} or C{scalar} radius in C{meter}).
894 @kwarg name: Optional name (C{str}).
895 @kwarg LatLon_and_kwds: Optional C{LatLon} class for the intersection
896 point plus C{LatLon} keyword arguments, include
897 B{C{datum}} if different from B{C{earth}}.
899 @return: The earth intersection (L{Vector3d}, C{Cartesian type} of
900 B{C{pov}} or B{C{LatLon}}).
902 @raise IntersectionError: Null B{C{pov}} or B{C{los}} vector, B{C{pov}}
903 is inside the earth or B{C{los}} points outside
904 the earth or points in an opposite direction.
906 @raise TypeError: Invalid B{C{pov}}, B{C{los}} or B{C{earth}}.
908 @see: Function L{pygeodesy.hartzell4}, L{pygeodesy.tyr3d} for B{C{los}},
909 method L{Ellipsoid.hartzell4} and U{I{Satellite Line-of-Sight
910 Intersection with Earth}<https://StephenHartzell.Medium.com/
911 satellite-line-of-sight-intersection-with-earth-d786b4a6a9b6>}.
912 '''
913 D = earth if isinstance(earth, Datum) else \
914 _spherical_datum(earth, name=hartzell.__name__)
915 try:
916 r, _ = _MODS.triaxials._hartzell3d2(pov, los, D.ellipsoid._triaxial)
917 except Exception as x:
918 raise IntersectionError(pov=pov, los=los, earth=earth, cause=x)
920# else:
921# E = D.ellipsoid
922# # Triaxial(a, b, c) == (E.a, E.a, E.b)
923#
924# def _Error(txt):
925# return IntersectionError(pov=pov, los=los, earth=earth, txt=txt)
926#
927# a2 = b2 = E.a2 # earth' x, y, ...
928# c2 = E.b2 # ... z semi-axis squared
929# q2 = E.b2_a2 # == c2 / a2
930# bc = E.a * E.b # == b * c
931#
932# V3 = _MODS.vector3d._otherV3d
933# p3 = V3(pov=pov)
934# u3 = V3(los=los) if los else p3.negate()
935# u3 = u3.unit() # unit vector, opposing signs
936#
937# x2, y2, z2 = p3.x2y2z2 # p3.times_(p3).xyz
938# ux, vy, wz = u3.times_(p3).xyz
939# u2, v2, w2 = u3.x2y2z2 # u3.times_(u3).xyz
940#
941# t = c2, c2, b2
942# m = fdot(t, u2, v2, w2) # a2 factored out
943# if m < EPS0: # zero or near-null LOS vector
944# raise _Error(_near_(_null_))
945#
946# # a2 and b2 factored out, b2 == a2 and b2 / a2 == 1
947# r = fsum_(b2 * w2, c2 * v2, -v2 * z2, vy * wz * 2,
948# c2 * u2, -u2 * z2, -w2 * x2, ux * wz * 2,
949# -w2 * y2, -u2 * y2 * q2, -v2 * x2 * q2, ux * vy * 2 * q2, floats=True)
950# if r > 0:
951# r = sqrt(r) * bc # == a * a * b * c / a2
952# elif r < 0: # LOS pointing away from or missing the earth
953# raise _Error(_opposite_ if max(ux, vy, wz) > 0 else _outside_)
954#
955# d = Fdot(t, ux, vy, wz).fadd_(r).fover(m) # -r for antipode, a2 factored out
956# if d > 0: # POV inside or LOS missing, outside the earth
957# s = fsum_(_1_0, x2 / a2, y2 / b2, z2 / c2, _N_2_0, floats=True) # like _sideOf
958# raise _Error(_outside_ if s > 0 else _inside_)
959# elif fsum_(x2, y2, z2) < d**2: # d past earth center
960# raise _Error(_too_(_distant_))
961#
962# r = p3.minus(u3.times(d))
963# # h = p3.minus(r).length # distance to ellipsoid
965 r = _xnamed(r, name or hartzell.__name__)
966 if LatLon_and_kwds:
967 c = _MODS.cartesianBase.CartesianBase(r, datum=D, name=r.name)
968 r = c.toLatLon(**LatLon_and_kwds)
969 return r
972def haversine(lat1, lon1, lat2, lon2, radius=R_M, wrap=False):
973 '''Compute the distance between two (spherical) points using the
974 U{Haversine<https://www.Movable-Type.co.UK/scripts/latlong.html>}
975 formula.
977 @arg lat1: Start latitude (C{degrees}).
978 @arg lon1: Start longitude (C{degrees}).
979 @arg lat2: End latitude (C{degrees}).
980 @arg lon2: End longitude (C{degrees}).
981 @kwarg radius: Mean earth radius, ellipsoid or datum (C{meter},
982 L{Ellipsoid}, L{Ellipsoid2}, L{Datum} or L{a_f2Tuple}).
983 @kwarg wrap: Wrap and L{pygeodesy.unroll180} longitudes (C{bool}).
985 @return: Distance (C{meter}, same units as B{C{radius}}).
987 @raise TypeError: Invalid B{C{radius}}.
989 @see: U{Distance between two (spherical) points
990 <https://www.EdWilliams.org/avform.htm#Dist>}, functions
991 L{cosineLaw}, L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert},
992 L{equirectangular}, L{euclidean}, L{flatLocal}/L{hubeny}, L{flatPolar},
993 L{thomas} and L{vincentys} and methods L{Ellipsoid.distance2},
994 C{LatLon.distanceTo*} and C{LatLon.equirectangularTo}.
996 @note: See note at function L{vincentys_}.
997 '''
998 return _distanceToS(haversine_, lat1, lat2, radius,
999 *unroll180(lon1, lon2, wrap=wrap))
1002def haversine_(phi2, phi1, lam21):
1003 '''Compute the I{angular} distance between two (spherical) points
1004 using the U{Haversine<https://www.Movable-Type.co.UK/scripts/latlong.html>}
1005 formula.
1007 @arg phi2: End latitude (C{radians}).
1008 @arg phi1: Start latitude (C{radians}).
1009 @arg lam21: Longitudinal delta, M{end-start} (C{radians}).
1011 @return: Angular distance (C{radians}).
1013 @see: Functions L{haversine}, L{cosineAndoyerLambert_},
1014 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_},
1015 L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_},
1016 L{flatPolar_}, L{thomas_} and L{vincentys_}.
1018 @note: See note at function L{vincentys_}.
1019 '''
1020 def _hsin(rad):
1021 return sin(rad * _0_5)**2
1023 h = _hsin(phi2 - phi1) + cos(phi1) * cos(phi2) * _hsin(lam21) # haversine
1024 return atan2(sqrt0(h), sqrt0(_1_0 - h)) * _2_0 # == asin(sqrt(h)) * 2
1027def heightOf(angle, distance, radius=R_M):
1028 '''Determine the height above the (spherical) earth' surface after
1029 traveling along a straight line at a given tilt.
1031 @arg angle: Tilt angle above horizontal (C{degrees}).
1032 @arg distance: Distance along the line (C{meter} or same units as
1033 B{C{radius}}).
1034 @kwarg radius: Optional mean earth radius (C{meter}).
1036 @return: Height (C{meter}, same units as B{C{distance}} and B{C{radius}}).
1038 @raise ValueError: Invalid B{C{angle}}, B{C{distance}} or B{C{radius}}.
1040 @see: U{MultiDop geog_lib.GeogBeamHt<https://GitHub.com/NASA/MultiDop>}
1041 (U{Shapiro et al. 2009, JTECH
1042 <https://Journals.AMetSoc.org/doi/abs/10.1175/2009JTECHA1256.1>}
1043 and U{Potvin et al. 2012, JTECH
1044 <https://Journals.AMetSoc.org/doi/abs/10.1175/JTECH-D-11-00019.1>}).
1045 '''
1046 r = h = Radius(radius)
1047 d = fabs(Distance(distance))
1048 if d > h:
1049 d, h = h, d
1051 if d > EPS0: # and h > EPS0
1052 d = d / h # /= h chokes PyChecker
1053 s = sin(Phi_(angle=angle, clip=_180_0))
1054 s = fsum_(_1_0, _2_0 * s * d, d**2)
1055 if s > 0:
1056 return h * sqrt(s) - r
1058 raise _ValueError(angle=angle, distance=distance, radius=radius)
1061def horizon(height, radius=R_M, refraction=False):
1062 '''Determine the distance to the horizon from a given altitude
1063 above the (spherical) earth.
1065 @arg height: Altitude (C{meter} or same units as B{C{radius}}).
1066 @kwarg radius: Optional mean earth radius (C{meter}).
1067 @kwarg refraction: Consider atmospheric refraction (C{bool}).
1069 @return: Distance (C{meter}, same units as B{C{height}} and B{C{radius}}).
1071 @raise ValueError: Invalid B{C{height}} or B{C{radius}}.
1073 @see: U{Distance to horizon<https://www.EdWilliams.org/avform.htm#Horizon>}.
1074 '''
1075 h, r = Height(height), Radius(radius)
1076 if min(h, r) < 0:
1077 raise _ValueError(height=height, radius=radius)
1079 if refraction:
1080 d2 = 2.415750694528 * h * r # 2.0 / 0.8279
1081 else:
1082 d2 = h * fsum_(r, r, h)
1083 return sqrt0(d2)
1086def _idlmn5(datum, lat1, lon1, lat2, lon2, small, wrap, s):
1087 '''(INTERNAL) Helper for C{intersection2} and C{intersections2}.
1088 '''
1089 m = small if small is _100km else Meter_(small=small)
1090 n = (intersections2 if s else intersection2).__name__
1091 if datum is None or euclidean(lat1, lon1, lat2, lon2, wrap=wrap) < m:
1092 d, m = None, _MODS.vector3d
1093 _i = m._intersects2 if s else m._intersect3d3
1094 _, lon2 = unroll180(lon1, lon2, wrap=wrap)
1095 elif isscalar(datum) and datum < 0 and not s:
1096 d = _spherical_datum(-datum, name=n)
1097 m = _MODS.sphericalNvector
1098 _i = m.intersection
1099 else:
1100 d = _spherical_datum(datum, name=n)
1101 if d.isSpherical:
1102 m = _MODS.sphericalTrigonometry
1103 _i = m._intersects2 if s else m._intersect
1104 elif d.isEllipsoidal:
1105 try:
1106 if d.ellipsoid.geodesic:
1107 pass
1108 m = _MODS.ellipsoidalKarney
1109 except ImportError:
1110 m = _MODS.ellipsoidalExact
1111 _i = m._intersections2 if s else m._intersection3 # ellispoidalBaseDi
1112 else:
1113 raise _AssertionError(datum=datum)
1114 return _i, d, lon2, m, n
1117def intersection2(lat1, lon1, bearing1,
1118 lat2, lon2, bearing2, datum=None, wrap=True, small=_100km):
1119 '''I{Conveniently} compute the intersection of two lines each defined
1120 by a (geodetic) point and a bearing from North, using either ...
1122 1) L{vector3d.intersection3d3} for B{C{small}} distances (below 100 Km
1123 or about 0.88 degrees) or if no B{C{datum}} is specified, or ...
1125 2) L{sphericalTrigonometry.intersection} for a spherical B{C{datum}}
1126 or if B{C{datum}} is a C{scalar} representing the earth
1127 radius, conventionally in C{meter} or ...
1129 3) L{sphericalNvector.intersection} if B{C{datum}} is a I{negative}
1130 C{scalar}, (negative) earth radius, conventionally in C{meter} or ...
1132 4) L{ellipsoidalKarney.intersection3} for an ellipsoidal B{C{datum}}
1133 and if I{Karney}'s U{geographiclib<https://PyPI.org/project/geographiclib>}
1134 is installed, otherwise ...
1136 5) L{ellipsoidalExact.intersection3}, provided B{C{datum}} is ellipsoidal.
1138 @arg lat1: Latitude of the first point (C{degrees}).
1139 @arg lon1: Longitude of the first point (C{degrees}).
1140 @arg bearing1: Bearing at the first point (compass C{degrees}).
1141 @arg lat2: Latitude of the second point (C{degrees}).
1142 @arg lon2: Longitude of the second point (C{degrees}).
1143 @arg bearing2: Bearing at the second point (compass C{degrees}).
1144 @kwarg datum: Optional ellipsoidal or spherical datum (L{Datum},
1145 L{Ellipsoid}, L{Ellipsoid2}, L{a_f2Tuple} or
1146 C{scalar} earth radius in C{meter}) or C{None}.
1147 @kwarg wrap: Wrap and unroll longitudes (C{bool}).
1148 @kwarg small: Upper limit for small distances (C{meter}).
1150 @return: A L{LatLon2Tuple}C{(lat, lon)} with the lat- and
1151 longitude of the intersection point.
1153 @raise IntersectionError: Ambiguous or infinite intersection
1154 or colinear, parallel or otherwise
1155 non-intersecting lines.
1157 @raise TypeError: Invalid B{C{datum}}.
1159 @raise UnitError: Invalid B{C{lat1}}, B{C{lon1}}, B{C{bearing1}},
1160 B{C{lat2}}, B{C{lon2}} or B{C{bearing2}}.
1162 @see: Method L{RhumbLine.intersection2}.
1164 @note: The returned intersections may be near-antipodal.
1165 '''
1166 b1, b2 = Bearing(bearing1=bearing1), Bearing(bearing2=bearing2)
1167 try:
1168 _i, d, l2, m, n = _idlmn5(datum, lat1, lon1, lat2, lon2, small, wrap,
1169 False)
1170 if d is None:
1171 t, _, _ = _i(m.Vector3d(lon1, lat1, 0), b1,
1172 m.Vector3d(l2, lat2, 0), b2, useZ=False)
1173 t = LatLon2Tuple(t.y, t.x, name=n)
1175 else:
1176 t = _i(m.LatLon(lat1, lon1, datum=d), b1,
1177 m.LatLon(lat2, lon2, datum=d), b2, height=0, wrap=wrap)
1178 if isinstance(t, Intersection3Tuple): # ellipsoidal
1179 t, _, _ = t
1180 t = LatLon2Tuple(t.lat, t.lon, name=n)
1182 except (TypeError, ValueError) as x:
1183 raise _xError(x, lat1=lat1, lon1=lon1, bearing1=bearing1,
1184 lat2=lat2, lon2=lon2, bearing2=bearing2,
1185 datum=datum, small=small, wrap=wrap)
1186 return t
1189def intersections2(lat1, lon1, radius1,
1190 lat2, lon2, radius2, datum=None, wrap=True, small=_100km):
1191 '''I{Conveniently} compute the intersections of two circles each defined
1192 by a (geodetic) center point and a radius, using either ...
1194 1) L{vector3d.intersections2} for B{C{small}} distances (below 100 Km
1195 or about 0.88 degrees) or if no B{C{datum}} is specified, or ...
1197 2) L{sphericalTrigonometry.intersections2} for a spherical B{C{datum}}
1198 or if B{C{datum}} is a C{scalar} representing the earth radius,
1199 conventionally in C{meter} or ...
1201 3) L{ellipsoidalKarney.intersections2} for an ellipsoidal B{C{datum}}
1202 and if I{Karney}'s U{geographiclib<https://PyPI.org/project/geographiclib>}
1203 is installed, otherwise ...
1205 4) L{ellipsoidalExact.intersections2}, provided B{C{datum}} is ellipsoidal.
1207 @arg lat1: Latitude of the first circle center (C{degrees}).
1208 @arg lon1: Longitude of the first circle center (C{degrees}).
1209 @arg radius1: Radius of the first circle (C{meter}, conventionally).
1210 @arg lat2: Latitude of the second circle center (C{degrees}).
1211 @arg lon2: Longitude of the second circle center (C{degrees}).
1212 @arg radius2: Radius of the second circle (C{meter}, same units as B{C{radius1}}).
1213 @kwarg datum: Optional ellipsoidal or spherical datum (L{Datum},
1214 L{Ellipsoid}, L{Ellipsoid2}, L{a_f2Tuple} or
1215 C{scalar} earth radius in C{meter}, same units as
1216 B{C{radius1}} and B{C{radius2}}) or C{None}.
1217 @kwarg wrap: Wrap and unroll longitudes (C{bool}).
1218 @kwarg small: Upper limit for small distances (C{meter}).
1220 @return: 2-Tuple of the intersection points, each a
1221 L{LatLon2Tuple}C{(lat, lon)}. For abutting circles, the
1222 points are the same instance, aka the I{radical center}.
1224 @raise IntersectionError: Concentric, antipodal, invalid or
1225 non-intersecting circles or no
1226 convergence.
1228 @raise TypeError: Invalid B{C{datum}}.
1230 @raise UnitError: Invalid B{C{lat1}}, B{C{lon1}}, B{C{radius1}},
1231 B{C{lat2}}, B{C{lon2}} or B{C{radius2}}.
1232 '''
1233 r1, r2 = Radius_(radius1=radius1), Radius_(radius2=radius2)
1234 try:
1235 _i, d, l2, m, n = _idlmn5(datum, lat1, lon1, lat2, lon2, small, wrap,
1236 True)
1237 if d is None:
1238 r1 = m2degrees(r1, radius=R_M, lat=lat1)
1239 r2 = m2degrees(r2, radius=R_M, lat=lat2)
1241 def _V2T(x, y, _, **unused): # _ == z unused
1242 return LatLon2Tuple(y, x, name=n)
1244 t = _i(m.Vector3d(lon1, lat1, 0), r1,
1245 m.Vector3d(l2, lat2, 0), r2, sphere=False,
1246 Vector=_V2T)
1247 else:
1249 def _LL2T(lat, lon, **unused):
1250 return LatLon2Tuple(lat, lon, name=n)
1252 t = _i(m.LatLon(lat1, lon1, datum=d), r1,
1253 m.LatLon(lat2, lon2, datum=d), r2,
1254 LatLon=_LL2T, height=0, wrap=wrap)
1256 except (TypeError, ValueError) as x:
1257 raise _xError(x, lat1=lat1, lon1=lon1, radius1=radius1,
1258 lat2=lat2, lon2=lon2, radius2=radius2,
1259 datum=datum, small=small, wrap=wrap)
1260 return t
1263def isantipode(lat1, lon1, lat2, lon2, eps=EPS):
1264 '''Check whether two points are antipodal, on diametrically
1265 opposite sides of the earth.
1267 @arg lat1: Latitude of one point (C{degrees}).
1268 @arg lon1: Longitude of one point (C{degrees}).
1269 @arg lat2: Latitude of the other point (C{degrees}).
1270 @arg lon2: Longitude of the other point (C{degrees}).
1271 @kwarg eps: Tolerance for near-equality (C{degrees}).
1273 @return: C{True} if points are antipodal within the
1274 B{C{eps}} tolerance, C{False} otherwise.
1276 @see: Functions L{isantipode_} and L{antipode}.
1277 '''
1278 return True if (fabs(lat1 + lat2) <= eps and
1279 fabs(lon1 + lon2) <= eps) else \
1280 _MODS.latlonBase._isequalTo(antipode(lat1, lon1),
1281 normal(lat2, lon2), eps=eps)
1284def isantipode_(phi1, lam1, phi2, lam2, eps=EPS):
1285 '''Check whether two points are antipodal, on diametrically
1286 opposite sides of the earth.
1288 @arg phi1: Latitude of one point (C{radians}).
1289 @arg lam1: Longitude of one point (C{radians}).
1290 @arg phi2: Latitude of the other point (C{radians}).
1291 @arg lam2: Longitude of the other point (C{radians}).
1292 @kwarg eps: Tolerance for near-equality (C{radians}).
1294 @return: C{True} if points are antipodal within the
1295 B{C{eps}} tolerance, C{False} otherwise.
1297 @see: Functions L{isantipode} and L{antipode_}.
1298 '''
1299 return True if (fabs(phi1 + phi2) <= eps and
1300 fabs(lam1 + lam2) <= eps) else \
1301 _MODS.latlonBase._isequalTo_(antipode_(phi1, lam1),
1302 normal_(phi2, lam2), eps=eps)
1305def isnormal(lat, lon, eps=0):
1306 '''Check whether B{C{lat}} I{and} B{C{lon}} are within the I{normal}
1307 range in C{degrees}.
1309 @arg lat: Latitude (C{degrees}).
1310 @arg lon: Longitude (C{degrees}).
1311 @kwarg eps: Optional tolerance below C{90} and C{180 degrees}.
1313 @return: C{True} if C{(abs(B{lat}) + B{eps}) <= 90} and
1314 C{(abs(B{lon}) + B{eps}) <= 180}, C{False} othwerwise.
1316 @see: Functions L{isnormal_} and L{normal}.
1317 '''
1318 return (_90_0 - fabs(lat)) >= eps and (_180_0 - fabs(lon)) >= eps
1321def isnormal_(phi, lam, eps=0):
1322 '''Check whether B{C{phi}} I{and} B{C{lam}} are within the I{normal}
1323 range in C{radians}.
1325 @arg phi: Latitude (C{radians}).
1326 @arg lam: Longitude (C{radians}).
1327 @kwarg eps: Optional tolerance below C{PI/2} and C{PI radians}.
1329 @return: C{True} if C{(abs(B{phi}) + B{eps}) <= PI/2} and
1330 C{(abs(B{lam}) + B{eps}) <= PI}, C{False} othwerwise.
1332 @see: Functions L{isnormal} and L{normal_}.
1333 '''
1334 return (PI_2 - fabs(phi)) >= eps and (PI - fabs(lam)) >= eps
1337def latlon2n_xyz(lat, lon, name=NN):
1338 '''Convert lat-, longitude to C{n-vector} (I{normal} to the
1339 earth's surface) X, Y and Z components.
1341 @arg lat: Latitude (C{degrees}).
1342 @arg lon: Longitude (C{degrees}).
1343 @kwarg name: Optional name (C{str}).
1345 @return: A L{Vector3Tuple}C{(x, y, z)}.
1347 @see: Function L{philam2n_xyz}.
1349 @note: These are C{n-vector} x, y and z components,
1350 I{NOT} geocentric ECEF x, y and z coordinates!
1351 '''
1352 return _2n_xyz(name, *sincos2d_(lat, lon))
1355def _normal2(a, b, n_2, n, n2):
1356 '''(INTERNAL) Helper for C{normal} and C{normal_}.
1357 '''
1358 if fabs(b) > n:
1359 b = remainder(b, n2)
1360 r = remainder(a, n) if fabs(a) > n_2 else a
1361 if r != a:
1362 r = -r
1363 b -= n if b > 0 else -n
1364 return float0_(r, b)
1367def normal(lat, lon, name=NN):
1368 '''Normalize a lat- I{and} longitude pair in C{degrees}.
1370 @arg lat: Latitude (C{degrees}).
1371 @arg lon: Longitude (C{degrees}).
1372 @kwarg name: Optional name (C{str}).
1374 @return: L{LatLon2Tuple}C{(lat, lon)} with C{abs(lat) <= 90}
1375 and C{abs(lon) <= 180}.
1377 @see: Functions L{normal_} and L{isnormal}.
1378 '''
1379 return LatLon2Tuple(*_normal2(lat, lon, _90_0, _180_0, _360_0), name=name)
1382def normal_(phi, lam, name=NN):
1383 '''Normalize a lat- I{and} longitude pair in C{radians}.
1385 @arg phi: Latitude (C{radians}).
1386 @arg lam: Longitude (C{radians}).
1387 @kwarg name: Optional name (C{str}).
1389 @return: L{PhiLam2Tuple}C{(phi, lam)} with C{abs(phi) <= PI/2}
1390 and C{abs(lam) <= PI}.
1392 @see: Functions L{normal} and L{isnormal_}.
1393 '''
1394 return PhiLam2Tuple(*_normal2(phi, lam, PI_2, PI, PI2), name=name)
1397def _2n_xyz(name, sa, ca, sb, cb):
1398 '''(INTERNAL) Helper for C{latlon2n_xyz} and C{philam2n_xyz}.
1399 '''
1400 # Kenneth Gade eqn 3, but using right-handed
1401 # vector x -> 0°E,0°N, y -> 90°E,0°N, z -> 90°N
1402 return Vector3Tuple(ca * cb, ca * sb, sa, name=name)
1405def n_xyz2latlon(x, y, z, name=NN):
1406 '''Convert C{n-vector} components to lat- and longitude in C{degrees}.
1408 @arg x: X component (C{scalar}).
1409 @arg y: Y component (C{scalar}).
1410 @arg z: Z component (C{scalar}).
1411 @kwarg name: Optional name (C{str}).
1413 @return: A L{LatLon2Tuple}C{(lat, lon)}.
1415 @see: Function L{n_xyz2philam}.
1416 '''
1417 return LatLon2Tuple(atan2d(z, hypot(x, y)), atan2d(y, x), name=name)
1420def n_xyz2philam(x, y, z, name=NN):
1421 '''Convert C{n-vector} components to lat- and longitude in C{radians}.
1423 @arg x: X component (C{scalar}).
1424 @arg y: Y component (C{scalar}).
1425 @arg z: Z component (C{scalar}).
1426 @kwarg name: Optional name (C{str}).
1428 @return: A L{PhiLam2Tuple}C{(phi, lam)}.
1430 @see: Function L{n_xyz2latlon}.
1431 '''
1432 return PhiLam2Tuple(atan2(z, hypot(x, y)), atan2(y, x), name=name)
1435def _opposes(d, m, n, n2):
1436 '''(INETNAL) Helper for C{opposing} and C{opposing_}.
1437 '''
1438 d = d % n2 # -20 % 360 == 340, -1 % PI2 == PI2 - 1
1439 return False if d < m or d > (n2 - m) else (
1440 True if (n - m) < d < (n + m) else None)
1443def opposing(bearing1, bearing2, margin=_90_0):
1444 '''Compare the direction of two bearings given in C{degrees}.
1446 @arg bearing1: First bearing (compass C{degrees}).
1447 @arg bearing2: Second bearing (compass C{degrees}).
1448 @kwarg margin: Optional, interior angle bracket (C{degrees}).
1450 @return: C{True} if both bearings point in opposite, C{False} if
1451 in similar or C{None} if in perpendicular directions.
1453 @see: Function L{opposing_}.
1454 '''
1455 m = Degrees_(margin=margin, low=EPS0, high=_90_0)
1456 return _opposes(bearing2 - bearing1, m,_180_0, _360_0)
1459def opposing_(radians1, radians2, margin=PI_2):
1460 '''Compare the direction of two bearings given in C{radians}.
1462 @arg radians1: First bearing (C{radians}).
1463 @arg radians2: Second bearing (C{radians}).
1464 @kwarg margin: Optional, interior angle bracket (C{radians}).
1466 @return: C{True} if both bearings point in opposite, C{False} if
1467 in similar or C{None} if in perpendicular directions.
1469 @see: Function L{opposing}.
1470 '''
1471 m = Radians_(margin=margin, low=EPS0, high=PI_2)
1472 return _opposes(radians2 - radians1, m, PI, PI2)
1475def philam2n_xyz(phi, lam, name=NN):
1476 '''Convert lat-, longitude to C{n-vector} (I{normal} to the
1477 earth's surface) X, Y and Z components.
1479 @arg phi: Latitude (C{radians}).
1480 @arg lam: Longitude (C{radians}).
1481 @kwarg name: Optional name (C{str}).
1483 @return: A L{Vector3Tuple}C{(x, y, z)}.
1485 @see: Function L{latlon2n_xyz}.
1487 @note: These are C{n-vector} x, y and z components,
1488 I{NOT} geocentric ECEF x, y and z coordinates!
1489 '''
1490 return _2n_xyz(name, *sincos2_(phi, lam))
1493def _radical2(d, r1, r2): # in .ellipsoidalBaseDI, .sphericalTrigonometry, .vector3d
1494 # (INTERNAL) See C{radical2} below
1495 # assert d > EPS0
1496 r = fsum_(_1_0, (r1 / d)**2, -(r2 / d)**2) * _0_5
1497 return Radical2Tuple(max(_0_0, min(_1_0, r)), r * d)
1500def radical2(distance, radius1, radius2):
1501 '''Compute the I{radical ratio} and I{radical line} of two
1502 U{intersecting circles<https://MathWorld.Wolfram.com/
1503 Circle-CircleIntersection.html>}.
1505 The I{radical line} is perpendicular to the axis thru the
1506 centers of the circles at C{(0, 0)} and C{(B{distance}, 0)}.
1508 @arg distance: Distance between the circle centers (C{scalar}).
1509 @arg radius1: Radius of the first circle (C{scalar}).
1510 @arg radius2: Radius of the second circle (C{scalar}).
1512 @return: A L{Radical2Tuple}C{(ratio, xline)} where C{0.0 <=
1513 ratio <= 1.0} and C{xline} is along the B{C{distance}}.
1515 @raise IntersectionError: The B{C{distance}} exceeds the sum
1516 of B{C{radius1}} and B{C{radius2}}.
1518 @raise UnitError: Invalid B{C{distance}}, B{C{radius1}} or
1519 B{C{radius2}}.
1521 @see: U{Circle-Circle Intersection
1522 <https://MathWorld.Wolfram.com/Circle-CircleIntersection.html>}.
1523 '''
1524 d = Distance_(distance, low=_0_0)
1525 r1 = Radius_(radius1=radius1)
1526 r2 = Radius_(radius2=radius2)
1527 if d > (r1 + r2):
1528 raise IntersectionError(distance=d, radius1=r1, radius2=r2,
1529 txt=_too_(_distant_))
1530 return _radical2(d, r1, r2) if d > EPS0 else \
1531 Radical2Tuple(_0_5, _0_0)
1534class Radical2Tuple(_NamedTuple):
1535 '''2-Tuple C{(ratio, xline)} of the I{radical} C{ratio} and
1536 I{radical} C{xline}, both C{scalar} and C{0.0 <= ratio <= 1.0}
1537 '''
1538 _Names_ = (_ratio_, _xline_)
1539 _Units_ = ( Scalar, Scalar)
1542def _scale_deg(lat1, lat2): # degrees
1543 # scale factor cos(mean of lats) for delta lon
1544 m = fabs(lat1 + lat2) * _0_5
1545 return cos(radians(m)) if m < 90 else _0_0
1548def _scale_rad(phi1, phi2): # radians, by .frechet, .hausdorff, .heights
1549 # scale factor cos(mean of phis) for delta lam
1550 m = fabs(phi1 + phi2) * _0_5
1551 return cos(m) if m < PI_2 else _0_0
1554def _sincosa6(phi2, phi1, lam21):
1555 '''(INTERNAL) C{sin}es, C{cos}ines and C{acos}ine.
1556 '''
1557 s2, c2, s1, c1, _, c21 = sincos2_(phi2, phi1, lam21)
1558 return s2, c2, s1, c1, acos1(s1 * s2 + c1 * c2 * c21), c21
1561def thomas(lat1, lon1, lat2, lon2, datum=_WGS84, wrap=False):
1562 '''Compute the distance between two (ellipsoidal) points using
1563 U{Thomas'<https://apps.DTIC.mil/dtic/tr/fulltext/u2/703541.pdf>}
1564 formula.
1566 @arg lat1: Start latitude (C{degrees}).
1567 @arg lon1: Start longitude (C{degrees}).
1568 @arg lat2: End latitude (C{degrees}).
1569 @arg lon2: End longitude (C{degrees}).
1570 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid},
1571 L{Ellipsoid2} or L{a_f2Tuple}).
1572 @kwarg wrap: Wrap and L{pygeodesy.unroll180} longitudes (C{bool}).
1574 @return: Distance (C{meter}, same units as the B{C{datum}}'s
1575 ellipsoid axes).
1577 @raise TypeError: Invalid B{C{datum}}.
1579 @see: Functions L{thomas_}, L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert},
1580 L{cosineLaw}, L{equirectangular}, L{euclidean}, L{flatLocal}/L{hubeny},
1581 L{flatPolar}, L{haversine}, L{vincentys} and method L{Ellipsoid.distance2}.
1582 '''
1583 return _distanceToE(thomas_, lat1, lat2, datum,
1584 *unroll180(lon1, lon2, wrap=wrap))
1587def thomas_(phi2, phi1, lam21, datum=_WGS84):
1588 '''Compute the I{angular} distance between two (ellipsoidal) points using
1589 U{Thomas'<https://apps.DTIC.mil/dtic/tr/fulltext/u2/703541.pdf>}
1590 formula.
1592 @arg phi2: End latitude (C{radians}).
1593 @arg phi1: Start latitude (C{radians}).
1594 @arg lam21: Longitudinal delta, M{end-start} (C{radians}).
1595 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid},
1596 L{Ellipsoid2} or L{a_f2Tuple}).
1598 @return: Angular distance (C{radians}).
1600 @raise TypeError: Invalid B{C{datum}}.
1602 @see: Functions L{thomas}, L{cosineAndoyerLambert_},
1603 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_},
1604 L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_},
1605 L{flatPolar_}, L{haversine_} and L{vincentys_} and U{Geodesy-PHP
1606 <https://GitHub.com/jtejido/geodesy-php/blob/master/src/Geodesy/
1607 Distance/ThomasFormula.php>}.
1608 '''
1609 s2, c2, s1, c1, r, _ = _sincosa6(phi2, phi1, lam21)
1610 if r and isnon0(c1) and isnon0(c2):
1611 E = _ellipsoidal(datum, thomas_)
1612 if E.f:
1613 r1 = atan2(E.b_a * s1, c1)
1614 r2 = atan2(E.b_a * s2, c2)
1616 j = (r2 + r1) * _0_5
1617 k = (r2 - r1) * _0_5
1618 sj, cj, sk, ck, h, _ = sincos2_(j, k, lam21 * _0_5)
1620 h = fsum_(sk**2, (ck * h)**2, -(sj * h)**2)
1621 u = _1_0 - h
1622 if isnon0(u) and isnon0(h):
1623 r = atan(sqrt0(h / u)) * _2_0 # == acos(1 - 2 * h)
1624 sr, cr = sincos2(r)
1625 if isnon0(sr):
1626 u = 2 * (sj * ck)**2 / u
1627 h = 2 * (sk * cj)**2 / h
1628 x = u + h
1629 y = u - h
1631 s = r / sr
1632 e = 4 * s**2
1633 d = 2 * cr
1634 a = e * d
1635 b = 2 * r
1636 c = s - (a - d) * _0_5
1637 f = E.f * _0_25
1639 t = fsum_(a * x, -b * y, c * x**2, -d * y**2, e * x * y)
1640 r -= fsum_(s * x, -y, -t * f * _0_25) * f * sr
1641 return r
1644def vincentys(lat1, lon1, lat2, lon2, radius=R_M, wrap=False):
1645 '''Compute the distance between two (spherical) points using
1646 U{Vincenty's<https://WikiPedia.org/wiki/Great-circle_distance>}
1647 spherical formula.
1649 @arg lat1: Start latitude (C{degrees}).
1650 @arg lon1: Start longitude (C{degrees}).
1651 @arg lat2: End latitude (C{degrees}).
1652 @arg lon2: End longitude (C{degrees}).
1653 @kwarg radius: Mean earth radius, ellipsoid or datum (C{meter},
1654 L{Ellipsoid}, L{Ellipsoid2}, L{Datum} or L{a_f2Tuple}).
1655 @kwarg wrap: Wrap and L{pygeodesy.unroll180} longitudes (C{bool}).
1657 @return: Distance (C{meter}, same units as B{C{radius}}).
1659 @raise UnitError: Invalid B{C{radius}}.
1661 @see: Functions L{vincentys_}, L{cosineAndoyerLambert},
1662 L{cosineForsytheAndoyerLambert},L{cosineLaw}, L{equirectangular},
1663 L{euclidean}, L{flatLocal}/L{hubeny}, L{flatPolar},
1664 L{haversine} and L{thomas} and methods L{Ellipsoid.distance2},
1665 C{LatLon.distanceTo*} and C{LatLon.equirectangularTo}.
1667 @note: See note at function L{vincentys_}.
1668 '''
1669 return _distanceToS(vincentys_, lat1, lat2, radius,
1670 *unroll180(lon1, lon2, wrap=wrap))
1673def vincentys_(phi2, phi1, lam21):
1674 '''Compute the I{angular} distance between two (spherical) points using
1675 U{Vincenty's<https://WikiPedia.org/wiki/Great-circle_distance>}
1676 spherical formula.
1678 @arg phi2: End latitude (C{radians}).
1679 @arg phi1: Start latitude (C{radians}).
1680 @arg lam21: Longitudinal delta, M{end-start} (C{radians}).
1682 @return: Angular distance (C{radians}).
1684 @see: Functions L{vincentys}, L{cosineAndoyerLambert_},
1685 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_},
1686 L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_},
1687 L{flatPolar_}, L{haversine_} and L{thomas_}.
1689 @note: Functions L{vincentys_}, L{haversine_} and L{cosineLaw_}
1690 produce equivalent results, but L{vincentys_} is suitable
1691 for antipodal points and slightly more expensive (M{3 cos,
1692 3 sin, 1 hypot, 1 atan2, 6 mul, 2 add}) than L{haversine_}
1693 (M{2 cos, 2 sin, 2 sqrt, 1 atan2, 5 mul, 1 add}) and
1694 L{cosineLaw_} (M{3 cos, 3 sin, 1 acos, 3 mul, 1 add}).
1695 '''
1696 s1, c1, s2, c2, s21, c21 = sincos2_(phi1, phi2, lam21)
1698 c = c2 * c21
1699 x = s1 * s2 + c1 * c
1700 y = c1 * s2 - s1 * c
1701 return atan2(hypot(c2 * s21, y), x)
1703# **) MIT License
1704#
1705# Copyright (C) 2016-2023 -- mrJean1 at Gmail -- All Rights Reserved.
1706#
1707# Permission is hereby granted, free of charge, to any person obtaining a
1708# copy of this software and associated documentation files (the "Software"),
1709# to deal in the Software without restriction, including without limitation
1710# the rights to use, copy, modify, merge, publish, distribute, sublicense,
1711# and/or sell copies of the Software, and to permit persons to whom the
1712# Software is furnished to do so, subject to the following conditions:
1713#
1714# The above copyright notice and this permission notice shall be included
1715# in all copies or substantial portions of the Software.
1716#
1717# THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
1718# OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
1719# FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
1720# THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR
1721# OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE,
1722# ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR
1723# OTHER DEALINGS IN THE SOFTWARE.