1 2# -*- coding: utf-8 -*- 3 4u'''Formulary of basic geodesy functions and approximations. 5''' 6# make sure int/int division yields float quotient, see .basics 7from __future__ import division 8 9from pygeodesy.basics import isnon0 as _non0 10from pygeodesy.datums import Datum, _ellipsoidal_datum, _mean_radius, \ 11 _spherical_datum, _WGS84 12from pygeodesy.ellipsoids import Ellipsoid 13from pygeodesy.errors import _AssertionError, IntersectionError, \ 14 LimitError, _limiterrors, _ValueError 15from pygeodesy.fmath import euclid, fdot, fsum_, hypot, hypot2, sqrt0 16from pygeodesy.interns import EPS, EPS0, EPS1, NN, PI, PI2, PI3, PI_2, R_M, \ 17 _distant_, _inside_, _near_, _null_, _outside_, \ 18 _too_, _0_0, _0_125, _0_25, _0_5, _1_0, \ 19 _2_0, _4_0, _32_0, _90_0, _180_0, _360_0 20from pygeodesy.lazily import _ALL_LAZY 21from pygeodesy.named import _NamedTuple, _xnamed 22from pygeodesy.namedTuples import Bearing2Tuple, Distance4Tuple, \ 23 LatLon2Tuple, PhiLam2Tuple, Vector3Tuple 24from pygeodesy.streprs import unstr 25from pygeodesy.units import Degrees_, Distance, Distance_, Height, Lam_, Lat, \ 26 Lon, Phi_, Radians, Radians_, Radius, Radius_, \ 27 Scalar, _100km 28from pygeodesy.utily import acos1, atan2b, degrees2m, degrees90, degrees180, \ 29 m2degrees, sincos2, sincos2_, tan_2, unroll180, \ 30 unrollPI, wrap90, wrap180, wrapPI, wrapPI_2 31 32from math import atan, atan2, cos, degrees, radians, sin, sqrt # pow 33 34__all__ = _ALL_LAZY.formy 35__version__ = '21.09.14' 36 37_opposite_ = 'opposite' 38_ratio_ = 'ratio' 39_xline_ = 'xline' 40 41 42def antipode(lat, lon): 43 '''Return the antipode, the point diametrically opposite 44 to a given point in C{degrees}. 45 46 @arg lat: Latitude (C{degrees}). 47 @arg lon: Longitude (C{degrees}). 48 49 @return: A L{LatLon2Tuple}C{(lat, lon)}. 50 51 @see: U{Geosphere<https://CRAN.R-Project.org/web/packages/geosphere/geosphere.pdf>}. 52 ''' 53 return LatLon2Tuple(-wrap90(lat), wrap180(lon + _180_0)) 54 55 56def antipode_(phi, lam): 57 '''Return the antipode, the point diametrically opposite 58 to a given point in C{radians}. 59 60 @arg phi: Latitude (C{radians}). 61 @arg lam: Longitude (C{radians}). 62 63 @return: A L{PhiLam2Tuple}C{(phi, lam)}. 64 65 @see: U{Geosphere<https://CRAN.R-Project.org/web/packages/geosphere/geosphere.pdf>}. 66 ''' 67 return PhiLam2Tuple(-wrapPI_2(phi), wrapPI(lam + PI)) 68 69 70def _area_or_(func_, lat1, lat2, radius, d_lon, unused): 71 '''(INTERNAL) Helper for area and spherical excess. 72 ''' 73 r = func_(Phi_(lat2=lat2), 74 Phi_(lat1=lat1), radians(d_lon)) 75 if radius: 76 r *= _mean_radius(radius, lat1, lat2)**2 77 return r 78 79 80def bearing(lat1, lon1, lat2, lon2, **options): 81 '''Compute the initial or final bearing (forward or reverse 82 azimuth) between a (spherical) start and end point. 83 84 @arg lat1: Start latitude (C{degrees}). 85 @arg lon1: Start longitude (C{degrees}). 86 @arg lat2: End latitude (C{degrees}). 87 @arg lon2: End longitude (C{degrees}). 88 @kwarg options: Optional keyword arguments for function 89 L{bearing_}. 90 91 @return: Initial or final bearing (compass C{degrees360}) or 92 zero if start and end point coincide. 93 ''' 94 r = bearing_(Phi_(lat1=lat1), Lam_(lon1=lon1), 95 Phi_(lat2=lat2), Lam_(lon2=lon2), **options) 96 return degrees(r) 97 98 99def bearing_(phi1, lam1, phi2, lam2, final=False, wrap=False): 100 '''Compute the initial or final bearing (forward or reverse 101 azimuth) between a (spherical) start and end point. 102 103 @arg phi1: Start latitude (C{radians}). 104 @arg lam1: Start longitude (C{radians}). 105 @arg phi2: End latitude (C{radians}). 106 @arg lam2: End longitude (C{radians}). 107 @kwarg final: Return final bearing if C{True}, initial 108 otherwise (C{bool}). 109 @kwarg wrap: Wrap and L{unrollPI} longitudes (C{bool}). 110 111 @return: Initial or final bearing (compass C{radiansPI2}) or 112 zero if start and end point coincide. 113 ''' 114 if final: # swap plus PI 115 phi1, lam1, phi2, lam2 = phi2, lam2, phi1, lam1 116 r = PI3 117 else: 118 r = PI2 119 120 db, _ = unrollPI(lam1, lam2, wrap=wrap) 121 sa1, ca1, sa2, ca2, sdb, cdb = sincos2_(phi1, phi2, db) 122 123 # see <https://MathForum.org/library/drmath/view/55417.html> 124 x = ca1 * sa2 - sa1 * ca2 * cdb 125 y = sdb * ca2 126 return (atan2(y, x) + r) % PI2 # wrapPI2 127 128 129def _bearingTo2(p1, p2, wrap=False): # for points.ispolar, sphericalTrigonometry.areaOf 130 '''(INTERNAL) Compute initial and final bearing. 131 ''' 132 try: # for LatLon_ and ellipsoidal LatLon 133 return p1.bearingTo2(p2, wrap=wrap) 134 except AttributeError: 135 pass 136 # XXX spherical version, OK for ellipsoidal ispolar? 137 a1, b1 = p1.philam 138 a2, b2 = p2.philam 139 return Bearing2Tuple(degrees(bearing_(a1, b1, a2, b2, final=False, wrap=wrap)), 140 degrees(bearing_(a1, b1, a2, b2, final=True, wrap=wrap)), 141 name=_bearingTo2.__name__) 142 143 144def compassAngle(lat1, lon1, lat2, lon2, adjust=True, wrap=False): 145 '''Return the angle from North for the direction vector 146 M{(lon2 - lon1, lat2 - lat1)} between two points. 147 148 Suitable only for short, not near-polar vectors up to a few 149 hundred Km or Miles. Use function L{bearing} for longer 150 vectors. 151 152 @arg lat1: From latitude (C{degrees}). 153 @arg lon1: From longitude (C{degrees}). 154 @arg lat2: To latitude (C{degrees}). 155 @arg lon2: To longitude (C{degrees}). 156 @kwarg adjust: Adjust the longitudinal delta by the 157 cosine of the mean latitude (C{bool}). 158 @kwarg wrap: Wrap and L{unroll180} longitudes (C{bool}). 159 160 @return: Compass angle from North (C{degrees360}). 161 162 @note: Courtesy Martin Schultz. 163 164 @see: U{Local, flat earth approximation 165 <https://www.EdWilliams.org/avform.htm#flat>}. 166 ''' 167 d_lon, _ = unroll180(lon1, lon2, wrap=wrap) 168 if adjust: # scale delta lon 169 d_lon *= _scale_deg(lat1, lat2) 170 return atan2b(d_lon, lat2 - lat1) 171 172 173def cosineAndoyerLambert(lat1, lon1, lat2, lon2, datum=_WGS84, wrap=False): 174 '''Compute the distance between two (ellipsoidal) points using the 175 U{Andoyer-Lambert correction<https://navlib.net/wp-content/uploads/ 176 2013/10/admiralty-manual-of-navigation-vol-1-1964-english501c.pdf>} of the 177 U{Law of Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>} 178 fromula. 179 180 @arg lat1: Start latitude (C{degrees}). 181 @arg lon1: Start longitude (C{degrees}). 182 @arg lat2: End latitude (C{degrees}). 183 @arg lon2: End longitude (C{degrees}). 184 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid}, 185 L{Ellipsoid2} or L{a_f2Tuple}). 186 @kwarg wrap: Wrap and L{unroll180} longitudes (C{bool}). 187 188 @return: Distance (C{meter}, same units as the B{C{datum}}'s 189 ellipsoid axes or C{radians} if B{C{datum}} is C{None}). 190 191 @raise TypeError: Invalid B{C{datum}}. 192 193 @see: Functions L{cosineAndoyerLambert_}, L{cosineForsytheAndoyerLambert}, 194 L{cosineLaw}, L{equirectangular}, L{euclidean}, L{flatLocal}/L{hubeny}, 195 L{flatPolar}, L{haversine}, L{thomas} and L{vincentys} and method 196 L{Ellipsoid.distance2}. 197 ''' 198 return _distanceToE(cosineAndoyerLambert_, lat1, lat2, datum, 199 *unroll180(lon1, lon2, wrap=wrap)) 200 201 202def cosineAndoyerLambert_(phi2, phi1, lam21, datum=_WGS84): 203 '''Compute the I{angular} distance between two (ellipsoidal) points using the 204 U{Andoyer-Lambert correction<https://navlib.net/wp-content/uploads/2013/10/ 205 admiralty-manual-of-navigation-vol-1-1964-english501c.pdf>} of the U{Law 206 of Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>} 207 fromula. 208 209 @arg phi2: End latitude (C{radians}). 210 @arg phi1: Start latitude (C{radians}). 211 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 212 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid}, 213 L{Ellipsoid2} or L{a_f2Tuple}). 214 215 @return: Angular distance (C{radians}). 216 217 @raise TypeError: Invalid B{C{datum}}. 218 219 @see: Functions L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert_}, 220 L{cosineLaw_}, L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_}, 221 L{flatPolar_}, L{haversine_}, L{thomas_} and L{vincentys_} and U{Geodesy-PHP 222 <https://GitHub.com/jtejido/geodesy-php/blob/master/src/Geodesy/ 223 Distance/AndoyerLambert.php>}. 224 ''' 225 s2, c2, s1, c1, r, c21 = _sincosa6(phi2, phi1, lam21) 226 if _non0(c1) and _non0(c2): 227 E = _ellipsoidal(datum, cosineAndoyerLambert_) 228 if E.f: # ellipsoidal 229 r2 = atan2(E.b_a * s2, c2) 230 r1 = atan2(E.b_a * s1, c1) 231 s2, c2, s1, c1 = sincos2_(r2, r1) 232 r = acos1(s1 * s2 + c1 * c2 * c21) 233 if r: 234 sr, _, sr_2, cr_2 = sincos2_(r, r * _0_5) 235 if _non0(sr_2) and _non0(cr_2): 236 s = (sr + r) * ((s1 - s2) / sr_2)**2 237 c = (sr - r) * ((s1 + s2) / cr_2)**2 238 r += (c - s) * E.f * _0_125 239 return r 240 241 242def cosineForsytheAndoyerLambert(lat1, lon1, lat2, lon2, datum=_WGS84, wrap=False): 243 '''Compute the distance between two (ellipsoidal) points using the 244 U{Forsythe-Andoyer-Lambert correction<https://www2.UNB.CA/gge/Pubs/TR77.pdf>} of 245 the U{Law of Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>} 246 formula. 247 248 @arg lat1: Start latitude (C{degrees}). 249 @arg lon1: Start longitude (C{degrees}). 250 @arg lat2: End latitude (C{degrees}). 251 @arg lon2: End longitude (C{degrees}). 252 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid}, 253 L{Ellipsoid2} or L{a_f2Tuple}). 254 @kwarg wrap: Wrap and L{unroll180} longitudes (C{bool}). 255 256 @return: Distance (C{meter}, same units as the B{C{datum}}'s 257 ellipsoid axes or C{radians} if B{C{datum}} is C{None}). 258 259 @raise TypeError: Invalid B{C{datum}}. 260 261 @see: Functions L{cosineForsytheAndoyerLambert_}, L{cosineAndoyerLambert}, 262 L{cosineLaw}, L{equirectangular}, L{euclidean}, L{flatLocal}/L{hubeny}, 263 L{flatPolar}, L{haversine}, L{thomas} and L{vincentys} and method 264 L{Ellipsoid.distance2}. 265 ''' 266 return _distanceToE(cosineForsytheAndoyerLambert_, lat1, lat2, datum, 267 *unroll180(lon1, lon2, wrap=wrap)) 268 269 270def cosineForsytheAndoyerLambert_(phi2, phi1, lam21, datum=_WGS84): 271 '''Compute the I{angular} distance between two (ellipsoidal) points using the 272 U{Forsythe-Andoyer-Lambert correction<https://www2.UNB.CA/gge/Pubs/TR77.pdf>} of 273 the U{Law of Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>} 274 formula. 275 276 @arg phi2: End latitude (C{radians}). 277 @arg phi1: Start latitude (C{radians}). 278 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 279 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid}, 280 L{Ellipsoid2} or L{a_f2Tuple}). 281 282 @return: Angular distance (C{radians}). 283 284 @raise TypeError: Invalid B{C{datum}}. 285 286 @see: Functions L{cosineForsytheAndoyerLambert}, L{cosineAndoyerLambert_}, 287 L{cosineLaw_}, L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_}, 288 L{flatPolar_}, L{haversine_}, L{thomas_} and L{vincentys_} and U{Geodesy-PHP 289 <https://GitHub.com/jtejido/geodesy-php/blob/master/src/Geodesy/ 290 Distance/ForsytheCorrection.php>}. 291 ''' 292 s2, c2, s1, c1, r, _ = _sincosa6(phi2, phi1, lam21) 293 if r and _non0(c1) and _non0(c2): 294 E = _ellipsoidal(datum, cosineForsytheAndoyerLambert_) 295 if E.f: # ellipsoidal 296 sr, cr, s2r, _ = sincos2_(r, r * _2_0) 297 if _non0(sr) and abs(cr) < EPS1: 298 s = (s1 + s2)**2 / (1 + cr) 299 t = (s1 - s2)**2 / (1 - cr) 300 x = s + t 301 y = s - t 302 303 s = 8 * r**2 / sr 304 a = 64 * r + _2_0 * s * cr # 16 * r**2 / tan(r) 305 d = 48 * sr + s # 8 * r**2 / tan(r) 306 b = -2 * d 307 e = 30 * s2r 308 c = fsum_(30 * r, e * _0_5, s * cr) # 8 * r**2 / tan(r) 309 310 t = fsum_( a * x, b * y, -c * x**2, d * x * y, e * y**2) 311 r += fsum_(-r * x, 3 * y * sr, t * E.f / _32_0) * E.f * _0_25 312 return r 313 314 315def cosineLaw(lat1, lon1, lat2, lon2, radius=R_M, wrap=False): 316 '''Compute the distance between two points using the 317 U{spherical Law of Cosines 318 <https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>} 319 formula. 320 321 @arg lat1: Start latitude (C{degrees}). 322 @arg lon1: Start longitude (C{degrees}). 323 @arg lat2: End latitude (C{degrees}). 324 @arg lon2: End longitude (C{degrees}). 325 @kwarg radius: Mean earth radius, ellipsoid or datum 326 (C{meter}, L{Ellipsoid}, L{Ellipsoid2}, 327 L{Datum} or L{a_f2Tuple}). 328 @kwarg wrap: Wrap and L{unroll180} longitudes (C{bool}). 329 330 @return: Distance (C{meter}, same units as B{C{radius}} or 331 the ellipsoid or datum axes). 332 333 @raise TypeError: Invalid B{C{radius}}. 334 335 @see: Functions L{cosineLaw_}, L{cosineAndoyerLambert}, 336 L{cosineForsytheAndoyerLambert}, L{equirectangular}, L{euclidean}, 337 L{flatLocal}/L{hubeny}, L{flatPolar}, L{haversine}, L{thomas} and 338 L{vincentys} and method L{Ellipsoid.distance2}. 339 340 @note: See note at function L{vincentys_}. 341 ''' 342 return _distanceToS(cosineLaw_, lat1, lat2, radius, 343 *unroll180(lon1, lon2, wrap=wrap)) 344 345 346def cosineLaw_(phi2, phi1, lam21): 347 '''Compute the I{angular} distance between two points using the 348 U{spherical Law of Cosines 349 <https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>} 350 formula. 351 352 @arg phi2: End latitude (C{radians}). 353 @arg phi1: Start latitude (C{radians}). 354 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 355 356 @return: Angular distance (C{radians}). 357 358 @see: Functions L{cosineLaw}, L{cosineAndoyerLambert_}, 359 L{cosineForsytheAndoyerLambert_}, L{equirectangular_}, 360 L{euclidean_}, L{flatLocal_}/L{hubeny_}, L{flatPolar_}, 361 L{haversine_}, L{thomas_} and L{vincentys_}. 362 363 @note: See note at function L{vincentys_}. 364 ''' 365 return _sincosa6(phi2, phi1, lam21)[4] 366 367 368def _distanceToE(func_, lat1, lat2, earth, d_lon, unused): 369 '''(INTERNAL) Helper for ellipsoidal distances. 370 ''' 371 E = _ellipsoidal(earth, func_) 372 r = func_(Phi_(lat2=lat2), 373 Phi_(lat1=lat1), radians(d_lon), datum=E) 374 return r * E.a 375 376 377def _distanceToS(func_, lat1, lat2, earth, d_lon, unused, **adjust): 378 '''(INTERNAL) Helper for spherical distances. 379 ''' 380 r = func_(Phi_(lat2=lat2), 381 Phi_(lat1=lat1), radians(d_lon), **adjust) 382 return r * _mean_radius(earth, lat1, lat2) 383 384 385def _ellipsoidal(earth, where): 386 '''(INTERNAL) Helper for distances. 387 ''' 388 return earth if isinstance(earth, Ellipsoid) else ( 389 earth if isinstance(earth, Datum) else 390 _ellipsoidal_datum(earth, name=where.__name__)).ellipsoid # PYCHOK indent 391 392 393def equirectangular(lat1, lon1, lat2, lon2, radius=R_M, **options): 394 '''Compute the distance between two points using 395 the U{Equirectangular Approximation / Projection 396 <https://www.Movable-Type.co.UK/scripts/latlong.html#equirectangular>}. 397 398 @arg lat1: Start latitude (C{degrees}). 399 @arg lon1: Start longitude (C{degrees}). 400 @arg lat2: End latitude (C{degrees}). 401 @arg lon2: End longitude (C{degrees}). 402 @kwarg radius: Mean earth radius, ellipsoid or datum 403 (C{meter}, L{Ellipsoid}, L{Ellipsoid2}, 404 L{Datum} or L{a_f2Tuple}). 405 @kwarg options: Optional keyword arguments for function 406 L{equirectangular_}. 407 408 @return: Distance (C{meter}, same units as B{C{radius}} or 409 the ellipsoid or datum axes). 410 411 @raise TypeError: Invalid B{C{radius}}. 412 413 @see: Function L{equirectangular_} for more details, the 414 available B{C{options}}, errors, restrictions and other, 415 approximate or accurate distance functions. 416 ''' 417 d = sqrt(equirectangular_(Lat(lat1=lat1), Lon(lon1=lon1), 418 Lat(lat2=lat2), Lon(lon2=lon2), 419 **options).distance2) # PYCHOK 4 vs 2-3 420 return degrees2m(d, radius=_mean_radius(radius, lat1, lat2)) 421 422 423def equirectangular_(lat1, lon1, lat2, lon2, 424 adjust=True, limit=45, wrap=False): 425 '''Compute the distance between two points using 426 the U{Equirectangular Approximation / Projection 427 <https://www.Movable-Type.co.UK/scripts/latlong.html#equirectangular>}. 428 429 This approximation is valid for short distance of several 430 hundred Km or Miles, see the B{C{limit}} keyword argument and 431 the L{LimitError}. 432 433 @arg lat1: Start latitude (C{degrees}). 434 @arg lon1: Start longitude (C{degrees}). 435 @arg lat2: End latitude (C{degrees}). 436 @arg lon2: End longitude (C{degrees}). 437 @kwarg adjust: Adjust the wrapped, unrolled longitudinal delta 438 by the cosine of the mean latitude (C{bool}). 439 @kwarg limit: Optional limit for lat- and longitudinal deltas 440 (C{degrees}) or C{None} or C{0} for unlimited. 441 @kwarg wrap: Wrap and L{unroll180} longitudes (C{bool}). 442 443 @return: A L{Distance4Tuple}C{(distance2, delta_lat, delta_lon, 444 unroll_lon2)}. 445 446 @raise LimitError: If the lat- and/or longitudinal delta exceeds 447 the B{C{-limit..+limit}} range and L{limiterrors} 448 set to C{True}. 449 450 @see: U{Local, flat earth approximation 451 <https://www.EdWilliams.org/avform.htm#flat>}, functions 452 L{equirectangular}, L{cosineAndoyerLambert}, 453 L{cosineForsytheAndoyerLambert}, L{cosineLaw}, L{euclidean}, 454 L{flatLocal}/L{hubeny}, L{flatPolar}, L{haversine}, L{thomas} 455 and L{vincentys} and methods L{Ellipsoid.distance2}, 456 C{LatLon.distanceTo*} and C{LatLon.equirectangularTo}. 457 ''' 458 d_lat = lat2 - lat1 459 d_lon, ulon2 = unroll180(lon1, lon2, wrap=wrap) 460 461 if limit and _limiterrors \ 462 and max(abs(d_lat), abs(d_lon)) > limit > 0: 463 t = unstr(equirectangular_.__name__, 464 lat1, lon1, lat2, lon2, limit=limit) 465 raise LimitError('delta exceeds limit', txt=t) 466 467 if adjust: # scale delta lon 468 d_lon *= _scale_deg(lat1, lat2) 469 470 d2 = hypot2(d_lat, d_lon) # degrees squared! 471 return Distance4Tuple(d2, d_lat, d_lon, ulon2 - lon2) 472 473 474def euclidean(lat1, lon1, lat2, lon2, radius=R_M, adjust=True, wrap=False): 475 '''Approximate the C{Euclidean} distance between two (spherical) points. 476 477 @arg lat1: Start latitude (C{degrees}). 478 @arg lon1: Start longitude (C{degrees}). 479 @arg lat2: End latitude (C{degrees}). 480 @arg lon2: End longitude (C{degrees}). 481 @kwarg radius: Mean earth radius, ellipsoid or datum 482 (C{meter}, L{Ellipsoid}, L{Ellipsoid2}, 483 L{Datum} or L{a_f2Tuple}). 484 @kwarg adjust: Adjust the longitudinal delta by the cosine 485 of the mean latitude (C{bool}). 486 @kwarg wrap: Wrap and L{unroll180} longitudes (C{bool}). 487 488 @return: Distance (C{meter}, same units as B{C{radius}} or 489 the ellipsoid or datum axes). 490 491 @raise TypeError: Invalid B{C{radius}}. 492 493 @see: U{Distance between two (spherical) points 494 <https://www.EdWilliams.org/avform.htm#Dist>}, functions L{euclid}, 495 L{euclidean_}, L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert}, 496 L{cosineLaw}, L{equirectangular}, L{flatLocal}/L{hubeny}, L{flatPolar}, 497 L{haversine}, L{thomas} and L{vincentys} and methods L{Ellipsoid.distance2}, 498 C{LatLon.distanceTo*} and C{LatLon.equirectangularTo}. 499 ''' 500 return _distanceToS(euclidean_, lat1, lat2, radius, 501 *unroll180(lon1, lon2, wrap=wrap), 502 adjust=adjust) 503 504 505def euclidean_(phi2, phi1, lam21, adjust=True): 506 '''Approximate the I{angular} C{Euclidean} distance between two 507 (spherical) points. 508 509 @arg phi2: End latitude (C{radians}). 510 @arg phi1: Start latitude (C{radians}). 511 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 512 @kwarg adjust: Adjust the longitudinal delta by the cosine 513 of the mean latitude (C{bool}). 514 515 @return: Angular distance (C{radians}). 516 517 @see: Functions L{euclid}, L{euclidean}, L{cosineAndoyerLambert_}, 518 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, L{equirectangular_}, 519 L{flatLocal_}/L{hubeny_}, L{flatPolar_}, L{haversine_}, L{thomas_} 520 and L{vincentys_}. 521 ''' 522 if adjust: 523 lam21 *= _scale_rad(phi2, phi1) 524 return euclid(phi2 - phi1, lam21) 525 526 527def excessAbc(A, b, c): 528 '''Compute the I{spherical excess} C{E} of a (spherical) triangle 529 from two sides and the included angle. 530 531 @arg A: An interior triangle angle (C{radians}). 532 @arg b: Frist adjacent triangle side (C{radians}). 533 @arg c: Second adjacent triangle side (C{radians}). 534 535 @return: Spherical excess ({radians}). 536 537 @raise UnitError: Invalid B{C{A}}, B{C{b}} or B{C{c}}. 538 539 @see: Function L{excessGirard}, L{excessLHuilier}, U{Spherical 540 trigonometry<https://WikiPedia.org/wiki/Spherical_trigonometry>}. 541 ''' 542 sA, cA, sb, cb, sc, cc = sincos2_(Radians_(A=A), Radians_(b=b) * _0_5, 543 Radians_(c=c) * _0_5) 544 return atan2(sA * sb * sc, cb * cc + cA * sb * sc) * _2_0 545 546 547def excessGirard(A, B, C): 548 '''Compute the I{spherical excess} C{E} of a (spherical) triangle using 549 U{Girard's<https://MathWorld.Wolfram.com/GirardsSphericalExcessFormula.html>} 550 formula. 551 552 @arg A: First interior triangle angle (C{radians}). 553 @arg B: Second interior triangle angle (C{radians}). 554 @arg C: Third interior triangle angle (C{radians}). 555 556 @return: Spherical excess ({radians}). 557 558 @raise UnitError: Invalid B{C{A}}, B{C{B}} or B{C{C}}. 559 560 @see: Function L{excessLHuilier}, U{Spherical trigonometry 561 <https://WikiPedia.org/wiki/Spherical_trigonometry>}. 562 ''' 563 return Radians(Girard=fsum_(Radians_(A=A), 564 Radians_(B=B), 565 Radians_(C=C), -PI)) 566 567 568def excessLHuilier(a, b, c): 569 '''Compute the I{spherical excess} C{E} of a (spherical) triangle using 570 U{L'Huilier's<https://MathWorld.Wolfram.com/LHuiliersTheorem.html>} 571 Theorem. 572 573 @arg a: First triangle side (C{radians}). 574 @arg b: Second triangle side (C{radians}). 575 @arg c: Third triangle side (C{radians}). 576 577 @return: Spherical excess ({radians}). 578 579 @raise UnitError: Invalid B{C{a}}, B{C{b}} or B{C{c}}. 580 581 @see: Function L{excessGirard}, U{Spherical trigonometry 582 <https://WikiPedia.org/wiki/Spherical_trigonometry>}. 583 ''' 584 a = Radians_(a=a) 585 b = Radians_(b=b) 586 c = Radians_(c=c) 587 588 s = fsum_(a, b, c) * _0_5 589 r = tan_2(s) * tan_2(s - a) * tan_2(s - b) * tan_2(s - c) 590 return Radians(LHuilier=atan(sqrt(r)) * _4_0) 591 592 593def excessKarney(lat1, lon1, lat2, lon2, radius=R_M, wrap=False): 594 '''Compute the surface area of a (spherical) quadrilateral bounded by 595 a segment of a great circle, two meridians and the equator using U{Karney's 596 <http://OSGeo-org.1560.x6.Nabble.com/Area-of-a-spherical-polygon-td3841625.html>} 597 method. 598 599 @arg lat1: Start latitude (C{degrees}). 600 @arg lon1: Start longitude (C{degrees}). 601 @arg lat2: End latitude (C{degrees}). 602 @arg lon2: End longitude (C{degrees}). 603 @kwarg radius: Mean earth radius, ellipsoid or datum 604 (C{meter}, L{Ellipsoid}, L{Ellipsoid2}, 605 L{Datum} or L{a_f2Tuple}) or C{None}. 606 @kwarg wrap: Wrap and L{unroll180} longitudes (C{bool}). 607 608 @return: Surface area, I{signed} (I{square} C{meter}, or units of 609 B{C{radius}} I{squared}) or I{spherical excess} (C{radians}) 610 if B{C{radius}} is C{None} or C{0}. 611 612 @raise TypeError: Invalid B{C{radius}}. 613 614 @raise UnitError: Invalid B{C{lat2}} or B{C{lat1}}. 615 616 @raise ValueError: Semi-circular longitudinal delta. 617 618 @see: Function L{excessKarney_} and L{excessQuad}. 619 ''' 620 return _area_or_(excessKarney_, lat1, lat2, radius, 621 *unroll180(lon1, lon2, wrap=wrap)) 622 623 624def excessKarney_(phi2, phi1, lam21): 625 '''Compute the I{spherical excess} C{E} of a (spherical) quadrilateral bounded 626 by a segment of a great circle, two meridians and the equator using U{Karney's 627 <http://OSGeo-org.1560.x6.Nabble.com/Area-of-a-spherical-polygon-td3841625.html>} 628 method. 629 630 @arg phi2: End latitude (C{radians}). 631 @arg phi1: Start latitude (C{radians}). 632 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 633 634 @return: Spherical excess, I{signed} (C{radians}). 635 636 @raise ValueError: Semi-circular longitudinal delta B{C{lam21}}. 637 638 @see: Function L{excessKarney}, U{Area of a spherical polygon 639 <http://OSGeo-org.1560.x6.Nabble.com/Area-of-a-spherical-polygon-td3841625.html>}. 640 ''' 641 # from: Veness <https://www.Movable-Type.co.UK/scripts/latlong.html> 642 # Area method due to Karney: for each edge of the polygon, 643 # 644 # tan(Δλ/2) · (tan(φ1/2) + tan(φ2/2)) 645 # tan(E/2) = ------------------------------------ 646 # 1 + tan(φ1/2) · tan(φ2/2) 647 # 648 # where E is the spherical excess of the trapezium obtained by 649 # extending the edge to the equator-circle vector for each edge. 650 t2 = tan_2(phi2) 651 t1 = tan_2(phi1) 652 t = tan_2(lam21, lam21=None) 653 return Radians(Karney=atan2(t * (t1 + t2), 654 _1_0 + (t1 * t2)) * _2_0) 655 656 657def excessQuad(lat1, lon1, lat2, lon2, radius=R_M, wrap=False): 658 '''Compute the surface area of a (spherical) quadrilateral bounded 659 by a segment of a great circle, two meridians and the equator. 660 661 @arg lat1: Start latitude (C{degrees}). 662 @arg lon1: Start longitude (C{degrees}). 663 @arg lat2: End latitude (C{degrees}). 664 @arg lon2: End longitude (C{degrees}). 665 @kwarg radius: Mean earth radius, ellipsoid or datum 666 (C{meter}, L{Ellipsoid}, L{Ellipsoid2}, 667 L{Datum} or L{a_f2Tuple}) or C{None}. 668 @kwarg wrap: Wrap and L{unroll180} longitudes (C{bool}). 669 670 @return: Surface area, I{signed} (I{square} C{meter}, or units of 671 B{C{radius}} I{squared}) or I{spherical excess} (C{radians}) 672 if B{C{radius}} is C{None} or C{0}. 673 674 @raise TypeError: Invalid B{C{radius}}. 675 676 @raise UnitError: Invalid B{C{lat2}} or B{C{lat1}}. 677 678 @see: Function L{excessQuad_} and L{excessKarney}. 679 ''' 680 return _area_or_(excessQuad_, lat1, lat2, radius, 681 *unroll180(lon1, lon2, wrap=wrap)) 682 683 684def excessQuad_(phi2, phi1, lam21): 685 '''Compute the I{spherical excess} C{E} of a (spherical) quadrilateral bounded 686 by a segment of a great circle, two meridians and the equator. 687 688 @arg phi2: End latitude (C{radians}). 689 @arg phi1: Start latitude (C{radians}). 690 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 691 692 @return: Spherical excess, I{signed} (C{radians}). 693 694 @see: Function L{excessQuad}, U{Spherical trigonometry 695 <https://WikiPedia.org/wiki/Spherical_trigonometry>}. 696 ''' 697 s = sin((phi2 + phi1) * _0_5) 698 c = cos((phi2 - phi1) * _0_5) 699 return Radians(Quad=atan2(tan_2(lam21) * s, c) * _2_0) 700 701 702def flatLocal(lat1, lon1, lat2, lon2, datum=_WGS84, wrap=False): 703 '''Compute the distance between two (ellipsoidal) points using 704 the U{ellipsoidal Earth to plane projection<https://WikiPedia.org/ 705 wiki/Geographical_distance#Ellipsoidal_Earth_projected_to_a_plane>} 706 aka U{Hubeny<https://www.OVG.AT/de/vgi/files/pdf/3781/>} formula. 707 708 @arg lat1: Start latitude (C{degrees}). 709 @arg lon1: Start longitude (C{degrees}). 710 @arg lat2: End latitude (C{degrees}). 711 @arg lon2: End longitude (C{degrees}). 712 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid}, 713 L{Ellipsoid2} or L{a_f2Tuple}). 714 @kwarg wrap: Wrap and L{unroll180} longitudes (C{bool}). 715 716 @return: Distance (C{meter}, same units as the B{C{datum}}'s 717 ellipsoid axes). 718 719 @raise TypeError: Invalid B{C{datum}}. 720 721 @note: The meridional and prime_vertical radii of curvature 722 are taken and scaled at the mean of both latitude. 723 724 @see: Functions L{flatLocal_}/L{hubeny_}, L{cosineLaw}, 725 L{flatPolar}, L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert}, 726 L{equirectangular}, L{euclidean}, L{haversine}, L{thomas}, L{vincentys}, 727 method L{Ellipsoid.distance2} and U{local, flat earth approximation 728 <https://www.EdWilliams.org/avform.htm#flat>}. 729 ''' 730 d, _ = unroll180(lon1, lon2, wrap=wrap) 731 return flatLocal_(Phi_(lat2=lat2), 732 Phi_(lat1=lat1), radians(d), datum=datum) 733 734 735hubeny = flatLocal # for Karl Hubeny 736 737 738def flatLocal_(phi2, phi1, lam21, datum=_WGS84): 739 '''Compute the I{angular} distance between two (ellipsoidal) points using 740 the U{ellipsoidal Earth to plane projection<https://WikiPedia.org/ 741 wiki/Geographical_distance#Ellipsoidal_Earth_projected_to_a_plane>} 742 aka U{Hubeny<https://www.OVG.AT/de/vgi/files/pdf/3781/>} formula. 743 744 @arg phi2: End latitude (C{radians}). 745 @arg phi1: Start latitude (C{radians}). 746 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 747 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid}, 748 L{Ellipsoid2} or L{a_f2Tuple}). 749 750 @return: Angular distance (C{radians}). 751 752 @raise TypeError: Invalid B{C{datum}}. 753 754 @note: The meridional and prime_vertical radii of curvature 755 are taken and scaled I{at the mean of both latitude}. 756 757 @see: Functions L{flatLocal}/L{hubeny}, L{cosineAndoyerLambert_}, 758 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, 759 L{flatPolar_}, L{equirectangular_}, L{euclidean_}, 760 L{haversine_}, L{thomas_} and L{vincentys_} and U{local, flat 761 earth approximation <https://www.EdWilliams.org/avform.htm#flat>}. 762 ''' 763 E = _ellipsoidal(datum, flatLocal_) 764 m, n = E.roc2_((phi2 + phi1) * _0_5, scaled=True) 765 return hypot(m * (phi2 - phi1), n * lam21) 766 767 768hubeny_ = flatLocal_ # for Karl Hubeny 769 770 771def flatPolar(lat1, lon1, lat2, lon2, radius=R_M, wrap=False): 772 '''Compute the distance between two (spherical) points using 773 the U{polar coordinate flat-Earth <https://WikiPedia.org/wiki/ 774 Geographical_distance#Polar_coordinate_flat-Earth_formula>} 775 formula. 776 777 @arg lat1: Start latitude (C{degrees}). 778 @arg lon1: Start longitude (C{degrees}). 779 @arg lat2: End latitude (C{degrees}). 780 @arg lon2: End longitude (C{degrees}). 781 @kwarg radius: Mean earth radius, ellipsoid or datum 782 (C{meter}, L{Ellipsoid}, L{Ellipsoid2}, 783 L{Datum} or L{a_f2Tuple}). 784 @kwarg wrap: Wrap and L{unroll180} longitudes (C{bool}). 785 786 @return: Distance (C{meter}, same units as B{C{radius}} or 787 the ellipsoid or datum axes). 788 789 @raise TypeError: Invalid B{C{radius}}. 790 791 @see: Functions L{flatPolar_}, L{cosineAndoyerLambert}, 792 L{cosineForsytheAndoyerLambert},L{cosineLaw}, 793 L{flatLocal}/L{hubeny}, L{equirectangular}, 794 L{euclidean}, L{haversine}, L{thomas} and 795 L{vincentys}. 796 ''' 797 return _distanceToS(flatPolar_, lat1, lat2, radius, 798 *unroll180(lon1, lon2, wrap=wrap)) 799 800 801def flatPolar_(phi2, phi1, lam21): 802 '''Compute the I{angular} distance between two (spherical) points 803 using the U{polar coordinate flat-Earth<https://WikiPedia.org/wiki/ 804 Geographical_distance#Polar_coordinate_flat-Earth_formula>} 805 formula. 806 807 @arg phi2: End latitude (C{radians}). 808 @arg phi1: Start latitude (C{radians}). 809 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 810 811 @return: Angular distance (C{radians}). 812 813 @see: Functions L{flatPolar}, L{cosineAndoyerLambert_}, 814 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, 815 L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_}, 816 L{haversine_}, L{thomas_} and L{vincentys_}. 817 ''' 818 a1 = PI_2 - phi1 # co-latitude 819 a2 = PI_2 - phi2 # co-latitude 820 ab = _2_0 * a1 * a2 * cos(lam21) 821 r2 = fsum_(a1**2, a2**2, -abs(ab)) 822 return sqrt0(r2) 823 824 825def hartzell(pov, los=None, earth=_WGS84, LatLon=None, **LatLon_kwds): 826 '''Compute the intersection of a Line-Of-Sight from a Point-Of-View in 827 space with the surface of the earth. 828 829 @arg pov: Point-Of-View outside the earth (C{Cartesian}, L{Ecef9Tuple} 830 or L{Vector3d}). 831 @kwarg los: Line-Of-Sight, I{direction} to earth (L{Vector3d}) or 832 C{None} to point to the earth' center. 833 @kwarg earth: The earth model (L{Datum}, L{Ellipsoid}, L{Ellipsoid2}, 834 L{a_f2Tuple} or C{scalar} radius in C{meter}. 835 @kwarg LatLon: Class to convert insection point (C{LatLon}, L{LatLon_}). 836 @kwarg LatLon_kwds: Optional, additional B{C{LatLon}} keyword 837 arguments, ignored if C{B{LatLon} is None}. 838 839 @return: The earth intersection (L{Vector3d}, C{Cartesian type} of 840 B{C{pov}} or B{C{LatLon}}). 841 842 @raise IntersectionError: Null B{C{pov}} or B{C{los}} vector, B{C{pov}} 843 is inside the earth or B{C{los}} points outside 844 the earth or points in an opposite direction. 845 846 @raise TypeError: Invalid B{C{pov}}, B{C{los}} or B{C{earth}}. 847 848 @see: U{Stephen Hartzell<https://StephenHartzell.medium.com/ 849 satellite-line-of-sight-intersection-with-earth-d786b4a6a9b6>}. 850 ''' 851 from pygeodesy.vector3d import _otherV3d 852 853 D = earth if isinstance(earth, Datum) else \ 854 _spherical_datum(earth, name=hartzell.__name__) 855 E = D.ellipsoid 856 857 a2 = b2 = E.a2 # earth x, y, ... 858 c2 = E.b2 # ... z half-axis squared 859 q2 = E.b2_a2 # == c2 / a2 860 bc = E.a * E.b # == b * c 861 862 p3 = _otherV3d(pov=pov) 863 u3 = _otherV3d(los=los) if los else p3.negate() 864 u3 = u3.unit() # unit vector, opposing signs 865 866 x2, y2, z2 = p3.times_(p3).xyz # == p3.x2y2z2 867 ux, vy, wz = u3.times_(p3).xyz 868 u2, v2, w2 = u3.times_(u3).xyz # == u3.x2y2z2 869 870 t = c2, c2, b2 # a2 factored out 871 m = fdot(t, u2, v2, w2) 872 if m < EPS0: # zero or near-null LOS vector 873 raise IntersectionError(pov=pov, los=los, earth=earth, txt=_near_(_null_)) 874 875 # a2 and b2 factored out, b2 == a2 and b2 / a2 == 1 876 r = fsum_(b2 * w2, c2 * v2, -v2 * z2, vy * wz * 2, 877 c2 * u2, -u2 * z2, -w2 * x2, ux * wz * 2, 878 -w2 * y2, -u2 * y2 * q2, -v2 * x2 * q2, ux * vy * 2 * q2) 879 if r > 0: 880 r = bc * sqrt(r) 881 elif r < 0: # LOS pointing away from or missing the earth 882 t = _opposite_ if max(ux, vy, wz) > 0 else _outside_ 883 raise IntersectionError(pov=pov, los=los, earth=earth, txt=t) 884 885 n = fdot(t, ux, vy, wz) 886 d = (n + r) / m # (n - r) / m for antipode 887 if d > 0: # POV inside or LOS missing the earth 888 t = _inside_ if min(x2 - a2, y2 - b2, z2 - c2) < EPS else _outside_ 889 raise IntersectionError(pov=pov, los=los, earth=earth, txt=t) 890 891 if fsum_(x2, y2, z2) < d**2: # d beyond earth center 892 raise IntersectionError(pov=pov, los=los, earth=earth, txt=_too_(_distant_)) 893 894 r = _xnamed(p3.minus(u3.times(d)), hartzell.__name__) 895 if LatLon is not None: 896 from pygeodesy.cartesianBase import CartesianBase as _CB 897 # earth datum is overidden in LatLon if datum is specified in LatLon_kwds 898 r = _CB(r, datum=D, name=r.name).toLatLon(LatLon=LatLon, **LatLon_kwds) 899 return r 900 901 902def haversine(lat1, lon1, lat2, lon2, radius=R_M, wrap=False): 903 '''Compute the distance between two (spherical) points using the 904 U{Haversine<https://www.Movable-Type.co.UK/scripts/latlong.html>} 905 formula. 906 907 @arg lat1: Start latitude (C{degrees}). 908 @arg lon1: Start longitude (C{degrees}). 909 @arg lat2: End latitude (C{degrees}). 910 @arg lon2: End longitude (C{degrees}). 911 @kwarg radius: Mean earth radius, ellipsoid or datum 912 (C{meter}, L{Ellipsoid}, L{Ellipsoid2}, 913 L{Datum} or L{a_f2Tuple}). 914 @kwarg wrap: Wrap and L{unroll180} longitudes (C{bool}). 915 916 @return: Distance (C{meter}, same units as B{C{radius}}). 917 918 @raise TypeError: Invalid B{C{radius}}. 919 920 @see: U{Distance between two (spherical) points 921 <https://www.EdWilliams.org/avform.htm#Dist>}, functions 922 L{cosineLaw}, L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert}, 923 L{equirectangular}, L{euclidean}, L{flatLocal}/L{hubeny}, L{flatPolar}, 924 L{thomas} and L{vincentys} and methods L{Ellipsoid.distance2}, 925 C{LatLon.distanceTo*} and C{LatLon.equirectangularTo}. 926 927 @note: See note at function L{vincentys_}. 928 ''' 929 return _distanceToS(haversine_, lat1, lat2, radius, 930 *unroll180(lon1, lon2, wrap=wrap)) 931 932 933def haversine_(phi2, phi1, lam21): 934 '''Compute the I{angular} distance between two (spherical) points 935 using the U{Haversine<https://www.Movable-Type.co.UK/scripts/latlong.html>} 936 formula. 937 938 @arg phi2: End latitude (C{radians}). 939 @arg phi1: Start latitude (C{radians}). 940 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 941 942 @return: Angular distance (C{radians}). 943 944 @see: Functions L{haversine}, L{cosineAndoyerLambert_}, 945 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, 946 L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_}, 947 L{flatPolar_}, L{thomas_} and L{vincentys_}. 948 949 @note: See note at function L{vincentys_}. 950 ''' 951 def _hsin(rad): 952 return sin(rad * _0_5)**2 953 954 h = _hsin(phi2 - phi1) + cos(phi1) * cos(phi2) * _hsin(lam21) # haversine 955 return atan2(sqrt0(h), sqrt0(_1_0 - h)) * _2_0 # == asin(sqrt(h)) * 2 956 957 958def heightOf(angle, distance, radius=R_M): 959 '''Determine the height above the (spherical) earth' surface after 960 traveling along a straight line at a given tilt. 961 962 @arg angle: Tilt angle above horizontal (C{degrees}). 963 @arg distance: Distance along the line (C{meter} or same units as 964 B{C{radius}}). 965 @kwarg radius: Optional mean earth radius (C{meter}). 966 967 @return: Height (C{meter}, same units as B{C{distance}} and B{C{radius}}). 968 969 @raise ValueError: Invalid B{C{angle}}, B{C{distance}} or B{C{radius}}. 970 971 @see: U{MultiDop geog_lib.GeogBeamHt<https://GitHub.com/NASA/MultiDop>} 972 (U{Shapiro et al. 2009, JTECH 973 <https://Journals.AMetSoc.org/doi/abs/10.1175/2009JTECHA1256.1>} 974 and U{Potvin et al. 2012, JTECH 975 <https://Journals.AMetSoc.org/doi/abs/10.1175/JTECH-D-11-00019.1>}). 976 ''' 977 r = h = Radius(radius) 978 d = abs(Distance(distance)) 979 if d > h: 980 d, h = h, d 981 982 if d > EPS0: 983 d = d / h # /= h chokes PyChecker 984 s = sin(Phi_(angle=angle, clip=_180_0)) 985 s = fsum_(_1_0, _2_0 * s * d, d**2) 986 if s > 0: 987 return h * sqrt(s) - r 988 989 raise _ValueError(angle=angle, distance=distance, radius=radius) 990 991 992def horizon(height, radius=R_M, refraction=False): 993 '''Determine the distance to the horizon from a given altitude 994 above the (spherical) earth. 995 996 @arg height: Altitude (C{meter} or same units as B{C{radius}}). 997 @kwarg radius: Optional mean earth radius (C{meter}). 998 @kwarg refraction: Consider atmospheric refraction (C{bool}). 999 1000 @return: Distance (C{meter}, same units as B{C{height}} and B{C{radius}}). 1001 1002 @raise ValueError: Invalid B{C{height}} or B{C{radius}}. 1003 1004 @see: U{Distance to horizon<https://www.EdWilliams.org/avform.htm#Horizon>}. 1005 ''' 1006 h, r = Height(height), Radius(radius) 1007 if min(h, r) < 0: 1008 raise _ValueError(height=height, radius=radius) 1009 1010 if refraction: 1011 d2 = 2.415750694528 * h * r # 2.0 / 0.8279 1012 else: 1013 d2 = h * fsum_(r, r, h) 1014 return sqrt0(d2) 1015 1016 1017def intersections2(lat1, lon1, radius1, 1018 lat2, lon2, radius2, datum=None, wrap=True): 1019 '''Conveniently compute the intersections of two circles each defined 1020 by a (geodetic) center point and a radius, using either ... 1021 1022 1) L{vector3d.intersections2} for small distances (below 100 KM or 1023 about 0.9 degrees) or if no B{C{datum}} is specified, or ... 1024 1025 2) L{sphericalTrigonometry.intersections2} for a spherical B{C{datum}} 1026 or if B{C{datum}} is a C{scalar} representing the earth radius, or ... 1027 1028 3) L{ellipsoidalKarney.intersections2} for an ellipsoidal B{C{datum}} 1029 and if I{Karney}'s U{geographiclib<https://PyPI.org/project/geographiclib>} 1030 is installed, or ... 1031 1032 4) L{ellipsoidalExact.intersections2} otherwise provided B{C{datum}} 1033 is ellipsoidal. 1034 1035 @arg lat1: Latitude of the first circle center (C{degrees}). 1036 @arg lon1: Longitude of the first circle center (C{degrees}). 1037 @arg radius1: Radius of the first circle (C{meter}, conventionally). 1038 @arg lat2: Latitude of the second circle center (C{degrees}). 1039 @arg lon2: Longitude of the second circle center (C{degrees}). 1040 @arg radius2: Radius of the second circle (C{meter}, same units as B{C{radius1}}). 1041 @kwarg datum: Optional ellipsoidal or spherical datum (L{Datum}, 1042 L{Ellipsoid}, L{Ellipsoid2}, L{a_f2Tuple} or 1043 C{scalar} earth radius in C{meter}, same units as 1044 B{C{radius1}} and B{C{radius2}}) or C{None}. 1045 @kwarg wrap: Wrap and unroll longitudes (C{bool}). 1046 1047 @return: 2-Tuple of the intersection points, each a 1048 L{LatLon2Tuple}C{(lat, lon)}. For abutting circles, 1049 the points are the same instance, aka I{radical center}. 1050 1051 @raise IntersectionError: Concentric, antipodal, invalid or 1052 non-intersecting circles or no 1053 convergence. 1054 1055 @raise TypeError: Invalid B{C{datum}}. 1056 1057 @raise UnitError: Invalid B{C{lat1}}, B{C{lon1}}, B{C{radius1}} 1058 B{C{lat2}}, B{C{lon2}} or B{C{radius2}}. 1059 ''' 1060 if datum is None or euclidean(lat1, lon1, lat1, lon2, radius=R_M, 1061 adjust=True, wrap=wrap) < _100km: 1062 import pygeodesy.vector3d as m 1063 1064 def _V2T(x, y, _, **unused): # _ == z unused 1065 return LatLon2Tuple(y, x, name=intersections2.__name__) 1066 1067 r1 = m2degrees(Radius_(radius1=radius1), radius=R_M, lat=lat1) 1068 r2 = m2degrees(Radius_(radius2=radius2), radius=R_M, lat=lat2) 1069 1070 _, lon2 = unroll180(lon1, lon2, wrap=wrap) 1071 t = m.intersections2(m.Vector3d(lon1, lat1, 0), r1, 1072 m.Vector3d(lon2, lat2, 0), r2, sphere=False, 1073 Vector=_V2T) 1074 1075 else: 1076 def _LL2T(lat, lon, **unused): 1077 return LatLon2Tuple(lat, lon, name=intersections2.__name__) 1078 1079 d = _spherical_datum(datum, name=intersections2.__name__) 1080 if d.isSpherical: 1081 import pygeodesy.sphericalTrigonometry as m 1082 elif d.isEllipsoidal: 1083 try: 1084 if d.ellipsoid.geodesic: 1085 pass 1086 import pygeodesy.ellipsoidalKarney as m 1087 except ImportError: 1088 import pygeodesy.ellipsoidalExact as m 1089 else: 1090 raise _AssertionError(datum=d) 1091 1092 t = m.intersections2(m.LatLon(lat1, lon1, datum=d), radius1, 1093 m.LatLon(lat2, lon2, datum=d), radius2, 1094 LatLon=_LL2T, height=0, wrap=wrap) 1095 return t 1096 1097 1098def isantipode(lat1, lon1, lat2, lon2, eps=EPS): 1099 '''Check whether two points are antipodal, on diametrically 1100 opposite sides of the earth. 1101 1102 @arg lat1: Latitude of one point (C{degrees}). 1103 @arg lon1: Longitude of one point (C{degrees}). 1104 @arg lat2: Latitude of the other point (C{degrees}). 1105 @arg lon2: Longitude of the other point (C{degrees}). 1106 @kwarg eps: Tolerance for near-equality (C{degrees}). 1107 1108 @return: C{True} if points are antipodal within the 1109 B{C{eps}} tolerance, C{False} otherwise. 1110 1111 @see: U{Geosphere<https://CRAN.R-Project.org/web/packages/geosphere/geosphere.pdf>}. 1112 ''' 1113 return abs(wrap90(lat1) + wrap90(lat2)) < eps and \ 1114 abs(abs(wrap180(lon1) - wrap180(lon2)) % _360_0 - _180_0) < eps 1115 1116 1117def isantipode_(phi1, lam1, phi2, lam2, eps=EPS): 1118 '''Check whether two points are antipodal, on diametrically 1119 opposite sides of the earth. 1120 1121 @arg phi1: Latitude of one point (C{radians}). 1122 @arg lam1: Longitude of one point (C{radians}). 1123 @arg phi2: Latitude of the other point (C{radians}). 1124 @arg lam2: Longitude of the other point (C{radians}). 1125 @kwarg eps: Tolerance for near-equality (C{radians}). 1126 1127 @return: C{True} if points are antipodal within the 1128 B{C{eps}} tolerance, C{False} otherwise. 1129 1130 @see: U{Geosphere<https://CRAN.R-Project.org/web/packages/geosphere/geosphere.pdf>}. 1131 ''' 1132 return abs(wrapPI_2(phi1) + wrapPI_2(phi2)) < eps and abs( 1133 abs(wrapPI(lam1) - wrapPI(lam2)) % PI2 - PI) < eps 1134 1135 1136def latlon2n_xyz(lat, lon, name=NN): 1137 '''Convert lat-, longitude to C{n-vector} (normal to the 1138 earth's surface) X, Y and Z components. 1139 1140 @arg lat: Latitude (C{degrees}). 1141 @arg lon: Longitude (C{degrees}). 1142 @kwarg name: Optional name (C{str}). 1143 1144 @return: A L{Vector3Tuple}C{(x, y, z)}. 1145 1146 @see: Function L{philam2n_xyz}. 1147 1148 @note: These are C{n-vector} x, y and z components, 1149 I{NOT} geocentric ECEF x, y and z coordinates! 1150 ''' 1151 return philam2n_xyz(radians(lat), radians(lon), name=name) 1152 1153 1154def n_xyz2latlon(x, y, z, name=NN): 1155 '''Convert C{n-vector} components to lat- and longitude in C{degrees}. 1156 1157 @arg x: X component (C{scalar}). 1158 @arg y: Y component (C{scalar}). 1159 @arg z: Z component (C{scalar}). 1160 @kwarg name: Optional name (C{str}). 1161 1162 @return: A L{LatLon2Tuple}C{(lat, lon)}. 1163 1164 @see: Function L{n_xyz2philam}. 1165 ''' 1166 a, b = n_xyz2philam(x, y, z) # PYCHOK PhiLam2Tuple 1167 return LatLon2Tuple(degrees90(a), degrees180(b), name=name) 1168 1169 1170def n_xyz2philam(x, y, z, name=NN): 1171 '''Convert C{n-vector} components to lat- and longitude in C{radians}. 1172 1173 @arg x: X component (C{scalar}). 1174 @arg y: Y component (C{scalar}). 1175 @arg z: Z component (C{scalar}). 1176 @kwarg name: Optional name (C{str}). 1177 1178 @return: A L{PhiLam2Tuple}C{(phi, lam)}. 1179 1180 @see: Function L{n_xyz2latlon}. 1181 ''' 1182 return PhiLam2Tuple(atan2(z, hypot(x, y)), atan2(y, x), name=name) 1183 1184 1185def opposing(bearing1, bearing2, margin=None): 1186 '''Compare the direction of two bearings given in C{degrees}. 1187 1188 @arg bearing1: First bearing (compass C{degrees}). 1189 @arg bearing2: Second bearing (compass C{degrees}). 1190 @kwarg margin: Optional, interior angle bracket (C{degrees}), 1191 default C{90}. 1192 1193 @return: C{True} if both bearings point in opposite, C{False} if 1194 in similar or C{None} if in perpendicular directions. 1195 1196 @see: Function L{opposing_}. 1197 ''' 1198 m = Degrees_(margin=margin, low=EPS0, high=_90_0) if margin else _90_0 1199 d = (bearing2 - bearing1) % _360_0 # note -20 % 360 == 340 1200 return False if d < m or d > (_360_0 - m) else ( 1201 True if (_180_0 - m) < d < (_180_0 + m) else None) 1202 1203 1204def opposing_(radians1, radians2, margin=None): 1205 '''Compare the direction of two bearings given in C{radians}. 1206 1207 @arg radians1: First bearing (C{radians}). 1208 @arg radians2: Second bearing (C{radians}). 1209 @kwarg margin: Optional, interior angle bracket (C{radians}), 1210 default C{PI_2}. 1211 1212 @return: C{True} if both bearings point in opposite, C{False} if 1213 in similar or C{None} if in perpendicular directions. 1214 1215 @see: Function L{opposing}. 1216 ''' 1217 m = Radians_(margin=margin, low=EPS0, high=PI_2) if margin else PI_2 1218 r = (radians2 - radians1) % PI2 # note -1 % PI2 == PI2 - 1 1219 return False if r < m or r > (PI2 - m) else ( 1220 True if (PI - m) < r < (PI + m) else None) 1221 1222 1223def philam2n_xyz(phi, lam, name=NN): 1224 '''Convert lat-, longitude to C{n-vector} (normal to the 1225 earth's surface) X, Y and Z components. 1226 1227 @arg phi: Latitude (C{radians}). 1228 @arg lam: Longitude (C{radians}). 1229 @kwarg name: Optional name (C{str}). 1230 1231 @return: A L{Vector3Tuple}C{(x, y, z)}. 1232 1233 @see: Function L{latlon2n_xyz}. 1234 1235 @note: These are C{n-vector} x, y and z components, 1236 I{NOT} geocentric ECEF x, y and z coordinates! 1237 ''' 1238 # Kenneth Gade eqn 3, but using right-handed 1239 # vector x -> 0°E,0°N, y -> 90°E,0°N, z -> 90°N 1240 sa, ca, sb, cb = sincos2_(phi, lam) 1241 return Vector3Tuple(ca * cb, ca * sb, sa, name=name) 1242 1243 1244def _radical2(d, r1, r2): # in .ellipsoidalBaseDI, .sphericalTrigonometry, .vector3d 1245 # (INTERNAL) See C{radical2} below 1246 # assert d > EPS0 1247 r = fsum_(_1_0, (r1 / d)**2, -(r2 / d)**2) * _0_5 1248 return Radical2Tuple(max(_0_0, min(_1_0, r)), r * d) 1249 1250 1251def radical2(distance, radius1, radius2): 1252 '''Compute the I{radical ratio} and I{radical line} of two 1253 U{intersecting circles<https://MathWorld.Wolfram.com/ 1254 Circle-CircleIntersection.html>}. 1255 1256 The I{radical line} is perpendicular to the axis thru the 1257 centers of the circles at C{(0, 0)} and C{(B{distance}, 0)}. 1258 1259 @arg distance: Distance between the circle centers (C{scalar}). 1260 @arg radius1: Radius of the first circle (C{scalar}). 1261 @arg radius2: Radius of the second circle (C{scalar}). 1262 1263 @return: A L{Radical2Tuple}C{(ratio, xline)} where C{0.0 <= 1264 ratio <= 1.0} and C{xline} is along the B{C{distance}}. 1265 1266 @raise IntersectionError: The B{C{distance}} exceeds the sum 1267 of B{C{radius1}} and B{C{radius2}}. 1268 1269 @raise UnitError: Invalid B{C{distance}}, B{C{radius1}} or 1270 B{C{radius2}}. 1271 1272 @see: U{Circle-Circle Intersection 1273 <https://MathWorld.Wolfram.com/Circle-CircleIntersection.html>}. 1274 ''' 1275 d = Distance_(distance, low=_0_0) 1276 r1 = Radius_(radius1=radius1) 1277 r2 = Radius_(radius2=radius2) 1278 if d > (r1 + r2): 1279 raise IntersectionError(distance=d, radius1=r1, radius2=r2, 1280 txt=_too_(_distant_)) 1281 return _radical2(d, r1, r2) if d > EPS0 else \ 1282 Radical2Tuple(_0_5, _0_0) 1283 1284 1285class Radical2Tuple(_NamedTuple): 1286 '''2-Tuple C{(ratio, xline)} of the I{radical} C{ratio} and 1287 I{radical} C{xline}, both C{scalar} and C{0.0 <= ratio <= 1.0} 1288 ''' 1289 _Names_ = (_ratio_, _xline_) 1290 _Units_ = ( Scalar, Scalar) 1291 1292 1293def _scale_deg(lat1, lat2): # degrees 1294 # scale factor cos(mean of lats) for delta lon 1295 m = abs(lat1 + lat2) * _0_5 1296 return cos(radians(m)) if m < _90_0 else _0_0 1297 1298 1299def _scale_rad(phi1, phi2): # radians, by .frechet, .hausdorff, .heights 1300 # scale factor cos(mean of phis) for delta lam 1301 m = abs(phi1 + phi2) * _0_5 1302 return cos(m) if m < PI_2 else _0_0 1303 1304 1305def _sincosa6(phi2, phi1, lam21): 1306 '''(INTERNAL) C{sin}es, C{cos}ines and C{acos}ine. 1307 ''' 1308 s2, c2, s1, c1, _, c21 = sincos2_(phi2, phi1, lam21) 1309 return s2, c2, s1, c1, acos1(s1 * s2 + c1 * c2 * c21), c21 1310 1311 1312def thomas(lat1, lon1, lat2, lon2, datum=_WGS84, wrap=False): 1313 '''Compute the distance between two (ellipsoidal) points using 1314 U{Thomas'<https://apps.DTIC.mil/dtic/tr/fulltext/u2/703541.pdf>} 1315 formula. 1316 1317 @arg lat1: Start latitude (C{degrees}). 1318 @arg lon1: Start longitude (C{degrees}). 1319 @arg lat2: End latitude (C{degrees}). 1320 @arg lon2: End longitude (C{degrees}). 1321 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid}, 1322 L{Ellipsoid2} or L{a_f2Tuple}). 1323 @kwarg wrap: Wrap and L{unroll180} longitudes (C{bool}). 1324 1325 @return: Distance (C{meter}, same units as the B{C{datum}}'s 1326 ellipsoid axes). 1327 1328 @raise TypeError: Invalid B{C{datum}}. 1329 1330 @see: Functions L{thomas_}, L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert}, 1331 L{cosineLaw}, L{equirectangular}, L{euclidean}, L{flatLocal}/L{hubeny}, 1332 L{flatPolar}, L{haversine}, L{vincentys} and method L{Ellipsoid.distance2}. 1333 ''' 1334 return _distanceToE(thomas_, lat1, lat2, datum, 1335 *unroll180(lon1, lon2, wrap=wrap)) 1336 1337 1338def thomas_(phi2, phi1, lam21, datum=_WGS84): 1339 '''Compute the I{angular} distance between two (ellipsoidal) points using 1340 U{Thomas'<https://apps.DTIC.mil/dtic/tr/fulltext/u2/703541.pdf>} 1341 formula. 1342 1343 @arg phi2: End latitude (C{radians}). 1344 @arg phi1: Start latitude (C{radians}). 1345 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 1346 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid}, 1347 L{Ellipsoid2} or L{a_f2Tuple}). 1348 1349 @return: Angular distance (C{radians}). 1350 1351 @raise TypeError: Invalid B{C{datum}}. 1352 1353 @see: Functions L{thomas}, L{cosineAndoyerLambert_}, 1354 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, 1355 L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_}, 1356 L{flatPolar_}, L{haversine_} and L{vincentys_} and U{Geodesy-PHP 1357 <https://GitHub.com/jtejido/geodesy-php/blob/master/src/Geodesy/ 1358 Distance/ThomasFormula.php>}. 1359 ''' 1360 s2, c2, s1, c1, r, _ = _sincosa6(phi2, phi1, lam21) 1361 if r and _non0(c1) and _non0(c2): 1362 E = _ellipsoidal(datum, thomas_) 1363 if E.f: 1364 r1 = atan2(E.b_a * s1, c1) 1365 r2 = atan2(E.b_a * s2, c2) 1366 1367 j = (r2 + r1) * _0_5 1368 k = (r2 - r1) * _0_5 1369 sj, cj, sk, ck, h, _ = sincos2_(j, k, lam21 * _0_5) 1370 1371 h = fsum_(sk**2, (ck * h)**2, -(sj * h)**2) 1372 u = _1_0 - h 1373 if _non0(u) and _non0(h): 1374 r = atan(sqrt0(h / u)) * _2_0 # == acos(1 - 2 * h) 1375 sr, cr = sincos2(r) 1376 if _non0(sr): 1377 u = 2 * (sj * ck)**2 / u 1378 h = 2 * (sk * cj)**2 / h 1379 x = u + h 1380 y = u - h 1381 1382 s = r / sr 1383 e = 4 * s**2 1384 d = 2 * cr 1385 a = e * d 1386 b = 2 * r 1387 c = s - (a - d) * _0_5 1388 f = E.f * _0_25 1389 1390 t = fsum_(a * x, -b * y, c * x**2, -d * y**2, e * x * y) 1391 r -= fsum_(s * x, -y, -t * f * _0_25) * f * sr 1392 return r 1393 1394 1395def vincentys(lat1, lon1, lat2, lon2, radius=R_M, wrap=False): 1396 '''Compute the distance between two (spherical) points using 1397 U{Vincenty's<https://WikiPedia.org/wiki/Great-circle_distance>} 1398 spherical formula. 1399 1400 @arg lat1: Start latitude (C{degrees}). 1401 @arg lon1: Start longitude (C{degrees}). 1402 @arg lat2: End latitude (C{degrees}). 1403 @arg lon2: End longitude (C{degrees}). 1404 @kwarg radius: Mean earth radius, ellipsoid or datum 1405 (C{meter}, L{Ellipsoid}, L{Ellipsoid2}, 1406 L{Datum} or L{a_f2Tuple}). 1407 @kwarg wrap: Wrap and L{unroll180} longitudes (C{bool}). 1408 1409 @return: Distance (C{meter}, same units as B{C{radius}}). 1410 1411 @raise UnitError: Invalid B{C{radius}}. 1412 1413 @see: Functions L{vincentys_}, L{cosineAndoyerLambert}, 1414 L{cosineForsytheAndoyerLambert},L{cosineLaw}, L{equirectangular}, 1415 L{euclidean}, L{flatLocal}/L{hubeny}, L{flatPolar}, 1416 L{haversine} and L{thomas} and methods L{Ellipsoid.distance2}, 1417 C{LatLon.distanceTo*} and C{LatLon.equirectangularTo}. 1418 1419 @note: See note at function L{vincentys_}. 1420 ''' 1421 return _distanceToS(vincentys_, lat1, lat2, radius, 1422 *unroll180(lon1, lon2, wrap=wrap)) 1423 1424 1425def vincentys_(phi2, phi1, lam21): 1426 '''Compute the I{angular} distance between two (spherical) points using 1427 U{Vincenty's<https://WikiPedia.org/wiki/Great-circle_distance>} 1428 spherical formula. 1429 1430 @arg phi2: End latitude (C{radians}). 1431 @arg phi1: Start latitude (C{radians}). 1432 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 1433 1434 @return: Angular distance (C{radians}). 1435 1436 @see: Functions L{vincentys}, L{cosineAndoyerLambert_}, 1437 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, 1438 L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_}, 1439 L{flatPolar_}, L{haversine_} and L{thomas_}. 1440 1441 @note: Functions L{vincentys_}, L{haversine_} and L{cosineLaw_} 1442 produce equivalent results, but L{vincentys_} is suitable 1443 for antipodal points and slightly more expensive (M{3 cos, 1444 3 sin, 1 hypot, 1 atan2, 6 mul, 2 add}) than L{haversine_} 1445 (M{2 cos, 2 sin, 2 sqrt, 1 atan2, 5 mul, 1 add}) and 1446 L{cosineLaw_} (M{3 cos, 3 sin, 1 acos, 3 mul, 1 add}). 1447 ''' 1448 s1, c1, s2, c2, s21, c21 = sincos2_(phi1, phi2, lam21) 1449 1450 c = c2 * c21 1451 x = s1 * s2 + c1 * c 1452 y = c1 * s2 - s1 * c 1453 return atan2(hypot(c2 * s21, y), x) 1454 1455# **) MIT License 1456# 1457# Copyright (C) 2016-2021 -- mrJean1 at Gmail -- All Rights Reserved. 1458# 1459# Permission is hereby granted, free of charge, to any person obtaining a 1460# copy of this software and associated documentation files (the "Software"), 1461# to deal in the Software without restriction, including without limitation 1462# the rights to use, copy, modify, merge, publish, distribute, sublicense, 1463# and/or sell copies of the Software, and to permit persons to whom the 1464# Software is furnished to do so, subject to the following conditions: 1465# 1466# The above copyright notice and this permission notice shall be included 1467# in all copies or substantial portions of the Software. 1468# 1469# THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS 1470# OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, 1471# FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL 1472# THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR 1473# OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, 1474# ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR 1475# OTHER DEALINGS IN THE SOFTWARE. 1476