Coverage for pygeodesy/formy.py: 99%

409 statements  

« prev     ^ index     » next       coverage.py v7.2.2, created at 2023-06-07 08:37 -0400

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 as _; del _ # PYCHOK semicolon 

8 

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 IntersectionError, LimitError, limiterrors, \ 

18 _TypeError, _ValueError, \ 

19 _xError, _xkwds, _xkwds_pop 

20from pygeodesy.fmath import euclid, hypot, hypot2, sqrt0 

21from pygeodesy.fsums import fsumf_, isscalar 

22from pygeodesy.interns import NN, _distant_, _SPACE_, _too_ 

23from pygeodesy.lazily import _ALL_LAZY, _ALL_MODS as _MODS 

24from pygeodesy.named import _NamedTuple, _xnamed, Fmt, unstr 

25from pygeodesy.namedTuples import Bearing2Tuple, Distance4Tuple, \ 

26 Intersection3Tuple, LatLon2Tuple, \ 

27 PhiLam2Tuple, Vector3Tuple 

28# from pygeodesy.streprs import Fmt, unstr # from .named 

29from pygeodesy.units import Bearing, Degrees_, Distance, Distance_, Height, \ 

30 Lam_, Lat, Lon, Meter_, Phi_, Radians, Radians_, \ 

31 Radius, Radius_, Scalar, _100km 

32from pygeodesy.utily import acos1, atan2b, atan2d, degrees2m, m2degrees, \ 

33 tan_2, sincos2, sincos2_, sincos2d_, _Wrap 

34 

35from contextlib import contextmanager 

36from math import atan, atan2, cos, degrees, fabs, radians, sin, sqrt # pow 

37 

38__all__ = _ALL_LAZY.formy 

39__version__ = '23.05.15' 

40 

41_delta_ = 'delta' 

42_D2_R2 = (PI / _180_0)**2 # degrees- to radians-squared 

43_EWGS84 = _WGS84.ellipsoid 

44_ratio_ = 'ratio' 

45_xline_ = 'xline' 

46 

47 

48def _anti2(a, b, n_2, n, n2): 

49 '''(INTERNAL) Helper for C{antipode} and C{antipode_}. 

50 ''' 

51 r = remainder(a, n) if fabs(a) > n_2 else a 

52 if r == a: 

53 r = -r 

54 b += n 

55 if fabs(b) > n: 

56 b = remainder(b, n2) 

57 return float0_(r, b) 

58 

59 

60def antipode(lat, lon, name=NN): 

61 '''Return the antipode, the point diametrically opposite 

62 to a given point in C{degrees}. 

63 

64 @arg lat: Latitude (C{degrees}). 

65 @arg lon: Longitude (C{degrees}). 

66 @kwarg name: Optional name (C{str}). 

67 

68 @return: A L{LatLon2Tuple}C{(lat, lon)}. 

69 

70 @see: Functions L{antipode_} and L{normal} and U{Geosphere 

71 <https://CRAN.R-Project.org/web/packages/geosphere/geosphere.pdf>}. 

72 ''' 

73 return LatLon2Tuple(*_anti2(lat, lon, _90_0, _180_0, _360_0), name=name) 

74 

75 

76def antipode_(phi, lam, name=NN): 

77 '''Return the antipode, the point diametrically opposite 

78 to a given point in C{radians}. 

79 

80 @arg phi: Latitude (C{radians}). 

81 @arg lam: Longitude (C{radians}). 

82 @kwarg name: Optional name (C{str}). 

83 

84 @return: A L{PhiLam2Tuple}C{(phi, lam)}. 

85 

86 @see: Functions L{antipode} and L{normal_} and U{Geosphere 

87 <https://CRAN.R-Project.org/web/packages/geosphere/geosphere.pdf>}. 

88 ''' 

89 return PhiLam2Tuple(*_anti2(phi, lam, PI_2, PI, PI2), name=name) 

90 

91 

92def bearing(lat1, lon1, lat2, lon2, **final_wrap): 

93 '''Compute the initial or final bearing (forward or reverse 

94 azimuth) between a (spherical) start and end point. 

95 

96 @arg lat1: Start latitude (C{degrees}). 

97 @arg lon1: Start longitude (C{degrees}). 

98 @arg lat2: End latitude (C{degrees}). 

99 @arg lon2: End longitude (C{degrees}). 

100 @kwarg final_wrap: Optional keyword arguments for function 

101 L{pygeodesy.bearing_}. 

102 

103 @return: Initial or final bearing (compass C{degrees360}) or 

104 zero if start and end point coincide. 

105 ''' 

106 r = bearing_(Phi_(lat1=lat1), Lam_(lon1=lon1), 

107 Phi_(lat2=lat2), Lam_(lon2=lon2), **final_wrap) 

108 return degrees(r) 

109 

110 

111def bearing_(phi1, lam1, phi2, lam2, final=False, wrap=False): 

112 '''Compute the initial or final bearing (forward or reverse azimuth) 

113 between a (spherical) start and end point. 

114 

115 @arg phi1: Start latitude (C{radians}). 

116 @arg lam1: Start longitude (C{radians}). 

117 @arg phi2: End latitude (C{radians}). 

118 @arg lam2: End longitude (C{radians}). 

119 @kwarg final: Return final bearing if C{True}, initial otherwise (C{bool}). 

120 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll B{C{phi2}} and 

121 B{C{lam2}} (C{bool}). 

122 

123 @return: Initial or final bearing (compass C{radiansPI2}) or zero if start 

124 and end point coincide. 

125 

126 @see: U{Bearing<https://www.Movable-Type.co.UK/scripts/latlong.html>}, U{Course 

127 between two points<https://www.EdWilliams.org/avform147.htm#Crs>} and 

128 U{Bearing Between Two Points<https://web.Archive.org/web/20020630205931/ 

129 https://MathForum.org/library/drmath/view/55417.html>}. 

130 ''' 

131 db, phi2, lam2 = _Wrap.philam3(lam1, phi2, lam2, wrap) 

132 if final: # swap plus PI 

133 phi1, lam1, phi2, lam2, db = phi2, lam2, phi1, lam1, -db 

134 r = PI3 

135 else: 

136 r = PI2 

137 sa1, ca1, sa2, ca2, sdb, cdb = sincos2_(phi1, phi2, db) 

138 

139 x = ca1 * sa2 - sa1 * ca2 * cdb 

140 y = sdb * ca2 

141 return _umod_PI2(atan2(y, x) + r) # .utily.wrapPI2 

142 

143 

144def _bearingTo2(p1, p2, wrap=False): # for points.ispolar, sphericalTrigonometry.areaOf 

145 '''(INTERNAL) Compute initial and final bearing. 

146 ''' 

147 try: # for LatLon_ and ellipsoidal LatLon 

148 return p1.bearingTo2(p2, wrap=wrap) 

149 except AttributeError: 

150 pass 

151 # XXX spherical version, OK for ellipsoidal ispolar? 

152 a1, b1 = p1.philam 

153 a2, b2 = p2.philam 

154 return Bearing2Tuple(degrees(bearing_(a1, b1, a2, b2, final=False, wrap=wrap)), 

155 degrees(bearing_(a1, b1, a2, b2, final=True, wrap=wrap)), 

156 name=_bearingTo2.__name__) 

157 

158 

159def compassAngle(lat1, lon1, lat2, lon2, adjust=True, wrap=False): 

160 '''Return the angle from North for the direction vector M{(lon2 - lon1, 

161 lat2 - lat1)} between two points. 

162 

163 Suitable only for short, not near-polar vectors up to a few hundred 

164 Km or Miles. Use function L{pygeodesy.bearing} for longer vectors. 

165 

166 @arg lat1: From latitude (C{degrees}). 

167 @arg lon1: From longitude (C{degrees}). 

168 @arg lat2: To latitude (C{degrees}). 

169 @arg lon2: To longitude (C{degrees}). 

170 @kwarg adjust: Adjust the longitudinal delta by the cosine of the 

171 mean latitude (C{bool}). 

172 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll B{C{lat2}} 

173 and B{C{lon2}} (C{bool}). 

174 

175 @return: Compass angle from North (C{degrees360}). 

176 

177 @note: Courtesy of Martin Schultz. 

178 

179 @see: U{Local, flat earth approximation 

180 <https://www.EdWilliams.org/avform.htm#flat>}. 

181 ''' 

182 d_lon, lat2, lon2 = _Wrap.latlon3(lon1, lat2, lon2, wrap) 

183 if adjust: # scale delta lon 

184 d_lon *= _scale_deg(lat1, lat2) 

185 return atan2b(d_lon, lat2 - lat1) 

186 

187 

188def cosineAndoyerLambert(lat1, lon1, lat2, lon2, datum=_WGS84, wrap=False): 

189 '''Compute the distance between two (ellipsoidal) points using the 

190 U{Andoyer-Lambert correction<https://NavLib.net/wp-content/uploads/2013/10/ 

191 admiralty-manual-of-navigation-vol-1-1964-english501c.pdf>} of the U{Law of 

192 Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>} formula. 

193 

194 @arg lat1: Start latitude (C{degrees}). 

195 @arg lon1: Start longitude (C{degrees}). 

196 @arg lat2: End latitude (C{degrees}). 

197 @arg lon2: End longitude (C{degrees}). 

198 @kwarg datum: Datum (L{Datum}) or ellipsoid (L{Ellipsoid}, 

199 L{Ellipsoid2} or L{a_f2Tuple}) to use. 

200 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll 

201 B{C{lat2}} and B{C{lon2}} (C{bool}). 

202 

203 @return: Distance (C{meter}, same units as the B{C{datum}}'s 

204 ellipsoid axes or C{radians} if B{C{datum}} is C{None}). 

205 

206 @raise TypeError: Invalid B{C{datum}}. 

207 

208 @see: Functions L{cosineAndoyerLambert_}, L{cosineForsytheAndoyerLambert}, 

209 L{cosineLaw}, L{equirectangular}, L{euclidean}, L{flatLocal}/L{hubeny}, 

210 L{flatPolar}, L{haversine}, L{thomas} and L{vincentys} and method 

211 L{Ellipsoid.distance2}. 

212 ''' 

213 return _dE(cosineAndoyerLambert_, datum, wrap, lat1, lon1, lat2, lon2) 

214 

215 

216def cosineAndoyerLambert_(phi2, phi1, lam21, datum=_WGS84): 

217 '''Compute the I{angular} distance between two (ellipsoidal) points using the 

218 U{Andoyer-Lambert correction<https://NavLib.net/wp-content/uploads/2013/10/ 

219 admiralty-manual-of-navigation-vol-1-1964-english501c.pdf>} of the U{Law of 

220 Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>} formula. 

221 

222 @arg phi2: End latitude (C{radians}). 

223 @arg phi1: Start latitude (C{radians}). 

224 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 

225 @kwarg datum: Datum (L{Datum}) or ellipsoid (L{Ellipsoid}, 

226 L{Ellipsoid2} or L{a_f2Tuple}) to use. 

227 

228 @return: Angular distance (C{radians}). 

229 

230 @raise TypeError: Invalid B{C{datum}}. 

231 

232 @see: Functions L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert_}, 

233 L{cosineLaw_}, L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_}, 

234 L{flatPolar_}, L{haversine_}, L{thomas_} and L{vincentys_} and U{Geodesy-PHP 

235 <https://GitHub.com/jtejido/geodesy-php/blob/master/src/Geodesy/Distance/ 

236 AndoyerLambert.php>}. 

237 ''' 

