Coverage for pygeodesy/geoids.py: 96%
689 statements
« prev ^ index » next coverage.py v7.6.1, created at 2025-04-25 13:15 -0400
« prev ^ index » next coverage.py v7.6.1, created at 2025-04-25 13:15 -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 _isin, len2, min2, isodd, _splituple, \
93 ub2str as _ub2str
94from pygeodesy.constants import EPS, _float as _F, _1_0, _N_90_0, _180_0, \
95 _N_180_0, _360_0
96from pygeodesy.datums import Datums, _ellipsoidal_datum, _WGS84
97# from pygeodesy.dms import parseDMS2 # _MODS
98from pygeodesy.errors import _incompatible, LenError, RangeError, _SciPyIssue, \
99 _xkwds_pop2
100from pygeodesy.fmath import favg, Fdot, fdot, Fhorner, frange
101# from pygoedesy.formy import heightOrthometric # _MODS
102from pygeodesy.heights import _as_llis2, _ascalar, _HeightBase, HeightError, \
103 _Wrap
104# from pygeodesy.internals import typename, _version2 # _MODS
105from pygeodesy.interns import NN, _COLONSPACE_, _COMMASPACE_, _DMAIN_, _E_, \
106 _height_, _in_, _kind_, _lat_, _lon_, _mean_, _N_, \
107 _n_a_, _numpy_, _on_, _outside_, _S_, _s_, _scipy_, \
108 _SPACE_, _stdev_, _tbd_, _W_, _width_, _4_
109from pygeodesy.lazily import _ALL_DOCS, _ALL_LAZY, _ALL_MODS as _MODS, _FOR_DOCS
110from pygeodesy.named import _name__, _Named, _NamedTuple
111# from pygeodesy.namedTuples import LatLon3Tuple # _MODS
112from pygeodesy.props import Property_RO, property_RO, property_ROver
113from pygeodesy.streprs import attrs, Fmt, fstr, pairs
114from pygeodesy.units import Height, Int_, Lat, Lon
115# from pygeodesy.utily import _Wrap # from .heights
117from math import floor as _floor
118# from os import SEEK_CUR, SEEK_SET # _MODS
119# import os.path # _MODS
120from struct import calcsize as _calcsize, unpack as _unpack
121try:
122 from StringIO import StringIO as _BytesIO # reads bytes
123 _ub2str = str # PYCHOK convert bytes to str for egm*.pgm text
125except ImportError: # Python 3+
126 from io import BytesIO as _BytesIO # PYCHOK expected
128__all__ = _ALL_LAZY.geoids
129__version__ = '25.04.14'
131_assert_ = 'assert'
132_bHASH_ = b'#'
133_endian_ = 'endian'
134_format_ = '%s %r'
135_header_ = 'header'
136_intCs = {} # cache int value, del below
137_lli_ = 'lli'
138_rb_ = 'rb'
139_supported_ = 'supported'
142class GeoidError(HeightError):
143 '''Geoid interpolator C{Geoid...} or interpolation issue.
144 '''
145 pass
148class _GeoidBase(_HeightBase):
149 '''(INTERNAL) Base class for C{Geoid...}s.
150 '''
151 _center = None # (lat, lon, height)
152 _cropped = None
153# _datum = _WGS84 # from _HeightBase
154 _egm = None # open C{egm*.pgm} geoid file
155 _endian = _tbd_
156 _Error = GeoidError # in ._HeightBase._as_lls, ...
157 _geoid = _n_a_
158 _hs_y_x = None # numpy 2darray, row-major order
159 _iscipy = True # scipy or Karney's interpolation
160 _kind = 3 # order for interp2d, RectBivariateSpline
161# _kmin = 2 # min number of knots
162 _knots = 0 # nlat * nlon
163 _mean = None # fixed in GeoidKarney
164 _nBytes = 0 # numpy size in bytes, float64
165 _pgm = None # PGM attributes, C{_PGM} or C{None}
166 _sizeB = 0 # geoid file size in bytes
167 _smooth = 0 # used only for RectBivariateSpline
168 _stdev = None # fixed in GeoidKarney
169 _u2B = 0 # np.itemsize or undefined
170 _yx_hits = None # cache hits, ala Karney's
172# _lat_d = _0_0 # increment, +tive
173# _lat_lo = _0_0 # lower lat, south
174# _lat_hi = _0_0 # upper lat, noth
175# _lon_d = _0_0 # increment, +tive
176# _lon_lo = _0_0 # left lon, west
177# _lon_hi = _0_0 # right lon, east
178# _lon_of = _0_0 # forward lon offset
179# _lon_og = _0_0 # reverse lon offset
181 def __init__(self, hs, p):
182 '''(INTERNAL) Set up the grid axes, the C{SciPy} interpolator and
183 several internal geoid attributes.
185 @arg hs: Grid knots with known height (C{numpy 2darray}).
186 @arg p: The C{slat, wlon, nlat, nlon, dlat, dlon} and other
187 geoid parameters (C{INTERNAL}).
188 '''
189 spi = self.scipy_interpolate
190 # for 2d scipy.interpolate.interp2d(xs, ys, hs, ...) and
191 # scipy.interpolate.RectBivariateSpline(ys, xs, hs, ...)
192 # require the shape of hs to be (len(ys), len(xs)), note
193 # the different (xs, ys, ...) and (ys, xs, ...) orders
194 if (p.nlat, p.nlon) != hs.shape:
195 raise GeoidError(shape=hs.shape, txt=_incompatible((p.nlat, p.nlon)))
197 # both axes and bounding box
198 ys, self._lat_d = self._gaxis2(p.slat, p.dlat, p.nlat, _lat_ + _s_)
199 xs, self._lon_d = self._gaxis2(p.wlon, p.dlon, p.nlon, _lon_ + _s_)
201 bb = ys[0], ys[-1], xs[0], xs[-1] + p.dlon # fudge lon_hi
202 # geoid grids are typically stored in row-major order, some
203 # with rows (90..-90) reversed and columns (0..360) wrapped
204 # to Easten longitude, 0 <= east < 180 and 180 <= west < 360
205 k = self.kind
206 if k in self._k2interp2d: # see _HeightBase
207 self._interp2d(xs, ys, hs, kind=k)
208 else: # XXX order ys and xs, see HeightLSQBiSpline
209 k = self._kxky(k)
210 self._ev = spi.RectBivariateSpline(ys, xs, hs, bbox=bb, ky=k, kx=k,
211 s=self._smooth).ev
212 self._hs_y_x = hs # numpy 2darray, row-major
213 self._nBytes = hs.nbytes # numpy size in bytes
214 self._knots = p.knots # grid knots == len(hs)
215 self._lon_of = float(p.flon) # forward offset
216 self._lon_og = g = float(p.glon) # reverse offset
217 # shrink the bounding box by 1 unit on every side:
218 # +self._lat_d, -self._lat_d, +self._lon_d, -self._lon_d
219 self._lat_lo, \
220 self._lat_hi, \
221 self._lon_lo, \
222 self._lon_hi = map(float, bb)
223 self._lon_lo -= g
224 self._lon_hi -= g
226 def __call__(self, *llis, **wrap_H):
227 '''Interpolate the geoid (or orthometric) height for one or more locations.
229 @arg llis: One or several locations (each C{LatLon}), all positional.
230 @kwarg wrap_H: Keyword arguments C{B{wrap}=False} (C{bool}) and
231 C{B{H}=False} (C{bool}). Use C{B{wrap}=True} to wrap
232 or I{normalize} all B{C{llis}} locations. If C{B{H}
233 is True}, return the I{orthometric} height instead of
234 the I{geoid} height at each location.
236 @return: A single geoid (or orthometric) height (C{float}) or
237 a list or tuple of geoid (or orthometric) heights (each
238 C{float}).
240 @raise GeoidError: Insufficient number of B{C{llis}}, an invalid
241 B{C{lli}} or the C{egm*.pgm} geoid file is closed.
243 @raise RangeError: An B{C{lli}} is outside this geoid's lat- or
244 longitude range.
246 @raise SciPyError: A C{scipy} issue.
248 @raise SciPyWarning: A C{scipy} warning as exception.
250 @note: To obtain I{orthometric} heights, each B{C{llis}} location
251 must have an ellipsoid C{height} or C{h} attribute, otherwise
252 C{height=0} is used.
254 @see: Function L{pygeodesy.heightOrthometric}.
255 '''
256 return self._called(llis, **wrap_H)
258 def __enter__(self):
259 '''Open context.
260 '''
261 return self
263 def __exit__(self, *unused): # PYCHOK exc_type, exc_value, exc_traceback)
264 '''Close context.
265 '''
266 self.close()
267 # return None # XXX False
269 def __repr__(self):
270 return self.toStr()
272 def __str__(self):
273 return Fmt.PAREN(self.classname, repr(self.name))
275 def _called(self, llis, wrap=False, H=False):
276 # handle __call__
277 _H = self._heightOrthometric if H else None
278 _as, llis = _as_llis2(llis, Error=GeoidError)
279 _w, hs = _Wrap._latlonop(wrap), []
280 _h, _a = self._hGeoid, hs.append
281 try:
282 for i, lli in enumerate(llis):
283 N = _h(*_w(lli.lat, lli.lon))
284 # orthometric or geoid height
285 _a(_H(lli, N) if _H else N)
286 return _as(hs)
287 except (GeoidError, RangeError) as x:
288 # XXX avoid str(LatLon()) degree symbols
289 n = _lli_ if _as is _ascalar else Fmt.INDEX(llis=i)
290 t = fstr((lli.lat, lli.lon), strepr=repr)
291 raise type(x)(n, t, wrap=wrap, H=H, cause=x)
292 except Exception as x:
293 if self._iscipy and self.scipy:
294 raise _SciPyIssue(x, self._ev_name)
295 else:
296 raise
298 @Property_RO
299 def _center(self):
300 ''' Cache for method L{center}.
301 '''
302 return self._llh3(favg(self._lat_lo, self._lat_hi),
303 favg(self._lon_lo, self._lon_hi))
305 def center(self, LatLon=None):
306 '''Return the center location and height of this geoid.
308 @kwarg LatLon: Optional class to return the location and height
309 (C{LatLon}) or C{None}.
311 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon,
312 height)} otherwise a B{C{LatLon}} instance with the lat-,
313 longitude and geoid height of the center grid location.
314 '''
315 return self._llh3LL(self._center, LatLon)
317 def close(self):
318 '''Close the C{egm*.pgm} geoid file if open (and applicable).
319 '''
320 if not self.closed:
321 self._egm.close()
322 self._egm = None
324 @property_RO
325 def closed(self):
326 '''Get the C{egm*.pgm} geoid file status.
327 '''
328 return self._egm is None
330 @Property_RO
331 def cropped(self):
332 '''Is geoid cropped (C{bool} or C{None} if crop not supported).
333 '''
334 return self._cropped
336 @Property_RO
337 def dtype(self):
338 '''Get the grid C{scipy} U{dtype<https://docs.SciPy.org/doc/numpy/
339 reference/generated/numpy.ndarray.dtype.html>} (C{numpy.dtype}).
340 '''
341 return self._hs_y_x.dtype
343 @Property_RO
344 def endian(self):
345 '''Get the geoid endianess and U{dtype<https://docs.SciPy.org/
346 doc/numpy/reference/generated/numpy.dtype.html>} (C{str}).
347 '''
348 return self._endian
350 def _ev(self, y, x): # PYCHOK overwritten with .RectBivariateSpline.ev
351 # see methods _HeightBase._ev and -._interp2d
352 return self._ev2d(x, y) # (y, x) flipped!
354 def _gaxis2(self, lo, d, n, name):
355 # build grid axis, hi = lo + (n - 1) * d
356 m, a = len2(frange(lo, n, d))
357 if m != n:
358 raise LenError(self.__class__, grid=m, **{name: n})
359 if d < 0:
360 d, a = -d, list(reversed(a))
361 a = self.numpy.array(a)
362 m, i = min2(*map(float, a[1:] - a[:-1]))
363 if m < EPS: # non-increasing axis
364 i = Fmt.INDEX(name, i + 1)
365 raise GeoidError(i, m, txt_not_='increasing')
366 return a, d
368 def _g2ll2(self, lat, lon): # PYCHOK no cover
369 '''(INTERNAL) I{Must be overloaded}.'''
370 self._notOverloaded(lat, lon)
372 def _gyx2g2(self, y, x):
373 # convert grid (y, x) indices to grid (lat, lon)
374 return ((self._lat_lo + self._lat_d * y),
375 (self._lon_lo + self._lon_of + self._lon_d * x))
377 def height(self, lats, lons, **wrap):
378 '''Interpolate the geoid height for one or several lat-/longitudes.
380 @arg lats: Latitude or latitudes (each C{degrees}).
381 @arg lons: Longitude or longitudes (each C{degrees}).
382 @kwarg wrap: Use C{B{wrap}=True} to wrap or I{normalize} all
383 B{C{lats}} and B{C{lons}}.
385 @return: A single geoid height (C{float}) or a list of geoid
386 heights (each C{float}).
388 @raise GeoidError: Insufficient or unequal number of B{C{lats}}
389 and B{C{lons}}.
391 @raise RangeError: A B{C{lat}} or B{C{lon}} is outside this geoid's
392 lat- or longitude range.
394 @raise SciPyError: A C{scipy} issue.
396 @raise SciPyWarning: A C{scipy} warning as exception.
397 '''
398 return _HeightBase.height(self, lats, lons, **wrap)
400 def height_(self, *latlons, **wrap):
401 '''Interpolate the geoid height for each M{(latlons[i], latlons[i+1])
402 pair for i in range(0, len(latlons), B{2})}.
404 @arg latlons: Alternating lat-/longitude pairs (each C{degrees}),
405 all positional.
407 @see: Method L{height} for further details.
409 @return: A tuple of geoid heights (each C{float}).
410 '''
411 lls = tuple(self._as_lls(latlons[0::2], *latlons[1::2]))
412 return self._called(lls, **wrap)
414 @property_ROver
415 def _heightOrthometric(self):
416 return _MODS.formy.heightOrthometric # overwrite property_ROver
418 def _hGeoid(self, lat, lon):
419 out = self.outside(lat, lon)
420 if out: # XXX avoid str(LatLon()) degree symbols
421 t = fstr((lat, lon), strepr=repr)
422 raise RangeError(lli=t, txt=_SPACE_(_outside_, _on_, out))
423 return float(self._ev(*self._ll2g2(lat, lon)))
425 @Property_RO
426 def _highest(self):
427 '''(INTERNAL) Cache for C{.highest}.
428 '''
429 return self._LL3T(self._llh3minmax(True), name__=self.highest)
431 def highest(self, LatLon=None, **unused):
432 '''Return the location and largest height of this geoid.
434 @kwarg LatLon: Optional class to return the location and height
435 (C{LatLon}) or C{None}.
437 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon,
438 height)} otherwise a B{C{LatLon}} instance with the lat-,
439 longitude and geoid height of the highest grid location.
440 '''
441 return self._llh3LL(self._highest, LatLon)
443 @Property_RO
444 def hits(self):
445 '''Get the number of cache hits (C{int} or C{None}).
446 '''
447 return self._yx_hits
449 @Property_RO
450 def kind(self):
451 '''Get the interpolator kind and order (C{int}).
452 '''
453 return self._kind
455 @Property_RO
456 def knots(self):
457 '''Get the number of grid knots (C{int}).
458 '''
459 return self._knots
461 def _ll2g2(self, lat, lon): # PYCHOK no cover
462 '''(INTERNAL) I{Must be overloaded}.'''
463 self._notOverloaded(lat, lon)
465 @property_ROver
466 def _LL3T(self):
467 '''(INTERNAL) Get L{LatLon3Tuple}, I{once}.
468 '''
469 return _MODS.namedTuples.LatLon3Tuple # overwrite property_ROver
471 def _llh3(self, lat, lon):
472 return self._LL3T(lat, lon, self._hGeoid(lat, lon), name=self.name)
474 def _llh3LL(self, llh, LatLon):
475 return llh if LatLon is None else self._xnamed(LatLon(*llh))
477 def _llh3minmax(self, highest, *unused):
478 hs, np = self._hs_y_x, self.numpy
479 # <https://docs.SciPy.org/doc/numpy/reference/generated/
480 # numpy.argmin.html#numpy.argmin>
481 arg = np.argmax if highest else np.argmin
482 y, x = np.unravel_index(arg(hs, axis=None), hs.shape)
483 return self._g2ll2(*self._gyx2g2(y, x)) + (float(hs[y, x]),)
485 def _load(self, g, dtype=float, n=-1, offset=0, **sep): # sep=NN
486 # numpy.fromfile, like .frombuffer
487 g.seek(offset, _MODS.os.SEEK_SET)
488 return self.numpy.fromfile(g, dtype, count=n, **sep)
490 @Property_RO
491 def _lowerleft(self):
492 '''(INTERNAL) Cache for C{.lowerleft}.
493 '''
494 return self._llh3(self._lat_lo, self._lon_lo)
496 def lowerleft(self, LatLon=None):
497 '''Return the lower-left location and height of this geoid.
499 @kwarg LatLon: Optional class to return the location
500 (C{LatLon}) and height or C{None}.
502 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon, height)}
503 otherwise a B{C{LatLon}} instance with the lat-, longitude and
504 geoid height of the lower-left, SW grid corner.
505 '''
506 return self._llh3LL(self._lowerleft, LatLon)
508 @Property_RO
509 def _loweright(self):
510 '''(INTERNAL) Cache for C{.loweright}.
511 '''
512 return self._llh3(self._lat_lo, self._lon_hi)
514 def loweright(self, LatLon=None):
515 '''Return the lower-right location and height of this geoid.
517 @kwarg LatLon: Optional class to return the location and height
518 (C{LatLon}) or C{None}.
520 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon, height)}
521 otherwise a B{C{LatLon}} instance with the lat-, longitude and
522 geoid height of the lower-right, SE grid corner.
523 '''
524 return self._llh3LL(self._loweright, LatLon)
526 lowerright = loweright # synonymous
528 @Property_RO
529 def _lowest(self):
530 '''(INTERNAL) Cache for C{.lowest}.
531 '''
532 return self._LL3T(self._llh3minmax(False), name__=self.lowest)
534 def lowest(self, LatLon=None, **unused):
535 '''Return the location and lowest height of this geoid.
537 @kwarg LatLon: Optional class to return the location and height
538 (C{LatLon}) or C{None}.
540 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon,
541 height)} otherwise a B{C{LatLon}} instance with the lat-,
542 longitude and geoid height of the lowest grid location.
543 '''
544 return self._llh3LL(self._lowest, LatLon)
546 @Property_RO
547 def mean(self):
548 '''Get the mean of this geoid's heights (C{float}).
549 '''
550 if self._mean is None: # see GeoidKarney
551 self._mean = float(self.numpy.mean(self._hs_y_x))
552 return self._mean
554 @property_RO
555 def name(self):
556 '''Get the name of this geoid (C{str}).
557 '''
558 return _HeightBase.name.fget(self) or self._geoid # recursion
560 @Property_RO
561 def nBytes(self):
562 '''Get the grid in-memory size in bytes (C{int}).
563 '''
564 return self._nBytes
566 def _open(self, geoid, datum, kind, name, smooth):
567 # open the geoid file
568 try:
569 self._geoid = _MODS.os.path.basename(geoid)
570 self._sizeB = _MODS.os.path.getsize(geoid)
571 g = open(geoid, _rb_)
572 except (IOError, OSError) as x:
573 raise GeoidError(geoid=geoid, cause=x)
575 if not _isin(datum, None, self._datum):
576 self._datum = _ellipsoidal_datum(datum, name=name)
577 self._kind = int(kind)
578 if name:
579 _HeightBase.name.fset(self, name) # rename
580 if smooth:
581 self._smooth = Int_(smooth=smooth, Error=GeoidError, low=0)
583 return g
585 def outside(self, lat, lon):
586 '''Check whether a location is outside this geoid's lat-/longitude
587 or crop range.
589 @arg lat: The latitude (C{degrees}).
590 @arg lon: The longitude (C{degrees}).
592 @return: A 1- or 2-character C{str} if outside, an empty C{str} otherwise.
593 '''
594 lat = _S_ if lat < self._lat_lo else (_N_ if lat > self._lat_hi else NN)
595 lon = _W_ if lon < self._lon_lo else (_E_ if lon > self._lon_hi else NN)
596 return NN(lat, lon) if lat and lon else (lat or lon)
598 @Property_RO
599 def pgm(self):
600 '''Get the PGM attributes (C{_PGM} or C{None} if not available/applicable).
601 '''
602 return self._pgm
604 @Property_RO
605 def sizeB(self):
606 '''Get the geoid grid file size in bytes (C{int}).
607 '''
608 return self._sizeB
610 @Property_RO
611 def smooth(self):
612 '''Get the C{RectBivariateSpline} smoothing (C{int}).
613 '''
614 return self._smooth
616 @Property_RO
617 def stdev(self):
618 '''Get the standard deviation of this geoid's heights (C{float}) or C{None}.
619 '''
620 if self._stdev is None: # see GeoidKarney
621 self._stdev = float(self.numpy.std(self._hs_y_x))
622 return self._stdev
624 def _swne(self, crop):
625 # crop box to 4-tuple (s, w, n, e)
626 try:
627 if len(crop) == 2:
628 try: # sw, ne LatLons
629 swne = (crop[0].lat, crop[0].lon,
630 crop[1].lat, crop[1].lon)
631 except AttributeError: # (s, w), (n, e)
632 swne = tuple(crop[0]) + tuple(crop[1])
633 else: # (s, w, n, e)
634 swne = crop
635 if len(swne) == 4:
636 s, w, n, e = map(float, swne)
637 if _N_90_0 <= s <= (n - _1_0) <= 89.0 and \
638 _N_180_0 <= w <= (e - _1_0) <= 179.0:
639 return s, w, n, e
640 except (IndexError, TypeError, ValueError):
641 pass
642 raise GeoidError(crop=crop)
644 def toStr(self, prec=3, sep=_COMMASPACE_): # PYCHOK signature
645 '''This geoid and all geoid attributes as a string.
647 @kwarg prec: Number of decimal digits (0..9 or C{None} for
648 default). Trailing zero decimals are stripped
649 for B{C{prec}} values of 1 and above, but kept
650 for negative B{C{prec}} values.
651 @kwarg sep: Separator to join (C{str}).
653 @return: Geoid name and attributes (C{str}).
654 '''
655 s = 1 if self.kind < 0 else 2
656 t = _MODS.internals.typename
657 t = tuple(Fmt.PAREN(t(m), fstr(m(), prec=prec)) for m in
658 (self.lowerleft, self.upperright,
659 self.center,
660 self.highest, self.lowest)) + \
661 attrs( _mean_, _stdev_, prec=prec, Nones=False) + \
662 attrs((_kind_, 'smooth')[:s], prec=prec, Nones=False) + \
663 attrs( 'cropped', 'dtype', _endian_, 'hits', 'knots', 'nBytes',
664 'sizeB', _scipy_, _numpy_, prec=prec, Nones=False)
665 return _COLONSPACE_(self, sep.join(t))
667 @Property_RO
668 def u2B(self):
669 '''Get the PGM itemsize in bytes (C{int}).
670 '''
671 return self._u2B
673 @Property_RO
674 def _upperleft(self):
675 '''(INTERNAL) Cache for C{.upperleft}.
676 '''
677 return self._llh3(self._lat_hi, self._lon_lo)
679 def upperleft(self, LatLon=None):
680 '''Return the upper-left location and height of this geoid.
682 @kwarg LatLon: Optional class to return the location and height
683 (C{LatLon}) or C{None}.
685 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon, height)}
686 otherwise a B{C{LatLon}} instance with the lat-, longitude and
687 geoid height of the upper-left, NW grid corner.
688 '''
689 return self._llh3LL(self._upperleft, LatLon)
691 @Property_RO
692 def _upperright(self):
693 '''(INTERNAL) Cache for C{.upperright}.
694 '''
695 return self._llh3(self._lat_hi, self._lon_hi)
697 def upperright(self, LatLon=None):
698 '''Return the upper-right location and height of this geoid.
700 @kwarg LatLon: Optional class to return the location and height
701 (C{LatLon}) or C{None}.
703 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon, height)}
704 otherwise a B{C{LatLon}} instance with the lat-, longitude and
705 geoid height of the upper-right, NE grid corner.
706 '''
707 return self._llh3LL(self._upperright, LatLon)
710class GeoidEGM96(_GeoidBase):
711 '''Geoid height interpolator for the EGM96 U{15 Minute Interpolation Grid<https://earth-info.NGA.mil>}
712 based on C{SciPy} interpolation U{RectBivariateSpline<https://docs.SciPy.org/doc/scipy/reference/
713 generated/scipy.interpolate.RectBivariateSpline.html>}, U{interp2d<https://docs.SciPy.org/doc/scipy/
714 reference/generated/scipy.interpolate.interp2d.html>} or U{bisplrep/-ev<https://docs.scipy.org/doc/
715 scipy/reference/generated/scipy.interpolate.bisplrep.html>}.
717 Use only the C{WW15MGH.GRD} file, unzipped from the EGM96 U{15 Minute Interpolation Grid
718 <https://earth-info.NGA.mil/index.php?dir=wgs84&action=wgs84>} download.
719 '''
720 def __init__(self, EGM96_grd, datum=_WGS84, kind=3, smooth=0, **name_crop):
721 '''New L{GeoidEGM96} interpolator.
723 @arg EGM96_grd: An C{EGM96_grd} grid file name (C{.GRD}).
724 @kwarg datum: Optional grid datum (L{Datum}, L{Ellipsoid}, L{Ellipsoid2} or L{a_f2Tuple}),
725 overriding C{WGS84}.
726 @kwarg kind: C{scipy.interpolate} order (C{int}), use 1..5 for U{RectBivariateSpline
727 <https://docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.
728 RectBivariateSpline.html>} or -1, -3 or -5 for U{bisplrep/-ev<https://
729 docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.bisplrep.html>}
730 or U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/scipy.
731 interpolate.interp2d.html>} C{linear}, C{cubic} respectively C{quintic},
732 see note for more details.
733 @kwarg smooth: Smoothing factor for C{B{kind}=1..5} only (C{int}).
734 @kwarg name_crop: Optional geoid C{B{name}=NN} (C{str}) and UNSUPPORTED keyword argument
735 C{B{crop}=None}.
737 @raise GeoidError: Invalid B{C{crop}}, B{C{kind}} or B{C{smooth}} or a ECM96 grid file
738 B{C{ECM96_grd}} issue.
740 @raise ImportError: Package C{numpy} or C{scipy} not found or not installed.
742 @raise LenError: Grid file B{C{EGM96_grd}} axis mismatch.
744 @raise SciPyError: A C{scipy} issue.
746 @raise SciPyWarning: A C{scipy} warning as exception.
748 @raise TypeError: Invalid B{C{datum}}.
750 @note: Specify C{B{kind}=-1, -3 or -5} to use C{scipy.interpolate.interp2d}
751 before or C{scipy.interpolate.bisplrep/-ev} since C{Scipy} version 1.14.
752 '''
753 crop, name = _xkwds_pop2(name_crop, crop=None)
754 if crop is not None:
755 raise GeoidError(crop=crop, txt_not_=_supported_)
757 g = self._open(EGM96_grd, datum, kind, _name__(**name), smooth)
758 _ = self.numpy # import numpy for .fromfile, .reshape
760 try:
761 p, hs = _Gpars(), self._load(g, sep=_SPACE_) # text
762 p.slat, n, p.wlon, e, p.dlat, p.dlon = hs[:6] # n-s, 0-E
763 p.nlat = int((n - p.slat) / p.dlat) + 1 # include S
764 p.nlon = int((e - p.wlon) / p.dlon) + 1 # include W
765 p.knots = p.nlat * p.nlon # inverted lats N downto S
766 p.glon = _180_0 # Eastern lons 0-360
767 hs = hs[6:].reshape(p.nlat, p.nlon)
768 _GeoidBase.__init__(self, hs, p)
770 except Exception as x:
771 raise _SciPyIssue(x, _in_, repr(EGM96_grd))
772 finally:
773 g.close()
775 def _g2ll2(self, lat, lon):
776 # convert grid (lat, lon) to earth (lat, lon)
777 while lon > _180_0: # Eastern
778 lon -= _360_0
779 return -lat, lon # invert lat
781 def _ll2g2(self, lat, lon):
782 # convert earth (lat, lon) to grid (lat, lon)
783 while lon < 0: # Eastern
784 lon += _360_0
785 return -lat, lon # invert lat
787 if _FOR_DOCS:
788 __call__ = _GeoidBase.__call__
789 height = _GeoidBase.height
792class GeoidG2012B(_GeoidBase):
793 '''Geoid height interpolator for U{GEOID12B Model
794 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/>} grids U{CONUS
795 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_CONUS.shtml>},
796 U{Alaska<https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_AK.shtml>},
797 U{Hawaii<https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_HI.shtml>},
798 U{Guam and Northern Mariana Islands
799 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_GMNI.shtml>},
800 U{Puerto Rico and U.S. Virgin Islands
801 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_PRVI.shtml>} and
802 U{American Samoa<https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_AS.shtml>}
803 based on C{SciPy} interpolation U{RectBivariateSpline<https://docs.SciPy.org/doc/
804 scipy/reference/generated/scipy.interpolate.RectBivariateSpline.html>}, U{interp2d
805 <https://docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.interp2d.html>}
806 or U{bisplrep/-ev<https://docs.scipy.org/doc/scipy/reference/generated/scipy.interpolate.
807 bisplrep.html>}.
809 Use any of the C{le} (little endian) or C{be} (big endian) C{g2012b*.bin} binary grid files.
810 '''
811 _datum = Datums.NAD83
813 def __init__(self, g2012b_bin, datum=Datums.NAD83, kind=3, smooth=0, **name_crop):
814 '''New L{GeoidG2012B} interpolator.
816 @arg g2012b_bin: A C{GEOID12B} grid file name (C{.bin}).
817 @kwarg datum: Optional grid datum (L{Datum}, L{Ellipsoid}, L{Ellipsoid2} or
818 L{a_f2Tuple}), overriding C{NAD83}.
819 @kwarg kind: C{scipy.interpolate} order (C{int}), use 1..5 for U{RectBivariateSpline
820 <https://docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.
821 RectBivariateSpline.html>} or -1, -3 or -5 for U{bisplrep/-ev<https://
822 docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.bisplrep.html>}
823 or U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/scipy.
824 interpolate.interp2d.html>} C{linear}, C{cubic} respectively C{quintic},
825 see note for more details.
826 @kwarg smooth: Smoothing factor for C{B{kind}=1..5} only (C{int}).
827 @kwarg name_crop: Optional geoid C{B{name}=NN} (C{str}) and UNSUPPORTED keyword argument
828 C{B{crop}=None}.
830 @raise GeoidError: Invalid B{C{crop}}, B{C{kind}} or B{C{smooth}} or a G2012B grid file
831 B{C{g2012b_bin}} issue.
833 @raise ImportError: Package C{numpy} or C{scipy} not found or not installed.
835 @raise LenError: Grid file B{C{g2012b_bin}} axis mismatch.
837 @raise SciPyError: A C{scipy} issue.
839 @raise SciPyWarning: A C{scipy} warning as exception.
841 @raise TypeError: Invalid B{C{datum}}.
843 @note: Specify C{B{kind}=-1, -3 or -5} to use C{scipy.interpolate.interp2d}
844 before or C{scipy.interpolate.bisplrep/-ev} since C{Scipy} version 1.14.
845 '''
846 crop, name = _xkwds_pop2(name_crop, crop=None)
847 if crop is not None:
848 raise GeoidError(crop=crop, txt_not_=_supported_)
850 g = self._open(g2012b_bin, datum, kind, _name__(**name), smooth)
851 _ = self.numpy # import numpy for ._load and
853 try:
854 p = _Gpars()
855 n = (self.sizeB // 4) - 11 # number of f4 heights
856 # U{numpy dtype formats are different from Python struct formats
857 # <https://docs.SciPy.org/doc/numpy-1.15.0/reference/arrays.dtypes.html>}
858 for en_ in ('<', '>'):
859 # skip 4xf8, get 3xi4
860 p.nlat, p.nlon, ien = map(int, self._load(g, en_+'i4', 3, 32))
861 if ien == 1: # correct endian
862 p.knots = p.nlat * p.nlon
863 if p.knots == n and 1 < p.nlat < n \
864 and 1 < p.nlon < n:
865 self._endian = en_+'f4'
866 break
867 else: # couldn't validate endian
868 raise GeoidError(_endian_)
870 # get the first 4xf8
871 p.slat, p.wlon, p.dlat, p.dlon = map(float, self._load(g, en_+'f8', 4))
872 # read all f4 heights, ignoring the first 4xf8 and 3xi4
873 hs = self._load(g, self._endian, n, 44).reshape(p.nlat, p.nlon)
874 p.wlon -= _360_0 # western-most East longitude to earth (..., lon)
875 _GeoidBase.__init__(self, hs, p)
877 except Exception as x:
878 raise _SciPyIssue(x, _in_, repr(g2012b_bin))
879 finally:
880 g.close()
882 def _g2ll2(self, lat, lon):
883 # convert grid (lat, lon) to earth (lat, lon)
884 return lat, lon
886 def _ll2g2(self, lat, lon):
887 # convert earth (lat, lon) to grid (lat, lon)
888 return lat, lon
890 if _FOR_DOCS:
891 __call__ = _GeoidBase.__call__
892 height = _GeoidBase.height
895class GeoidHeight5Tuple(_NamedTuple): # .geoids.py
896 '''5-Tuple C{(lat, lon, egm84, egm96, egm2008)} for U{GeoidHeights.dat
897 <https://SourceForge.net/projects/geographiclib/files/testdata/>}
898 tests with the heights for 3 different EGM grids at C{degrees90}
899 and C{degrees180} degrees (after converting C{lon} from original
900 C{0 <= EasterLon <= 360}).
901 '''
902 _Names_ = (_lat_, _lon_, 'egm84', 'egm96', 'egm2008')
903 _Units_ = ( Lat, Lon, Height, Height, Height)
906def _I(i):
907 '''(INTERNAL) Cache a single C{int} constant.
908 '''
909 i = int(i)
910 return _intCs.setdefault(i, i) # PYCHOK undefined due to del _intCs
913def _T(cs):
914 '''(INTERNAL) Cache a tuple of single C{int} constants.
915 '''
916 return tuple(map(_I, _splituple(cs)))
918_T0s12 = (_I(0),) * 12 # PYCHOK _T('0, 0, ..., 0')
921class GeoidKarney(_GeoidBase):
922 '''Geoid height interpolator for I{Karney}'s U{GeographicLib Earth
923 Gravitational Model (EGM)<https://GeographicLib.SourceForge.io/C++/doc/
924 geoid.html>} geoid U{egm*.pgm<https://GeographicLib.SourceForge.io/
925 C++/doc/geoid.html#geoidinst>} datasets using bilinear or U{cubic
926 <https://dl.ACM.org/citation.cfm?id=368443>} interpolation and U{caching
927 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidcache>}
928 in pure Python, transcoded from I{Karney}'s U{C++ class Geoid
929 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinterp>}.
931 Use any of the geoid U{egm84-, egm96- or egm2008-*.pgm
932 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinst>}
933 datasets.
934 '''
935 _C0 = _F(372), _F(240), _F(372) # n, _ and s common denominators
936 # matrices c3n_, c3, c3s_, transposed from GeographicLib/Geoid.cpp
937 _C3 = ((_T('0, 0, 62, 124, 124, 62, 0, 0, 0, 0, 0, 0'),
938 _T0s12,
939 _T('-131, 7, -31, -62, -62, -31, 45, 216, 156, -45, -55, -7'),
940 _T0s12,
941 _T('138, -138, 0, 0, 0, 0, -183, 33, 153, -3, 48, -48'), # PYCHOK indent
942 _T('144, 42, -62, -124, -124, -62, -9, 87, 99, 9, 42, -42'),
943 _T0s12,
944 _T('0, 0, 0, 0, 0, 0, 93, -93, -93, 93, 0, 0'),
945 _T('-102, 102, 0, 0, 0, 0, 18, 12, -12, -18, -84, 84'),
946 _T('-31, -31, 31, 62, 62, 31, 0, -93, -93, 0, 31, 31')), # PYCHOK indent
948 (_T('9, -9, 9, 186, 54, -9, -9, 54, -54, 9, -9, 9'),
949 _T('-18, 18, -88, -42, 162, -32, 8, -78, 78, -8, 18, -18'),
950 _T('-88, 8, -18, -42, -78, 18, 18, 162, 78, -18, -32, -8'),
951 _T('0, 0, 90, -150, 30, 30, 30, -90, 90, -30, 0, 0'),
952 _T('96, -96, 96, -96, -24, 24, -96, -24, 144, -24, 24, -24'), # PYCHOK indent
953 _T('90, 30, 0, -150, -90, 0, 0, 30, 90, 0, 30, -30'),
954 _T('0, 0, -20, 60, -60, 20, -20, 60, -60, 20, 0, 0'),
955 _T('0, 0, -60, 60, 60, -60, 60, -60, -60, 60, 0, 0'),
956 _T('-60, 60, 0, 60, -60, 0, 0, 60, -60, 0, -60, 60'),
957 _T('-20, -20, 0, 60, 60, 0, 0, -60, -60, 0, 20, 20')),
959 (_T('18, -18, 36, 210, 162, -36, 0, 0, 0, 0, -18, 18'), # PYCHOK indent
960 _T('-36, 36, -165, 45, 141, -21, 0, 0, 0, 0, 36, -36'),
961 _T('-122, -2, -27, -111, -75, 27, 62, 124, 124, 62, -64, 2'),
962 _T('0, 0, 93, -93, -93, 93, 0, 0, 0, 0, 0, 0'),
963 _T('120, -120, 147, -57, -129, 39, 0, 0, 0, 0, 66, -66'), # PYCHOK indent
964 _T('135, 51, -9, -192, -180, 9, 31, 62, 62, 31, 51, -51'),
965 _T0s12,
966 _T('0, 0, -93, 93, 93, -93, 0, 0, 0, 0, 0, 0'),
967 _T('-84, 84, 18, 12, -12, -18, 0, 0, 0, 0, -102, 102'),
968 _T('-31, -31, 0, 93, 93, 0, -31, -62, -62, -31, 31, 31')))
970 _BT = (_T('0, 0'), # bilinear 4-tuple [i, j] indices
971 _T('1, 0'),
972 _T('0, 1'),
973 _T('1, 1'))
975 _CM = (_T(' 0, -1'), # 10x12 cubic matrix [i, j] indices
976 _T(' 1, -1'),
977 _T('-1, 0'),
978 _T(' 0, 0'),
979 _T(' 1, 0'),
980 _T(' 2, 0'),
981 _T('-1, 1'),
982 _T(' 0, 1'),
983 _T(' 1, 1'),
984 _T(' 2, 1'),
985 _T(' 0, 2'),
986 _T(' 1, 2'))
988# _cropped = None
989 _endian = '>H' # struct.unpack 1 ushort (big endian, unsigned short)
990 _4endian = '>4H' # struct.unpack 4 ushorts
991 _Rendian = NN # struct.unpack a row of ushorts
992# _highest = (-8.4, 147.367, 85.839) if egm2008-1.pgm else (
993# (-8.167, 147.25, 85.422) if egm96-5.pgm else
994# (-4.5, 148.75, 81.33)) # egm84-15.pgm
995 _iscipy = False
996# _lowest = (4.7, 78.767, -106.911) if egm2008-1.pgm else (
997# (4.667, 78.833, -107.043) if egm96-5.pgm else
998# (4.75, 79.25, -107.34)) # egm84-15.pgm
999 _mean = _F(-1.317) # from egm2008-1, -1.438 egm96-5, -0.855 egm84-15
1000 _nBytes = None # not applicable
1001 _nterms = len(_C3[0]) # columns length, number of rows
1002 _smooth = None # not applicable
1003 _stdev = _F(29.244) # from egm2008-1, 29.227 egm96-5, 29.183 egm84-15
1004 _u2B = _calcsize(_endian) # pixelsize_ in bytes
1005 _4u2B = _calcsize(_4endian) # 4 pixelsize_s in bytes
1006 _Ru2B = 0 # row of pixelsize_s in bytes
1007 _yx_hits = 0 # cache hits
1008 _yx_i = () # cached (y, x) indices
1009 _yx_t = () # cached 4- or 10-tuple for _ev2k resp. _ev3k
1011 def __init__(self, egm_pgm, crop=None, datum=_WGS84, kind=3, **name_smooth):
1012 '''New L{GeoidKarney} interpolator.
1014 @arg egm_pgm: An U{EGM geoid dataset<https://GeographicLib.SourceForge.io/
1015 C++/doc/geoid.html#geoidinst>} file name (C{egm*.pgm}), see
1016 note below.
1017 @kwarg crop: Optional box to limit geoid locations, a 4-tuple (C{south,
1018 west, north, east}), 2-tuple (C{(south, west), (north, east)})
1019 with 2 C{degrees90} lat- and C{degrees180} longitudes or as
1020 2-tuple (C{LatLonSW, LatLonNE}) of C{LatLon} instances.
1021 @kwarg datum: Optional grid datum (C{Datum}, L{Ellipsoid}, L{Ellipsoid2} or
1022 L{a_f2Tuple}), overriding C{WGS84}.
1023 @kwarg kind: Interpolation order (C{int}), 2 for C{bilinear} or 3 for C{cubic}.
1024 @kwarg name_smooth: Optional geoid C{B{name}=NN} (C{str}) and UNSUPPORTED
1025 keyword argument C{B{smooth}}, use C{B{smooth}=None} to ignore.
1027 @raise GeoidError: EGM dataset B{C{egm_pgm}} issue or invalid B{C{crop}},
1028 B{C{kind}} or B{C{smooth}}.
1030 @raise TypeError: Invalid B{C{datum}}.
1032 @see: Class L{GeoidPGM} and function L{egmGeoidHeights}.
1034 @note: Geoid file B{C{egm_pgm}} remains open and I{must be closed} by calling
1035 method C{close} or by using C{with B{GeoidKarney}(...) as ...:} context.
1036 '''
1037 smooth, name = _xkwds_pop2(name_smooth, smooth=None)
1038 if smooth is not None:
1039 raise GeoidError(smooth=smooth, txt_not_=_supported_)
1041 if _isin(kind, 2):
1042 self._ev2d = self._ev2k # see ._ev_name
1043 elif not _isin(kind, 3):
1044 raise GeoidError(kind=kind)
1046 self._egm = g = self._open(egm_pgm, datum, kind, _name__(**name), None)
1047 self._pgm = p = _PGM(g, pgm=egm_pgm, itemsize=self.u2B, sizeB=self.sizeB)
1049 self._Rendian = self._4endian.replace(_4_, str(p.nlon))
1050 self._Ru2B = _calcsize(self._Rendian)
1052 self._knots = p.knots # grid knots
1053 self._lon_of = float(p.flon) # forward offset
1054 self._lon_og = float(p.glon) # reverse offset
1055 # set earth (lat, lon) limits (s, w, n, e)
1056 self._lat_lo, self._lon_lo, \
1057 self._lat_hi, self._lon_hi = self._swne(crop if crop else p.crop4)
1058 self._cropped = bool(crop)
1060 def _c0c3v(self, y, x):
1061 # get the common denominator, the 10x12 cubic matrix and
1062 # the 12 cubic v-coefficients around geoid index (y, x)
1063 p = self._pgm
1064 if 0 < x < (p.nlon - 2) and 0 < y < (p.nlat - 2):
1065 # read 4x4 ushorts, drop the 4 corners
1066 S = _MODS.os.SEEK_SET
1067 e = self._4endian
1068 g = self._egm
1069 n = self._4u2B
1070 R = self._Ru2B
1071 b = self._seek(y - 1, x - 1)
1072 v = _unpack(e, g.read(n))[1:3]
1073 b += R
1074 g.seek(b, S)
1075 v += _unpack(e, g.read(n))
1076 b += R
1077 g.seek(b, S)
1078 v += _unpack(e, g.read(n))
1079 b += R
1080 g.seek(b, S)
1081 v += _unpack(e, g.read(n))[1:3]
1082 j = 1
1084 else: # likely some wrapped y and/or x's
1085 v = self._raws(y, x, GeoidKarney._CM)
1086 j = 0 if y < 1 else (1 if y < (p.nlat - 2) else 2)
1088 return GeoidKarney._C0[j], GeoidKarney._C3[j], v
1090 @Property_RO
1091 def dtype(self):
1092 '''Get the geoid's grid data type (C{str}).
1093 '''
1094 return 'ushort'
1096 def _ev(self, lat, lon): # PYCHOK expected
1097 # interpolate the geoid height at grid (lat, lon)
1098 fy, fx = self._g2yx2(lat, lon)
1099 y, x = int(_floor(fy)), int(_floor(fx))
1100 fy -= y
1101 fx -= x
1102 H = self._ev2d(fy, fx, y, x) # PYCHOK ._ev2k or ._ev3k
1103 H *= self._pgm.Scale # H.fmul(self._pgm.Scale)
1104 H += self._pgm.Offset # H.fadd(self._pgm.Offset)
1105 return H.fsum() # float(H)
1107 def _ev2k(self, fy, fx, *yx):
1108 # compute the bilinear 4-tuple and interpolate raw H
1109 if self._yx_i == yx:
1110 self._yx_hits += 1
1111 else:
1112 y, x = self._yx_i = yx
1113 self._yx_t = self._raws(y, x, GeoidKarney._BT)
1114 t = self._yx_t
1115 v = _1_0, -fx, fx
1116 H = Fdot(v, t[0], t[0], t[1]).fmul(_1_0 - fy) # c = a * (1 - fy)
1117 H += Fdot(v, t[2], t[2], t[3]).fmul(fy) # c += b * fy
1118 return H
1120 def _ev3k(self, fy, fx, *yx):
1121 # compute the cubic 10-tuple and interpolate raw H
1122 if self._yx_i == yx:
1123 self._yx_hits += 1
1124 else:
1125 c0, c3, v = self._c0c3v(*yx)
1126 # assert len(c3) == self._nterms
1127 self._yx_t = tuple(fdot(v, *r3) / c0 for r3 in c3)
1128 self._yx_i = yx
1129 # GeographicLib/Geoid.cpp Geoid::height(lat, lon) ...
1130 # real h = t[0] + fx * (t[1] + fx * (t[3] + fx * t[6])) +
1131 # fy * (t[2] + fx * (t[4] + fx * t[7]) +
1132 # fy * (t[5] + fx * t[8] + fy * t[9]));
1133 t = self._yx_t
1134 v = _1_0, fx, fy
1135 H = Fdot(v, t[5], t[8], t[9])
1136 H *= fy
1137 H += Fhorner(fx, t[2], t[4], t[7])
1138 H *= fy
1139 H += Fhorner(fx, t[0], t[1], t[3], t[6])
1140 return H
1142 _ev2d = _ev3k # overriden for kind=2, see ._ev_name
1144 def _g2ll2(self, lat, lon):
1145 # convert grid (lat, lon) to earth (lat, lon), uncropped
1146 while lon > _180_0:
1147 lon -= _360_0
1148 return lat, lon
1150 def _g2yx2(self, lat, lon):
1151 # convert grid (lat, lon) to grid (y, x) indices
1152 p = self._pgm
1153 # note, slat = +90, rlat < 0 makes y >=0
1154 return ((lat - p.slat) * p.rlat), ((lon - p.wlon) * p.rlon)
1156 def _gyx2g2(self, y, x):
1157 # convert grid (y, x) indices to grid (lat, lon)
1158 p = self._pgm
1159 return (p.slat + p.dlat * y), (p.wlon + p.dlon * x)
1161 @Property_RO
1162 def _highest_ltd(self):
1163 '''(INTERNAL) Cache for C{.highest}.
1164 '''
1165 return self._LL3T(self._llh3minmax(True, -12, -4), name__=self.highest)
1167 def highest(self, LatLon=None, full=False): # PYCHOK full
1168 '''Return the location and largest height of this geoid.
1170 @kwarg LatLon: Optional class to return the location and height
1171 (C{LatLon}) or C{None}.
1172 @kwarg full: Search the full or limited latitude range (C{bool}).
1174 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon,
1175 height)} otherwise a B{C{LatLon}} instance with the lat-,
1176 longitude and geoid height of the highest grid location.
1177 '''
1178 llh = self._highest if full or self.cropped else self._highest_ltd
1179 return self._llh3LL(llh, LatLon)
1181 def _lat2y2(self, lat2):
1182 # convert earth lat(s) to min and max grid y indices
1183 ys, m = [], self._pgm.nlat - 1
1184 for lat in lat2:
1185 y, _ = self._g2yx2(*self._ll2g2(lat, 0))
1186 ys.append(max(min(int(y), m), 0))
1187 return min(ys), max(ys) + 1
1189 def _ll2g2(self, lat, lon):
1190 # convert earth (lat, lon) to grid (lat, lon), uncropped
1191 while lon < 0:
1192 lon += _360_0
1193 return lat, lon
1195 def _llh3minmax(self, highest, *lat2):
1196 # find highest or lowest, takes 10+ secs for egm2008-1.pgm geoid
1197 # (Python 2.7.16, macOS 10.13.6 High Sierra, iMac 3 GHz Core i3)
1198 if highest:
1199 def _mt(r, h):
1200 m = max(r)
1201 return m, (m > h)
1203 else: # lowest
1204 def _mt(r, h): # PYCHOK redef
1205 m = min(r)
1206 return m, (m < h)
1208 y = x = 0
1209 h = self._raw(y, x)
1210 for j, r in self._raw2(*lat2):
1211 m, t = _mt(r, h)
1212 if t:
1213 h, y, x = m, j, r.index(m)
1214 h *= self._pgm.Scale
1215 h += self._pgm.Offset
1216 return self._g2ll2(*self._gyx2g2(y, x)) + (h,)
1218 @Property_RO
1219 def _lowest_ltd(self):
1220 '''(INTERNAL) Cache for C{.lowest}.
1221 '''
1222 return self._LL3T(self._llh3minmax(False, 0, 8), name__=self.lowest)
1224 def lowest(self, LatLon=None, full=False): # PYCHOK full
1225 '''Return the location and lowest height of this geoid.
1227 @kwarg LatLon: Optional class to return the location and height
1228 (C{LatLon}) or C{None}.
1229 @kwarg full: Search the full or limited latitude range (C{bool}).
1231 @return: If C{B{LatLon} is None}, a L{LatLon3Tuple}C{(lat, lon,
1232 height)} otherwise a B{C{LatLon}} instance with the lat-,
1233 longitude and geoid height of the lowest grid location.
1234 '''
1235 llh = self._lowest if full or self.cropped else self._lowest_ltd
1236 return self._llh3LL(llh, LatLon)
1238 def _raw(self, y, x):
1239 # get the ushort geoid height at geoid index (y, x),
1240 # like GeographicLib/Geoid.hpp real rawval(is, iy)
1241 p = self._pgm
1242 if x < 0:
1243 x += p.nlon
1244 elif x >= p.nlon:
1245 x -= p.nlon
1246 h = p.nlon // 2
1247 if y < 0:
1248 y = -y
1249 elif y >= p.nlat:
1250 y = (p.nlat - 1) * 2 - y
1251 else:
1252 h = 0
1253 x += h if x < h else -h
1254 self._seek(y, x)
1255 h = _unpack(self._endian, self._egm.read(self._u2B))
1256 return h[0]
1258 def _raws(self, y, x, ijs):
1259 # get bilinear 4-tuple or 10x12 cubic matrix
1260 return tuple(self._raw(y + j, x + i) for i, j in ijs)
1262 def _raw2(self, *lat2):
1263 # yield a 2-tuple (y, ushorts) for each row or for
1264 # the rows between two (or more) earth lat values
1265 p = self._pgm
1266 g = self._egm
1267 e = self._Rendian
1268 n = self._Ru2B
1269 # min(lat2) <= lat <= max(lat2) or 0 <= y < p.nlat
1270 s, t = self._lat2y2(lat2) if lat2 else (0, p.nlat)
1271 self._seek(s, 0) # to start of row s
1272 for y in range(s, t):
1273 yield y, _unpack(e, g.read(n))
1275 def _seek(self, y, x):
1276 # position geoid to grid index (y, x)
1277 p, g = self._pgm, self._egm
1278 if g:
1279 b = p.skip + (y * p.nlon + x) * self._u2B
1280 g.seek(b, _MODS.os.SEEK_SET)
1281 return b # position
1282 raise GeoidError('closed file', txt=repr(p.egm)) # IOError
1285class GeoidPGM(_GeoidBase):
1286 '''Geoid height interpolator for I{Karney}'s U{GeographicLib Earth
1287 Gravitational Model (EGM)<https://GeographicLib.SourceForge.io/C++/doc/geoid.html>}
1288 geoid U{egm*.pgm<https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinst>}
1289 datasets but based on C{SciPy} U{RectBivariateSpline<https://docs.SciPy.org/doc/scipy/
1290 reference/generated/scipy.interpolate.RectBivariateSpline.html>}, U{bisplrep/-ev
1291 <https://docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.bisplrep.html>}
1292 or U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.
1293 interp2d.html>} interpolation.
1295 Use any of the U{egm84-, egm96- or egm2008-*.pgm <https://GeographicLib.SourceForge.io/
1296 C++/doc/geoid.html#geoidinst>} datasets. However, unless cropped, an entire C{egm*.pgm}
1297 dataset is loaded into the C{SciPy} interpolator and converted from 2-byte C{int} to
1298 8-byte C{dtype float64}. Therefore, internal memory usage is 4x the U{egm*.pgm
1299 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinst>} file size and may
1300 exceed the available memory, especially with 32-bit Python, see properties C{.nBytes}
1301 and C{.sizeB}.
1302 '''
1303 _cropped = False
1304 _endian = '>u2'
1306 def __init__(self, egm_pgm, crop=None, datum=_WGS84, kind=3, smooth=0, **name):
1307 '''New L{GeoidPGM} interpolator.
1309 @arg egm_pgm: An U{EGM geoid dataset<https://GeographicLib.SourceForge.io/
1310 C++/doc/geoid.html#geoidinst>} file name (C{egm*.pgm}).
1311 @kwarg crop: Optional box to crop B{C{egm_pgm}}, a 4-tuple (C{south, west,
1312 north, east}) or 2-tuple (C{(south, west), (north, east)}),
1313 in C{degrees90} lat- and C{degrees180} longitudes or a 2-tuple
1314 (C{LatLonSW, LatLonNE}) of C{LatLon} instances.
1315 @kwarg datum: Optional grid datum (L{Datum}, L{Ellipsoid}, L{Ellipsoid2} or
1316 L{a_f2Tuple}), overriding C{WGS84}.
1317 @kwarg kind: C{scipy.interpolate} order (C{int}), use 1..5 for U{RectBivariateSpline
1318 <https://docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.
1319 RectBivariateSpline.html>} or -1, -3 or -5 for U{bisplrep/-ev<https://
1320 docs.SciPy.org/doc/scipy/reference/generated/scipy.interpolate.bisplrep.html>}
1321 or U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/scipy.
1322 interpolate.interp2d.html>} C{linear}, C{cubic} respectively C{quintic},
1323 see note for more details.
1324 @kwarg smooth: Smoothing factor for C{B{kind}=1..5} only (C{int}).
1325 @kwarg name: Optional geoid C{B{name}=NN} (C{str}).
1327 @raise GeoidError: EGM dataset B{C{egm_pgm}} issue or invalid B{C{crop}}, B{C{kind}}
1328 or B{C{smooth}}.
1330 @raise ImportError: Package C{numpy} or C{scipy} not found or not installed.
1332 @raise LenError: EGM dataset B{C{egm_pgm}} axis mismatch.
1334 @raise SciPyError: A C{scipy} issue.
1336 @raise SciPyWarning: A C{scipy} warning as exception.
1338 @raise TypeError: Invalid B{C{datum}} or unexpected argument.
1340 @note: Specify C{B{kind}=-1, -3 or -5} to use C{scipy.interpolate.interp2d}
1341 before or C{scipy.interpolate.bisplrep/-ev} since C{Scipy} version 1.14.
1343 @note: The U{GeographicLib egm*.pgm<https://GeographicLib.SourceForge.io/C++/doc/
1344 geoid.html#geoidinst>} file sizes are based on a 2-byte C{int} height
1345 converted to 8-byte C{dtype float64} for C{scipy} interpolators. Therefore,
1346 internal memory usage is 4 times the C{egm*.pgm} file size and may exceed
1347 the available memory, especially with 32-bit Python. To reduce memory
1348 usage, set keyword argument B{C{crop}} to the region of interest. For
1349 example C{B{crop}=(20, -125, 50, -65)} covers the U{conterminous US
1350 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/maps/GEOID12B_CONUS_grids.png>}
1351 (CONUS), less than 3% of the entire C{egm2008-1.pgm} dataset.
1353 @see: Class L{GeoidKarney} and function L{egmGeoidHeights}.
1354 '''
1355 np = self.numpy
1356 self._u2B = np.dtype(self.endian).itemsize
1358 g = self._open(egm_pgm, datum, kind, _name__(**name), smooth)
1359 self._pgm = p = _PGM(g, pgm=egm_pgm, itemsize=self.u2B, sizeB=self.sizeB)
1360 if crop:
1361 g = p._cropped(g, abs(kind) + 1, *self._swne(crop))
1362 if _MODS.internals._version2(np.__version__) < (1, 9):
1363 g = open(g.name, _rb_) # reopen tempfile for numpy 1.8.0-
1364 self._cropped = True
1365 try:
1366 # U{numpy dtype formats are different from Python struct formats
1367 # <https://docs.SciPy.org/doc/numpy-1.15.0/reference/arrays.dtypes.html>}
1368 # read all heights, skipping the PGM header lines, converted to float
1369 hs = self._load(g, self.endian, p.knots, p.skip).reshape(p.nlat, p.nlon) * p.Scale
1370 if p.Offset: # offset
1371 hs = p.Offset + hs
1372 if p.dlat < 0: # flip the rows
1373 hs = np.flipud(hs)
1374 _GeoidBase.__init__(self, hs, p)
1375 except Exception as x:
1376 raise _SciPyIssue(x, _in_, repr(egm_pgm))
1377 finally:
1378 g.close()
1380 def _g2ll2(self, lat, lon):
1381 # convert grid (lat, lon) to earth (lat, lon), un-/cropped
1382 if self._cropped:
1383 lon -= self._lon_of
1384 else:
1385 while lon > _180_0:
1386 lon -= _360_0
1387 return lat, lon
1389 def _ll2g2(self, lat, lon):
1390 # convert earth (lat, lon) to grid (lat, lon), un-/cropped
1391 if self._cropped:
1392 lon += self._lon_of
1393 else:
1394 while lon < 0:
1395 lon += _360_0
1396 return lat, lon
1398 if _FOR_DOCS:
1399 __call__ = _GeoidBase.__call__
1400 height = _GeoidBase.height
1403class _Gpars(_Named):
1404 '''(INTERNAL) Basic geoid parameters.
1405 '''
1406 # interpolator parameters
1407 dlat = 0 # +/- latitude resolution in C{degrees}
1408 dlon = 0 # longitude resolution in C{degrees}
1409 nlat = 1 # number of latitude knots (C{int})
1410 nlon = 0 # number of longitude knots (C{int})
1411 rlat = 0 # +/- latitude resolution in C{float}, 1 / .dlat
1412 rlon = 0 # longitude resolution in C{float}, 1 / .dlon
1413 slat = 0 # nothern- or southern most latitude (C{degrees90})
1414 wlon = 0 # western-most longitude in Eastern lon (C{degrees360})
1416 flon = 0 # forward, earth to grid longitude offset
1417 glon = 0 # reverse, grid to earth longitude offset
1419 knots = 0 # number of knots, nlat * nlon (C{int})
1420 skip = 0 # header bytes to skip (C{int})
1422 def __repr__(self):
1423 t = _COMMASPACE_.join(pairs((a, getattr(self, a)) for
1424 a in dir(self.__class__)
1425 if a[:1].isupper()))
1426 return _COLONSPACE_(self, t)
1428 def __str__(self):
1429 return Fmt.PAREN(self.classname, repr(self.name))
1432class _PGM(_Gpars):
1433 '''(INTERNAL) Parse an C{egm*.pgm} geoid dataset file.
1435 # Geoid file in PGM format for the GeographicLib::Geoid class
1436 # Description WGS84 EGM96, 5-minute grid
1437 # URL https://Earth-Info.NGA.mil/GandG/wgs84/gravitymod/egm96/egm96.html
1438 # DateTime 2009-08-29 18:45:03
1439 # MaxBilinearError 0.140
1440 # RMSBilinearError 0.005
1441 # MaxCubicError 0.003
1442 # RMSCubicError 0.001
1443 # Offset -108
1444 # Scale 0.003
1445 # Origin 90N 0E
1446 # AREA_OR_POINT Point
1447 # Vertical_Datum WGS84
1448 <width> <height>
1449 <pixel>
1450 ...
1451 '''
1452 crop4 = () # 4-tuple (C{south, west, north, east}).
1453 egm = None
1454 glon = 180 # reverse offset, uncropped
1455# pgm = NN # name
1456 sizeB = 0
1457 u2B = 2 # item size of grid height (C{int}).
1459 @staticmethod
1460 def _llstr2floats(latlon):
1461 # llstr to (lat, lon) floats
1462 lat, lon = latlon.split()
1463 return _MODS.dms.parseDMS2(lat, lon)
1465 # PGM file attributes, CamelCase but not .istitle()
1466 AREA_OR_POINT = str
1467 DateTime = str
1468 Description = str # 'WGS84 EGM96, 5-minute grid'
1469 Geoid = str # 'file in PGM format for the GeographicLib::Geoid class'
1470 MaxBilinearError = float
1471 MaxCubicError = float
1472 Offset = float
1473 Origin = _llstr2floats
1474 Pixel = 0
1475 RMSBilinearError = float
1476 RMSCubicError = float
1477 Scale = float
1478 URL = str # 'https://Earth-Info.NGA.mil/GandG/wgs84/...'
1479 Vertical_Datum = str
1481 def __init__(self, g, pgm=NN, itemsize=0, sizeB=0): # MCCABE 22
1482 '''(INTERNAL) New C{_PGM} parsed C{egm*.pgm} geoid dataset.
1483 '''
1484 self.name = pgm # geoid file name
1485 if itemsize:
1486 self._u2B = itemsize
1487 if sizeB:
1488 self.sizeB = sizeB
1490 t = g.readline() # make sure newline == '\n'
1491 if t != b'P5\n' and t.strip() != b'P5':
1492 raise self._Errorf(_format_, _header_, t)
1494 while True: # read all # Attr ... lines,
1495 try: # ignore empty ones or comments
1496 t = g.readline().strip()
1497 if t.startswith(_bHASH_):
1498 t = t.lstrip(_bHASH_).lstrip()
1499 a, v = map(_ub2str, t.split(None, 1))
1500 f = getattr(_PGM, a, None)
1501 if callable(f) and a[:1].isupper():
1502 setattr(self, a, f(v))
1503 elif t:
1504 break
1505 except (TypeError, ValueError):
1506 raise self._Errorf(_format_, 'Attr', t)
1507 else: # should never get here
1508 raise self._Errorf(_format_, _header_, g.tell())
1510 try: # must be (even) width and (odd) height
1511 nlon, nlat = map(int, t.split())
1512 if nlon < 2 or nlon > (360 * 60) or isodd(nlon) or \
1513 nlat < 2 or nlat > (181 * 60) or not isodd(nlat):
1514 raise ValueError
1515 except (TypeError, ValueError):
1516 raise self._Errorf(_format_, _SPACE_(_width_, _height_), t)
1518 try: # must be 16 bit pixel height
1519 t = g.readline().strip()
1520 self.Pixel = int(t)
1521 if not 255 < self.Pixel < 65536: # >u2 or >H only
1522 raise ValueError
1523 except (TypeError, ValueError):
1524 raise self._Errorf(_format_, 'pixel', t)
1526 for a in dir(_PGM): # set undefined # Attr ... to None
1527 if a[:1].isupper() and callable(getattr(self, a)):
1528 setattr(self, a, None)
1530 if self.Origin is None:
1531 raise self._Errorf(_format_, 'Origin', self.Origin)
1532 if self.Offset is None or self.Offset > 0:
1533 raise self._Errorf(_format_, 'Offset', self.Offset)
1534 if self.Scale is None or self.Scale < EPS:
1535 raise self._Errorf(_format_, 'Scale', self.Scale)
1537 self.skip = g.tell()
1538 self.knots = nlat * nlon
1540 self.nlat, self.nlon = nlat, nlon
1541 self.slat, self.wlon = self.Origin
1542 # note, negative .dlat and .rlat since rows
1543 # are from .slat 90N down in decreasing lat
1544 self.dlat, self.dlon = _180_0 / (1 - nlat), _360_0 / nlon
1545 self.rlat, self.rlon = (1 - nlat) / _180_0, nlon / _360_0
1547 # grid corners in earth (lat, lon), .slat = 90, .dlat < 0
1548 n = float(self.slat)
1549 s = n + self.dlat * (nlat - 1)
1550 w = self.wlon - self.glon
1551 e = w + self.dlon * nlon
1552 self.crop4 = s, w, n, e
1554 n = self.sizeB - self.skip
1555 if n > 0 and n != (self.knots * self.u2B):
1556 raise self._Errorf('%s(%s x %s != %s)', _assert_, nlat, nlon, n)
1558 def _cropped(self, g, k1, south, west, north, east): # MCCABE 15
1559 '''Crop the geoid to (south, west, north, east) box.
1560 '''
1561 # flon offset for both west and east
1562 f = 360 if west < 0 else 0
1563 # earth (lat, lon) to grid indices (y, x),
1564 # note y is decreasing, i.e. n < s
1565 s, w = self._lle2yx2(south, west, f)
1566 n, e = self._lle2yx2(north, east, f)
1567 s += 1 # s > n
1568 e += 1 # e > w
1570 hi, wi = self.nlat, self.nlon
1571 # handle special cases
1572 if (s - n) > hi:
1573 n, s = 0, hi # entire lat range
1574 if (e - w) > wi:
1575 w, e, f = 0, wi, 180 # entire lon range
1576 if s == hi and w == n == 0 and e == wi:
1577 return g # use entire geoid as-is
1579 if (e - w) < k1 or (s - n) < (k1 + 1):
1580 raise self._Errorf(_format_, 'swne', (north - south, east - west))
1582 if e > wi > w: # wrap around
1583 # read w..wi and 0..e
1584 r, p = (wi - w), (e - wi)
1585 elif e > w:
1586 r, p = (e - w), 0
1587 else:
1588 raise self._Errorf('%s(%s < %s)', _assert_, w, e)
1590 # convert to bytes
1591 r *= self.u2B
1592 p *= self.u2B
1593 q = wi * self.u2B # stride
1594 # number of rows and cols to skip from
1595 # the original (.slat, .wlon) origin
1596 z = self.skip + (n * wi + w) * self.u2B
1597 # sanity check
1598 if r < 2 or p < 0 or q < 2 or z < self.skip \
1599 or z > self.sizeB:
1600 raise self._Errorf(_format_, _assert_, (r, p, q, z))
1602 # can't use _BytesIO since numpy
1603 # needs .fileno attr in .fromfile
1604 t, c = 0, self._tmpfile()
1605 # reading (s - n) rows, forward
1606 for y in range(n, s): # PYCHOK y unused
1607 g.seek(z, _MODS.os.SEEK_SET)
1608 # Python 2 tmpfile.write returns None
1609 t += c.write(g.read(r)) or r
1610 if p: # wrap around to start of row
1611 g.seek(-q, _MODS.os.SEEK_CUR)
1612 # assert(g.tell() == (z - w * self.u2B))
1613 # Python 2 tmpfile.write returns None
1614 t += c.write(g.read(p)) or p
1615 z += q
1616 c.flush()
1617 g.close()
1619 s -= n # nlat
1620 e -= w # nlon
1621 k = s * e # knots
1622 z = k * self.u2B
1623 if t != z:
1624 raise self._Errorf('%s(%s != %s) %s', _assert_, t, z, self)
1626 # update the _Gpars accordingly, note attributes
1627 # .dlat, .dlon, .rlat and .rlon remain unchanged
1628 self.slat += n * self.dlat
1629 self.wlon += w * self.dlon
1630 self.nlat = s
1631 self.nlon = e
1632 self.flon = self.glon = f
1634 self.crop4 = south, west, north, east
1635 self.knots = k
1636 self.skip = 0 # no header lines in c
1638 c.seek(0, _MODS.os.SEEK_SET)
1639 # c = open(c.name, _rb_) # reopen for numpy 1.8.0-
1640 return c
1642 def _Errorf(self, fmt, *args): # PYCHOK no cover
1643 t = fmt % args
1644 e = self.pgm or NN
1645 if e:
1646 t = _SPACE_(t, _in_, repr(e))
1647 return PGMError(t)
1649 def _lle2yx2(self, lat, lon, flon):
1650 # earth (lat, lon) to grid indices (y, x)
1651 # with .dlat decreasing from 90N .slat
1652 lat -= self.slat
1653 lon += flon - self.wlon
1654 return (min(self.nlat - 1, max(0, int(lat * self.rlat))),
1655 max(0, int(lon * self.rlon)))
1657 def _tmpfile(self):
1658 # create a tmpfile to hold the cropped geoid grid
1659 try:
1660 from tempfile import NamedTemporaryFile as tmpfile
1661 except ImportError: # Python 2.7.16-
1662 from _MODS.os import tmpfile
1663 t = _MODS.os.path.basename(self.pgm)
1664 t = _MODS.os.path.splitext(t)[0]
1665 f = tmpfile(mode='w+b', prefix=t or 'egm')
1666 f.seek(0, _MODS.os.SEEK_SET) # force overwrite
1667 return f
1669 @Property_RO
1670 def pgm(self):
1671 '''Get the geoid file name (C{str}).
1672 '''
1673 return self.name
1676class PGMError(GeoidError):
1677 '''An issue while parsing or cropping an C{egm*.pgm} geoid dataset.
1678 '''
1679 pass
1682def egmGeoidHeights(GeoidHeights_dat):
1683 '''Generate geoid U{egm*.pgm<https://GeographicLib.SourceForge.io/
1684 C++/doc/geoid.html#geoidinst>} height tests from U{GeoidHeights.dat
1685 <https://SourceForge.net/projects/geographiclib/files/testdata/>}
1686 U{Test data for Geoids<https://GeographicLib.SourceForge.io/C++/doc/
1687 geoid.html#testgeoid>}.
1689 @arg GeoidHeights_dat: The un-gz-ed C{GeoidHeights.dat} file
1690 (C{str} or C{file} handle).
1692 @return: For each test, yield a L{GeoidHeight5Tuple}C{(lat, lon,
1693 egm84, egm96, egm2008)}.
1695 @raise GeoidError: Invalid B{C{GeoidHeights_dat}}.
1697 @note: Function L{egmGeoidHeights} is used to test the geoids
1698 L{GeoidKarney} and L{GeoidPGM}, see PyGeodesy module
1699 C{test/testGeoids.py}.
1700 '''
1701 dat = GeoidHeights_dat
1702 if isinstance(dat, bytes):
1703 dat = _BytesIO(dat)
1705 try:
1706 dat.seek(0, _MODS.os.SEEK_SET) # reset
1707 except AttributeError as x:
1708 raise GeoidError(GeoidHeights_dat=type(dat), cause=x)
1710 for t in dat.readlines():
1711 t = t.strip()
1712 if t and not t.startswith(_bHASH_):
1713 lat, lon, egm84, egm96, egm2008 = map(float, t.split())
1714 while lon > _180_0: # EasternLon to earth lon
1715 lon -= _360_0
1716 yield GeoidHeight5Tuple(lat, lon, egm84, egm96, egm2008)
1719__all__ += _ALL_DOCS(_GeoidBase)
1721if __name__ == _DMAIN_: # MCCABE 14
1723 from pygeodesy.internals import printf, _secs2str, _versions, _sys
1724 from time import time
1726 _crop = {}
1727 _GeoidEGM = GeoidKarney
1728 _kind = 3
1730 geoids = _sys.argv[1:]
1731 while geoids:
1732 G = geoids.pop(0)
1733 g = G.lower()
1735 if '-crop'.startswith(g):
1736 _crop = dict(crp=(20, -125, 50, -65)) # CONUS
1738 elif '-egm96'.startswith(g):
1739 _GeoidEGM = GeoidEGM96
1741 elif '-karney'.startswith(g):
1742 _GeoidEGM = GeoidKarney
1744 elif '-kind'.startswith(g):
1745 _kind = int(geoids.pop(0))
1747 elif '-pgm'.startswith(g):
1748 _GeoidEGM = GeoidPGM
1750 elif _isin(g[-4:], '.pgm', '.grd'):
1751 g = _GeoidEGM(G, kind=_kind, **_crop)
1752 t = time()
1753 _ = g.highest()
1754 t = _secs2str(time() - t)
1755 printf('%s: %s (%s)', g.toStr(), t, _versions(), nl=1, nt=1)
1756 t = g.pgm
1757 if t:
1758 printf(repr(t), nt=1)
1759 # <https://GeographicLib.SourceForge.io/cgi-bin/GeoidEval>:
1760 # The height of the EGM96 geoid at Timbuktu
1761 # echo 16:46:33N 3:00:34W | GeoidEval
1762 # => 28.7068 -0.02e-6 -1.73e-6
1763 # The 1st number is the height of the geoid, the 2nd and
1764 # 3rd are its slopes in northerly and easterly direction
1765 t = 'Timbuktu %s' % (g,)
1766 k = {'egm84-15.pgm': '31.2979',
1767 'egm96-5.pgm': '28.7067',
1768 'egm2008-1.pgm': '28.7880'}.get(g.name.lower(), '28.7880')
1769 ll = _MODS.dms.parseDMS2('16:46:33N', '3:00:34W', sep=':')
1770 for ll in (ll, (16.776, -3.009),):
1771 try:
1772 h, ll = g.height(*ll), fstr(ll, prec=6)
1773 printf('%s.height(%s): %.4F vs %s', t, ll, h, k)
1774 except (GeoidError, RangeError) as x:
1775 printf(_COLONSPACE_(t, str(x)))
1777 elif _isin(g[-4:], '.bin'):
1778 g = GeoidG2012B(G, kind=_kind)
1779 printf(g.toStr())
1781 else:
1782 raise GeoidError(grid=repr(G))
1784_I = int # PYCHOK unused _I
1785del _intCs, _T, _T0s12 # trash ints cache and map
1788# <https://GeographicLib.SourceForge.io/cgi-bin/GeoidEval>
1789# _lowerleft = -90, -179, -30.1500 # egm2008-1.pgm
1790# _lowerleft = -90, -179, -29.5350 # egm96-5.pgm
1791# _lowerleft = -90, -179, -29.7120 # egm84-15.pgm
1793# _center = 0, 0, 17.2260 # egm2008-1.pgm
1794# _center = 0, 0, 17.1630 # egm96-5.pgm
1795# _center = 0, 0, 18.3296 # egm84-15.pgm
1797# _upperright = 90, 180, 14.8980 # egm2008-1.pgm
1798# _upperright = 90, 180, 13.6050 # egm96-5.pgm
1799# _upperright = 90, 180, 13.0980 # egm84-15.pgm
1802# % python3.12 -m pygeodesy.geoids -egm96 ../testGeoids/WW15MGH.GRD
1803#
1804# 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)
1805#
1806# Timbuktu GeoidEGM96('WW15MGH.GRD').height(16.775833, -3.009444): 28.7073 vs 28.7880
1807# Timbuktu GeoidEGM96('WW15MGH.GRD').height(16.776, -3.009): 28.7072 vs 28.7880
1810# % python3.12 -m pygeodesy.geoids -Karney ../testGeoids/egm*.pgm
1811#
1812# 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)
1813#
1814# _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'
1815#
1816# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1817# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1818#
1819# 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)
1820#
1821# _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'
1822#
1823# Timbuktu GeoidKarney('egm84-15.pgm').height(16.775833, -3.009444): 31.2983 vs 31.2979
1824# Timbuktu GeoidKarney('egm84-15.pgm').height(16.776, -3.009): 31.2979 vs 31.2979
1825#
1826# 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)
1827#
1828# _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'
1829#
1830# Timbuktu GeoidKarney('egm96-5.pgm').height(16.775833, -3.009444): 28.7068 vs 28.7067
1831# Timbuktu GeoidKarney('egm96-5.pgm').height(16.776, -3.009): 28.7067 vs 28.7067
1834# % python3.8 -m pygeodesy.geoids -egm96 ../testGeoids/WW15MGH.GRD
1835#
1836# 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)
1837#
1838# Timbuktu GeoidEGM96('WW15MGH.GRD').height(16.775833, -3.009444): 28.7073 vs 28.7880
1839# Timbuktu GeoidEGM96('WW15MGH.GRD').height(16.776, -3.009): 28.7072 vs 28.7880
1842# % python3.8 -m pygeodesy.geoids -Karney ../testGeoids/egm*.pgm
1843#
1844# 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)
1845#
1846# _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'
1847#
1848# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1849# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1850#
1851# 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)
1852#
1853# _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'
1854#
1855# Timbuktu GeoidKarney('egm84-15.pgm').height(16.775833, -3.009444): 31.2983 vs 31.2979
1856# Timbuktu GeoidKarney('egm84-15.pgm').height(16.776, -3.009): 31.2979 vs 31.2979
1857#
1858# 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)
1859#
1860# _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'
1861#
1862# Timbuktu GeoidKarney('egm96-5.pgm').height(16.775833, -3.009444): 28.7068 vs 28.7067
1863# Timbuktu GeoidKarney('egm96-5.pgm').height(16.776, -3.009): 28.7067 vs 28.7067
1866# % python2 -m pygeodesy.geoids -Karney ../testGeoids/egm*.pgm
1867#
1868# 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)
1869#
1870# _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'
1871#
1872# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1873# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1874#
1875# 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)
1876#
1877# _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'
1878#
1879# Timbuktu GeoidKarney('egm84-15.pgm').height(16.775833, -3.009444): 31.2983 vs 31.2979
1880# Timbuktu GeoidKarney('egm84-15.pgm').height(16.776, -3.009): 31.2979 vs 31.2979
1881#
1882# 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)
1883#
1884# _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'
1885#
1886# Timbuktu GeoidKarney('egm96-5.pgm').height(16.775833, -3.009444): 28.7068 vs 28.7067
1887# Timbuktu GeoidKarney('egm96-5.pgm').height(16.776, -3.009): 28.7067 vs 28.7067
1890# % python3.12 -m pygeodesy.geoids -PGM ../testGeoids/egm*.pgm
1891#
1892# 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)
1893#
1894# _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'
1895#
1896# Timbuktu GeoidPGM('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1897# Timbuktu GeoidPGM('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1898#
1899# 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)
1900#
1901# _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'
1902#
1903# Timbuktu GeoidPGM('egm84-15.pgm').height(16.775833, -3.009444): 31.2979 vs 31.2979
1904# Timbuktu GeoidPGM('egm84-15.pgm').height(16.776, -3.009): 31.2975 vs 31.2979
1905#
1906# 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)
1907#
1908# _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'
1909#
1910# Timbuktu GeoidPGM('egm96-5.pgm').height(16.775833, -3.009444): 28.7065 vs 28.7067
1911# Timbuktu GeoidPGM('egm96-5.pgm').height(16.776, -3.009): 28.7064 vs 28.7067
1913# **) MIT License
1914#
1915# Copyright (C) 2016-2025 -- mrJean1 at Gmail -- All Rights Reserved.
1916#
1917# Permission is hereby granted, free of charge, to any person obtaining a
1918# copy of this software and associated documentation files (the "Software"),
1919# to deal in the Software without restriction, including without limitation
1920# the rights to use, copy, modify, merge, publish, distribute, sublicense,
1921# and/or sell copies of the Software, and to permit persons to whom the
1922# Software is furnished to do so, subject to the following conditions:
1923#
1924# The above copyright notice and this permission notice shall be included
1925# in all copies or substantial portions of the Software.
1926#
1927# THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
1928# OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
1929# FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
1930# THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR
1931# OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE,
1932# ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR
1933# OTHER DEALINGS IN THE SOFTWARE.