Coverage for pygeodesy/geoids.py: 96%
688 statements
« prev ^ index » next coverage.py v7.6.1, created at 2025-04-09 11:05 -0400
« prev ^ index » next coverage.py v7.6.1, created at 2025-04-09 11:05 -0400
2# -*- coding: utf-8 -*-
4u'''Geoid models and geoid height interpolations.
6Classes L{GeoidEGM96}, L{GeoidG2012B}, L{GeoidKarney} and L{GeoidPGM} to
7interpolate the height of various U{geoid<https://WikiPedia.org/wiki/Geoid>}s
8at C{LatLon} locations or separate lat-/longitudes using various interpolation
9methods and C{geoid} model files.
11L{GeoidKarney} is a transcoding of I{Charles Karney}'s C++ class U{Geoid
12<https://GeographicLib.SourceForge.io/C++/doc/geoid.html>} to pure Python.
14The L{GeoidEGM96}, L{GeoidG2012B} and L{GeoidPGM} interpolators both depend on
15U{scipy<https://SciPy.org>} and U{numpy<https://PyPI.org/project/numpy>} and
16require those packages to be installed.
18In addition, each geoid interpolator needs C{grid knots} (down)loaded from
19a C{geoid} model file, I{specific to the interpolator}. More details below
20and in the documentation of the interpolator class. For each interpolator,
21there are several interpolation choices, like I{linear}, I{cubic}, etc.
23Typical usage
24=============
261. Choose an interpolator class L{GeoidEGM96}, L{GeoidG2012B}, L{GeoidKarney}
27or L{GeoidPGM} and download a C{geoid} model file, containing locations with
28known heights also referred to as the C{grid knots}. See the documentation of
29the interpolator class for references to available C{grid} models.
31C{>>> from pygeodesy import GeoidEGM96 # or -G2012B, -Karney or -PGM as GeoidXyz}
332. Instantiate an interpolator with the C{geoid} model file and use keyword
34arguments to select different interpolation options
36C{>>> ginterpolator = GeoidXyz(geoid_model_file, **options)}
383. Get the interpolated geoid height of C{LatLon} location(s) with
40C{>>> ll = LatLon(1, 2, ...)}
42C{>>> h = ginterpolator(ll)}
44or
46C{>>> h1, h2, h3, ... = ginterpolator(ll1, ll2, ll3, ...)}
48or a list, tuple, generator, etc. of C{LatLon}s
50C{>>> hs = ginterpolator(lls)}
524. For separate lat- and longitudes invoke the C{height} method as
54C{>>> h = ginterpolator.height(lat, lon)}
56or as 2 lists, 2 tuples, etc.
58C{>>> hs = ginterpolator.height(lats, lons)}
60or for several positionals use the C{height_} method
62C{>>> h1, h2, ... = ginterpolator.height_(lat1, lon1, lat2, lon2, ...)}
645. An example is in U{issue #64<https://GitHub.com/mrJean1/PyGeodesy/issues/64>},
65courtesy of SBFRF.
67@note: Classes L{GeoidEGM96}, L{GeoidG2012B} and L{GeoidPGM} require both U{numpy
68 <https://PyPI.org/project/numpy>} and U{scipy<https://PyPI.org/project/scipy>}
69 to be installed.
71@note: Errors from C{scipy} are raised as L{SciPyError}s. Warnings issued by C{scipy} can
72 be thrown as L{SciPyWarning} exceptions, provided Python C{warnings} are filtered
73 accordingly, see L{SciPyWarning}.
75@see: I{Karney}'s U{GeographicLib<https://GeographicLib.SourceForge.io/C++/doc/index.html>},
76 U{Geoid height<https://GeographicLib.SourceForge.io/C++/doc/geoid.html>} and U{Installing
77 the Geoid datasets<https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinst>},
78 World Geodetic System 1984 (WG84) and U{Earth Gravitational Model 96 (EGM96) Data and
79 Apps<https://earth-info.NGA.mil/index.php?dir=wgs84&action=wgs84>},
80 U{SciPy<https://docs.SciPy.org/doc/scipy/reference/interpolate.html>} interpolation
81 U{RectBivariateSpline<https://docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.
82 RectBivariateSpline.html>}, U{bisplrep/-ev<https://docs.scipy.org/doc/scipy/reference/generated/
83 scipy.interpolate.bisplrep.html>} and U{interp2d<https://docs.SciPy.org/doc/scipy/reference/
84 generated/scipy.interpolate.interp2d.html>}, functions L{elevations.elevation2} and
85 L{elevations.geoidHeight2}, U{I{Ellispoid vs Orthometric Elevations}<https://www.YouTube.com/
86 watch?v=dX6a6kCk3Po>} and U{6.22.1 Avoiding Pitfalls Related to Ellipsoid Height and Height
87 Above Mean Sea Level<https://Wiki.ROS.org/mavros>}.
88'''
89# make sure int/int division yields float quotient, see .basics
90from __future__ import division as _; del _ # PYCHOK semicolon
92from pygeodesy.basics import len2, min2, isodd, _splituple, ub2str as _ub2str
93from pygeodesy.constants import EPS, _float as _F, _1_0, _N_90_0, _180_0, \
94 _N_180_0, _360_0
95from pygeodesy.datums import Datums, _ellipsoidal_datum, _WGS84
96# from pygeodesy.dms import parseDMS2 # _MODS
97from pygeodesy.errors import _incompatible, LenError, RangeError, _SciPyIssue, \
98 _xkwds_pop2
99from pygeodesy.fmath import favg, Fdot, fdot, Fhorner, frange
100# from pygoedesy.formy import heightOrthometric # _MODS
101from pygeodesy.heights import _as_llis2, _ascalar, _HeightBase, HeightError, \
102 _Wrap
103# from pygeodesy.internals import _version2 # _MODS
104from pygeodesy.interns import NN, _COLONSPACE_, _COMMASPACE_, _E_, _height_, \
105 _in_, _kind_, _lat_, _lon_, _mean_, _N_, _n_a_, \
106 _numpy_, _on_, _outside_, _S_, _s_, _scipy_, \
107 _SPACE_, _stdev_, _tbd_, _W_, _width_, _4_
108from pygeodesy.lazily import _ALL_DOCS, _ALL_LAZY, _ALL_MODS as _MODS, _FOR_DOCS
109from pygeodesy.named import _name__, _Named, _NamedTuple
110# from pygeodesy.namedTuples import LatLon3Tuple # _MODS
111from pygeodesy.props import Property_RO, property_RO, property_ROver
112from pygeodesy.streprs import attrs, Fmt, fstr, pairs
113from pygeodesy.units import Height, Int_, Lat, Lon
114# from pygeodesy.utily import _Wrap # from .heights
116from math import floor as _floor
117# from os import SEEK_CUR, SEEK_SET # _MODS
118# import os.path # _MODS
119from struct import calcsize as _calcsize, unpack as _unpack
120try:
121 from StringIO import StringIO as _BytesIO # reads bytes
122 _ub2str = str # PYCHOK convert bytes to str for egm*.pgm text
124except ImportError: # Python 3+
125 from io import BytesIO as _BytesIO # PYCHOK expected
127__all__ = _ALL_LAZY.geoids
128__version__ = '25.01.15'
130_assert_ = 'assert'
131_bHASH_ = b'#'
132_endian_ = 'endian'
133_format_ = '%s %r'
134_header_ = 'header'
135_intCs = {} # cache int value, del below
136_lli_ = 'lli'
137_rb_ = 'rb'
138_supported_ = 'supported'
141class GeoidError(HeightError):
142 '''Geoid interpolator C{Geoid...} or interpolation issue.
143 '''
144 pass
147class _GeoidBase(_HeightBase):
148 '''(INTERNAL) Base class for C{Geoid...}s.
149 '''
150 _center = None # (lat, lon, height)
151 _cropped = None
152# _datum = _WGS84 # from _HeightBase
153 _egm = None # open C{egm*.pgm} geoid file
154 _endian = _tbd_
155 _Error = GeoidError # in ._HeightBase._as_lls, ...
156 _geoid = _n_a_
157 _hs_y_x = None # numpy 2darray, row-major order
158 _iscipy = True # scipy or Karney's interpolation
159 _kind = 3 # order for interp2d, RectBivariateSpline
160# _kmin = 2 # min number of knots
161 _knots = 0 # nlat * nlon
162 _mean = None # fixed in GeoidKarney
163 _nBytes = 0 # numpy size in bytes, float64
164 _pgm = None # PGM attributes, C{_PGM} or C{None}
165 _sizeB = 0 # geoid file size in bytes
166 _smooth = 0 # used only for RectBivariateSpline
167 _stdev = None # fixed in GeoidKarney
168 _u2B = 0 # np.itemsize or undefined
169 _yx_hits = None # cache hits, ala Karney's
171# _lat_d = _0_0 # increment, +tive
172# _lat_lo = _0_0 # lower lat, south
173# _lat_hi = _0_0 # upper lat, noth
174# _lon_d = _0_0 # increment, +tive
175# _lon_lo = _0_0 # left lon, west
176# _lon_hi = _0_0 # right lon, east
177# _lon_of = _0_0 # forward lon offset
178# _lon_og = _0_0 # reverse lon offset
180 def __init__(self, hs, p):
181 '''(INTERNAL) Set up the grid axes, the C{SciPy} interpolator and
182 several internal geoid attributes.
184 @arg hs: Grid knots with known height (C{numpy 2darray}).
185 @arg p: The C{slat, wlon, nlat, nlon, dlat, dlon} and other
186 geoid parameters (C{INTERNAL}).
187 '''
188 spi = self.scipy_interpolate
189 # for 2d scipy.interpolate.interp2d(xs, ys, hs, ...) and
190 # scipy.interpolate.RectBivariateSpline(ys, xs, hs, ...)
191 # require the shape of hs to be (len(ys), len(xs)), note
192 # the different (xs, ys, ...) and (ys, xs, ...) orders
193 if (p.nlat, p.nlon) != hs.shape:
194 raise GeoidError(shape=hs.shape, txt=_incompatible((p.nlat, p.nlon)))
196 # both axes and bounding box
197 ys, self._lat_d = self._gaxis2(p.slat, p.dlat, p.nlat, _lat_ + _s_)
198 xs, self._lon_d = self._gaxis2(p.wlon, p.dlon, p.nlon, _lon_ + _s_)
200 bb = ys[0], ys[-1], xs[0], xs[-1] + p.dlon # fudge lon_hi
201 # geoid grids are typically stored in row-major order, some
202 # with rows (90..-90) reversed and columns (0..360) wrapped
203 # to Easten longitude, 0 <= east < 180 and 180 <= west < 360
204 k = self.kind
205 if k in self._k2interp2d: # see _HeightBase
206 self._interp2d(xs, ys, hs, kind=k)
207 else: # XXX order ys and xs, see HeightLSQBiSpline
208 k = self._kxky(k)
209 self._ev = spi.RectBivariateSpline(ys, xs, hs, bbox=bb, ky=k, kx=k,
210 s=self._smooth).ev
211 self._hs_y_x = hs # numpy 2darray, row-major
212 self._nBytes = hs.nbytes # numpy size in bytes
213 self._knots = p.knots # grid knots == len(hs)
214 self._lon_of = float(p.flon) # forward offset
215 self._lon_og = g = float(p.glon) # reverse offset
216 # shrink the bounding box by 1 unit on every side:
217 # +self._lat_d, -self._lat_d, +self._lon_d, -self._lon_d
218 self._lat_lo, \
219 self._lat_hi, \
220 self._lon_lo, \
221 self._lon_hi = map(float, bb)
222 self._lon_lo -= g
223 self._lon_hi -= g
225 def __call__(self, *llis, **wrap_H):
226 '''Interpolate the geoid (or orthometric) height for one or more locations.
228 @arg llis: One or several locations (each C{LatLon}), all positional.
229 @kwarg wrap_H: Keyword arguments C{B{wrap}=False} (C{bool}) and
230 C{B{H}=False} (C{bool}). Use C{B{wrap}=True} to wrap
231 or I{normalize} all B{C{llis}} locations. If C{B{H}
232 is True}, return the I{orthometric} height instead of
233 the I{geoid} height at each location.
235 @return: A single geoid (or orthometric) height (C{float}) or
236 a list or tuple of geoid (or orthometric) heights (each
237 C{float}).
239 @raise GeoidError: Insufficient number of B{C{llis}}, an invalid
240 B{C{lli}} or the C{egm*.pgm} geoid file is closed.
242 @raise RangeError: An B{C{lli}} is outside this geoid's lat- or
243 longitude range.
245 @raise SciPyError: A C{scipy} issue.
247 @raise SciPyWarning: A C{scipy} warning as exception.
249 @note: To obtain I{orthometric} heights, each B{C{llis}} location
250 must have an ellipsoid C{height} or C{h} attribute, otherwise
251 C{height=0} is used.
253 @see: Function L{pygeodesy.heightOrthometric}.
254 '''
255 return self._called(llis, **wrap_H)
257 def __enter__(self):
258 '''Open context.
259 '''
260 return self
262 def __exit__(self, *unused): # PYCHOK exc_type, exc_value, exc_traceback)
263 '''Close context.
264 '''
265 self.close()
266 # return None # XXX False
268 def __repr__(self):
269 return self.toStr()
271 def __str__(self):
272 return Fmt.PAREN(self.classname, repr(self.name))
274 def _called(self, llis, wrap=False, H=False):
275 # handle __call__
276 _H = self._heightOrthometric if H else None
277 _as, llis = _as_llis2(llis, Error=GeoidError)
278 _w, hs = _Wrap._latlonop(wrap), []
279 _h, _a = self._hGeoid, hs.append
280 try:
281 for i, lli in enumerate(llis):
282 N = _h(*_w(lli.lat, lli.lon))
283 # orthometric or geoid height
284 _a(_H(lli, N) if _H else N)
285 return _as(hs)
286 except (GeoidError, RangeError) as x:
287 # XXX avoid str(LatLon()) degree symbols
288 n = _lli_ if _as is _ascalar else Fmt.INDEX(llis=i)
289 t = fstr((lli.lat, lli.lon), strepr=repr)
290 raise type(x)(n, t, wrap=wrap, H=H, cause=x)
291 except Exception as x:
292 if self._iscipy and self.scipy:
293 raise _SciPyIssue(x, self._ev_name)
294 else:
295 raise
297 @Property_RO
298 def _center(self):
299 ''' Cache for method L{center}.
300 '''
301 return self._llh3(favg(self._lat_lo, self._lat_hi),
302 favg(self._lon_lo, self._lon_hi))
304 def center(self, LatLon=None):
305 '''Return the center location and height of this geoid.
307 @kwarg LatLon: Optional class to return the location and height
308 (C{LatLon}) or C{None}.
310 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon,
311 height)} otherwise a B{C{LatLon}} instance with the lat-,
312 longitude and geoid height of the center grid location.
313 '''
314 return self._llh3LL(self._center, LatLon)
316 def close(self):
317 '''Close the C{egm*.pgm} geoid file if open (and applicable).
318 '''
319 if not self.closed:
320 self._egm.close()
321 self._egm = None
323 @property_RO
324 def closed(self):
325 '''Get the C{egm*.pgm} geoid file status.
326 '''
327 return self._egm is None
329 @Property_RO
330 def cropped(self):
331 '''Is geoid cropped (C{bool} or C{None} if crop not supported).
332 '''
333 return self._cropped
335 @Property_RO
336 def dtype(self):
337 '''Get the grid C{scipy} U{dtype<https://docs.SciPy.org/doc/numpy/
338 reference/generated/numpy.ndarray.dtype.html>} (C{numpy.dtype}).
339 '''
340 return self._hs_y_x.dtype
342 @Property_RO
343 def endian(self):
344 '''Get the geoid endianess and U{dtype<https://docs.SciPy.org/
345 doc/numpy/reference/generated/numpy.dtype.html>} (C{str}).
346 '''
347 return self._endian
349 def _ev(self, y, x): # PYCHOK overwritten with .RectBivariateSpline.ev
350 # see methods _HeightBase._ev and -._interp2d
351 return self._ev2d(x, y) # (y, x) flipped!
353 def _gaxis2(self, lo, d, n, name):
354 # build grid axis, hi = lo + (n - 1) * d
355 m, a = len2(frange(lo, n, d))
356 if m != n:
357 raise LenError(self.__class__, grid=m, **{name: n})
358 if d < 0:
359 d, a = -d, list(reversed(a))
360 a = self.numpy.array(a)
361 m, i = min2(*map(float, a[1:] - a[:-1]))
362 if m < EPS: # non-increasing axis
363 i = Fmt.INDEX(name, i + 1)
364 raise GeoidError(i, m, txt_not_='increasing')
365 return a, d
367 def _g2ll2(self, lat, lon): # PYCHOK no cover
368 '''(INTERNAL) I{Must be overloaded}.'''
369 self._notOverloaded(lat, lon)
371 def _gyx2g2(self, y, x):
372 # convert grid (y, x) indices to grid (lat, lon)
373 return ((self._lat_lo + self._lat_d * y),
374 (self._lon_lo + self._lon_of + self._lon_d * x))
376 def height(self, lats, lons, **wrap):
377 '''Interpolate the geoid height for one or several lat-/longitudes.
379 @arg lats: Latitude or latitudes (each C{degrees}).
380 @arg lons: Longitude or longitudes (each C{degrees}).
381 @kwarg wrap: Use C{B{wrap}=True} to wrap or I{normalize} all
382 B{C{lats}} and B{C{lons}}.
384 @return: A single geoid height (C{float}) or a list of geoid
385 heights (each C{float}).
387 @raise GeoidError: Insufficient or unequal number of B{C{lats}}
388 and B{C{lons}}.
390 @raise RangeError: A B{C{lat}} or B{C{lon}} is outside this geoid's
391 lat- or longitude range.
393 @raise SciPyError: A C{scipy} issue.
395 @raise SciPyWarning: A C{scipy} warning as exception.
396 '''
397 return _HeightBase.height(self, lats, lons, **wrap)
399 def height_(self, *latlons, **wrap):
400 '''Interpolate the geoid height for each M{(latlons[i], latlons[i+1])
401 pair for i in range(0, len(latlons), B{2})}.
403 @arg latlons: Alternating lat-/longitude pairs (each C{degrees}),
404 all positional.
406 @see: Method L{height} for further details.
408 @return: A tuple of geoid heights (each C{float}).
409 '''
410 lls = tuple(self._as_lls(latlons[0::2], *latlons[1::2]))
411 return self._called(lls, **wrap)
413 @property_ROver
414 def _heightOrthometric(self):
415 return _MODS.formy.heightOrthometric # overwrite property_ROver
417 def _hGeoid(self, lat, lon):
418 out = self.outside(lat, lon)
419 if out: # XXX avoid str(LatLon()) degree symbols
420 t = fstr((lat, lon), strepr=repr)
421 raise RangeError(lli=t, txt=_SPACE_(_outside_, _on_, out))
422 return float(self._ev(*self._ll2g2(lat, lon)))
424 @Property_RO
425 def _highest(self):
426 '''(INTERNAL) Cache for C{.highest}.
427 '''
428 return self._LL3T(self._llh3minmax(True), name__=self.highest)
430 def highest(self, LatLon=None, **unused):
431 '''Return the location and largest height of this geoid.
433 @kwarg LatLon: Optional class to return the location and height
434 (C{LatLon}) or C{None}.
436 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon,
437 height)} otherwise a B{C{LatLon}} instance with the lat-,
438 longitude and geoid height of the highest grid location.
439 '''
440 return self._llh3LL(self._highest, LatLon)
442 @Property_RO
443 def hits(self):
444 '''Get the number of cache hits (C{int} or C{None}).
445 '''
446 return self._yx_hits
448 @Property_RO
449 def kind(self):
450 '''Get the interpolator kind and order (C{int}).
451 '''
452 return self._kind
454 @Property_RO
455 def knots(self):
456 '''Get the number of grid knots (C{int}).
457 '''
458 return self._knots
460 def _ll2g2(self, lat, lon): # PYCHOK no cover
461 '''(INTERNAL) I{Must be overloaded}.'''
462 self._notOverloaded(lat, lon)
464 @property_ROver
465 def _LL3T(self):
466 '''(INTERNAL) Get L{LatLon3Tuple}, I{once}.
467 '''
468 return _MODS.namedTuples.LatLon3Tuple # overwrite property_ROver
470 def _llh3(self, lat, lon):
471 return self._LL3T(lat, lon, self._hGeoid(lat, lon), name=self.name)
473 def _llh3LL(self, llh, LatLon):
474 return llh if LatLon is None else self._xnamed(LatLon(*llh))
476 def _llh3minmax(self, highest, *unused):
477 hs, np = self._hs_y_x, self.numpy
478 # <https://docs.SciPy.org/doc/numpy/reference/generated/
479 # numpy.argmin.html#numpy.argmin>
480 arg = np.argmax if highest else np.argmin
481 y, x = np.unravel_index(arg(hs, axis=None), hs.shape)
482 return self._g2ll2(*self._gyx2g2(y, x)) + (float(hs[y, x]),)
484 def _load(self, g, dtype=float, n=-1, offset=0, **sep): # sep=NN
485 # numpy.fromfile, like .frombuffer
486 g.seek(offset, _MODS.os.SEEK_SET)
487 return self.numpy.fromfile(g, dtype, count=n, **sep)
489 @Property_RO
490 def _lowerleft(self):
491 '''(INTERNAL) Cache for C{.lowerleft}.
492 '''
493 return self._llh3(self._lat_lo, self._lon_lo)
495 def lowerleft(self, LatLon=None):
496 '''Return the lower-left location and height of this geoid.
498 @kwarg LatLon: Optional class to return the location
499 (C{LatLon}) and height or C{None}.
501 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon, height)}
502 otherwise a B{C{LatLon}} instance with the lat-, longitude and
503 geoid height of the lower-left, SW grid corner.
504 '''
505 return self._llh3LL(self._lowerleft, LatLon)
507 @Property_RO
508 def _loweright(self):
509 '''(INTERNAL) Cache for C{.loweright}.
510 '''
511 return self._llh3(self._lat_lo, self._lon_hi)
513 def loweright(self, LatLon=None):
514 '''Return the lower-right location and height of this geoid.
516 @kwarg LatLon: Optional class to return the location and height
517 (C{LatLon}) or C{None}.
519 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon, height)}
520 otherwise a B{C{LatLon}} instance with the lat-, longitude and
521 geoid height of the lower-right, SE grid corner.
522 '''
523 return self._llh3LL(self._loweright, LatLon)
525 lowerright = loweright # synonymous
527 @Property_RO
528 def _lowest(self):
529 '''(INTERNAL) Cache for C{.lowest}.
530 '''
531 return self._LL3T(self._llh3minmax(False), name__=self.lowest)
533 def lowest(self, LatLon=None, **unused):
534 '''Return the location and lowest height of this geoid.
536 @kwarg LatLon: Optional class to return the location and height
537 (C{LatLon}) or C{None}.
539 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon,
540 height)} otherwise a B{C{LatLon}} instance with the lat-,
541 longitude and geoid height of the lowest grid location.
542 '''
543 return self._llh3LL(self._lowest, LatLon)
545 @Property_RO
546 def mean(self):
547 '''Get the mean of this geoid's heights (C{float}).
548 '''
549 if self._mean is None: # see GeoidKarney
550 self._mean = float(self.numpy.mean(self._hs_y_x))
551 return self._mean
553 @property_RO
554 def name(self):
555 '''Get the name of this geoid (C{str}).
556 '''
557 return _HeightBase.name.fget(self) or self._geoid # recursion
559 @Property_RO
560 def nBytes(self):
561 '''Get the grid in-memory size in bytes (C{int}).
562 '''
563 return self._nBytes
565 def _open(self, geoid, datum, kind, name, smooth):
566 # open the geoid file
567 try:
568 self._geoid = _MODS.os.path.basename(geoid)
569 self._sizeB = _MODS.os.path.getsize(geoid)
570 g = open(geoid, _rb_)
571 except (IOError, OSError) as x:
572 raise GeoidError(geoid=geoid, cause=x)
574 if datum not in (None, self._datum):
575 self._datum = _ellipsoidal_datum(datum, name=name)
576 self._kind = int(kind)
577 if name:
578 _HeightBase.name.fset(self, name) # rename
579 if smooth:
580 self._smooth = Int_(smooth=smooth, Error=GeoidError, low=0)
582 return g
584 def outside(self, lat, lon):
585 '''Check whether a location is outside this geoid's lat-/longitude
586 or crop range.
588 @arg lat: The latitude (C{degrees}).
589 @arg lon: The longitude (C{degrees}).
591 @return: A 1- or 2-character C{str} if outside, an empty C{str} otherwise.
592 '''
593 lat = _S_ if lat < self._lat_lo else (_N_ if lat > self._lat_hi else NN)
594 lon = _W_ if lon < self._lon_lo else (_E_ if lon > self._lon_hi else NN)
595 return NN(lat, lon) if lat and lon else (lat or lon)
597 @Property_RO
598 def pgm(self):
599 '''Get the PGM attributes (C{_PGM} or C{None} if not available/applicable).
600 '''
601 return self._pgm
603 @Property_RO
604 def sizeB(self):
605 '''Get the geoid grid file size in bytes (C{int}).
606 '''
607 return self._sizeB
609 @Property_RO
610 def smooth(self):
611 '''Get the C{RectBivariateSpline} smoothing (C{int}).
612 '''
613 return self._smooth
615 @Property_RO
616 def stdev(self):
617 '''Get the standard deviation of this geoid's heights (C{float}) or C{None}.
618 '''
619 if self._stdev is None: # see GeoidKarney
620 self._stdev = float(self.numpy.std(self._hs_y_x))
621 return self._stdev
623 def _swne(self, crop):
624 # crop box to 4-tuple (s, w, n, e)
625 try:
626 if len(crop) == 2:
627 try: # sw, ne LatLons
628 swne = (crop[0].lat, crop[0].lon,
629 crop[1].lat, crop[1].lon)
630 except AttributeError: # (s, w), (n, e)
631 swne = tuple(crop[0]) + tuple(crop[1])
632 else: # (s, w, n, e)
633 swne = crop
634 if len(swne) == 4:
635 s, w, n, e = map(float, swne)
636 if _N_90_0 <= s <= (n - _1_0) <= 89.0 and \
637 _N_180_0 <= w <= (e - _1_0) <= 179.0:
638 return s, w, n, e
639 except (IndexError, TypeError, ValueError):
640 pass
641 raise GeoidError(crop=crop)
643 def toStr(self, prec=3, sep=_COMMASPACE_): # PYCHOK signature
644 '''This geoid and all geoid attributes as a string.
646 @kwarg prec: Number of decimal digits (0..9 or C{None} for
647 default). Trailing zero decimals are stripped
648 for B{C{prec}} values of 1 and above, but kept
649 for negative B{C{prec}} values.
650 @kwarg sep: Separator to join (C{str}).
652 @return: Geoid name and attributes (C{str}).
653 '''
654 s = 1 if self.kind < 0 else 2
655 t = tuple(Fmt.PAREN(m.__name__, fstr(m(), prec=prec)) for m in
656 (self.lowerleft, self.upperright,
657 self.center,
658 self.highest, self.lowest)) + \
659 attrs( _mean_, _stdev_, prec=prec, Nones=False) + \
660 attrs((_kind_, 'smooth')[:s], prec=prec, Nones=False) + \
661 attrs( 'cropped', 'dtype', _endian_, 'hits', 'knots', 'nBytes',
662 'sizeB', _scipy_, _numpy_, prec=prec, Nones=False)
663 return _COLONSPACE_(self, sep.join(t))
665 @Property_RO
666 def u2B(self):
667 '''Get the PGM itemsize in bytes (C{int}).
668 '''
669 return self._u2B
671 @Property_RO
672 def _upperleft(self):
673 '''(INTERNAL) Cache for C{.upperleft}.
674 '''
675 return self._llh3(self._lat_hi, self._lon_lo)
677 def upperleft(self, LatLon=None):
678 '''Return the upper-left location and height of this geoid.
680 @kwarg LatLon: Optional class to return the location and height
681 (C{LatLon}) or C{None}.
683 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon, height)}
684 otherwise a B{C{LatLon}} instance with the lat-, longitude and
685 geoid height of the upper-left, NW grid corner.
686 '''
687 return self._llh3LL(self._upperleft, LatLon)
689 @Property_RO
690 def _upperright(self):
691 '''(INTERNAL) Cache for C{.upperright}.
692 '''
693 return self._llh3(self._lat_hi, self._lon_hi)
695 def upperright(self, LatLon=None):
696 '''Return the upper-right location and height of this geoid.
698 @kwarg LatLon: Optional class to return the location and height
699 (C{LatLon}) or C{None}.
701 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon, height)}
702 otherwise a B{C{LatLon}} instance with the lat-, longitude and
703 geoid height of the upper-right, NE grid corner.
704 '''
705 return self._llh3LL(self._upperright, LatLon)
708class GeoidEGM96(_GeoidBase):
709 '''Geoid height interpolator for the EGM96 U{15 Minute Interpolation Grid<https://earth-info.NGA.mil>}
710 based on C{SciPy} interpolation U{RectBivariateSpline<https://docs.SciPy.org/doc/scipy/reference/
711 generated/scipy.interpolate.RectBivariateSpline.html>}, U{interp2d<https://docs.SciPy.org/doc/scipy/
712 reference/generated/scipy.interpolate.interp2d.html>} or U{bisplrep/-ev<https://docs.scipy.org/doc/
713 scipy/reference/generated/scipy.interpolate.bisplrep.html>}.
715 Use only the C{WW15MGH.GRD} file, unzipped from the EGM96 U{15 Minute Interpolation Grid
716 <https://earth-info.NGA.mil/index.php?dir=wgs84&action=wgs84>} download.
717 '''
718 def __init__(self, EGM96_grd, datum=_WGS84, kind=3, smooth=0, **name_crop):
719 '''New L{GeoidEGM96} interpolator.
721 @arg EGM96_grd: An C{EGM96_grd} grid file name (C{.GRD}).
722 @kwarg datum: Optional grid datum (L{Datum}, L{Ellipsoid}, L{Ellipsoid2} or L{a_f2Tuple}),
723 overriding C{WGS84}.
724 @kwarg kind: C{scipy.interpolate} order (C{int}), use 1..5 for U{RectBivariateSpline
725 <https://docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.
726 RectBivariateSpline.html>} or -1, -3 or -5 for U{bisplrep/-ev<https://
727 docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.bisplrep.html>}
728 or U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/scipy.
729 interpolate.interp2d.html>} C{linear}, C{cubic} respectively C{quintic},
730 see note for more details.
731 @kwarg smooth: Smoothing factor for C{B{kind}=1..5} only (C{int}).
732 @kwarg name_crop: Optional geoid C{B{name}=NN} (C{str}) and UNSUPPORTED keyword argument
733 C{B{crop}=None}.
735 @raise GeoidError: Invalid B{C{crop}}, B{C{kind}} or B{C{smooth}} or a ECM96 grid file
736 B{C{ECM96_grd}} issue.
738 @raise ImportError: Package C{numpy} or C{scipy} not found or not installed.
740 @raise LenError: Grid file B{C{EGM96_grd}} axis mismatch.
742 @raise SciPyError: A C{scipy} issue.
744 @raise SciPyWarning: A C{scipy} warning as exception.
746 @raise TypeError: Invalid B{C{datum}}.
748 @note: Specify C{B{kind}=-1, -3 or -5} to use C{scipy.interpolate.interp2d}
749 before or C{scipy.interpolate.bisplrep/-ev} since C{Scipy} version 1.14.
750 '''
751 crop, name = _xkwds_pop2(name_crop, crop=None)
752 if crop is not None:
753 raise GeoidError(crop=crop, txt_not_=_supported_)
755 g = self._open(EGM96_grd, datum, kind, _name__(**name), smooth)
756 _ = self.numpy # import numpy for .fromfile, .reshape
758 try:
759 p, hs = _Gpars(), self._load(g, sep=_SPACE_) # text
760 p.slat, n, p.wlon, e, p.dlat, p.dlon = hs[:6] # n-s, 0-E
761 p.nlat = int((n - p.slat) / p.dlat) + 1 # include S
762 p.nlon = int((e - p.wlon) / p.dlon) + 1 # include W
763 p.knots = p.nlat * p.nlon # inverted lats N downto S
764 p.glon = _180_0 # Eastern lons 0-360
765 hs = hs[6:].reshape(p.nlat, p.nlon)
766 _GeoidBase.__init__(self, hs, p)
768 except Exception as x:
769 raise _SciPyIssue(x, _in_, repr(EGM96_grd))
770 finally:
771 g.close()
773 def _g2ll2(self, lat, lon):
774 # convert grid (lat, lon) to earth (lat, lon)
775 while lon > _180_0: # Eastern
776 lon -= _360_0
777 return -lat, lon # invert lat
779 def _ll2g2(self, lat, lon):
780 # convert earth (lat, lon) to grid (lat, lon)
781 while lon < 0: # Eastern
782 lon += _360_0
783 return -lat, lon # invert lat
785 if _FOR_DOCS:
786 __call__ = _GeoidBase.__call__
787 height = _GeoidBase.height
790class GeoidG2012B(_GeoidBase):
791 '''Geoid height interpolator for U{GEOID12B Model
792 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/>} grids U{CONUS
793 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_CONUS.shtml>},
794 U{Alaska<https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_AK.shtml>},
795 U{Hawaii<https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_HI.shtml>},
796 U{Guam and Northern Mariana Islands
797 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_GMNI.shtml>},
798 U{Puerto Rico and U.S. Virgin Islands
799 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_PRVI.shtml>} and
800 U{American Samoa<https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_AS.shtml>}
801 based on C{SciPy} interpolation U{RectBivariateSpline<https://docs.SciPy.org/doc/
802 scipy/reference/generated/scipy.interpolate.RectBivariateSpline.html>}, U{interp2d
803 <https://docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.interp2d.html>}
804 or U{bisplrep/-ev<https://docs.scipy.org/doc/scipy/reference/generated/scipy.interpolate.
805 bisplrep.html>}.
807 Use any of the C{le} (little endian) or C{be} (big endian) C{g2012b*.bin} binary grid files.
808 '''
809 _datum = Datums.NAD83
811 def __init__(self, g2012b_bin, datum=Datums.NAD83, kind=3, smooth=0, **name_crop):
812 '''New L{GeoidG2012B} interpolator.
814 @arg g2012b_bin: A C{GEOID12B} grid file name (C{.bin}).
815 @kwarg datum: Optional grid datum (L{Datum}, L{Ellipsoid}, L{Ellipsoid2} or
816 L{a_f2Tuple}), overriding C{NAD83}.
817 @kwarg kind: C{scipy.interpolate} order (C{int}), use 1..5 for U{RectBivariateSpline
818 <https://docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.
819 RectBivariateSpline.html>} or -1, -3 or -5 for U{bisplrep/-ev<https://
820 docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.bisplrep.html>}
821 or U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/scipy.
822 interpolate.interp2d.html>} C{linear}, C{cubic} respectively C{quintic},
823 see note for more details.
824 @kwarg smooth: Smoothing factor for C{B{kind}=1..5} only (C{int}).
825 @kwarg name_crop: Optional geoid C{B{name}=NN} (C{str}) and UNSUPPORTED keyword argument
826 C{B{crop}=None}.
828 @raise GeoidError: Invalid B{C{crop}}, B{C{kind}} or B{C{smooth}} or a G2012B grid file
829 B{C{g2012b_bin}} issue.
831 @raise ImportError: Package C{numpy} or C{scipy} not found or not installed.
833 @raise LenError: Grid file B{C{g2012b_bin}} axis mismatch.
835 @raise SciPyError: A C{scipy} issue.
837 @raise SciPyWarning: A C{scipy} warning as exception.
839 @raise TypeError: Invalid B{C{datum}}.
841 @note: Specify C{B{kind}=-1, -3 or -5} to use C{scipy.interpolate.interp2d}
842 before or C{scipy.interpolate.bisplrep/-ev} since C{Scipy} version 1.14.
843 '''
844 crop, name = _xkwds_pop2(name_crop, crop=None)
845 if crop is not None:
846 raise GeoidError(crop=crop, txt_not_=_supported_)
848 g = self._open(g2012b_bin, datum, kind, _name__(**name), smooth)
849 _ = self.numpy # import numpy for ._load and
851 try:
852 p = _Gpars()
853 n = (self.sizeB // 4) - 11 # number of f4 heights
854 # U{numpy dtype formats are different from Python struct formats
855 # <https://docs.SciPy.org/doc/numpy-1.15.0/reference/arrays.dtypes.html>}
856 for en_ in ('<', '>'):
857 # skip 4xf8, get 3xi4
858 p.nlat, p.nlon, ien = map(int, self._load(g, en_+'i4', 3, 32))
859 if ien == 1: # correct endian
860 p.knots = p.nlat * p.nlon
861 if p.knots == n and 1 < p.nlat < n \
862 and 1 < p.nlon < n:
863 self._endian = en_+'f4'
864 break
865 else: # couldn't validate endian
866 raise GeoidError(_endian_)
868 # get the first 4xf8
869 p.slat, p.wlon, p.dlat, p.dlon = map(float, self._load(g, en_+'f8', 4))
870 # read all f4 heights, ignoring the first 4xf8 and 3xi4
871 hs = self._load(g, self._endian, n, 44).reshape(p.nlat, p.nlon)
872 p.wlon -= _360_0 # western-most East longitude to earth (..., lon)
873 _GeoidBase.__init__(self, hs, p)
875 except Exception as x:
876 raise _SciPyIssue(x, _in_, repr(g2012b_bin))
877 finally:
878 g.close()
880 def _g2ll2(self, lat, lon):
881 # convert grid (lat, lon) to earth (lat, lon)
882 return lat, lon
884 def _ll2g2(self, lat, lon):
885 # convert earth (lat, lon) to grid (lat, lon)
886 return lat, lon
888 if _FOR_DOCS:
889 __call__ = _GeoidBase.__call__
890 height = _GeoidBase.height
893class GeoidHeight5Tuple(_NamedTuple): # .geoids.py
894 '''5-Tuple C{(lat, lon, egm84, egm96, egm2008)} for U{GeoidHeights.dat
895 <https://SourceForge.net/projects/geographiclib/files/testdata/>}
896 tests with the heights for 3 different EGM grids at C{degrees90}
897 and C{degrees180} degrees (after converting C{lon} from original
898 C{0 <= EasterLon <= 360}).
899 '''
900 _Names_ = (_lat_, _lon_, 'egm84', 'egm96', 'egm2008')
901 _Units_ = ( Lat, Lon, Height, Height, Height)
904def _I(i):
905 '''(INTERNAL) Cache a single C{int} constant.
906 '''
907 i = int(i)
908 return _intCs.setdefault(i, i) # PYCHOK undefined due to del _intCs
911def _T(cs):
912 '''(INTERNAL) Cache a tuple of single C{int} constants.
913 '''
914 return tuple(map(_I, _splituple(cs)))
916_T0s12 = (_I(0),) * 12 # PYCHOK _T('0, 0, ..., 0')
919class GeoidKarney(_GeoidBase):
920 '''Geoid height interpolator for I{Karney}'s U{GeographicLib Earth
921 Gravitational Model (EGM)<https://GeographicLib.SourceForge.io/C++/doc/
922 geoid.html>} geoid U{egm*.pgm<https://GeographicLib.SourceForge.io/
923 C++/doc/geoid.html#geoidinst>} datasets using bilinear or U{cubic
924 <https://dl.ACM.org/citation.cfm?id=368443>} interpolation and U{caching
925 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidcache>}
926 in pure Python, transcoded from I{Karney}'s U{C++ class Geoid
927 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinterp>}.
929 Use any of the geoid U{egm84-, egm96- or egm2008-*.pgm
930 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinst>}
931 datasets.
932 '''
933 _C0 = _F(372), _F(240), _F(372) # n, _ and s common denominators
934 # matrices c3n_, c3, c3s_, transposed from GeographicLib/Geoid.cpp
935 _C3 = ((_T('0, 0, 62, 124, 124, 62, 0, 0, 0, 0, 0, 0'),
936 _T0s12,
937 _T('-131, 7, -31, -62, -62, -31, 45, 216, 156, -45, -55, -7'),
938 _T0s12,
939 _T('138, -138, 0, 0, 0, 0, -183, 33, 153, -3, 48, -48'), # PYCHOK indent
940 _T('144, 42, -62, -124, -124, -62, -9, 87, 99, 9, 42, -42'),
941 _T0s12,
942 _T('0, 0, 0, 0, 0, 0, 93, -93, -93, 93, 0, 0'),
943 _T('-102, 102, 0, 0, 0, 0, 18, 12, -12, -18, -84, 84'),
944 _T('-31, -31, 31, 62, 62, 31, 0, -93, -93, 0, 31, 31')), # PYCHOK indent
946 (_T('9, -9, 9, 186, 54, -9, -9, 54, -54, 9, -9, 9'),
947 _T('-18, 18, -88, -42, 162, -32, 8, -78, 78, -8, 18, -18'),
948 _T('-88, 8, -18, -42, -78, 18, 18, 162, 78, -18, -32, -8'),
949 _T('0, 0, 90, -150, 30, 30, 30, -90, 90, -30, 0, 0'),
950 _T('96, -96, 96, -96, -24, 24, -96, -24, 144, -24, 24, -24'), # PYCHOK indent
951 _T('90, 30, 0, -150, -90, 0, 0, 30, 90, 0, 30, -30'),
952 _T('0, 0, -20, 60, -60, 20, -20, 60, -60, 20, 0, 0'),
953 _T('0, 0, -60, 60, 60, -60, 60, -60, -60, 60, 0, 0'),
954 _T('-60, 60, 0, 60, -60, 0, 0, 60, -60, 0, -60, 60'),
955 _T('-20, -20, 0, 60, 60, 0, 0, -60, -60, 0, 20, 20')),
957 (_T('18, -18, 36, 210, 162, -36, 0, 0, 0, 0, -18, 18'), # PYCHOK indent
958 _T('-36, 36, -165, 45, 141, -21, 0, 0, 0, 0, 36, -36'),
959 _T('-122, -2, -27, -111, -75, 27, 62, 124, 124, 62, -64, 2'),
960 _T('0, 0, 93, -93, -93, 93, 0, 0, 0, 0, 0, 0'),
961 _T('120, -120, 147, -57, -129, 39, 0, 0, 0, 0, 66, -66'), # PYCHOK indent
962 _T('135, 51, -9, -192, -180, 9, 31, 62, 62, 31, 51, -51'),
963 _T0s12,
964 _T('0, 0, -93, 93, 93, -93, 0, 0, 0, 0, 0, 0'),
965 _T('-84, 84, 18, 12, -12, -18, 0, 0, 0, 0, -102, 102'),
966 _T('-31, -31, 0, 93, 93, 0, -31, -62, -62, -31, 31, 31')))
968 _BT = (_T('0, 0'), # bilinear 4-tuple [i, j] indices
969 _T('1, 0'),
970 _T('0, 1'),
971 _T('1, 1'))
973 _CM = (_T(' 0, -1'), # 10x12 cubic matrix [i, j] indices
974 _T(' 1, -1'),
975 _T('-1, 0'),
976 _T(' 0, 0'),
977 _T(' 1, 0'),
978 _T(' 2, 0'),
979 _T('-1, 1'),
980 _T(' 0, 1'),
981 _T(' 1, 1'),
982 _T(' 2, 1'),
983 _T(' 0, 2'),
984 _T(' 1, 2'))
986# _cropped = None
987 _endian = '>H' # struct.unpack 1 ushort (big endian, unsigned short)
988 _4endian = '>4H' # struct.unpack 4 ushorts
989 _Rendian = NN # struct.unpack a row of ushorts
990# _highest = (-8.4, 147.367, 85.839) if egm2008-1.pgm else (
991# (-8.167, 147.25, 85.422) if egm96-5.pgm else
992# (-4.5, 148.75, 81.33)) # egm84-15.pgm
993 _iscipy = False
994# _lowest = (4.7, 78.767, -106.911) if egm2008-1.pgm else (
995# (4.667, 78.833, -107.043) if egm96-5.pgm else
996# (4.75, 79.25, -107.34)) # egm84-15.pgm
997 _mean = _F(-1.317) # from egm2008-1, -1.438 egm96-5, -0.855 egm84-15
998 _nBytes = None # not applicable
999 _nterms = len(_C3[0]) # columns length, number of rows
1000 _smooth = None # not applicable
1001 _stdev = _F(29.244) # from egm2008-1, 29.227 egm96-5, 29.183 egm84-15
1002 _u2B = _calcsize(_endian) # pixelsize_ in bytes
1003 _4u2B = _calcsize(_4endian) # 4 pixelsize_s in bytes
1004 _Ru2B = 0 # row of pixelsize_s in bytes
1005 _yx_hits = 0 # cache hits
1006 _yx_i = () # cached (y, x) indices
1007 _yx_t = () # cached 4- or 10-tuple for _ev2k resp. _ev3k
1009 def __init__(self, egm_pgm, crop=None, datum=_WGS84, kind=3, **name_smooth):
1010 '''New L{GeoidKarney} interpolator.
1012 @arg egm_pgm: An U{EGM geoid dataset<https://GeographicLib.SourceForge.io/
1013 C++/doc/geoid.html#geoidinst>} file name (C{egm*.pgm}), see
1014 note below.
1015 @kwarg crop: Optional box to limit geoid locations, a 4-tuple (C{south,
1016 west, north, east}), 2-tuple (C{(south, west), (north, east)})
1017 with 2 C{degrees90} lat- and C{degrees180} longitudes or as
1018 2-tuple (C{LatLonSW, LatLonNE}) of C{LatLon} instances.
1019 @kwarg datum: Optional grid datum (C{Datum}, L{Ellipsoid}, L{Ellipsoid2} or
1020 L{a_f2Tuple}), overriding C{WGS84}.
1021 @kwarg kind: Interpolation order (C{int}), 2 for C{bilinear} or 3 for C{cubic}.
1022 @kwarg name_smooth: Optional geoid C{B{name}=NN} (C{str}) and UNSUPPORTED
1023 keyword argument C{B{smooth}}, use C{B{smooth}=None} to ignore.
1025 @raise GeoidError: EGM dataset B{C{egm_pgm}} issue or invalid B{C{crop}},
1026 B{C{kind}} or B{C{smooth}}.
1028 @raise TypeError: Invalid B{C{datum}}.
1030 @see: Class L{GeoidPGM} and function L{egmGeoidHeights}.
1032 @note: Geoid file B{C{egm_pgm}} remains open and I{must be closed} by calling
1033 method C{close} or by using C{with B{GeoidKarney}(...) as ...:} context.
1034 '''
1035 smooth, name = _xkwds_pop2(name_smooth, smooth=None)
1036 if smooth is not None:
1037 raise GeoidError(smooth=smooth, txt_not_=_supported_)
1039 if kind in (2,):
1040 self._ev2d = self._ev2k # see ._ev_name
1041 elif kind not in (3,):
1042 raise GeoidError(kind=kind)
1044 self._egm = g = self._open(egm_pgm, datum, kind, _name__(**name), None)
1045 self._pgm = p = _PGM(g, pgm=egm_pgm, itemsize=self.u2B, sizeB=self.sizeB)
1047 self._Rendian = self._4endian.replace(_4_, str(p.nlon))
1048 self._Ru2B = _calcsize(self._Rendian)
1050 self._knots = p.knots # grid knots
1051 self._lon_of = float(p.flon) # forward offset
1052 self._lon_og = float(p.glon) # reverse offset
1053 # set earth (lat, lon) limits (s, w, n, e)
1054 self._lat_lo, self._lon_lo, \
1055 self._lat_hi, self._lon_hi = self._swne(crop if crop else p.crop4)
1056 self._cropped = bool(crop)
1058 def _c0c3v(self, y, x):
1059 # get the common denominator, the 10x12 cubic matrix and
1060 # the 12 cubic v-coefficients around geoid index (y, x)
1061 p = self._pgm
1062 if 0 < x < (p.nlon - 2) and 0 < y < (p.nlat - 2):
1063 # read 4x4 ushorts, drop the 4 corners
1064 S = _MODS.os.SEEK_SET
1065 e = self._4endian
1066 g = self._egm
1067 n = self._4u2B
1068 R = self._Ru2B
1069 b = self._seek(y - 1, x - 1)
1070 v = _unpack(e, g.read(n))[1:3]
1071 b += R
1072 g.seek(b, S)
1073 v += _unpack(e, g.read(n))
1074 b += R
1075 g.seek(b, S)
1076 v += _unpack(e, g.read(n))
1077 b += R
1078 g.seek(b, S)
1079 v += _unpack(e, g.read(n))[1:3]
1080 j = 1
1082 else: # likely some wrapped y and/or x's
1083 v = self._raws(y, x, GeoidKarney._CM)
1084 j = 0 if y < 1 else (1 if y < (p.nlat - 2) else 2)
1086 return GeoidKarney._C0[j], GeoidKarney._C3[j], v
1088 @Property_RO
1089 def dtype(self):
1090 '''Get the geoid's grid data type (C{str}).
1091 '''
1092 return 'ushort'
1094 def _ev(self, lat, lon): # PYCHOK expected
1095 # interpolate the geoid height at grid (lat, lon)
1096 fy, fx = self._g2yx2(lat, lon)
1097 y, x = int(_floor(fy)), int(_floor(fx))
1098 fy -= y
1099 fx -= x
1100 H = self._ev2d(fy, fx, y, x) # PYCHOK ._ev2k or ._ev3k
1101 H *= self._pgm.Scale # H.fmul(self._pgm.Scale)
1102 H += self._pgm.Offset # H.fadd(self._pgm.Offset)
1103 return H.fsum() # float(H)
1105 def _ev2k(self, fy, fx, *yx):
1106 # compute the bilinear 4-tuple and interpolate raw H
1107 if self._yx_i == yx:
1108 self._yx_hits += 1
1109 else:
1110 y, x = self._yx_i = yx
1111 self._yx_t = self._raws(y, x, GeoidKarney._BT)
1112 t = self._yx_t
1113 v = _1_0, -fx, fx
1114 H = Fdot(v, t[0], t[0], t[1]).fmul(_1_0 - fy) # c = a * (1 - fy)
1115 H += Fdot(v, t[2], t[2], t[3]).fmul(fy) # c += b * fy
1116 return H
1118 def _ev3k(self, fy, fx, *yx):
1119 # compute the cubic 10-tuple and interpolate raw H
1120 if self._yx_i == yx:
1121 self._yx_hits += 1
1122 else:
1123 c0, c3, v = self._c0c3v(*yx)
1124 # assert len(c3) == self._nterms
1125 self._yx_t = tuple(fdot(v, *r3) / c0 for r3 in c3)
1126 self._yx_i = yx
1127 # GeographicLib/Geoid.cpp Geoid::height(lat, lon) ...
1128 # real h = t[0] + fx * (t[1] + fx * (t[3] + fx * t[6])) +
1129 # fy * (t[2] + fx * (t[4] + fx * t[7]) +
1130 # fy * (t[5] + fx * t[8] + fy * t[9]));
1131 t = self._yx_t
1132 v = _1_0, fx, fy
1133 H = Fdot(v, t[5], t[8], t[9])
1134 H *= fy
1135 H += Fhorner(fx, t[2], t[4], t[7])
1136 H *= fy
1137 H += Fhorner(fx, t[0], t[1], t[3], t[6])
1138 return H
1140 _ev2d = _ev3k # overriden for kind=2, see ._ev_name
1142 def _g2ll2(self, lat, lon):
1143 # convert grid (lat, lon) to earth (lat, lon), uncropped
1144 while lon > _180_0:
1145 lon -= _360_0
1146 return lat, lon
1148 def _g2yx2(self, lat, lon):
1149 # convert grid (lat, lon) to grid (y, x) indices
1150 p = self._pgm
1151 # note, slat = +90, rlat < 0 makes y >=0
1152 return ((lat - p.slat) * p.rlat), ((lon - p.wlon) * p.rlon)
1154 def _gyx2g2(self, y, x):
1155 # convert grid (y, x) indices to grid (lat, lon)
1156 p = self._pgm
1157 return (p.slat + p.dlat * y), (p.wlon + p.dlon * x)
1159 @Property_RO
1160 def _highest_ltd(self):
1161 '''(INTERNAL) Cache for C{.highest}.
1162 '''
1163 return self._LL3T(self._llh3minmax(True, -12, -4), name__=self.highest)
1165 def highest(self, LatLon=None, full=False): # PYCHOK full
1166 '''Return the location and largest height of this geoid.
1168 @kwarg LatLon: Optional class to return the location and height
1169 (C{LatLon}) or C{None}.
1170 @kwarg full: Search the full or limited latitude range (C{bool}).
1172 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon,
1173 height)} otherwise a B{C{LatLon}} instance with the lat-,
1174 longitude and geoid height of the highest grid location.
1175 '''
1176 llh = self._highest if full or self.cropped else self._highest_ltd
1177 return self._llh3LL(llh, LatLon)
1179 def _lat2y2(self, lat2):
1180 # convert earth lat(s) to min and max grid y indices
1181 ys, m = [], self._pgm.nlat - 1
1182 for lat in lat2:
1183 y, _ = self._g2yx2(*self._ll2g2(lat, 0))
1184 ys.append(max(min(int(y), m), 0))
1185 return min(ys), max(ys) + 1
1187 def _ll2g2(self, lat, lon):
1188 # convert earth (lat, lon) to grid (lat, lon), uncropped
1189 while lon < 0:
1190 lon += _360_0
1191 return lat, lon
1193 def _llh3minmax(self, highest, *lat2):
1194 # find highest or lowest, takes 10+ secs for egm2008-1.pgm geoid
1195 # (Python 2.7.16, macOS 10.13.6 High Sierra, iMac 3 GHz Core i3)
1196 if highest:
1197 def _mt(r, h):
1198 m = max(r)
1199 return m, (m > h)
1201 else: # lowest
1202 def _mt(r, h): # PYCHOK redef
1203 m = min(r)
1204 return m, (m < h)
1206 y = x = 0
1207 h = self._raw(y, x)
1208 for j, r in self._raw2(*lat2):
1209 m, t = _mt(r, h)
1210 if t:
1211 h, y, x = m, j, r.index(m)
1212 h *= self._pgm.Scale
1213 h += self._pgm.Offset
1214 return self._g2ll2(*self._gyx2g2(y, x)) + (h,)
1216 @Property_RO
1217 def _lowest_ltd(self):
1218 '''(INTERNAL) Cache for C{.lowest}.
1219 '''
1220 return self._LL3T(self._llh3minmax(False, 0, 8), name__=self.lowest)
1222 def lowest(self, LatLon=None, full=False): # PYCHOK full
1223 '''Return the location and lowest height of this geoid.
1225 @kwarg LatLon: Optional class to return the location and height
1226 (C{LatLon}) or C{None}.
1227 @kwarg full: Search the full or limited latitude range (C{bool}).
1229 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon,
1230 height)} otherwise a B{C{LatLon}} instance with the lat-,
1231 longitude and geoid height of the lowest grid location.
1232 '''
1233 llh = self._lowest if full or self.cropped else self._lowest_ltd
1234 return self._llh3LL(llh, LatLon)
1236 def _raw(self, y, x):
1237 # get the ushort geoid height at geoid index (y, x),
1238 # like GeographicLib/Geoid.hpp real rawval(is, iy)
1239 p = self._pgm
1240 if x < 0:
1241 x += p.nlon
1242 elif x >= p.nlon:
1243 x -= p.nlon
1244 h = p.nlon // 2
1245 if y < 0:
1246 y = -y
1247 elif y >= p.nlat:
1248 y = (p.nlat - 1) * 2 - y
1249 else:
1250 h = 0
1251 x += h if x < h else -h
1252 self._seek(y, x)
1253 h = _unpack(self._endian, self._egm.read(self._u2B))
1254 return h[0]
1256 def _raws(self, y, x, ijs):
1257 # get bilinear 4-tuple or 10x12 cubic matrix
1258 return tuple(self._raw(y + j, x + i) for i, j in ijs)
1260 def _raw2(self, *lat2):
1261 # yield a 2-tuple (y, ushorts) for each row or for
1262 # the rows between two (or more) earth lat values
1263 p = self._pgm
1264 g = self._egm
1265 e = self._Rendian
1266 n = self._Ru2B
1267 # min(lat2) <= lat <= max(lat2) or 0 <= y < p.nlat
1268 s, t = self._lat2y2(lat2) if lat2 else (0, p.nlat)
1269 self._seek(s, 0) # to start of row s
1270 for y in range(s, t):
1271 yield y, _unpack(e, g.read(n))
1273 def _seek(self, y, x):
1274 # position geoid to grid index (y, x)
1275 p, g = self._pgm, self._egm
1276 if g:
1277 b = p.skip + (y * p.nlon + x) * self._u2B
1278 g.seek(b, _MODS.os.SEEK_SET)
1279 return b # position
1280 raise GeoidError('closed file', txt=repr(p.egm)) # IOError
1283class GeoidPGM(_GeoidBase):
1284 '''Geoid height interpolator for I{Karney}'s U{GeographicLib Earth
1285 Gravitational Model (EGM)<https://GeographicLib.SourceForge.io/C++/doc/geoid.html>}
1286 geoid U{egm*.pgm<https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinst>}
1287 datasets but based on C{SciPy} U{RectBivariateSpline<https://docs.SciPy.org/doc/scipy/
1288 reference/generated/scipy.interpolate.RectBivariateSpline.html>}, U{bisplrep/-ev
1289 <https://docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.bisplrep.html>}
1290 or U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.
1291 interp2d.html>} interpolation.
1293 Use any of the U{egm84-, egm96- or egm2008-*.pgm <https://GeographicLib.SourceForge.io/
1294 C++/doc/geoid.html#geoidinst>} datasets. However, unless cropped, an entire C{egm*.pgm}
1295 dataset is loaded into the C{SciPy} interpolator and converted from 2-byte C{int} to
1296 8-byte C{dtype float64}. Therefore, internal memory usage is 4x the U{egm*.pgm
1297 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinst>} file size and may
1298 exceed the available memory, especially with 32-bit Python, see properties C{.nBytes}
1299 and C{.sizeB}.
1300 '''
1301 _cropped = False
1302 _endian = '>u2'
1304 def __init__(self, egm_pgm, crop=None, datum=_WGS84, kind=3, smooth=0, **name):
1305 '''New L{GeoidPGM} interpolator.
1307 @arg egm_pgm: An U{EGM geoid dataset<https://GeographicLib.SourceForge.io/
1308 C++/doc/geoid.html#geoidinst>} file name (C{egm*.pgm}).
1309 @kwarg crop: Optional box to crop B{C{egm_pgm}}, a 4-tuple (C{south, west,
1310 north, east}) or 2-tuple (C{(south, west), (north, east)}),
1311 in C{degrees90} lat- and C{degrees180} longitudes or a 2-tuple
1312 (C{LatLonSW, LatLonNE}) of C{LatLon} instances.
1313 @kwarg datum: Optional grid datum (L{Datum}, L{Ellipsoid}, L{Ellipsoid2} or
1314 L{a_f2Tuple}), overriding C{WGS84}.
1315 @kwarg kind: C{scipy.interpolate} order (C{int}), use 1..5 for U{RectBivariateSpline
1316 <https://docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.
1317 RectBivariateSpline.html>} or -1, -3 or -5 for U{bisplrep/-ev<https://
1318 docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.bisplrep.html>}
1319 or U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/scipy.
1320 interpolate.interp2d.html>} C{linear}, C{cubic} respectively C{quintic},
1321 see note for more details.
1322 @kwarg smooth: Smoothing factor for C{B{kind}=1..5} only (C{int}).
1323 @kwarg name: Optional geoid C{B{name}=NN} (C{str}).
1325 @raise GeoidError: EGM dataset B{C{egm_pgm}} issue or invalid B{C{crop}}, B{C{kind}}
1326 or B{C{smooth}}.
1328 @raise ImportError: Package C{numpy} or C{scipy} not found or not installed.
1330 @raise LenError: EGM dataset B{C{egm_pgm}} axis mismatch.
1332 @raise SciPyError: A C{scipy} issue.
1334 @raise SciPyWarning: A C{scipy} warning as exception.
1336 @raise TypeError: Invalid B{C{datum}} or unexpected argument.
1338 @note: Specify C{B{kind}=-1, -3 or -5} to use C{scipy.interpolate.interp2d}
1339 before or C{scipy.interpolate.bisplrep/-ev} since C{Scipy} version 1.14.
1341 @note: The U{GeographicLib egm*.pgm<https://GeographicLib.SourceForge.io/C++/doc/
1342 geoid.html#geoidinst>} file sizes are based on a 2-byte C{int} height
1343 converted to 8-byte C{dtype float64} for C{scipy} interpolators. Therefore,
1344 internal memory usage is 4 times the C{egm*.pgm} file size and may exceed
1345 the available memory, especially with 32-bit Python. To reduce memory
1346 usage, set keyword argument B{C{crop}} to the region of interest. For
1347 example C{B{crop}=(20, -125, 50, -65)} covers the U{conterminous US
1348 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/maps/GEOID12B_CONUS_grids.png>}
1349 (CONUS), less than 3% of the entire C{egm2008-1.pgm} dataset.
1351 @see: Class L{GeoidKarney} and function L{egmGeoidHeights}.
1352 '''
1353 np = self.numpy
1354 self._u2B = np.dtype(self.endian).itemsize
1356 g = self._open(egm_pgm, datum, kind, _name__(**name), smooth)
1357 self._pgm = p = _PGM(g, pgm=egm_pgm, itemsize=self.u2B, sizeB=self.sizeB)
1358 if crop:
1359 g = p._cropped(g, abs(kind) + 1, *self._swne(crop))
1360 if _MODS.internals._version2(np.__version__) < (1, 9):
1361 g = open(g.name, _rb_) # reopen tempfile for numpy 1.8.0-
1362 self._cropped = True
1363 try:
1364 # U{numpy dtype formats are different from Python struct formats
1365 # <https://docs.SciPy.org/doc/numpy-1.15.0/reference/arrays.dtypes.html>}
1366 # read all heights, skipping the PGM header lines, converted to float
1367 hs = self._load(g, self.endian, p.knots, p.skip).reshape(p.nlat, p.nlon) * p.Scale
1368 if p.Offset: # offset
1369 hs = p.Offset + hs
1370 if p.dlat < 0: # flip the rows
1371 hs = np.flipud(hs)
1372 _GeoidBase.__init__(self, hs, p)
1373 except Exception as x:
1374 raise _SciPyIssue(x, _in_, repr(egm_pgm))
1375 finally:
1376 g.close()
1378 def _g2ll2(self, lat, lon):
1379 # convert grid (lat, lon) to earth (lat, lon), un-/cropped
1380 if self._cropped:
1381 lon -= self._lon_of
1382 else:
1383 while lon > _180_0:
1384 lon -= _360_0
1385 return lat, lon
1387 def _ll2g2(self, lat, lon):
1388 # convert earth (lat, lon) to grid (lat, lon), un-/cropped
1389 if self._cropped:
1390 lon += self._lon_of
1391 else:
1392 while lon < 0:
1393 lon += _360_0
1394 return lat, lon
1396 if _FOR_DOCS:
1397 __call__ = _GeoidBase.__call__
1398 height = _GeoidBase.height
1401class _Gpars(_Named):
1402 '''(INTERNAL) Basic geoid parameters.
1403 '''
1404 # interpolator parameters
1405 dlat = 0 # +/- latitude resolution in C{degrees}
1406 dlon = 0 # longitude resolution in C{degrees}
1407 nlat = 1 # number of latitude knots (C{int})
1408 nlon = 0 # number of longitude knots (C{int})
1409 rlat = 0 # +/- latitude resolution in C{float}, 1 / .dlat
1410 rlon = 0 # longitude resolution in C{float}, 1 / .dlon
1411 slat = 0 # nothern- or southern most latitude (C{degrees90})
1412 wlon = 0 # western-most longitude in Eastern lon (C{degrees360})
1414 flon = 0 # forward, earth to grid longitude offset
1415 glon = 0 # reverse, grid to earth longitude offset
1417 knots = 0 # number of knots, nlat * nlon (C{int})
1418 skip = 0 # header bytes to skip (C{int})
1420 def __repr__(self):
1421 t = _COMMASPACE_.join(pairs((a, getattr(self, a)) for
1422 a in dir(self.__class__)
1423 if a[:1].isupper()))
1424 return _COLONSPACE_(self, t)
1426 def __str__(self):
1427 return Fmt.PAREN(self.classname, repr(self.name))
1430class _PGM(_Gpars):
1431 '''(INTERNAL) Parse an C{egm*.pgm} geoid dataset file.
1433 # Geoid file in PGM format for the GeographicLib::Geoid class
1434 # Description WGS84 EGM96, 5-minute grid
1435 # URL https://Earth-Info.NGA.mil/GandG/wgs84/gravitymod/egm96/egm96.html
1436 # DateTime 2009-08-29 18:45:03
1437 # MaxBilinearError 0.140
1438 # RMSBilinearError 0.005
1439 # MaxCubicError 0.003
1440 # RMSCubicError 0.001
1441 # Offset -108
1442 # Scale 0.003
1443 # Origin 90N 0E
1444 # AREA_OR_POINT Point
1445 # Vertical_Datum WGS84
1446 <width> <height>
1447 <pixel>
1448 ...
1449 '''
1450 crop4 = () # 4-tuple (C{south, west, north, east}).
1451 egm = None
1452 glon = 180 # reverse offset, uncropped
1453# pgm = NN # name
1454 sizeB = 0
1455 u2B = 2 # item size of grid height (C{int}).
1457 @staticmethod
1458 def _llstr2floats(latlon):
1459 # llstr to (lat, lon) floats
1460 lat, lon = latlon.split()
1461 return _MODS.dms.parseDMS2(lat, lon)
1463 # PGM file attributes, CamelCase but not .istitle()
1464 AREA_OR_POINT = str
1465 DateTime = str
1466 Description = str # 'WGS84 EGM96, 5-minute grid'
1467 Geoid = str # 'file in PGM format for the GeographicLib::Geoid class'
1468 MaxBilinearError = float
1469 MaxCubicError = float
1470 Offset = float
1471 Origin = _llstr2floats
1472 Pixel = 0
1473 RMSBilinearError = float
1474 RMSCubicError = float
1475 Scale = float
1476 URL = str # 'https://Earth-Info.NGA.mil/GandG/wgs84/...'
1477 Vertical_Datum = str
1479 def __init__(self, g, pgm=NN, itemsize=0, sizeB=0): # MCCABE 22
1480 '''(INTERNAL) New C{_PGM} parsed C{egm*.pgm} geoid dataset.
1481 '''
1482 self.name = pgm # geoid file name
1483 if itemsize:
1484 self._u2B = itemsize
1485 if sizeB:
1486 self.sizeB = sizeB
1488 t = g.readline() # make sure newline == '\n'
1489 if t != b'P5\n' and t.strip() != b'P5':
1490 raise self._Errorf(_format_, _header_, t)
1492 while True: # read all # Attr ... lines,
1493 try: # ignore empty ones or comments
1494 t = g.readline().strip()
1495 if t.startswith(_bHASH_):
1496 t = t.lstrip(_bHASH_).lstrip()
1497 a, v = map(_ub2str, t.split(None, 1))
1498 f = getattr(_PGM, a, None)
1499 if callable(f) and a[:1].isupper():
1500 setattr(self, a, f(v))
1501 elif t:
1502 break
1503 except (TypeError, ValueError):
1504 raise self._Errorf(_format_, 'Attr', t)
1505 else: # should never get here
1506 raise self._Errorf(_format_, _header_, g.tell())
1508 try: # must be (even) width and (odd) height
1509 nlon, nlat = map(int, t.split())
1510 if nlon < 2 or nlon > (360 * 60) or isodd(nlon) or \
1511 nlat < 2 or nlat > (181 * 60) or not isodd(nlat):
1512 raise ValueError
1513 except (TypeError, ValueError):
1514 raise self._Errorf(_format_, _SPACE_(_width_, _height_), t)
1516 try: # must be 16 bit pixel height
1517 t = g.readline().strip()
1518 self.Pixel = int(t)
1519 if not 255 < self.Pixel < 65536: # >u2 or >H only
1520 raise ValueError
1521 except (TypeError, ValueError):
1522 raise self._Errorf(_format_, 'pixel', t)
1524 for a in dir(_PGM): # set undefined # Attr ... to None
1525 if a[:1].isupper() and callable(getattr(self, a)):
1526 setattr(self, a, None)
1528 if self.Origin is None:
1529 raise self._Errorf(_format_, 'Origin', self.Origin)
1530 if self.Offset is None or self.Offset > 0:
1531 raise self._Errorf(_format_, 'Offset', self.Offset)
1532 if self.Scale is None or self.Scale < EPS:
1533 raise self._Errorf(_format_, 'Scale', self.Scale)
1535 self.skip = g.tell()
1536 self.knots = nlat * nlon
1538 self.nlat, self.nlon = nlat, nlon
1539 self.slat, self.wlon = self.Origin
1540 # note, negative .dlat and .rlat since rows
1541 # are from .slat 90N down in decreasing lat
1542 self.dlat, self.dlon = _180_0 / (1 - nlat), _360_0 / nlon
1543 self.rlat, self.rlon = (1 - nlat) / _180_0, nlon / _360_0
1545 # grid corners in earth (lat, lon), .slat = 90, .dlat < 0
1546 n = float(self.slat)
1547 s = n + self.dlat * (nlat - 1)
1548 w = self.wlon - self.glon
1549 e = w + self.dlon * nlon
1550 self.crop4 = s, w, n, e
1552 n = self.sizeB - self.skip
1553 if n > 0 and n != (self.knots * self.u2B):
1554 raise self._Errorf('%s(%s x %s != %s)', _assert_, nlat, nlon, n)
1556 def _cropped(self, g, k1, south, west, north, east): # MCCABE 15
1557 '''Crop the geoid to (south, west, north, east) box.
1558 '''
1559 # flon offset for both west and east
1560 f = 360 if west < 0 else 0
1561 # earth (lat, lon) to grid indices (y, x),
1562 # note y is decreasing, i.e. n < s
1563 s, w = self._lle2yx2(south, west, f)
1564 n, e = self._lle2yx2(north, east, f)
1565 s += 1 # s > n
1566 e += 1 # e > w
1568 hi, wi = self.nlat, self.nlon
1569 # handle special cases
1570 if (s - n) > hi:
1571 n, s = 0, hi # entire lat range
1572 if (e - w) > wi:
1573 w, e, f = 0, wi, 180 # entire lon range
1574 if s == hi and w == n == 0 and e == wi:
1575 return g # use entire geoid as-is
1577 if (e - w) < k1 or (s - n) < (k1 + 1):
1578 raise self._Errorf(_format_, 'swne', (north - south, east - west))
1580 if e > wi > w: # wrap around
1581 # read w..wi and 0..e
1582 r, p = (wi - w), (e - wi)
1583 elif e > w:
1584 r, p = (e - w), 0
1585 else:
1586 raise self._Errorf('%s(%s < %s)', _assert_, w, e)
1588 # convert to bytes
1589 r *= self.u2B
1590 p *= self.u2B
1591 q = wi * self.u2B # stride
1592 # number of rows and cols to skip from
1593 # the original (.slat, .wlon) origin
1594 z = self.skip + (n * wi + w) * self.u2B
1595 # sanity check
1596 if r < 2 or p < 0 or q < 2 or z < self.skip \
1597 or z > self.sizeB:
1598 raise self._Errorf(_format_, _assert_, (r, p, q, z))
1600 # can't use _BytesIO since numpy
1601 # needs .fileno attr in .fromfile
1602 t, c = 0, self._tmpfile()
1603 # reading (s - n) rows, forward
1604 for y in range(n, s): # PYCHOK y unused
1605 g.seek(z, _MODS.os.SEEK_SET)
1606 # Python 2 tmpfile.write returns None
1607 t += c.write(g.read(r)) or r
1608 if p: # wrap around to start of row
1609 g.seek(-q, _MODS.os.SEEK_CUR)
1610 # assert(g.tell() == (z - w * self.u2B))
1611 # Python 2 tmpfile.write returns None
1612 t += c.write(g.read(p)) or p
1613 z += q
1614 c.flush()
1615 g.close()
1617 s -= n # nlat
1618 e -= w # nlon
1619 k = s * e # knots
1620 z = k * self.u2B
1621 if t != z:
1622 raise self._Errorf('%s(%s != %s) %s', _assert_, t, z, self)
1624 # update the _Gpars accordingly, note attributes
1625 # .dlat, .dlon, .rlat and .rlon remain unchanged
1626 self.slat += n * self.dlat
1627 self.wlon += w * self.dlon
1628 self.nlat = s
1629 self.nlon = e
1630 self.flon = self.glon = f
1632 self.crop4 = south, west, north, east
1633 self.knots = k
1634 self.skip = 0 # no header lines in c
1636 c.seek(0, _MODS.os.SEEK_SET)
1637 # c = open(c.name, _rb_) # reopen for numpy 1.8.0-
1638 return c
1640 def _Errorf(self, fmt, *args): # PYCHOK no cover
1641 t = fmt % args
1642 e = self.pgm or NN
1643 if e:
1644 t = _SPACE_(t, _in_, repr(e))
1645 return PGMError(t)
1647 def _lle2yx2(self, lat, lon, flon):
1648 # earth (lat, lon) to grid indices (y, x)
1649 # with .dlat decreasing from 90N .slat
1650 lat -= self.slat
1651 lon += flon - self.wlon
1652 return (min(self.nlat - 1, max(0, int(lat * self.rlat))),
1653 max(0, int(lon * self.rlon)))
1655 def _tmpfile(self):
1656 # create a tmpfile to hold the cropped geoid grid
1657 try:
1658 from tempfile import NamedTemporaryFile as tmpfile
1659 except ImportError: # Python 2.7.16-
1660 from _MODS.os import tmpfile
1661 t = _MODS.os.path.basename(self.pgm)
1662 t = _MODS.os.path.splitext(t)[0]
1663 f = tmpfile(mode='w+b', prefix=t or 'egm')
1664 f.seek(0, _MODS.os.SEEK_SET) # force overwrite
1665 return f
1667 @Property_RO
1668 def pgm(self):
1669 '''Get the geoid file name (C{str}).
1670 '''
1671 return self.name
1674class PGMError(GeoidError):
1675 '''An issue while parsing or cropping an C{egm*.pgm} geoid dataset.
1676 '''
1677 pass
1680def egmGeoidHeights(GeoidHeights_dat):
1681 '''Generate geoid U{egm*.pgm<https://GeographicLib.SourceForge.io/
1682 C++/doc/geoid.html#geoidinst>} height tests from U{GeoidHeights.dat
1683 <https://SourceForge.net/projects/geographiclib/files/testdata/>}
1684 U{Test data for Geoids<https://GeographicLib.SourceForge.io/C++/doc/
1685 geoid.html#testgeoid>}.
1687 @arg GeoidHeights_dat: The un-gz-ed C{GeoidHeights.dat} file
1688 (C{str} or C{file} handle).
1690 @return: For each test, yield a L{GeoidHeight5Tuple}C{(lat, lon,
1691 egm84, egm96, egm2008)}.
1693 @raise GeoidError: Invalid B{C{GeoidHeights_dat}}.
1695 @note: Function L{egmGeoidHeights} is used to test the geoids
1696 L{GeoidKarney} and L{GeoidPGM}, see PyGeodesy module
1697 C{test/testGeoids.py}.
1698 '''
1699 dat = GeoidHeights_dat
1700 if isinstance(dat, bytes):
1701 dat = _BytesIO(dat)
1703 try:
1704 dat.seek(0, _MODS.os.SEEK_SET) # reset
1705 except AttributeError as x:
1706 raise GeoidError(GeoidHeights_dat=type(dat), cause=x)
1708 for t in dat.readlines():
1709 t = t.strip()
1710 if t and not t.startswith(_bHASH_):
1711 lat, lon, egm84, egm96, egm2008 = map(float, t.split())
1712 while lon > _180_0: # EasternLon to earth lon
1713 lon -= _360_0
1714 yield GeoidHeight5Tuple(lat, lon, egm84, egm96, egm2008)
1717__all__ += _ALL_DOCS(_GeoidBase)
1719if __name__ == '__main__': # MCCABE 14
1721 from pygeodesy.internals import printf, _secs2str, _versions, _sys
1722 from time import time
1724 _crop = {}
1725 _GeoidEGM = GeoidKarney
1726 _kind = 3
1728 geoids = _sys.argv[1:]
1729 while geoids:
1730 geoid = geoids.pop(0)
1732 if '-crop'.startswith(geoid.lower()):
1733 _crop = dict(crp=(20, -125, 50, -65)) # CONUS
1735 elif '-egm96'.startswith(geoid.lower()):
1736 _GeoidEGM = GeoidEGM96
1738 elif '-karney'.startswith(geoid.lower()):
1739 _GeoidEGM = GeoidKarney
1741 elif '-kind'.startswith(geoid.lower()):
1742 _kind = int(geoids.pop(0))
1744 elif '-pgm'.startswith(geoid.lower()):
1745 _GeoidEGM = GeoidPGM
1747 elif geoid[-4:].lower() in ('.pgm', '.grd'):
1748 g = _GeoidEGM(geoid, kind=_kind, **_crop)
1749 t = time()
1750 _ = g.highest()
1751 t = _secs2str(time() - t)
1752 printf('%s: %s (%s)', g.toStr(), t, _versions(), nl=1, nt=1)
1753 t = g.pgm
1754 if t:
1755 printf(repr(t), nt=1)
1756 # <https://GeographicLib.SourceForge.io/cgi-bin/GeoidEval>:
1757 # The height of the EGM96 geoid at Timbuktu
1758 # echo 16:46:33N 3:00:34W | GeoidEval
1759 # => 28.7068 -0.02e-6 -1.73e-6
1760 # The 1st number is the height of the geoid, the 2nd and
1761 # 3rd are its slopes in northerly and easterly direction
1762 t = 'Timbuktu %s' % (g,)
1763 k = {'egm84-15.pgm': '31.2979',
1764 'egm96-5.pgm': '28.7067',
1765 'egm2008-1.pgm': '28.7880'}.get(g.name.lower(), '28.7880')
1766 ll = _MODS.dms.parseDMS2('16:46:33N', '3:00:34W', sep=':')
1767 for ll in (ll, (16.776, -3.009),):
1768 try:
1769 h, ll = g.height(*ll), fstr(ll, prec=6)
1770 printf('%s.height(%s): %.4F vs %s', t, ll, h, k)
1771 except (GeoidError, RangeError) as x:
1772 printf(_COLONSPACE_(t, str(x)))
1774 elif geoid[-4:].lower() in ('.bin',):
1775 g = GeoidG2012B(geoid, kind=_kind)
1776 printf(g.toStr())
1778 else:
1779 raise GeoidError(grid=repr(geoid))
1781_I = int # PYCHOK unused _I
1782del _intCs, _T, _T0s12 # trash ints cache and map
1785# <https://GeographicLib.SourceForge.io/cgi-bin/GeoidEval>
1786# _lowerleft = -90, -179, -30.1500 # egm2008-1.pgm
1787# _lowerleft = -90, -179, -29.5350 # egm96-5.pgm
1788# _lowerleft = -90, -179, -29.7120 # egm84-15.pgm
1790# _center = 0, 0, 17.2260 # egm2008-1.pgm
1791# _center = 0, 0, 17.1630 # egm96-5.pgm
1792# _center = 0, 0, 18.3296 # egm84-15.pgm
1794# _upperright = 90, 180, 14.8980 # egm2008-1.pgm
1795# _upperright = 90, 180, 13.6050 # egm96-5.pgm
1796# _upperright = 90, 180, 13.0980 # egm84-15.pgm
1799# % python3.12 -m pygeodesy.geoids -egm96 ../testGeoids/WW15MGH.GRD
1800#
1801# GeoidEGM96('WW15MGH.GRD'): lowerleft(-90.0, -180.0, -29.534), upperright(90.0, 180.25, 13.606), center(0.0, 0.125, 17.125), highest(-8.25, -32.75, 85.391), lowest(4.75, -101.25, -106.991): 1.267 ms (pygeodesy 24.12.24 Python 3.12.7 64bit arm64 macOS 14.6.1)
1802#
1803# Timbuktu GeoidEGM96('WW15MGH.GRD').height(16.775833, -3.009444): 28.7073 vs 28.7880
1804# Timbuktu GeoidEGM96('WW15MGH.GRD').height(16.776, -3.009): 28.7072 vs 28.7880
1807# % python3.12 -m pygeodesy.geoids -Karney ../testGeoids/egm*.pgm
1808#
1809# GeoidKarney('egm2008-1.pgm'): lowerleft(-90.0, -180.0, -30.15), upperright(90.0, 180.0, 14.898), center(0.0, 0.0, 17.226), highest(-8.4, 147.367, 85.839), lowest(4.7, 78.767, -106.911): 204.334 ms (pygeodesy 24.8.24 Python 3.12.5 64bit arm64 macOS 14.6.1)
1810#
1811# _PGM('../testGeoids/egm2008-1.pgm'): AREA_OR_POINT='Point', DateTime='2009-08-31 06:54:00', Description='WGS84 EGM2008, 1-minute grid', Geoid='file in PGM format for the GeographicLib::Geoid class', MaxBilinearError=0.025, MaxCubicError=0.003, Offset=-108.0, Origin=LatLon2Tuple(lat=90.0, lon=0.0), Pixel=65535, RMSBilinearError=0.001, RMSCubicError=0.001, Scale=0.003, URL='http://earth-info.nga.mil/GandG/wgs84/gravitymod/egm2008', Vertical_Datum='WGS84'
1812#
1813# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1814# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1815#
1816# GeoidKarney('egm84-15.pgm'): lowerleft(-90.0, -180.0, -29.712), upperright(90.0, 180.0, 13.098), center(0.0, 0.0, 18.33), highest(-4.5, 148.75, 81.33), lowest(4.75, 79.25, -107.34): 1.007 ms (pygeodesy 24.8.24 Python 3.12.5 64bit arm64 macOS 14.6.1)
1817#
1818# _PGM('../testGeoids/egm84-15.pgm'): AREA_OR_POINT='Point', DateTime='2009-08-29 18:45:02', Description='WGS84 EGM84, 15-minute grid', Geoid='file in PGM format for the GeographicLib::Geoid class', MaxBilinearError=0.413, MaxCubicError=0.02, Offset=-108.0, Origin=LatLon2Tuple(lat=90.0, lon=0.0), Pixel=65535, RMSBilinearError=0.018, RMSCubicError=0.001, Scale=0.003, URL='http://earth-info.nga.mil/GandG/wgs84/gravitymod/wgs84_180/wgs84_180.html', Vertical_Datum='WGS84'
1819#
1820# Timbuktu GeoidKarney('egm84-15.pgm').height(16.775833, -3.009444): 31.2983 vs 31.2979
1821# Timbuktu GeoidKarney('egm84-15.pgm').height(16.776, -3.009): 31.2979 vs 31.2979
1822#
1823# GeoidKarney('egm96-5.pgm'): lowerleft(-90.0, -180.0, -29.535), upperright(90.0, 180.0, 13.605), center(0.0, 0.0, 17.163), highest(-8.167, 147.25, 85.422), lowest(4.667, 78.833, -107.043): 8.509 ms (pygeodesy 24.8.24 Python 3.12.5 64bit arm64 macOS 14.6.1)
1824#
1825# _PGM('../testGeoids/egm96-5.pgm'): AREA_OR_POINT='Point', DateTime='2009-08-29 18:45:03', Description='WGS84 EGM96, 5-minute grid', Geoid='file in PGM format for the GeographicLib::Geoid class', MaxBilinearError=0.14, MaxCubicError=0.003, Offset=-108.0, Origin=LatLon2Tuple(lat=90.0, lon=0.0), Pixel=65535, RMSBilinearError=0.005, RMSCubicError=0.001, Scale=0.003, URL='http://earth-info.nga.mil/GandG/wgs84/gravitymod/egm96/egm96.html', Vertical_Datum='WGS84'
1826#
1827# Timbuktu GeoidKarney('egm96-5.pgm').height(16.775833, -3.009444): 28.7068 vs 28.7067
1828# Timbuktu GeoidKarney('egm96-5.pgm').height(16.776, -3.009): 28.7067 vs 28.7067
1831# % python3.8 -m pygeodesy.geoids -egm96 ../testGeoids/WW15MGH.GRD
1832#
1833# GeoidEGM96('WW15MGH.GRD'): lowerleft(-90.0, -180.0, -29.534), upperright(90.0, 180.25, 13.606), center(0.0, 0.125, 17.125), highest(-8.25, -32.75, 85.391), lowest(4.75, -101.25, -106.991): 1.267 ms (pygeodesy 24.12.24 Python 3.8.10 64bit arm64_x86_64 macOS 10.16)
1834#
1835# Timbuktu GeoidEGM96('WW15MGH.GRD').height(16.775833, -3.009444): 28.7073 vs 28.7880
1836# Timbuktu GeoidEGM96('WW15MGH.GRD').height(16.776, -3.009): 28.7072 vs 28.7880
1839# % python3.8 -m pygeodesy.geoids -Karney ../testGeoids/egm*.pgm
1840#
1841# GeoidKarney('egm2008-1.pgm'): lowerleft(-90.0, -180.0, -30.15), upperright(90.0, 180.0, 14.898), center(0.0, 0.0, 17.226), highest(-8.4, 147.367, 85.839), lowest(4.7, 78.767, -106.911): 353.050 ms (pygeodesy 24.8.24 Python 3.8.10 64bit arm64_x86_64 macOS 10.16)
1842#
1843# _PGM('../testGeoids/egm2008-1.pgm'): AREA_OR_POINT='Point', DateTime='2009-08-31 06:54:00', Description='WGS84 EGM2008, 1-minute grid', Geoid='file in PGM format for the GeographicLib::Geoid class', MaxBilinearError=0.025, MaxCubicError=0.003, Offset=-108.0, Origin=LatLon2Tuple(lat=90.0, lon=0.0), Pixel=65535, RMSBilinearError=0.001, RMSCubicError=0.001, Scale=0.003, URL='http://earth-info.nga.mil/GandG/wgs84/gravitymod/egm2008', Vertical_Datum='WGS84'
1844#
1845# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1846# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1847#
1848# GeoidKarney('egm84-15.pgm'): lowerleft(-90.0, -180.0, -29.712), upperright(90.0, 180.0, 13.098), center(0.0, 0.0, 18.33), highest(-4.5, 148.75, 81.33), lowest(4.75, 79.25, -107.34): 1.727 ms (pygeodesy 24.8.24 Python 3.8.10 64bit arm64_x86_64 macOS 10.16)
1849#
1850# _PGM('../testGeoids/egm84-15.pgm'): AREA_OR_POINT='Point', DateTime='2009-08-29 18:45:02', Description='WGS84 EGM84, 15-minute grid', Geoid='file in PGM format for the GeographicLib::Geoid class', MaxBilinearError=0.413, MaxCubicError=0.02, Offset=-108.0, Origin=LatLon2Tuple(lat=90.0, lon=0.0), Pixel=65535, RMSBilinearError=0.018, RMSCubicError=0.001, Scale=0.003, URL='http://earth-info.nga.mil/GandG/wgs84/gravitymod/wgs84_180/wgs84_180.html', Vertical_Datum='WGS84'
1851#
1852# Timbuktu GeoidKarney('egm84-15.pgm').height(16.775833, -3.009444): 31.2983 vs 31.2979
1853# Timbuktu GeoidKarney('egm84-15.pgm').height(16.776, -3.009): 31.2979 vs 31.2979
1854#
1855# GeoidKarney('egm96-5.pgm'): lowerleft(-90.0, -180.0, -29.535), upperright(90.0, 180.0, 13.605), center(0.0, 0.0, 17.163), highest(-8.167, 147.25, 85.422), lowest(4.667, 78.833, -107.043): 14.807 ms (pygeodesy 24.8.24 Python 3.8.10 64bit arm64_x86_64 macOS 10.16)
1856#
1857# _PGM('../testGeoids/egm96-5.pgm'): AREA_OR_POINT='Point', DateTime='2009-08-29 18:45:03', Description='WGS84 EGM96, 5-minute grid', Geoid='file in PGM format for the GeographicLib::Geoid class', MaxBilinearError=0.14, MaxCubicError=0.003, Offset=-108.0, Origin=LatLon2Tuple(lat=90.0, lon=0.0), Pixel=65535, RMSBilinearError=0.005, RMSCubicError=0.001, Scale=0.003, URL='http://earth-info.nga.mil/GandG/wgs84/gravitymod/egm96/egm96.html', Vertical_Datum='WGS84'
1858#
1859# Timbuktu GeoidKarney('egm96-5.pgm').height(16.775833, -3.009444): 28.7068 vs 28.7067
1860# Timbuktu GeoidKarney('egm96-5.pgm').height(16.776, -3.009): 28.7067 vs 28.7067
1863# % python2 -m pygeodesy.geoids -Karney ../testGeoids/egm*.pgm
1864#
1865# GeoidKarney('egm2008-1.pgm'): lowerleft(-90.0, -180.0, -30.15), upperright(90.0, 180.0, 14.898), center(0.0, 0.0, 17.226), highest(-8.4, 147.367, 85.839), lowest(4.7, 78.767, -106.911): 283.362 ms (pygeodesy 24.8.24 Python 2.7.18 64bit arm64_x86_64 macOS 10.16)
1866#
1867# _PGM('../testGeoids/egm2008-1.pgm'): AREA_OR_POINT='Point', DateTime='2009-08-31 06:54:00', Description='WGS84 EGM2008, 1-minute grid', Geoid='file in PGM format for the GeographicLib::Geoid class', MaxBilinearError=0.025, MaxCubicError=0.003, Offset=-108.0, Origin=LatLon2Tuple(lat=90.0, lon=0.0), Pixel=65535, RMSBilinearError=0.001, RMSCubicError=0.001, Scale=0.003, URL='http://earth-info.nga.mil/GandG/wgs84/gravitymod/egm2008', Vertical_Datum='WGS84'
1868#
1869# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1870# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1871#
1872# GeoidKarney('egm84-15.pgm'): lowerleft(-90.0, -180.0, -29.712), upperright(90.0, 180.0, 13.098), center(0.0, 0.0, 18.33), highest(-4.5, 148.75, 81.33), lowest(4.75, 79.25, -107.34): 1.378 ms (pygeodesy 24.8.24 Python 2.7.18 64bit arm64_x86_64 macOS 10.16)
1873#
1874# _PGM('../testGeoids/egm84-15.pgm'): AREA_OR_POINT='Point', DateTime='2009-08-29 18:45:02', Description='WGS84 EGM84, 15-minute grid', Geoid='file in PGM format for the GeographicLib::Geoid class', MaxBilinearError=0.413, MaxCubicError=0.02, Offset=-108.0, Origin=LatLon2Tuple(lat=90.0, lon=0.0), Pixel=65535, RMSBilinearError=0.018, RMSCubicError=0.001, Scale=0.003, URL='http://earth-info.nga.mil/GandG/wgs84/gravitymod/wgs84_180/wgs84_180.html', Vertical_Datum='WGS84'
1875#
1876# Timbuktu GeoidKarney('egm84-15.pgm').height(16.775833, -3.009444): 31.2983 vs 31.2979
1877# Timbuktu GeoidKarney('egm84-15.pgm').height(16.776, -3.009): 31.2979 vs 31.2979
1878#
1879# GeoidKarney('egm96-5.pgm'): lowerleft(-90.0, -180.0, -29.535), upperright(90.0, 180.0, 13.605), center(0.0, 0.0, 17.163), highest(-8.167, 147.25, 85.422), lowest(4.667, 78.833, -107.043): 11.612 ms (pygeodesy 24.8.24 Python 2.7.18 64bit arm64_x86_64 macOS 10.16)
1880#
1881# _PGM('../testGeoids/egm96-5.pgm'): AREA_OR_POINT='Point', DateTime='2009-08-29 18:45:03', Description='WGS84 EGM96, 5-minute grid', Geoid='file in PGM format for the GeographicLib::Geoid class', MaxBilinearError=0.14, MaxCubicError=0.003, Offset=-108.0, Origin=LatLon2Tuple(lat=90.0, lon=0.0), Pixel=65535, RMSBilinearError=0.005, RMSCubicError=0.001, Scale=0.003, URL='http://earth-info.nga.mil/GandG/wgs84/gravitymod/egm96/egm96.html', Vertical_Datum='WGS84'
1882#
1883# Timbuktu GeoidKarney('egm96-5.pgm').height(16.775833, -3.009444): 28.7068 vs 28.7067
1884# Timbuktu GeoidKarney('egm96-5.pgm').height(16.776, -3.009): 28.7067 vs 28.7067
1887# % python3.12 -m pygeodesy.geoids -PGM ../testGeoids/egm*.pgm
1888#
1889# GeoidPGM('egm2008-1.pgm'): lowerleft(-90.0, -180.0, -30.15), upperright(90.0, 180.0, 14.898), center(0.0, 0.0, 17.226), highest(-8.4, -32.633, 85.839), lowest(4.683, -101.25, -106.911): 543.148 ms (pygeodesy 24.8.24 Python 3.12.5 64bit arm64 macOS 14.6.1)
1890#
1891# _PGM('../testGeoids/egm2008-1.pgm'): AREA_OR_POINT='Point', DateTime='2009-08-31 06:54:00', Description='WGS84 EGM2008, 1-minute grid', Geoid='file in PGM format for the GeographicLib::Geoid class', MaxBilinearError=0.025, MaxCubicError=0.003, Offset=-108.0, Origin=LatLon2Tuple(lat=90.0, lon=0.0), Pixel=65535, RMSBilinearError=0.001, RMSCubicError=0.001, Scale=0.003, URL='http://earth-info.nga.mil/GandG/wgs84/gravitymod/egm2008', Vertical_Datum='WGS84'
1892#
1893# Timbuktu GeoidPGM('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1894# Timbuktu GeoidPGM('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1895#
1896# GeoidPGM('egm84-15.pgm'): lowerleft(-90.0, -180.0, -29.712), upperright(90.0, 180.0, 13.098), center(0.0, 0.0, 18.33), highest(-4.5, -31.25, 81.33), lowest(4.75, -100.75, -107.34): 1.762 ms (pygeodesy 24.8.24 Python 3.12.5 64bit arm64 macOS 14.6.1)
1897#
1898# _PGM('../testGeoids/egm84-15.pgm'): AREA_OR_POINT='Point', DateTime='2009-08-29 18:45:02', Description='WGS84 EGM84, 15-minute grid', Geoid='file in PGM format for the GeographicLib::Geoid class', MaxBilinearError=0.413, MaxCubicError=0.02, Offset=-108.0, Origin=LatLon2Tuple(lat=90.0, lon=0.0), Pixel=65535, RMSBilinearError=0.018, RMSCubicError=0.001, Scale=0.003, URL='http://earth-info.nga.mil/GandG/wgs84/gravitymod/wgs84_180/wgs84_180.html', Vertical_Datum='WGS84'
1899#
1900# Timbuktu GeoidPGM('egm84-15.pgm').height(16.775833, -3.009444): 31.2979 vs 31.2979
1901# Timbuktu GeoidPGM('egm84-15.pgm').height(16.776, -3.009): 31.2975 vs 31.2979
1902#
1903# GeoidPGM('egm96-5.pgm'): lowerleft(-90.0, -180.0, -29.535), upperright(90.0, 180.0, 13.605), center(0.0, -0.0, 17.179), highest(-8.167, -32.75, 85.422), lowest(4.667, -101.167, -107.043): 12.594 ms (pygeodesy 24.8.24 Python 3.12.5 64bit arm64 macOS 14.6.1)
1904#
1905# _PGM('../testGeoids/egm96-5.pgm'): AREA_OR_POINT='Point', DateTime='2009-08-29 18:45:03', Description='WGS84 EGM96, 5-minute grid', Geoid='file in PGM format for the GeographicLib::Geoid class', MaxBilinearError=0.14, MaxCubicError=0.003, Offset=-108.0, Origin=LatLon2Tuple(lat=90.0, lon=0.0), Pixel=65535, RMSBilinearError=0.005, RMSCubicError=0.001, Scale=0.003, URL='http://earth-info.nga.mil/GandG/wgs84/gravitymod/egm96/egm96.html', Vertical_Datum='WGS84'
1906#
1907# Timbuktu GeoidPGM('egm96-5.pgm').height(16.775833, -3.009444): 28.7065 vs 28.7067
1908# Timbuktu GeoidPGM('egm96-5.pgm').height(16.776, -3.009): 28.7064 vs 28.7067
1910# **) MIT License
1911#
1912# Copyright (C) 2016-2025 -- mrJean1 at Gmail -- All Rights Reserved.
1913#
1914# Permission is hereby granted, free of charge, to any person obtaining a
1915# copy of this software and associated documentation files (the "Software"),
1916# to deal in the Software without restriction, including without limitation
1917# the rights to use, copy, modify, merge, publish, distribute, sublicense,
1918# and/or sell copies of the Software, and to permit persons to whom the
1919# Software is furnished to do so, subject to the following conditions:
1920#
1921# The above copyright notice and this permission notice shall be included
1922# in all copies or substantial portions of the Software.
1923#
1924# THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
1925# OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
1926# FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
1927# THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR
1928# OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE,
1929# ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR
1930# OTHER DEALINGS IN THE SOFTWARE.