238 s2, c2, s1, c1, r, c21 = _sincosa6(phi2, phi1, lam21) 

239 if isnon0(c1) and isnon0(c2): 

240 E = _ellipsoidal(datum, cosineAndoyerLambert_) 

241 if E.f: # ellipsoidal 

242 r2 = atan2(E.b_a * s2, c2) 

243 r1 = atan2(E.b_a * s1, c1) 

244 s2, c2, s1, c1 = sincos2_(r2, r1) 

245 r = acos1(s1 * s2 + c1 * c2 * c21) 

246 if r: 

247 sr, _, sr_2, cr_2 = sincos2_(r, r * _0_5) 

248 if isnon0(sr_2) and isnon0(cr_2): 

249 s = (sr + r) * ((s1 - s2) / sr_2)**2 

250 c = (sr - r) * ((s1 + s2) / cr_2)**2 

251 r += (c - s) * E.f * _0_125 

252 return r 

253 

254 

255def cosineForsytheAndoyerLambert(lat1, lon1, lat2, lon2, datum=_WGS84, wrap=False): 

256 '''Compute the distance between two (ellipsoidal) points using the 

257 U{Forsythe-Andoyer-Lambert correction<https://www2.UNB.Ca/gge/Pubs/TR77.pdf>} of 

258 the U{Law of Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>} 

259 formula. 

260 

261 @arg lat1: Start latitude (C{degrees}). 

262 @arg lon1: Start longitude (C{degrees}). 

263 @arg lat2: End latitude (C{degrees}). 

264 @arg lon2: End longitude (C{degrees}). 

265 @kwarg datum: Datum (L{Datum}) or ellipsoid (L{Ellipsoid}, 

266 L{Ellipsoid2} or L{a_f2Tuple}) to use. 

267 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll 

268 B{C{lat2}} and B{C{lon2}} (C{bool}). 

269 

270 @return: Distance (C{meter}, same units as the B{C{datum}}'s 

271 ellipsoid axes or C{radians} if B{C{datum}} is C{None}). 

272 

273 @raise TypeError: Invalid B{C{datum}}. 

274 

275 @see: Functions L{cosineForsytheAndoyerLambert_}, L{cosineAndoyerLambert}, 

276 L{cosineLaw}, L{equirectangular}, L{euclidean}, L{flatLocal}/L{hubeny}, 

277 L{flatPolar}, L{haversine}, L{thomas} and L{vincentys} and method 

278 L{Ellipsoid.distance2}. 

279 ''' 

280 return _dE(cosineForsytheAndoyerLambert_, datum, wrap, lat1, lon1, lat2, lon2) 

281 

282 

283def cosineForsytheAndoyerLambert_(phi2, phi1, lam21, datum=_WGS84): 

284 '''Compute the I{angular} distance between two (ellipsoidal) points using the 

285 U{Forsythe-Andoyer-Lambert correction<https://www2.UNB.Ca/gge/Pubs/TR77.pdf>} of 

286 the U{Law of Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>} 

287 formula. 

288 

289 @arg phi2: End latitude (C{radians}). 

290 @arg phi1: Start latitude (C{radians}). 

291 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 

292 @kwarg datum: Datum (L{Datum}) or ellipsoid to use (L{Ellipsoid}, 

293 L{Ellipsoid2} or L{a_f2Tuple}). 

294 

295 @return: Angular distance (C{radians}). 

296 

297 @raise TypeError: Invalid B{C{datum}}. 

298 

299 @see: Functions L{cosineForsytheAndoyerLambert}, L{cosineAndoyerLambert_}, 

300 L{cosineLaw_}, L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_}, 

301 L{flatPolar_}, L{haversine_}, L{thomas_} and L{vincentys_} and U{Geodesy-PHP 

302 <https://GitHub.com/jtejido/geodesy-php/blob/master/src/Geodesy/ 

303 Distance/ForsytheCorrection.php>}. 

304 ''' 

305 s2, c2, s1, c1, r, _ = _sincosa6(phi2, phi1, lam21) 

306 if r and isnon0(c1) and isnon0(c2): 

307 E = _ellipsoidal(datum, cosineForsytheAndoyerLambert_) 

308 if E.f: # ellipsoidal 

309 sr, cr, s2r, _ = sincos2_(r, r * 2) 

310 if isnon0(sr) and fabs(cr) < EPS1: 

311 s = (s1 + s2)**2 / (1 + cr) 

312 t = (s1 - s2)**2 / (1 - cr) 

313 x = s + t 

314 y = s - t 

315 

316 s = 8 * r**2 / sr 

317 a = 64 * r + s * cr * 2 # 16 * r**2 / tan(r) 

318 d = 48 * sr + s # 8 * r**2 / tan(r) 

319 b = -2 * d 

320 e = 30 * s2r 

321 c = fsumf_(30 * r, e * _0_5, s * cr) # 8 * r**2 / tan(r) 

322 

323 t = fsumf_( a * x, b * y, -c * x**2, d * x * y, e * y**2) 

324 r += fsumf_(-r * x, 3 * y * sr, t * E.f / _32_0) * E.f * _0_25 

325 return r 

326 

327 

328def cosineLaw(lat1, lon1, lat2, lon2, radius=R_M, wrap=False): 

329 '''Compute the distance between two points using the U{spherical Law of 

330 Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>} 

331 formula. 

332 

333 @arg lat1: Start latitude (C{degrees}). 

334 @arg lon1: Start longitude (C{degrees}). 

335 @arg lat2: End latitude (C{degrees}). 

336 @arg lon2: End longitude (C{degrees}). 

337 @kwarg radius: Mean earth radius (C{meter}), datum (L{Datum}) 

338 or ellipsoid (L{Ellipsoid}, L{Ellipsoid2} or 

339 L{a_f2Tuple}) to use. 

340 @kwarg wrap: If C{True}, wrap or I{normalize} and B{C{lat2}} 

341 and B{C{lon2}} (C{bool}). 

342 

343 @return: Distance (C{meter}, same units as B{C{radius}} or the 

344 ellipsoid or datum axes). 

345 

346 @raise TypeError: Invalid B{C{radius}}. 

347 

348 @see: Functions L{cosineLaw_}, L{cosineAndoyerLambert}, 

349 L{cosineForsytheAndoyerLambert}, L{equirectangular}, L{euclidean}, 

350 L{flatLocal}/L{hubeny}, L{flatPolar}, L{haversine}, L{thomas} and 

351 L{vincentys} and method L{Ellipsoid.distance2}. 

352 

353 @note: See note at function L{vincentys_}. 

354 ''' 

355 return _dS(cosineLaw_, radius, wrap, lat1, lon1, lat2, lon2) 

356 

357 

358def cosineLaw_(phi2, phi1, lam21): 

359 '''Compute the I{angular} distance between two points using the U{spherical 

360 Law of Cosines<https://www.Movable-Type.co.UK/scripts/latlong.html#cosine-law>} 

361 formula. 

362 

363 @arg phi2: End latitude (C{radians}). 

364 @arg phi1: Start latitude (C{radians}). 

365 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 

366 

367 @return: Angular distance (C{radians}). 

368 

369 @see: Functions L{cosineLaw}, L{cosineAndoyerLambert_}, 

370 L{cosineForsytheAndoyerLambert_}, L{equirectangular_}, 

371 L{euclidean_}, L{flatLocal_}/L{hubeny_}, L{flatPolar_}, 

372 L{haversine_}, L{thomas_} and L{vincentys_}. 

373 

374 @note: See note at function L{vincentys_}. 

375 ''' 

376 return _sincosa6(phi2, phi1, lam21)[4] 

377 

378 

379def _d3(wrap, lat1, lon1, lat2, lon2): 

380 '''(INTERNAL) Helper for _dE, _dS and _eA. 

381 ''' 

382 if wrap: 

383 d_lon, lat2, _ = _Wrap.latlon3(lon1, lat2, lon2, wrap) 

384 return radians(lat2), Phi_(lat1=lat1), radians(d_lon) 

385 else: # for backward compaibility 

386 return Phi_(lat2=lat2), Phi_(lat1=lat1), Phi_(d_lon=lon2 - lon1) 

387 

388 

389def _dE(func_, earth, *wrap_lls): 

390 '''(INTERNAL) Helper for ellipsoidal distances. 

391 ''' 

392 E = _ellipsoidal(earth, func_) 

393 r = func_(*_d3(*wrap_lls), datum=E) 

394 return r * E.a 

395 

396 

397def _dS(func_, radius, *wrap_lls, **adjust): 

398 '''(INTERNAL) Helper for spherical distances. 

399 ''' 

400 r = func_(*_d3(*wrap_lls), **adjust) 

401 if radius is not R_M: 

402 _, lat1, _, lat2, _ = wrap_lls 

403 radius = _mean_radius(radius, lat1, lat2) 

404 return r * radius 

405 

406 

407def _eA(excess_, radius, *wrap_lls): 

408 '''(INTERNAL) Helper for spherical excess or area. 

409 ''' 

410 r = excess_(*_d3(*wrap_lls)) 

411 if radius: 

412 _, lat1, _, lat2, _ = wrap_lls 

413 r *= _mean_radius(radius, lat1, lat2)**2 

414 return r 

415 

416 

417def _ellipsoidal(earth, where): 

418 '''(INTERNAL) Helper for distances. 

419 ''' 

420 return _EWGS84 if earth in (_WGS84, _EWGS84) else ( 

421 earth if isinstance(earth, Ellipsoid) else 

422 (earth if isinstance(earth, Datum) else # PYCHOK indent 

423 _ellipsoidal_datum(earth, name=where.__name__)).ellipsoid) 

424 

425 

426def equirectangular(lat1, lon1, lat2, lon2, radius=R_M, **adjust_limit_wrap): 

427 '''Compute the distance between two points using 

428 the U{Equirectangular Approximation / Projection 

429 <https://www.Movable-Type.co.UK/scripts/latlong.html#equirectangular>}. 

430 

431 @arg lat1: Start latitude (C{degrees}). 

432 @arg lon1: Start longitude (C{degrees}). 

433 @arg lat2: End latitude (C{degrees}). 

434 @arg lon2: End longitude (C{degrees}). 

435 @kwarg radius: Mean earth radius (C{meter}), datum (L{Datum}) 

436 or ellipsoid (L{Ellipsoid}, L{Ellipsoid2} or 

437 L{a_f2Tuple}). 

438 @kwarg adjust_limit_wrap: Optional keyword arguments for 

439 function L{equirectangular_}. 

440 

441 @return: Distance (C{meter}, same units as B{C{radius}} or 

442 the ellipsoid or datum axes). 

443 

444 @raise TypeError: Invalid B{C{radius}}. 

445 

446 @see: Function L{equirectangular_} for more details, the 

447 available B{C{options}}, errors, restrictions and other, 

448 approximate or accurate distance functions. 

449 ''' 

450 d = sqrt(equirectangular_(Lat(lat1=lat1), Lon(lon1=lon1), 

451 Lat(lat2=lat2), Lon(lon2=lon2), 

452 **adjust_limit_wrap).distance2) # PYCHOK 4 vs 2-3 

453 return degrees2m(d, radius=_mean_radius(radius, lat1, lat2)) 

454 

455 

456def _equirectangular(lat1, lon1, lat2, lon2, **adjust_limit_wrap): 

457 '''(INTERNAL) Helper for the L{frechet._FrecherMeterRadians} 

458 and L{hausdorff._HausdorffMeterRedians} classes. 

459 ''' 

460 return equirectangular_(lat1, lon1, lat2, lon2, **adjust_limit_wrap).distance2 * _D2_R2 

461 

462 

463def equirectangular_(lat1, lon1, lat2, lon2, adjust=True, limit=45, wrap=False): 

464 '''Compute the distance between two points using the U{Equirectangular 

465 Approximation / Projection 

466 <https://www.Movable-Type.co.UK/scripts/latlong.html#equirectangular>}. 

467 

468 This approximation is valid for short distance of several hundred Km 

469 or Miles, see the B{C{limit}} keyword argument and L{LimitError}. 

470 

471 @arg lat1: Start latitude (C{degrees}). 

472 @arg lon1: Start longitude (C{degrees}). 

473 @arg lat2: End latitude (C{degrees}). 

474 @arg lon2: End longitude (C{degrees}). 

475 @kwarg adjust: Adjust the wrapped, unrolled longitudinal delta 

476 by the cosine of the mean latitude (C{bool}). 

477 @kwarg limit: Optional limit for lat- and longitudinal deltas 

478 (C{degrees}) or C{None} or C{0} for unlimited. 

479 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll B{C{lat2}} 

480 and B{C{lon2}} (C{bool}). 

481 

482 @return: A L{Distance4Tuple}C{(distance2, delta_lat, delta_lon, 

483 unroll_lon2)} in C{degrees squared}. 

484 

485 @raise LimitError: If the lat- and/or longitudinal delta exceeds the 

486 B{C{-limit..limit}} range and L{pygeodesy.limiterrors} 

487 set to C{True}. 

488 

489 @see: U{Local, flat earth approximation 

490 <https://www.EdWilliams.org/avform.htm#flat>}, functions 

491 L{equirectangular}, L{cosineAndoyerLambert}, 

492 L{cosineForsytheAndoyerLambert}, L{cosineLaw}, L{euclidean}, 

493 L{flatLocal}/L{hubeny}, L{flatPolar}, L{haversine}, L{thomas} 

494 and L{vincentys} and methods L{Ellipsoid.distance2}, 

495 C{LatLon.distanceTo*} and C{LatLon.equirectangularTo}. 

496 ''' 

497 d_lon, lat2, ulon2 = _Wrap.latlon3(lon1, lat2, lon2, wrap) 

498 d_lat = lat2 - lat1 

499 

500 if limit and limit > 0 and limiterrors(): 

501 d = max(fabs(d_lat), fabs(d_lon)) 

502 if d > limit: 

503 t = _SPACE_(_delta_, Fmt.PAREN_g(d), Fmt.exceeds_limit(limit)) 

504 s = unstr(equirectangular_, lat1, lon1, lat2, lon2, 

505 limit=limit, wrap=wrap) 

506 raise LimitError(s, txt=t) 

507 

508 if adjust: # scale delta lon 

509 d_lon *= _scale_deg(lat1, lat2) 

510 

511 d2 = hypot2(d_lat, d_lon) # degrees squared! 

512 return Distance4Tuple(d2, d_lat, d_lon, ulon2 - lon2) 

513 

514 

515def euclidean(lat1, lon1, lat2, lon2, radius=R_M, adjust=True, wrap=False): 

516 '''Approximate the C{Euclidean} distance between two (spherical) points. 

517 

518 @arg lat1: Start latitude (C{degrees}). 

519 @arg lon1: Start longitude (C{degrees}). 

520 @arg lat2: End latitude (C{degrees}). 

521 @arg lon2: End longitude (C{degrees}). 

522 @kwarg radius: Mean earth radius (C{meter}), datum (L{Datum}) 

523 or ellipsoid (L{Ellipsoid}, L{Ellipsoid2} or 

524 L{a_f2Tuple}) to use. 

525 @kwarg adjust: Adjust the longitudinal delta by the cosine of 

526 the mean latitude (C{bool}). 

527 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll B{C{lat2}} 

528 and B{C{lon2}} (C{bool}). 

529 

530 @return: Distance (C{meter}, same units as B{C{radius}} or the 

531 ellipsoid or datum axes). 

532 

533 @raise TypeError: Invalid B{C{radius}}. 

534 

535 @see: U{Distance between two (spherical) points 

536 <https://www.EdWilliams.org/avform.htm#Dist>}, functions L{euclid}, 

537 L{euclidean_}, L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert}, 

538 L{cosineLaw}, L{equirectangular}, L{flatLocal}/L{hubeny}, L{flatPolar}, 

539 L{haversine}, L{thomas} and L{vincentys} and methods L{Ellipsoid.distance2}, 

540 C{LatLon.distanceTo*} and C{LatLon.equirectangularTo}. 

541 ''' 

542 return _dS(euclidean_, radius, wrap, lat1, lon1, lat2, lon2, adjust=adjust) 

543 

544 

545def euclidean_(phi2, phi1, lam21, adjust=True): 

546 '''Approximate the I{angular} C{Euclidean} distance between two 

547 (spherical) points. 

548 

549 @arg phi2: End latitude (C{radians}). 

550 @arg phi1: Start latitude (C{radians}). 

551 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 

552 @kwarg adjust: Adjust the longitudinal delta by the cosine 

553 of the mean latitude (C{bool}). 

554 

555 @return: Angular distance (C{radians}). 

556 

557 @see: Functions L{euclid}, L{euclidean}, L{cosineAndoyerLambert_}, 

558 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, L{equirectangular_}, 

559 L{flatLocal_}/L{hubeny_}, L{flatPolar_}, L{haversine_}, L{thomas_} 

560 and L{vincentys_}. 

561 ''' 

562 if adjust: 

563 lam21 *= _scale_rad(phi2, phi1) 

564 return euclid(phi2 - phi1, lam21) 

565 

566 

567def excessAbc_(A, b, c): 

568 '''Compute the I{spherical excess} C{E} of a (spherical) triangle 

569 from two sides and the included angle. 

570 

571 @arg A: An interior triangle angle (C{radians}). 

572 @arg b: Frist adjacent triangle side (C{radians}). 

573 @arg c: Second adjacent triangle side (C{radians}). 

574 

575 @return: Spherical excess (C{radians}). 

576 

577 @raise UnitError: Invalid B{C{A}}, B{C{b}} or B{C{c}}. 

578 

579 @see: Functions L{excessGirard_}, L{excessLHuilier_} and U{Spherical 

580 trigonometry<https://WikiPedia.org/wiki/Spherical_trigonometry>}. 

581 ''' 

582 sA, cA, sb, cb, sc, cc = sincos2_(Radians_(A=A), Radians_(b=b) * _0_5, 

583 Radians_(c=c) * _0_5) 

584 return atan2(sA * sb * sc, cb * cc + cA * sb * sc) * _2_0 

585 

586 

587def excessGirard_(A, B, C): 

588 '''Compute the I{spherical excess} C{E} of a (spherical) triangle using 

589 U{Girard's<https://MathWorld.Wolfram.com/GirardsSphericalExcessFormula.html>} 

590 formula. 

591 

592 @arg A: First interior triangle angle (C{radians}). 

593 @arg B: Second interior triangle angle (C{radians}). 

594 @arg C: Third interior triangle angle (C{radians}). 

595 

596 @return: Spherical excess (C{radians}). 

597 

598 @raise UnitError: Invalid B{C{A}}, B{C{B}} or B{C{C}}. 

599 

600 @see: Function L{excessLHuilier_} and U{Spherical trigonometry 

601 <https://WikiPedia.org/wiki/Spherical_trigonometry>}. 

602 ''' 

603 return Radians(Girard=fsumf_(Radians_(A=A), 

604 Radians_(B=B), 

605 Radians_(C=C), -PI)) 

606 

607 

608def excessLHuilier_(a, b, c): 

609 '''Compute the I{spherical excess} C{E} of a (spherical) triangle using 

610 U{L'Huilier's<https://MathWorld.Wolfram.com/LHuiliersTheorem.html>} 

611 Theorem. 

612 

613 @arg a: First triangle side (C{radians}). 

614 @arg b: Second triangle side (C{radians}). 

615 @arg c: Third triangle side (C{radians}). 

616 

617 @return: Spherical excess (C{radians}). 

618 

619 @raise UnitError: Invalid B{C{a}}, B{C{b}} or B{C{c}}. 

620 

621 @see: Function L{excessGirard_} and U{Spherical trigonometry 

622 <https://WikiPedia.org/wiki/Spherical_trigonometry>}. 

623 ''' 

624 a = Radians_(a=a) 

625 b = Radians_(b=b) 

626 c = Radians_(c=c) 

627 

628 s = fsumf_(a, b, c) * _0_5 

629 r = tan_2(s) * tan_2(s - a) * tan_2(s - b) * tan_2(s - c) 

630 r = atan(sqrt(r)) if r > 0 else _0_0 

631 return Radians(LHuilier=r * _4_0) 

632 

633 

634def excessKarney(lat1, lon1, lat2, lon2, radius=R_M, wrap=False): 

635 '''Compute the surface area of a (spherical) quadrilateral bounded by a 

636 segment of a great circle, two meridians and the equator using U{Karney's 

637 <https://MathOverflow.net/questions/97711/the-area-of-spherical-polygons>} 

638 method. 

639 

640 @arg lat1: Start latitude (C{degrees}). 

641 @arg lon1: Start longitude (C{degrees}). 

642 @arg lat2: End latitude (C{degrees}). 

643 @arg lon2: End longitude (C{degrees}). 

644 @kwarg radius: Mean earth radius (C{meter}), datum (L{Datum}) 

645 or ellipsoid (L{Ellipsoid}, L{Ellipsoid2} or 

646 L{a_f2Tuple}) or C{None}. 

647 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll 

648 B{C{lat2}} and B{C{lon2}} (C{bool}). 

649 

650 @return: Surface area, I{signed} (I{square} C{meter} or the same units as 

651 B{C{radius}} I{squared}) or the I{spherical excess} (C{radians}) 

652 if C{B{radius}=0} or C{None}. 

653 

654 @raise TypeError: Invalid B{C{radius}}. 

655 

656 @raise UnitError: Invalid B{C{lat2}} or B{C{lat1}}. 

657 

658 @raise ValueError: Semi-circular longitudinal delta. 

659 

660 @see: Functions L{excessKarney_} and L{excessQuad}. 

661 ''' 

662 return _eA(excessKarney_, radius, wrap, lat1, lon1, lat2, lon2) 

663 

664 

665def excessKarney_(phi2, phi1, lam21): 

666 '''Compute the I{spherical excess} C{E} of a (spherical) quadrilateral bounded 

667 by a segment of a great circle, two meridians and the equator using U{Karney's 

668 <https://MathOverflow.net/questions/97711/the-area-of-spherical-polygons>} 

669 method. 

670 

671 @arg phi2: End latitude (C{radians}). 

672 @arg phi1: Start latitude (C{radians}). 

673 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 

674 

675 @return: Spherical excess, I{signed} (C{radians}). 

676 

677 @raise ValueError: Semi-circular longitudinal delta B{C{lam21}}. 

678 

679 @see: Function L{excessKarney} and U{Area of a spherical polygon 

680 <https://MathOverflow.net/questions/97711/the-area-of-spherical-polygons>}. 

681 ''' 

682 # from: Veness <https://www.Movable-Type.co.UK/scripts/latlong.html> Area 

683 # method due to Karney: for each edge of the polygon, 

684 # 

685 # tan(Δλ / 2) · (tan(φ1 / 2) + tan(φ2 / 2)) 

686 # tan(E / 2) = ----------------------------------------- 

687 # 1 + tan(φ1 / 2) · tan(φ2 / 2) 

688 # 

689 # where E is the spherical excess of the trapezium obtained by extending 

690 # the edge to the equator-circle vector for each edge (see also ***). 

691 t2 = tan_2(phi2) 

692 t1 = tan_2(phi1) 

693 t = tan_2(lam21, lam21=None) 

694 return Radians(Karney=atan2(t * (t1 + t2), 

695 _1_0 + (t1 * t2)) * _2_0) 

696 

697 

698# ***) Original post no longer available, following is a copy of the main part 

699# <http://OSGeo-org.1560.x6.Nabble.com/Area-of-a-spherical-polygon-td3841625.html> 

700# 

701# The area of a polygon on a (unit) sphere is given by the spherical excess 

702# 

703# A = 2 * pi - sum(exterior angles) 

704# 

705# However this is badly conditioned if the polygon is small. In this case, use 

706# 

707# A = sum(S12{i, i+1}) over the edges of the polygon 

708# 

709# where S12 is the area of the quadrilateral bounded by an edge of the polygon, 

710# two meridians and the equator, i.e. with vertices (phi1, lambda1), (phi2, 

711# lambda2), (0, lambda1) and (0, lambda2). S12 is given by 

712# 

713# tan(S12 / 2) = tan(lambda21 / 2) * (tan(phi1 / 2) + tan(phi2 / 2)) / 

714# (tan(phi1 / 2) * tan(phi2 / 2) + 1) 

715# 

716# = tan(lambda21 / 2) * tanh((Lambertian(phi1) + 

717# Lambertian(phi2)) / 2) 

718# 

719# where lambda21 = lambda2 - lambda1 and lamb(x) is the Lambertian (or 

720# inverse Gudermannian) function 

721# 

722# Lambertian(x) = asinh(tan(x)) = atanh(sin(x)) = 2 * atanh(tan(x / 2)) 

723# 

724# Notes: The formula for S12 is exact, except that... 

725# - it is indeterminate if an edge is a semi-circle 

726# - the formula for A applies only if the polygon does not include a pole 

727# (if it does, then add +/- 2 * pi to the result) 

728# - in the limit of small phi and lambda, S12 reduces to the trapezoidal 

729# formula, S12 = (lambda2 - lambda1) * (phi1 + phi2) / 2 

730# - I derived this result from the equation for the area of a spherical 

731# triangle in terms of two edges and the included angle given by, e.g. 

732# U{Todhunter, I. - Spherical Trigonometry (1871), Sec. 103, Eq. (2) 

733# <http://Books.Google.com/books?id=3uBHAAAAIAAJ&pg=PA71>} 

734# - I would be interested to know if this formula for S12 is already known 

735# - Charles Karney 

736 

737 

738def excessQuad(lat1, lon1, lat2, lon2, radius=R_M, wrap=False): 

739 '''Compute the surface area of a (spherical) quadrilateral bounded by a segment 

740 of a great circle, two meridians and the equator. 

741 

742 @arg lat1: Start latitude (C{degrees}). 

743 @arg lon1: Start longitude (C{degrees}). 

744 @arg lat2: End latitude (C{degrees}). 

745 @arg lon2: End longitude (C{degrees}). 

746 @kwarg radius: Mean earth radius (C{meter}), datum (L{Datum}) 

747 or ellipsoid (L{Ellipsoid}, L{Ellipsoid2} or 

748 L{a_f2Tuple}) or C{None}. 

749 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll 

750 B{C{lat2}} and B{C{lon2}} (C{bool}). 

751 

752 @return: Surface area, I{signed} (I{square} C{meter} or the same units as 

753 B{C{radius}} I{squared}) or the I{spherical excess} (C{radians}) 

754 if C{B{radius}=0} or C{None}. 

755 

756 @raise TypeError: Invalid B{C{radius}}. 

757 

758 @raise UnitError: Invalid B{C{lat2}} or B{C{lat1}}. 

759 

760 @see: Function L{excessQuad_} and L{excessKarney}. 

761 ''' 

762 return _eA(excessQuad_, radius, wrap, lat1, lon1, lat2, lon2) 

763 

764 

765def excessQuad_(phi2, phi1, lam21): 

766 '''Compute the I{spherical excess} C{E} of a (spherical) quadrilateral bounded 

767 by a segment of a great circle, two meridians and the equator. 

768 

769 @arg phi2: End latitude (C{radians}). 

770 @arg phi1: Start latitude (C{radians}). 

771 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 

772 

773 @return: Spherical excess, I{signed} (C{radians}). 

774 

775 @see: Function L{excessQuad}, U{Spherical trigonometry 

776 <https://WikiPedia.org/wiki/Spherical_trigonometry>}. 

777 ''' 

778 s = sin((phi2 + phi1) * _0_5) 

779 c = cos((phi2 - phi1) * _0_5) 

780 return Radians(Quad=atan2(tan_2(lam21) * s, c) * _2_0) 

781 

782 

783def flatLocal(lat1, lon1, lat2, lon2, datum=_WGS84, scaled=True, wrap=False): 

784 '''Compute the distance between two (ellipsoidal) points using 

785 the U{ellipsoidal Earth to plane projection<https://WikiPedia.org/ 

786 wiki/Geographical_distance#Ellipsoidal_Earth_projected_to_a_plane>} 

787 aka U{Hubeny<https://www.OVG.AT/de/vgi/files/pdf/3781/>} formula. 

788 

789 @arg lat1: Start latitude (C{degrees}). 

790 @arg lon1: Start longitude (C{degrees}). 

791 @arg lat2: End latitude (C{degrees}). 

792 @arg lon2: End longitude (C{degrees}). 

793 @kwarg datum: Datum (L{Datum}) or ellipsoid (L{Ellipsoid}, 

794 L{Ellipsoid2} or L{a_f2Tuple}) to use. 

795 @kwarg scaled: Scale prime_vertical by C{cos(B{phi})} (C{bool}), 

796 see method L{pygeodesy.Ellipsoid.roc2_}. 

797 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll 

798 B{C{lat2}} and B{C{lon2}} (C{bool}). 

799 

800 @return: Distance (C{meter}, same units as the B{C{datum}}'s 

801 ellipsoid axes). 

802 

803 @raise TypeError: Invalid B{C{datum}}. 

804 

805 @note: The meridional and prime_vertical radii of curvature 

806 are taken and scaled at the mean of both latitude. 

807 

808 @see: Functions L{flatLocal_} or L{hubeny_}, L{cosineLaw}, L{flatPolar}, 

809 L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert}, 

810 L{equirectangular}, L{euclidean}, L{haversine}, L{thomas}, 

811 L{vincentys}, method L{Ellipsoid.distance2} and U{local, flat 

812 earth approximation<https://www.EdWilliams.org/avform.htm#flat>}. 

813 ''' 

814 E = _ellipsoidal(datum, flatLocal) 

815 return E._hubeny_2(*_d3(wrap, lat1, lon1, lat2, lon2), 

816 scaled=scaled, squared=False) * E.a 

817 

818hubeny = flatLocal # PYCHOK for Karl Hubeny 

819 

820 

821def flatLocal_(phi2, phi1, lam21, datum=_WGS84, scaled=True): 

822 '''Compute the I{angular} distance between two (ellipsoidal) points using 

823 the U{ellipsoidal Earth to plane projection<https://WikiPedia.org/ 

824 wiki/Geographical_distance#Ellipsoidal_Earth_projected_to_a_plane>} 

825 aka U{Hubeny<https://www.OVG.AT/de/vgi/files/pdf/3781/>} formula. 

826 

827 @arg phi2: End latitude (C{radians}). 

828 @arg phi1: Start latitude (C{radians}). 

829 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 

830 @kwarg datum: Datum (L{Datum}) or ellipsoid (L{Ellipsoid}, 

831 L{Ellipsoid2} or L{a_f2Tuple}) to use. 

832 @kwarg scaled: Scale prime_vertical by C{cos(B{phi})} (C{bool}), 

833 see method L{pygeodesy.Ellipsoid.roc2_}. 

834 

835 @return: Angular distance (C{radians}). 

836 

837 @raise TypeError: Invalid B{C{datum}}. 

838 

839 @note: The meridional and prime_vertical radii of curvature 

840 are taken and scaled I{at the mean of both latitude}. 

841 

842 @see: Functions L{flatLocal} or L{hubeny}, L{cosineAndoyerLambert_}, 

843 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, L{flatPolar_}, 

844 L{equirectangular_}, L{euclidean_}, L{haversine_}, L{thomas_} 

845 and L{vincentys_} and U{local, flat earth approximation 

846 <https://www.EdWilliams.org/avform.htm#flat>}. 

847 ''' 

848 E = _ellipsoidal(datum, flatLocal_) 

849 return E._hubeny_2(phi2, phi1, lam21, scaled=scaled, squared=False) 

850 

851hubeny_ = flatLocal_ # PYCHOK for Karl Hubeny 

852 

853 

854def flatPolar(lat1, lon1, lat2, lon2, radius=R_M, wrap=False): 

855 '''Compute the distance between two (spherical) points using 

856 the U{polar coordinate flat-Earth <https://WikiPedia.org/wiki/ 

857 Geographical_distance#Polar_coordinate_flat-Earth_formula>} 

858 formula. 

859 

860 @arg lat1: Start latitude (C{degrees}). 

861 @arg lon1: Start longitude (C{degrees}). 

862 @arg lat2: End latitude (C{degrees}). 

863 @arg lon2: End longitude (C{degrees}). 

864 @kwarg radius: Mean earth radius (C{meter}), datum (L{Datum}) 

865 or ellipsoid (L{Ellipsoid}, L{Ellipsoid2} or 

866 L{a_f2Tuple}) to use. 

867 @kwarg wrap: If C{True}, wrap or I{normalize} and B{C{lat2}} 

868 and B{C{lon2}} (C{bool}). 

869 

870 @return: Distance (C{meter}, same units as B{C{radius}} or the 

871 ellipsoid or datum axes). 

872 

873 @raise TypeError: Invalid B{C{radius}}. 

874 

875 @see: Functions L{flatPolar_}, L{cosineAndoyerLambert}, 

876 L{cosineForsytheAndoyerLambert},L{cosineLaw}, 

877 L{flatLocal}/L{hubeny}, L{equirectangular}, 

878 L{euclidean}, L{haversine}, L{thomas} and 

879 L{vincentys}. 

880 ''' 

881 return _dS(flatPolar_, radius, wrap, lat1, lon1, lat2, lon2) 

882 

883 

884def flatPolar_(phi2, phi1, lam21): 

885 '''Compute the I{angular} distance between two (spherical) points 

886 using the U{polar coordinate flat-Earth<https://WikiPedia.org/wiki/ 

887 Geographical_distance#Polar_coordinate_flat-Earth_formula>} 

888 formula. 

889 

890 @arg phi2: End latitude (C{radians}). 

891 @arg phi1: Start latitude (C{radians}). 

892 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 

893 

894 @return: Angular distance (C{radians}). 

895 

896 @see: Functions L{flatPolar}, L{cosineAndoyerLambert_}, 

897 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, 

898 L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_}, 

899 L{haversine_}, L{thomas_} and L{vincentys_}. 

900 ''' 

901 a = fabs(PI_2 - phi1) # co-latitude 

902 b = fabs(PI_2 - phi2) # co-latitude 

903 if a < b: 

904 a, b = b, a 

905 if a < EPS0: 

906 a = _0_0 

907 elif b > 0: 

908 b = b / a # /= chokes PyChecker 

909 c = b * cos(lam21) * _2_0 

910 c = fsumf_(_1_0, b**2, -fabs(c)) 

911 a *= sqrt0(c) 

912 return a 

913 

914 

915def hartzell(pov, los=None, earth=_WGS84, name=NN, **LatLon_and_kwds): 

916 '''Compute the intersection of the earth's surface and a Line-Of-Sight 

917 from a Point-Of-View in space. 

918 

919 @arg pov: Point-Of-View outside the earth (C{Cartesian}, L{Ecef9Tuple} 

920 or L{Vector3d}). 

921 @kwarg los: Line-Of-Sight, I{direction} to earth (L{Vector3d}) or 

922 C{None} to point to the earth' center. 

923 @kwarg earth: The earth model (L{Datum}, L{Ellipsoid}, L{Ellipsoid2}, 

924 L{a_f2Tuple} or C{scalar} radius in C{meter}). 

925 @kwarg name: Optional name (C{str}). 

926 @kwarg LatLon_and_kwds: Optional C{LatLon} class for the intersection 

927 point plus C{LatLon} keyword arguments, include 

928 B{C{datum}} if different from B{C{earth}}. 

929 

930 @return: The earth intersection (L{Vector3d}, C{Cartesian type} of 

931 B{C{pov}} or B{C{LatLon}}). 

932 

933 @raise IntersectionError: Null B{C{pov}} or B{C{los}} vector, B{C{pov}} 

934 is inside the earth or B{C{los}} points outside 

935 the earth or points in an opposite direction. 

936 

937 @raise TypeError: Invalid B{C{pov}}, B{C{los}} or B{C{earth}}. 

938 

939 @see: Function L{pygeodesy.hartzell4}, L{pygeodesy.tyr3d} for B{C{los}}, 

940 method L{Ellipsoid.hartzell4} and U{I{Satellite Line-of-Sight 

941 Intersection with Earth}<https://StephenHartzell.Medium.com/ 

942 satellite-line-of-sight-intersection-with-earth-d786b4a6a9b6>}. 

943 ''' 

944 D = earth if isinstance(earth, Datum) else \ 

945 _spherical_datum(earth, name=hartzell.__name__) 

946 try: 

947 r, _ = _MODS.triaxials._hartzell3d2(pov, los, D.ellipsoid._triaxial) 

948 except Exception as x: 

949 raise IntersectionError(pov=pov, los=los, earth=earth, cause=x) 

950 

951# else: 

952# E = D.ellipsoid 

953# # Triaxial(a, b, c) == (E.a, E.a, E.b) 

954# 

955# def _Error(txt): 

956# return IntersectionError(pov=pov, los=los, earth=earth, txt=txt) 

957# 

958# a2 = b2 = E.a2 # earth' x, y, ... 

959# c2 = E.b2 # ... z semi-axis squared 

960# q2 = E.b2_a2 # == c2 / a2 

961# bc = E.a * E.b # == b * c 

962# 

963# V3 = _MODS.vector3d._otherV3d 

964# p3 = V3(pov=pov) 

965# u3 = V3(los=los) if los else p3.negate() 

966# u3 = u3.unit() # unit vector, opposing signs 

967# 

968# x2, y2, z2 = p3.x2y2z2 # p3.times_(p3).xyz 

969# ux, vy, wz = u3.times_(p3).xyz 

970# u2, v2, w2 = u3.x2y2z2 # u3.times_(u3).xyz 

971# 

972# t = c2, c2, b2 

973# m = fdot(t, u2, v2, w2) # a2 factored out 

974# if m < EPS0: # zero or near-null LOS vector 

975# raise _Error(_near_(_null_)) 

976# 

977# # a2 and b2 factored out, b2 == a2 and b2 / a2 == 1 

978# r = fsumf_(b2 * w2, c2 * v2, -v2 * z2, vy * wz * 2, 

979# c2 * u2, -u2 * z2, -w2 * x2, ux * wz * 2, 

980# -w2 * y2, -u2 * y2 * q2, -v2 * x2 * q2, ux * vy * 2 * q2) 

981# if r > 0: 

982# r = sqrt(r) * bc # == a * a * b * c / a2 

983# elif r < 0: # LOS pointing away from or missing the earth 

984# raise _Error(_opposite_ if max(ux, vy, wz) > 0 else _outside_) 

985# 

986# d = Fdot(t, ux, vy, wz).fadd_(r).fover(m) # -r for antipode, a2 factored out 

987# if d > 0: # POV inside or LOS missing, outside the earth 

988# s = fsumf_(_1_0, x2 / a2, y2 / b2, z2 / c2, _N_2_0) # like _sideOf 

989# raise _Error(_outside_ if s > 0 else _inside_) 

990# elif fsumf_(x2, y2, z2) < d**2: # d past earth center 

991# raise _Error(_too_(_distant_)) 

992# 

993# r = p3.minus(u3.times(d)) 

994# # h = p3.minus(r).length # distance to ellipsoid 

995 

996 r = _xnamed(r, name or hartzell.__name__) 

997 if LatLon_and_kwds: 

998 c = _MODS.cartesianBase.CartesianBase(r, datum=D, name=r.name) 

999 r = c.toLatLon(**LatLon_and_kwds) 

1000 return r 

1001 

1002 

1003def haversine(lat1, lon1, lat2, lon2, radius=R_M, wrap=False): 

1004 '''Compute the distance between two (spherical) points using the 

1005 U{Haversine<https://www.Movable-Type.co.UK/scripts/latlong.html>} 

1006 formula. 

1007 

1008 @arg lat1: Start latitude (C{degrees}). 

1009 @arg lon1: Start longitude (C{degrees}). 

1010 @arg lat2: End latitude (C{degrees}). 

1011 @arg lon2: End longitude (C{degrees}). 

1012 @kwarg radius: Mean earth radius (C{meter}), datum (L{Datum}) 

1013 or ellipsoid (L{Ellipsoid}, L{Ellipsoid2} or 

1014 L{a_f2Tuple}) to use. 

1015 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll 

1016 B{C{lat2}} and B{C{lon2}} (C{bool}). 

1017 

1018 @return: Distance (C{meter}, same units as B{C{radius}}). 

1019 

1020 @raise TypeError: Invalid B{C{radius}}. 

1021 

1022 @see: U{Distance between two (spherical) points 

1023 <https://www.EdWilliams.org/avform.htm#Dist>}, functions 

1024 L{cosineLaw}, L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert}, 

1025 L{equirectangular}, L{euclidean}, L{flatLocal}/L{hubeny}, L{flatPolar}, 

1026 L{thomas} and L{vincentys} and methods L{Ellipsoid.distance2}, 

1027 C{LatLon.distanceTo*} and C{LatLon.equirectangularTo}. 

1028 

1029 @note: See note at function L{vincentys_}. 

1030 ''' 

1031 return _dS(haversine_, radius, wrap, lat1, lon1, lat2, lon2) 

1032 

1033 

1034def haversine_(phi2, phi1, lam21): 

1035 '''Compute the I{angular} distance between two (spherical) points 

1036 using the U{Haversine<https://www.Movable-Type.co.UK/scripts/latlong.html>} 

1037 formula. 

1038 

1039 @arg phi2: End latitude (C{radians}). 

1040 @arg phi1: Start latitude (C{radians}). 

1041 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 

1042 

1043 @return: Angular distance (C{radians}). 

1044 

1045 @see: Functions L{haversine}, L{cosineAndoyerLambert_}, 

1046 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, 

1047 L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_}, 

1048 L{flatPolar_}, L{thomas_} and L{vincentys_}. 

1049 

1050 @note: See note at function L{vincentys_}. 

1051 ''' 

1052 def _hsin(rad): 

1053 return sin(rad * _0_5)**2 

1054 

1055 h = _hsin(phi2 - phi1) + cos(phi1) * cos(phi2) * _hsin(lam21) # haversine 

1056 return atan2(sqrt0(h), sqrt0(_1_0 - h)) * _2_0 # == asin(sqrt(h)) * 2 

1057 

1058 

1059def heightOf(angle, distance, radius=R_M): 

1060 '''Determine the height above the (spherical) earth' surface after 

1061 traveling along a straight line at a given tilt. 

1062 

1063 @arg angle: Tilt angle above horizontal (C{degrees}). 

1064 @arg distance: Distance along the line (C{meter} or same units as 

1065 B{C{radius}}). 

1066 @kwarg radius: Optional mean earth radius (C{meter}). 

1067 

1068 @return: Height (C{meter}, same units as B{C{distance}} and B{C{radius}}). 

1069 

1070 @raise ValueError: Invalid B{C{angle}}, B{C{distance}} or B{C{radius}}. 

1071 

1072 @see: U{MultiDop geog_lib.GeogBeamHt<https://GitHub.com/NASA/MultiDop>} 

1073 (U{Shapiro et al. 2009, JTECH 

1074 <https://Journals.AMetSoc.org/doi/abs/10.1175/2009JTECHA1256.1>} 

1075 and U{Potvin et al. 2012, JTECH 

1076 <https://Journals.AMetSoc.org/doi/abs/10.1175/JTECH-D-11-00019.1>}). 

1077 ''' 

1078 r = h = Radius(radius) 

1079 d = fabs(Distance(distance)) 

1080 if d > h: 

1081 d, h = h, d 

1082 

1083 if d > EPS0: # and h > EPS0 

1084 d = d / h # /= h chokes PyChecker 

1085 s = sin(Phi_(angle=angle, clip=_180_0)) 

1086 s = fsumf_(_1_0, _2_0 * s * d, d**2) 

1087 if s > 0: 

1088 return h * sqrt(s) - r 

1089 

1090 raise _ValueError(angle=angle, distance=distance, radius=radius) 

1091 

1092 

1093def horizon(height, radius=R_M, refraction=False): 

1094 '''Determine the distance to the horizon from a given altitude 

1095 above the (spherical) earth. 

1096 

1097 @arg height: Altitude (C{meter} or same units as B{C{radius}}). 

1098 @kwarg radius: Optional mean earth radius (C{meter}). 

1099 @kwarg refraction: Consider atmospheric refraction (C{bool}). 

1100 

1101 @return: Distance (C{meter}, same units as B{C{height}} and B{C{radius}}). 

1102 

1103 @raise ValueError: Invalid B{C{height}} or B{C{radius}}. 

1104 

1105 @see: U{Distance to horizon<https://www.EdWilliams.org/avform.htm#Horizon>}. 

1106 ''' 

1107 h, r = Height(height), Radius(radius) 

1108 if min(h, r) < 0: 

1109 raise _ValueError(height=height, radius=radius) 

1110 

1111 if refraction: 

1112 d2 = 2.415750694528 * h * r # 2.0 / 0.8279 

1113 else: 

1114 d2 = h * fsumf_(r, r, h) 

1115 return sqrt0(d2) 

1116 

1117 

1118class _idllmn6(object): # see also .geodesicw._wargs, .vector2d._numpy 

1119 '''(INTERNAL) Helper for C{intersection2} and C{intersections2}. 

1120 ''' 

1121 @contextmanager # <https://www.python.org/dev/peps/pep-0343/> Examples 

1122 def __call__(self, datum, lat1, lon1, lat2, lon2, small, wrap, s, **kwds): 

1123 try: 

1124 if wrap: 

1125 _, lat2, lon2 = _Wrap.latlon3(lon1, lat2, lon2, wrap) 

1126 kwds = _xkwds(kwds, wrap=wrap) # for _xError 

1127 m = small if small is _100km else Meter_(small=small) 

1128 n = (intersections2 if s else intersection2).__name__ 

1129 if datum is None or euclidean(lat1, lon1, lat2, lon2) < m: 

1130 d, m = None, _MODS.vector3d 

1131 _i = m._intersects2 if s else m._intersect3d3 

1132 elif isscalar(datum) and datum < 0 and not s: 

1133 d = _spherical_datum(-datum, name=n) 

1134 m = _MODS.sphericalNvector 

1135 _i = m.intersection 

1136 else: 

1137 d = _spherical_datum(datum, name=n) 

1138 if d.isSpherical: 

1139 m = _MODS.sphericalTrigonometry 

1140 _i = m._intersects2 if s else m._intersect 

1141 elif d.isEllipsoidal: 

1142 try: 

1143 if d.ellipsoid.geodesic: 

1144 pass 

1145 m = _MODS.ellipsoidalKarney 

1146 except ImportError: 

1147 m = _MODS.ellipsoidalExact 

1148 _i = m._intersections2 if s else m._intersection3 # ellispoidalBaseDI 

1149 else: 

1150 raise _TypeError(datum=datum) 

1151 yield _i, d, lat2, lon2, m, n 

1152 

1153 except (TypeError, ValueError) as x: 

1154 raise _xError(x, lat1=lat1, lon1=lon1, datum=datum, 

1155 lat2=lat2, lon2=lon2, small=small, **kwds) 

1156 

1157_idllmn6 = _idllmn6() # PYCHOK singleton 

1158 

1159 

1160def intersection2(lat1, lon1, bearing1, 

1161 lat2, lon2, bearing2, datum=None, wrap=False, small=_100km): # was=True 

1162 '''I{Conveniently} compute the intersection of two lines each defined 

1163 by a (geodetic) point and a bearing from North, using either ... 

1164 

1165 1) L{vector3d.intersection3d3} for B{C{small}} distances (below 100 Km 

1166 or about 0.88 degrees) or if I{no} B{C{datum}} is specified, or ... 

1167 

1168 2) L{sphericalTrigonometry.intersection} for a spherical B{C{datum}} 

1169 or a C{scalar B{datum}} representing the earth radius, conventionally 

1170 in C{meter} or ... 

1171 

1172 3) L{sphericalNvector.intersection} if B{C{datum}} is a I{negative} 

1173 C{scalar}, (negative) earth radius, conventionally in C{meter} or ... 

1174 

1175 4) L{ellipsoidalKarney.intersection3} for an ellipsoidal B{C{datum}} 

1176 and if I{Karney}'s U{geographiclib<https://PyPI.org/project/geographiclib>} 

1177 is installed, otherwise ... 

1178 

1179 5) L{ellipsoidalExact.intersection3}, provided B{C{datum}} is ellipsoidal. 

1180 

1181 @arg lat1: Latitude of the first point (C{degrees}). 

1182 @arg lon1: Longitude of the first point (C{degrees}). 

1183 @arg bearing1: Bearing at the first point (compass C{degrees}). 

1184 @arg lat2: Latitude of the second point (C{degrees}). 

1185 @arg lon2: Longitude of the second point (C{degrees}). 

1186 @arg bearing2: Bearing at the second point (compass C{degrees}). 

1187 @kwarg datum: Optional datum (L{Datum}) or ellipsoid (L{Ellipsoid}, 

1188 L{Ellipsoid2} or L{a_f2Tuple}) or C{scalar} earth 

1189 radius (C{meter}, same units as B{C{radius1}} and 

1190 B{C{radius2}}) or C{None}. 

1191 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll B{C{lat2}} 

1192 and B{C{lon2}} (C{bool}). 

1193 @kwarg small: Upper limit for small distances (C{meter}). 

1194 

1195 @return: A L{LatLon2Tuple}C{(lat, lon)} with the lat- and 

1196 longitude of the intersection point. 

1197 

1198 @raise IntersectionError: Ambiguous or infinite intersection 

1199 or colinear, parallel or otherwise 

1200 non-intersecting lines. 

1201 

1202 @raise TypeError: Invalid B{C{datum}}. 

1203 

1204 @raise UnitError: Invalid B{C{lat1}}, B{C{lon1}}, B{C{bearing1}}, 

1205 B{C{lat2}}, B{C{lon2}} or B{C{bearing2}}. 

1206 

1207 @see: Method L{RhumbLine.intersection2}. 

1208 

1209 @note: The returned intersections may be near-antipodal. 

1210 ''' 

1211 b1 = Bearing(bearing1=bearing1) 

1212 b2 = Bearing(bearing2=bearing2) 

1213 with _idllmn6(datum, lat1, lon1, lat2, lon2, 

1214 small, wrap, False, bearing1=b1, bearing2=b2) as t: 

1215 _i, d, lat2, lon2, m, n = t 

1216 if d is None: 

1217 t, _, _ = _i(m.Vector3d(lon1, lat1, 0), b1, 

1218 m.Vector3d(lon2, lat2, 0), b2, useZ=False) 

1219 t = LatLon2Tuple(t.y, t.x, name=n) 

1220 

1221 else: 

1222 t = _i(m.LatLon(lat1, lon1, datum=d), b1, 

1223 m.LatLon(lat2, lon2, datum=d), b2, height=0, wrap=False) 

1224 if isinstance(t, Intersection3Tuple): # ellipsoidal 

1225 t, _, _ = t 

1226 t = LatLon2Tuple(t.lat, t.lon, name=n) 

1227 return t 

1228 

1229 

1230def intersections2(lat1, lon1, radius1, 

1231 lat2, lon2, radius2, datum=None, wrap=False, small=_100km): # was=True 

1232 '''I{Conveniently} compute the intersections of two circles each defined 

1233 by a (geodetic) center point and a radius, using either ... 

1234 

1235 1) L{vector3d.intersections2} for B{C{small}} distances (below 100 Km 

1236 or about 0.88 degrees) or if I{no} B{C{datum}} is specified, or ... 

1237 

1238 2) L{sphericalTrigonometry.intersections2} for a spherical B{C{datum}} 

1239 or a C{scalar B{datum}} representing the earth radius, conventionally 

1240 in C{meter} or ... 

1241 

1242 3) L{ellipsoidalKarney.intersections2} for an ellipsoidal B{C{datum}} 

1243 and if I{Karney}'s U{geographiclib<https://PyPI.org/project/geographiclib>} 

1244 is installed, otherwise ... 

1245 

1246 4) L{ellipsoidalExact.intersections2}, provided B{C{datum}} is ellipsoidal. 

1247 

1248 @arg lat1: Latitude of the first circle center (C{degrees}). 

1249 @arg lon1: Longitude of the first circle center (C{degrees}). 

1250 @arg radius1: Radius of the first circle (C{meter}, conventionally). 

1251 @arg lat2: Latitude of the second circle center (C{degrees}). 

1252 @arg lon2: Longitude of the second circle center (C{degrees}). 

1253 @arg radius2: Radius of the second circle (C{meter}, same units as B{C{radius1}}). 

1254 @kwarg datum: Optional datum (L{Datum}) or ellipsoid (L{Ellipsoid}, 

1255 L{Ellipsoid2} or L{a_f2Tuple}) or C{scalar} earth 

1256 radius (C{meter}, same units as B{C{radius1}} and 

1257 B{C{radius2}}) or C{None}. 

1258 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll B{C{lat2}} 

1259 and B{C{lon2}} (C{bool}). 

1260 @kwarg small: Upper limit for small distances (C{meter}). 

1261 

1262 @return: 2-Tuple of the intersection points, each a 

1263 L{LatLon2Tuple}C{(lat, lon)}. For abutting circles, the 

1264 points are the same instance, aka the I{radical center}. 

1265 

1266 @raise IntersectionError: Concentric, antipodal, invalid or 

1267 non-intersecting circles or no 

1268 convergence. 

1269 

1270 @raise TypeError: Invalid B{C{datum}}. 

1271 

1272 @raise UnitError: Invalid B{C{lat1}}, B{C{lon1}}, B{C{radius1}}, 

1273 B{C{lat2}}, B{C{lon2}} or B{C{radius2}}. 

1274 ''' 

1275 r1 = Radius_(radius1=radius1) 

1276 r2 = Radius_(radius2=radius2) 

1277 with _idllmn6(datum, lat1, lon1, lat2, lon2, 

1278 small, wrap, True, radius1=r1, radius2=r2) as t: 

1279 _i, d, lat2, lon2, m, n = t 

1280 if d is None: 

1281 r1 = m2degrees(r1, radius=R_M, lat=lat1) 

1282 r2 = m2degrees(r2, radius=R_M, lat=lat2) 

1283 

1284 def _V2T(x, y, _, **unused): # _ == z unused 

1285 return LatLon2Tuple(y, x, name=n) 

1286 

1287 t = _i(m.Vector3d(lon1, lat1, 0), r1, 

1288 m.Vector3d(lon2, lat2, 0), r2, sphere=False, 

1289 Vector=_V2T) 

1290 else: 

1291 def _LL2T(lat, lon, **unused): 

1292 return LatLon2Tuple(lat, lon, name=n) 

1293 

1294 t = _i(m.LatLon(lat1, lon1, datum=d), r1, 

1295 m.LatLon(lat2, lon2, datum=d), r2, 

1296 LatLon=_LL2T, height=0, wrap=False) 

1297 return t 

1298 

1299 

1300def isantipode(lat1, lon1, lat2, lon2, eps=EPS): 

1301 '''Check whether two points are I{antipodal}, on diametrically 

1302 opposite sides of the earth. 

1303 

1304 @arg lat1: Latitude of one point (C{degrees}). 

1305 @arg lon1: Longitude of one point (C{degrees}). 

1306 @arg lat2: Latitude of the other point (C{degrees}). 

1307 @arg lon2: Longitude of the other point (C{degrees}). 

1308 @kwarg eps: Tolerance for near-equality (C{degrees}). 

1309 

1310 @return: C{True} if points are antipodal within the 

1311 B{C{eps}} tolerance, C{False} otherwise. 

1312 

1313 @see: Functions L{isantipode_} and L{antipode}. 

1314 ''' 

1315 return (fabs(lat1 + lat2) <= eps and 

1316 fabs(lon1 + lon2) <= eps) or _isequalTo( 

1317 normal(lat1, lon1), antipode(lat2, lon2), eps) 

1318 

1319 

1320def isantipode_(phi1, lam1, phi2, lam2, eps=EPS): 

1321 '''Check whether two points are I{antipodal}, on diametrically 

1322 opposite sides of the earth. 

1323 

1324 @arg phi1: Latitude of one point (C{radians}). 

1325 @arg lam1: Longitude of one point (C{radians}). 

1326 @arg phi2: Latitude of the other point (C{radians}). 

1327 @arg lam2: Longitude of the other point (C{radians}). 

1328 @kwarg eps: Tolerance for near-equality (C{radians}). 

1329 

1330 @return: C{True} if points are antipodal within the 

1331 B{C{eps}} tolerance, C{False} otherwise. 

1332 

1333 @see: Functions L{isantipode} and L{antipode_}. 

1334 ''' 

1335 return (fabs(phi1 + phi2) <= eps and 

1336 fabs(lam1 + lam2) <= eps) or _isequalTo_( 

1337 normal_(phi1, lam1), antipode_(phi2, lam2), eps) 

1338 

1339 

1340def _isequalTo(p1, p2, eps=EPS): 

1341 '''Compare 2 point lat-/lons ignoring C{class}. 

1342 ''' 

1343 return (fabs(p1.lat - p2.lat) <= eps and 

1344 fabs(p1.lon - p2.lon) <= eps) if eps else (p1.latlon == p2.latlon) 

1345 

1346 

1347def _isequalTo_(p1, p2, eps=EPS): 

1348 '''(INTERNAL) Compare 2 point phi-/lams ignoring C{class}. 

1349 ''' 

1350 return (fabs(p1.phi - p2.phi) <= eps and 

1351 fabs(p1.lam - p2.lam) <= eps) if eps else (p1.philam == p2.philam) 

1352 

1353 

1354def isnormal(lat, lon, eps=0): 

1355 '''Check whether B{C{lat}} I{and} B{C{lon}} are within their 

1356 respective I{normal} range in C{degrees}. 

1357 

1358 @arg lat: Latitude (C{degrees}). 

1359 @arg lon: Longitude (C{degrees}). 

1360 @kwarg eps: Optional tolerance C{degrees}). 

1361 

1362 @return: C{True} if C{(abs(B{lat}) + B{eps}) <= 90} and 

1363 C{(abs(B{lon}) + B{eps}) <= 180}, C{False} othwerwise. 

1364 

1365 @see: Functions L{isnormal_} and L{normal}. 

1366 ''' 

1367 return (_90_0 - fabs(lat)) >= eps and (_180_0 - fabs(lon)) >= eps 

1368 

1369 

1370def isnormal_(phi, lam, eps=0): 

1371 '''Check whether B{C{phi}} I{and} B{C{lam}} are within their 

1372 respective I{normal} range in C{radians}. 

1373 

1374 @arg phi: Latitude (C{radians}). 

1375 @arg lam: Longitude (C{radians}). 

1376 @kwarg eps: Optional tolerance C{radians}). 

1377 

1378 @return: C{True} if C{(abs(B{phi}) + B{eps}) <= PI/2} and 

1379 C{(abs(B{lam}) + B{eps}) <= PI}, C{False} othwerwise. 

1380 

1381 @see: Functions L{isnormal} and L{normal_}. 

1382 ''' 

1383 return (PI_2 - fabs(phi)) >= eps and (PI - fabs(lam)) >= eps 

1384 

1385 

1386def latlon2n_xyz(lat, lon, name=NN): 

1387 '''Convert lat-, longitude to C{n-vector} (I{normal} to the 

1388 earth's surface) X, Y and Z components. 

1389 

1390 @arg lat: Latitude (C{degrees}). 

1391 @arg lon: Longitude (C{degrees}). 

1392 @kwarg name: Optional name (C{str}). 

1393 

1394 @return: A L{Vector3Tuple}C{(x, y, z)}. 

1395 

1396 @see: Function L{philam2n_xyz}. 

1397 

1398 @note: These are C{n-vector} x, y and z components, 

1399 I{NOT} geocentric ECEF x, y and z coordinates! 

1400 ''' 

1401 return _2n_xyz(name, *sincos2d_(lat, lon)) 

1402 

1403 

1404def _normal2(a, b, n_2, n, n2): 

1405 '''(INTERNAL) Helper for C{normal} and C{normal_}. 

1406 ''' 

1407 if fabs(b) > n: 

1408 b = remainder(b, n2) 

1409 if fabs(a) > n_2: 

1410 r = remainder(a, n) 

1411 if r != a: 

1412 a = -r 

1413 b -= n if b > 0 else -n 

1414 return float0_(a, b) 

1415 

1416 

1417def normal(lat, lon, name=NN): 

1418 '''Normalize a lat- I{and} longitude pair in C{degrees}. 

1419 

1420 @arg lat: Latitude (C{degrees}). 

1421 @arg lon: Longitude (C{degrees}). 

1422 @kwarg name: Optional name (C{str}). 

1423 

1424 @return: L{LatLon2Tuple}C{(lat, lon)} with C{abs(lat) <= 90} 

1425 and C{abs(lon) <= 180}. 

1426 

1427 @see: Functions L{normal_} and L{isnormal}. 

1428 ''' 

1429 return LatLon2Tuple(*_normal2(lat, lon, _90_0, _180_0, _360_0), 

1430 name=name or normal.__name__) 

1431 

1432 

1433def normal_(phi, lam, name=NN): 

1434 '''Normalize a lat- I{and} longitude pair in C{radians}. 

1435 

1436 @arg phi: Latitude (C{radians}). 

1437 @arg lam: Longitude (C{radians}). 

1438 @kwarg name: Optional name (C{str}). 

1439 

1440 @return: L{PhiLam2Tuple}C{(phi, lam)} with C{abs(phi) <= PI/2} 

1441 and C{abs(lam) <= PI}. 

1442 

1443 @see: Functions L{normal} and L{isnormal_}. 

1444 ''' 

1445 return PhiLam2Tuple(*_normal2(phi, lam, PI_2, PI, PI2), 

1446 name=name or normal_.__name__) 

1447 

1448 

1449def _2n_xyz(name, sa, ca, sb, cb): 

1450 '''(INTERNAL) Helper for C{latlon2n_xyz} and C{philam2n_xyz}. 

1451 ''' 

1452 # Kenneth Gade eqn 3, but using right-handed 

1453 # vector x -> 0°E,0°N, y -> 90°E,0°N, z -> 90°N 

1454 return Vector3Tuple(ca * cb, ca * sb, sa, name=name) 

1455 

1456 

1457def n_xyz2latlon(x, y, z, name=NN): 

1458 '''Convert C{n-vector} components to lat- and longitude in C{degrees}. 

1459 

1460 @arg x: X component (C{scalar}). 

1461 @arg y: Y component (C{scalar}). 

1462 @arg z: Z component (C{scalar}). 

1463 @kwarg name: Optional name (C{str}). 

1464 

1465 @return: A L{LatLon2Tuple}C{(lat, lon)}. 

1466 

1467 @see: Function L{n_xyz2philam}. 

1468 ''' 

1469 return LatLon2Tuple(atan2d(z, hypot(x, y)), atan2d(y, x), name=name) 

1470 

1471 

1472def n_xyz2philam(x, y, z, name=NN): 

1473 '''Convert C{n-vector} components to lat- and longitude in C{radians}. 

1474 

1475 @arg x: X component (C{scalar}). 

1476 @arg y: Y component (C{scalar}). 

1477 @arg z: Z component (C{scalar}). 

1478 @kwarg name: Optional name (C{str}). 

1479 

1480 @return: A L{PhiLam2Tuple}C{(phi, lam)}. 

1481 

1482 @see: Function L{n_xyz2latlon}. 

1483 ''' 

1484 return PhiLam2Tuple(atan2(z, hypot(x, y)), atan2(y, x), name=name) 

1485 

1486 

1487def _opposes(d, m, n, n2): 

1488 '''(INETNAL) Helper for C{opposing} and C{opposing_}. 

1489 ''' 

1490 d = d % n2 # -20 % 360 == 340, -1 % PI2 == PI2 - 1 

1491 return False if d < m or d > (n2 - m) else ( 

1492 True if (n - m) < d < (n + m) else None) 

1493 

1494 

1495def opposing(bearing1, bearing2, margin=_90_0): 

1496 '''Compare the direction of two bearings given in C{degrees}. 

1497 

1498 @arg bearing1: First bearing (compass C{degrees}). 

1499 @arg bearing2: Second bearing (compass C{degrees}). 

1500 @kwarg margin: Optional, interior angle bracket (C{degrees}). 

1501 

1502 @return: C{True} if both bearings point in opposite, C{False} if 

1503 in similar or C{None} if in perpendicular directions. 

1504 

1505 @see: Function L{opposing_}. 

1506 ''' 

1507 m = Degrees_(margin=margin, low=EPS0, high=_90_0) 

1508 return _opposes(bearing2 - bearing1, m,_180_0, _360_0) 

1509 

1510 

1511def opposing_(radians1, radians2, margin=PI_2): 

1512 '''Compare the direction of two bearings given in C{radians}. 

1513 

1514 @arg radians1: First bearing (C{radians}). 

1515 @arg radians2: Second bearing (C{radians}). 

1516 @kwarg margin: Optional, interior angle bracket (C{radians}). 

1517 

1518 @return: C{True} if both bearings point in opposite, C{False} if 

1519 in similar or C{None} if in perpendicular directions. 

1520 

1521 @see: Function L{opposing}. 

1522 ''' 

1523 m = Radians_(margin=margin, low=EPS0, high=PI_2) 

1524 return _opposes(radians2 - radians1, m, PI, PI2) 

1525 

1526 

1527def philam2n_xyz(phi, lam, name=NN): 

1528 '''Convert lat-, longitude to C{n-vector} (I{normal} to the 

1529 earth's surface) X, Y and Z components. 

1530 

1531 @arg phi: Latitude (C{radians}). 

1532 @arg lam: Longitude (C{radians}). 

1533 @kwarg name: Optional name (C{str}). 

1534 

1535 @return: A L{Vector3Tuple}C{(x, y, z)}. 

1536 

1537 @see: Function L{latlon2n_xyz}. 

1538 

1539 @note: These are C{n-vector} x, y and z components, 

1540 I{NOT} geocentric ECEF x, y and z coordinates! 

1541 ''' 

1542 return _2n_xyz(name, *sincos2_(phi, lam)) 

1543 

1544 

1545def _radical2(d, r1, r2): # in .ellipsoidalBaseDI, .sphericalTrigonometry, .vector3d 

1546 # (INTERNAL) See C{radical2} below 

1547 # assert d > EPS0 

1548 r = fsumf_(_1_0, (r1 / d)**2, -(r2 / d)**2) * _0_5 

1549 return Radical2Tuple(max(_0_0, min(_1_0, r)), r * d) 

1550 

1551 

1552def radical2(distance, radius1, radius2): 

1553 '''Compute the I{radical ratio} and I{radical line} of two 

1554 U{intersecting circles<https://MathWorld.Wolfram.com/ 

1555 Circle-CircleIntersection.html>}. 

1556 

1557 The I{radical line} is perpendicular to the axis thru the 

1558 centers of the circles at C{(0, 0)} and C{(B{distance}, 0)}. 

1559 

1560 @arg distance: Distance between the circle centers (C{scalar}). 

1561 @arg radius1: Radius of the first circle (C{scalar}). 

1562 @arg radius2: Radius of the second circle (C{scalar}). 

1563 

1564 @return: A L{Radical2Tuple}C{(ratio, xline)} where C{0.0 <= 

1565 ratio <= 1.0} and C{xline} is along the B{C{distance}}. 

1566 

1567 @raise IntersectionError: The B{C{distance}} exceeds the sum 

1568 of B{C{radius1}} and B{C{radius2}}. 

1569 

1570 @raise UnitError: Invalid B{C{distance}}, B{C{radius1}} or 

1571 B{C{radius2}}. 

1572 

1573 @see: U{Circle-Circle Intersection 

1574 <https://MathWorld.Wolfram.com/Circle-CircleIntersection.html>}. 

1575 ''' 

1576 d = Distance_(distance, low=_0_0) 

1577 r1 = Radius_(radius1=radius1) 

1578 r2 = Radius_(radius2=radius2) 

1579 if d > (r1 + r2): 

1580 raise IntersectionError(distance=d, radius1=r1, radius2=r2, 

1581 txt=_too_(_distant_)) 

1582 return _radical2(d, r1, r2) if d > EPS0 else \ 

1583 Radical2Tuple(_0_5, _0_0) 

1584 

1585 

1586class Radical2Tuple(_NamedTuple): 

1587 '''2-Tuple C{(ratio, xline)} of the I{radical} C{ratio} and 

1588 I{radical} C{xline}, both C{scalar} and C{0.0 <= ratio <= 1.0} 

1589 ''' 

1590 _Names_ = (_ratio_, _xline_) 

1591 _Units_ = ( Scalar, Scalar) 

1592 

1593 

1594def _radistance(inst): 

1595 '''(INTERNAL) Helper for the L{frechet._FrecherMeterRadians} 

1596 and L{hausdorff._HausdorffMeterRedians} classes. 

1597 ''' 

1598 kwds_ = _xkwds(inst._kwds, wrap=False) 

1599 wrap_ = _xkwds_pop(kwds_, wrap=False) 

1600 func_ = inst._func_ 

1601 try: # calling lower-overhead C{func_} 

1602 func_(0, _0_25, _0_5, **kwds_) 

1603 wrap_ = _Wrap._philamop(wrap_) 

1604 except TypeError: 

1605 return inst.distance 

1606 

1607 _rad = radians 

1608 

1609 def _philam(p): 

1610 return _rad(p.lat), _rad(p.lon) 

1611 

1612 def _radistance(point1, point2): 

1613 try: 

1614 phi1, lam1 = wrap_(point1.phi, point1.lam) 

1615 phi2, lam2 = wrap_(point2.phi, point2.lam) 

1616 except AttributeError: # no 'phi or .lam 

1617 phi1, lam1 = wrap_(*_philam(point1)) 

1618 phi2, lam2 = wrap_(*_philam(point2)) 

1619 return func_(phi2, phi1, lam2 - lam1, **kwds_) 

1620 

1621 inst._units = inst._units_ 

1622 return _radistance 

1623 

1624 

1625def _scale_deg(lat1, lat2): # degrees 

1626 # scale factor cos(mean of lats) for delta lon 

1627 m = fabs(lat1 + lat2) * _0_5 

1628 return cos(radians(m)) if m < 90 else _0_0 

1629 

1630 

1631def _scale_rad(phi1, phi2): # radians, by .frechet, .hausdorff, .heights 

1632 # scale factor cos(mean of phis) for delta lam 

1633 m = fabs(phi1 + phi2) * _0_5 

1634 return cos(m) if m < PI_2 else _0_0 

1635 

1636 

1637def _sincosa6(phi2, phi1, lam21): # [4] in cosineLaw 

1638 '''(INTERNAL) C{sin}es, C{cos}ines and C{acos}ine. 

1639 ''' 

1640 s2, c2, s1, c1, _, c21 = sincos2_(phi2, phi1, lam21) 

1641 return s2, c2, s1, c1, acos1(s1 * s2 + c1 * c2 * c21), c21 

1642 

1643 

1644def thomas(lat1, lon1, lat2, lon2, datum=_WGS84, wrap=False): 

1645 '''Compute the distance between two (ellipsoidal) points using 

1646 U{Thomas'<https://apps.DTIC.mil/dtic/tr/fulltext/u2/703541.pdf>} 

1647 formula. 

1648 

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 datum: Datum (L{Datum}) or ellipsoid (L{Ellipsoid}, 

1654 L{Ellipsoid2} or L{a_f2Tuple}) to use. 

1655 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll 

1656 B{C{lat2}} and B{C{lon2}} (C{bool}). 

1657 

1658 @return: Distance (C{meter}, same units as the B{C{datum}}'s 

1659 ellipsoid axes). 

1660 

1661 @raise TypeError: Invalid B{C{datum}}. 

1662 

1663 @see: Functions L{thomas_}, L{cosineAndoyerLambert}, L{cosineForsytheAndoyerLambert}, 

1664 L{cosineLaw}, L{equirectangular}, L{euclidean}, L{flatLocal}/L{hubeny}, 

1665 L{flatPolar}, L{haversine}, L{vincentys} and method L{Ellipsoid.distance2}. 

1666 ''' 

1667 return _dE(thomas_, datum, wrap, lat1, lon1, lat2, lon2) 

1668 

1669 

1670def thomas_(phi2, phi1, lam21, datum=_WGS84): 

1671 '''Compute the I{angular} distance between two (ellipsoidal) points using 

1672 U{Thomas'<https://apps.DTIC.mil/dtic/tr/fulltext/u2/703541.pdf>} 

1673 formula. 

1674 

1675 @arg phi2: End latitude (C{radians}). 

1676 @arg phi1: Start latitude (C{radians}). 

1677 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 

1678 @kwarg datum: Datum or ellipsoid to use (L{Datum}, L{Ellipsoid}, 

1679 L{Ellipsoid2} or L{a_f2Tuple}). 

1680 

1681 @return: Angular distance (C{radians}). 

1682 

1683 @raise TypeError: Invalid B{C{datum}}. 

1684 

1685 @see: Functions L{thomas}, L{cosineAndoyerLambert_}, 

1686 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, 

1687 L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_}, 

1688 L{flatPolar_}, L{haversine_} and L{vincentys_} and U{Geodesy-PHP 

1689 <https://GitHub.com/jtejido/geodesy-php/blob/master/src/Geodesy/ 

1690 Distance/ThomasFormula.php>}. 

1691 ''' 

1692 s2, c2, s1, c1, r, _ = _sincosa6(phi2, phi1, lam21) 

1693 if r and isnon0(c1) and isnon0(c2): 

1694 E = _ellipsoidal(datum, thomas_) 

1695 if E.f: 

1696 r1 = atan2(E.b_a * s1, c1) 

1697 r2 = atan2(E.b_a * s2, c2) 

1698 

1699 j = (r2 + r1) * _0_5 

1700 k = (r2 - r1) * _0_5 

1701 sj, cj, sk, ck, h, _ = sincos2_(j, k, lam21 * _0_5) 

1702 

1703 h = fsumf_(sk**2, (ck * h)**2, -(sj * h)**2) 

1704 u = _1_0 - h 

1705 if isnon0(u) and isnon0(h): 

1706 r = atan(sqrt0(h / u)) * 2 # == acos(1 - 2 * h) 

1707 sr, cr = sincos2(r) 

1708 if isnon0(sr): 

1709 u = 2 * (sj * ck)**2 / u 

1710 h = 2 * (sk * cj)**2 / h 

1711 x = u + h 

1712 y = u - h 

1713 

1714 s = r / sr 

1715 e = 4 * s**2 

1716 d = 2 * cr 

1717 a = e * d 

1718 b = 2 * r 

1719 c = s - (a - d) * _0_5 

1720 f = E.f * _0_25 

1721 

1722 t = fsumf_(a * x, -b * y, c * x**2, -d * y**2, e * x * y) 

1723 r -= fsumf_(s * x, -y, -t * f * _0_25) * f * sr 

1724 return r 

1725 

1726 

1727def vincentys(lat1, lon1, lat2, lon2, radius=R_M, wrap=False): 

1728 '''Compute the distance between two (spherical) points using 

1729 U{Vincenty's<https://WikiPedia.org/wiki/Great-circle_distance>} 

1730 spherical formula. 

1731 

1732 @arg lat1: Start latitude (C{degrees}). 

1733 @arg lon1: Start longitude (C{degrees}). 

1734 @arg lat2: End latitude (C{degrees}). 

1735 @arg lon2: End longitude (C{degrees}). 

1736 @kwarg radius: Mean earth radius (C{meter}), datum (L{Datum}) 

1737 or ellipsoid (L{Ellipsoid}, L{Ellipsoid2} or 

1738 L{a_f2Tuple}) to use. 

1739 @kwarg wrap: If C{True}, wrap or I{normalize} and unroll 

1740 B{C{lat2}} and B{C{lon2}} (C{bool}). 

1741 

1742 @return: Distance (C{meter}, same units as B{C{radius}}). 

1743 

1744 @raise UnitError: Invalid B{C{radius}}. 

1745 

1746 @see: Functions L{vincentys_}, L{cosineAndoyerLambert}, 

1747 L{cosineForsytheAndoyerLambert},L{cosineLaw}, L{equirectangular}, 

1748 L{euclidean}, L{flatLocal}/L{hubeny}, L{flatPolar}, 

1749 L{haversine} and L{thomas} and methods L{Ellipsoid.distance2}, 

1750 C{LatLon.distanceTo*} and C{LatLon.equirectangularTo}. 

1751 

1752 @note: See note at function L{vincentys_}. 

1753 ''' 

1754 return _dS(vincentys_, radius, wrap, lat1, lon1, lat2, lon2) 

1755 

1756 

1757def vincentys_(phi2, phi1, lam21): 

1758 '''Compute the I{angular} distance between two (spherical) points using 

1759 U{Vincenty's<https://WikiPedia.org/wiki/Great-circle_distance>} 

1760 spherical formula. 

1761 

1762 @arg phi2: End latitude (C{radians}). 

1763 @arg phi1: Start latitude (C{radians}). 

1764 @arg lam21: Longitudinal delta, M{end-start} (C{radians}). 

1765 

1766 @return: Angular distance (C{radians}). 

1767 

1768 @see: Functions L{vincentys}, L{cosineAndoyerLambert_}, 

1769 L{cosineForsytheAndoyerLambert_}, L{cosineLaw_}, 

1770 L{equirectangular_}, L{euclidean_}, L{flatLocal_}/L{hubeny_}, 

1771 L{flatPolar_}, L{haversine_} and L{thomas_}. 

1772 

1773 @note: Functions L{vincentys_}, L{haversine_} and L{cosineLaw_} 

1774 produce equivalent results, but L{vincentys_} is suitable 

1775 for antipodal points and slightly more expensive (M{3 cos, 

1776 3 sin, 1 hypot, 1 atan2, 6 mul, 2 add}) than L{haversine_} 

1777 (M{2 cos, 2 sin, 2 sqrt, 1 atan2, 5 mul, 1 add}) and 

1778 L{cosineLaw_} (M{3 cos, 3 sin, 1 acos, 3 mul, 1 add}). 

1779 ''' 

1780 s1, c1, s2, c2, s21, c21 = sincos2_(phi1, phi2, lam21) 

1781 

1782 c = c2 * c21 

1783 x = s1 * s2 + c1 * c 

1784 y = c1 * s2 - s1 * c 

1785 return atan2(hypot(c2 * s21, y), x) 

1786 

1787# **) MIT License 

1788# 

1789# Copyright (C) 2016-2023 -- mrJean1 at Gmail -- All Rights Reserved. 

1790# 

1791# Permission is hereby granted, free of charge, to any person obtaining a 

1792# copy of this software and associated documentation files (the "Software"), 

1793# to deal in the Software without restriction, including without limitation 

1794# the rights to use, copy, modify, merge, publish, distribute, sublicense, 

1795# and/or sell copies of the Software, and to permit persons to whom the 

1796# Software is furnished to do so, subject to the following conditions: 

1797# 

1798# The above copyright notice and this permission notice shall be included 

1799# in all copies or substantial portions of the Software. 

1800# 

1801# THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS 

1802# OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, 

1803# FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL 

1804# THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR 

1805# OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, 

1806# ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR 

1807# OTHER DEALINGS IN THE SOFTWARE.