Coverage for pygeodesy/geoids.py: 96%
666 statements
« prev ^ index » next coverage.py v7.2.2, created at 2023-04-05 15:46 -0400
« prev ^ index » next coverage.py v7.2.2, created at 2023-04-05 15:46 -0400
2# -*- coding: utf-8 -*-
4u'''Geoid models and geoid height interpolations.
6Classes L{GeoidG2012B}, L{GeoidKarney} and L{GeoidPGM} to interpolate the
7height of various U{geoid<https://WikiPedia.org/wiki/Geoid>}s at C{LatLon}
8locations or separate lat-/longitudes using different interpolation methods
9and 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.
13The L{GeoidG2012B} and L{GeoidPGM} interpolators both depend on U{scipy
14<https://SciPy.org>} and U{numpy<https://PyPI.org/project/numpy>} and
15require those packages to be installed.
17In addition, each geoid interpolator needs C{grid knots} (down)loaded from
18a C{geoid} model file, I{specific to the interpolator}, more details below
19and in the documentation of the interpolator class. For each interpolator,
20there are several interpolation choices, like I{linear}, I{cubic}, etc.
22Typical usage
23=============
251. Choose one of the interpolator classes L{GeoidG2012B}, L{GeoidKarney}
26or L{GeoidPGM} and download a C{geoid} model file, containing locations with
27known heights also referred to as the C{grid knots}. See the documentation
28of interpolator class for references to available C{grid} models.
30C{>>> from pygeodesy import GeoidG2012B # or -Karney or -PGM as GeoidXyz}
322. Instantiate an interpolator with the C{geoid} model file and use keyword
33arguments to select different interpolation options
35C{>>> ginterpolator = GeoidXyz(geoid_model_file, **options)}
373. Get the interpolated geoid height of other C{LatLon} location(s) with
39C{>>> ll = LatLon(1, 2, ...)}
40C{>>> h = ginterpolator(ll)}
42or
44C{>>> h0, h1, h2, ... = ginterpolator(ll0, ll1, ll2, ...)}
46or a list, tuple, generator, etc. of C{LatLon}s
48C{>>> hs = ginterpolator(lls)}
504. For separate lat- and longitudes invoke the C{.height} method as
52C{>>> h = ginterpolator.height(lat, lon)}
54or as 2 lists, 2 tuples, etc.
56C{>>> hs = ginterpolator.height(lats, lons)}
585. An example is in U{issue #64<https://GitHub.com/mrJean1/PyGeodesy/issues/64>},
59courtesy of SBFRF.
61@note: Classes L{GeoidG2012B} and L{GeoidPGM} require both U{numpy
62 <https://PyPI.org/project/numpy>} and U{scipy<https://PyPI.org/project/scipy>}
63 to be installed.
65@note: Errors from C{scipy} are raised as L{SciPyError}s. Warnings issued by
66 C{scipy} can be thrown as L{SciPyWarning} exceptions, provided Python
67 C{warnings} are filtered accordingly, see L{SciPyWarning}.
69@see: I{Karney}'s U{GeographicLib<https://GeographicLib.SourceForge.io/html/index.html>},
70 U{Geoid height<https://GeographicLib.SourceForge.io/html/geoid.html>} and U{Installing
71 the Geoid datasets<https://GeographicLib.SourceForge.io/html/geoid.html#geoidinst>},
72 U{SciPy<https://docs.SciPy.org/doc/scipy/reference/interpolate.html>} interpolation
73 U{RectBivariateSpline<https://docs.SciPy.org/doc/scipy/reference/generated/scipy.
74 interpolate.RectBivariateSpline.html>} and U{interp2d<https://docs.SciPy.org/doc/scipy/
75 reference/generated/scipy.interpolate.interp2d.html>}, functions L{elevations.elevation2}
76 and L{elevations.geoidHeight2}, U{I{Ellispoid vs Orthometric Elevations}<https://
77 www.YouTube.com/watch?v=dX6a6kCk3Po>} and U{I{Pitfalls Related to Ellipsoid Height
78 and Height Above Mean Sea Level (AMSL)}<https://Wiki.ROS.org/mavros#mavros.2FPlugins.
79 Avoiding_Pitfalls_Related_to_Ellipsoid_Height_and_Height_Above_Mean_Sea_Level>}.
80'''
81# make sure int/int division yields float quotient, see .basics
82from __future__ import division as _; del _ # PYCHOK semicolon
84from pygeodesy.basics import len2, map1, map2, isodd, ub2str as _ub2str
85from pygeodesy.constants import EPS, _float as _F, _0_0, _1_0, _180_0, _360_0
86from pygeodesy.datums import _ellipsoidal_datum, _WGS84
87from pygeodesy.dms import parseDMS2
88from pygeodesy.errors import _incompatible, LenError, RangeError, _SciPyIssue
89from pygeodesy.fmath import favg, Fdot, fdot, Fhorner, frange
90from pygeodesy.heights import _allis2, _ascalar, _HeightBase, HeightError
91from pygeodesy.interns import NN, _4_, _COLONSPACE_, _COMMASPACE_, _cubic_, \
92 _DOT_, _E_, _height_, _in_, _kind_, _knots_, \
93 _lat_, _linear_, _lon_, _mean_, _N_, _n_a_, _not_, \
94 _numpy_, _on_, _outside_, _S_, _s_, _scipy_, \
95 _SPACE_, _stdev_, _supported_, _tbd_, _W_, _width_
96from pygeodesy.lazily import _ALL_DOCS, _ALL_LAZY, _ALL_MODS as _MODS, _FOR_DOCS
97from pygeodesy.named import _Named, _NamedTuple, notOverloaded
98# from pygeodesy.namedTuples import LatLon3Tuple # _MODS
99from pygeodesy.props import Property_RO, property_RO
100from pygeodesy.streprs import attrs, Fmt, fstr, pairs
101from pygeodesy.units import Height, Int_, Lat, Lon
103from math import floor
104import os.path as _os_path
105from os import SEEK_CUR as _SEEK_CUR, SEEK_SET as _SEEK_SET
106from struct import calcsize as _calcsize, unpack as _unpack
108try:
109 from StringIO import StringIO as _BytesIO # reads bytes
110 _ub2str = str # PYCHOK convert bytes to str for egm*.pgm text
112except ImportError: # Python 3+
113 from io import BytesIO as _BytesIO # PYCHOK expected
115__all__ = _ALL_LAZY.geoids
116__version__ = '22.10.23'
118_assert_ = 'assert'
119_bHASH_ = b'#'
120_endian_ = 'endian'
121_format_ = '%s %r'
122_header_ = 'header'
123# temporarily hold a single instance for each int value
124_intCs = {}
125_interp2d_ks = {-2: _linear_,
126 -3: _cubic_,
127 -5: 'quintic'}
128_lli_ = 'lli'
129_non_increasing_ = 'non-increasing'
130_rb_ = 'rb'
133class _GeoidBase(_HeightBase):
134 '''(INTERNAL) Base class for C{Geoid...}s.
135 '''
136 _cropped = None
137 _datum = _WGS84
138 _egm = None # open C{egm*.pgm} geoid file
139 _endian = _tbd_
140 _geoid = _n_a_
141 _hs_y_x = None # numpy 2darray, row-major order
142 _interp2d = None # interp2d interpolation
143 _kind = 3 # order for interp2d, RectBivariateSpline
144 _knots = 0 # nlat * nlon
145 _mean = None # fixed in GeoidKarney
146# _name = NN # _Named
147 _nBytes = 0 # numpy size in bytes, float64
148 _pgm = None # PGM attributes, C{_PGM} or C{None}
149 _sizeB = 0 # geoid file size in bytes
150 _smooth = 0 # used only for RectBivariateSpline
151 _stdev = None # fixed in GeoidKarney
152 _u2B = 0 # np.itemsize or undefined
154 _lat_d = _0_0 # increment, +tive
155 _lat_lo = _0_0 # lower lat, south
156 _lat_hi = _0_0 # upper lat, noth
157 _lon_d = _0_0 # increment, +tive
158 _lon_lo = _0_0 # left lon, west
159 _lon_hi = _0_0 # right lon, east
160 _lon_of = _0_0 # forward lon offset
161 _lon_og = _0_0 # reverse lon offset
163 _center = None # (lat, lon, height)
164 _yx_hits = None # cache hits, ala Karney
166 def __init__(self, hs, p):
167 '''(INTERNAL) Set up the grid axes, the C{SciPy} interpolator
168 and several internal geoid attributes.
170 @arg hs: Grid knots with known height (C{numpy 2darray}).
171 @arg p: The C{slat, wlon, nlat, nlon, dlat, dlon} and
172 other geoid parameters (C{INTERNAL}).
174 @raise GeoidError: Incompatible grid B{C{hs}} shape or
175 invalid B{C{kind}}.
177 @raise LenError: Mismatch grid B{C{hs}} axis.
179 @raise SciPyError: A C{scipy.interpolate.inter2d} or
180 C{-.RectBivariateSpline} issue.
182 @raise SciPyWarning: A C{scipy.interpolate.inter2d} or
183 C{-.RectBivariateSpline} warning as
184 exception.
185 '''
186 spi = self.scipy_interpolate
187 # for 2d scipy.interpolate.interp2d(xs, ys, hs, ...) and
188 # scipy.interpolate.RectBivariateSpline(ys, xs, hs, ...)
189 # require the shape of hs to be (len(ys), len(xs)), note
190 # the different (xs, ys, ...) and (ys, xs, ...) orders
191 if (p.nlat, p.nlon) != hs.shape:
192 raise GeoidError(shape=hs.shape, txt=_incompatible((p.nlat, p.nlon)))
194 # both axes and bounding box
195 ys, self._lat_d = self._gaxis2(p.slat, p.dlat, p.nlat, _lat_ + _s_)
196 xs, self._lon_d = self._gaxis2(p.wlon, p.dlon, p.nlon, _lon_ + _s_)
198 bb = ys[0], ys[-1], xs[0], xs[-1] + p.dlon # fudge lon_hi
199 # geoid grids are typically stored in row-major order, some
200 # with rows (90..-90) reversed and columns (0..360) wrapped
201 # to Easten longitude, 0 <= east < 180 and 180 <= west < 360
202 k = self.kind
203 if k in _interp2d_ks:
204 self._interp2d = spi.interp2d(xs, ys, hs, kind=_interp2d_ks[k])
205 elif 1 <= k <= 5:
206 self._ev = spi.RectBivariateSpline(ys, xs, hs, bbox=bb, ky=k, kx=k,
207 s=self._smooth).ev
208 else:
209 raise GeoidError(kind=k)
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
214 self._lon_of = float(p.flon) # forward offset
215 self._lon_og = float(p.glon) # reverse offset
216 # shrink the box by 1 unit on every side
217 # bb += self._lat_d, -self._lat_d, self._lon_d, -self._lon_d
218 self._lat_lo = float(bb[0])
219 self._lat_hi = float(bb[1])
220 self._lon_lo = float(bb[2] - p.glon)
221 self._lon_hi = float(bb[3] - p.glon)
223 def __call__(self, *llis):
224 '''Interpolate the geoid height for one or several locations.
226 @arg llis: The location or locations (C{LatLon}, ... or
227 C{LatLon}s).
229 @return: A single interpolated geoid height (C{float}) or
230 a list or tuple of interpolated geoid heights
231 (C{float}s).
233 @raise GeoidError: Insufficient number of B{C{llis}}, an
234 invalid B{C{lli}} or the C{egm*.pgm}
235 geoid file is closed.
237 @raise RangeError: An B{C{lli}} is outside this geoid's lat-
238 or longitude range.
240 @raise SciPyError: A C{scipy.interpolate.inter2d} or
241 C{-.RectBivariateSpline} issue.
243 @raise SciPyWarning: A C{scipy.interpolate.inter2d} or
244 C{-.RectBivariateSpline} warning as
245 exception.
246 '''
247 return self._called(llis, True)
249 def __enter__(self):
250 '''Open context.
251 '''
252 return self
254 def __exit__(self, *unused): # PYCHOK exc_type, exc_value, exc_traceback)
255 '''Close context.
256 '''
257 self.close()
258 # return None # XXX False
260 def __repr__(self):
261 return self.toStr()
263 def __str__(self):
264 return Fmt.PAREN(self.classname, repr(self.name))
266 def _called(self, llis, scipy):
267 # handle __call__
268 _as, llis = _allis2(llis, Error=GeoidError)
269 try:
270 hs = []
271 for i, lli in enumerate(llis):
272 hs.append(self._hGeoid(lli.lat, lli.lon))
273 return _as(hs)
275 except (GeoidError, RangeError) as x:
276 # XXX avoid str(LatLon()) degree symbols
277 t = _lli_ if _as is _ascalar else Fmt.SQUARE(llis=i)
278 lli = fstr((lli.lat, lli.lon), strepr=repr)
279 raise type(x)(t, lli, cause=x)
280 except Exception as x:
281 if scipy and self.scipy:
282 raise _SciPyIssue(x)
283 else:
284 raise
286 @Property_RO
287 def _center(self):
288 ''' Cache for method L{center}.
289 '''
290 return self._llh3(favg(self._lat_lo, self._lat_hi),
291 favg(self._lon_lo, self._lon_hi))
293 def center(self, LatLon=None):
294 '''Return the center location and height of this geoid.
296 @kwarg LatLon: Optional class to return the location and height
297 (C{LatLon}) or C{None}.
299 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
300 lon, height)} otherwise a B{C{LatLon}} instance
301 with the lat-, longitude and height of the grid
302 center location.
303 '''
304 return self._llh3LL(self._center, LatLon)
306 def close(self):
307 '''Close the C{egm*.pgm} geoid file if open (and applicable).
308 '''
309 if not self.closed:
310 self._egm.close()
311 self._egm = None
313 @property_RO
314 def closed(self):
315 '''Get the C{egm*.pgm} geoid file status.
316 '''
317 return self._egm is None
319 @Property_RO
320 def cropped(self):
321 '''Is geoid cropped (C{bool} or C{None} if crop not supported).
322 '''
323 return self._cropped
325 @Property_RO
326 def dtype(self):
327 '''Get the grid C{scipy} U{dtype<https://docs.SciPy.org/doc/numpy/
328 reference/generated/numpy.ndarray.dtype.html>} (C{numpy.dtype}).
329 '''
330 return self._hs_y_x.dtype
332 @Property_RO
333 def endian(self):
334 '''Get the geoid endianess and U{dtype<https://docs.SciPy.org/
335 doc/numpy/reference/generated/numpy.dtype.html>} (C{str}).
336 '''
337 return self._endian
339 def _ev(self, y, x): # PYCHOK expected
340 # only used for .interpolate.interp2d, but
341 # overwritten for .RectBivariateSpline,
342 # note (y, x) must be flipped!
343 return self._interp2d(x, y)
345 def _gaxis2(self, lo, d, n, name):
346 # build grid axis, hi = lo + (n - 1) * d
347 m, a = len2(frange(lo, n, d))
348 if m != n:
349 raise LenError(self.__class__, grid=m, **{name: n})
350 if d < 0:
351 d, a = -d, list(reversed(a))
352 for i in range(1, m):
353 e = a[i] - a[i-1]
354 if e < EPS: # non-increasing axis
355 i = Fmt.SQUARE(name, i)
356 raise GeoidError(i, e, txt=_non_increasing_)
357 return self.numpy.array(a), d
359 def _g2ll2(self, lat, lon): # PYCHOK no cover
360 notOverloaded(self, lat, lon)
362 def _gyx2g2(self, y, x):
363 # convert grid (y, x) indices to grid (lat, lon)
364 return ((self._lat_lo + self._lat_d * y),
365 (self._lon_lo + self._lon_of + self._lon_d * x))
367 def height(self, lats, lons):
368 '''Interpolate the geoid height for one or several lat-/longitudes.
370 @arg lats: Latitude or latitudes (C{degrees} or C{degrees}s).
371 @arg lons: Longitude or longitudes (C{degrees} or C{degrees}s).
373 @return: A single interpolated geoid height (C{float}) or a
374 list of interpolated geoid heights (C{float}s).
376 @raise GeoidError: Insufficient or non-matching number of
377 B{C{lats}} and B{C{lons}}.
379 @raise RangeError: A B{C{lat}} or B{C{lon}} is outside this
380 geoid's lat- or longitude range.
382 @raise SciPyError: A C{scipy.interpolate.inter2d} or
383 C{-.RectBivariateSpline} issue.
385 @raise SciPyWarning: A C{scipy.interpolate.inter2d} or
386 C{-.RectBivariateSpline} warning as
387 exception.
388 '''
389 return _HeightBase._height(self, lats, lons, Error=GeoidError)
391 def _hGeoid(self, lat, lon):
392 out = self.outside(lat, lon)
393 if out:
394 lli = fstr((lat, lon), strepr=repr)
395 raise RangeError(lli=lli, txt=_SPACE_(_outside_, _on_, out))
396 return float(self._ev(*self._ll2g2(lat, lon)))
398 @Property_RO
399 def _highest(self):
400 '''(INTERNAL) Cache for L{highest} method.
401 '''
402 return self._llh3minmax(True)
404 def highest(self, LatLon=None, **unused):
405 '''Return the location and largest height of this geoid.
407 @kwarg LatLon: Optional class to return the location and height
408 (C{LatLon}) or C{None}.
410 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
411 lon, height)} otherwise a B{C{LatLon}} instance
412 with the lat-, longitude and height of the highest
413 grid location.
414 '''
415 return self._llh3LL(self._highest, LatLon)
417 @Property_RO
418 def hits(self):
419 '''Get the number of cache hits (C{int} or C{None}).
420 '''
421 return self._yx_hits
423 @Property_RO
424 def kind(self):
425 '''Get the interpolator kind and order (C{int}).
426 '''
427 return self._kind
429 @Property_RO
430 def knots(self):
431 '''Get the number of grid knots (C{int}).
432 '''
433 return self._knots
435 def _ll2g2(self, lat, lon): # PYCHOK no cover
436 notOverloaded(self, lat, lon)
438 @property_RO
439 def _LL3T(self):
440 '''(INTERNAL) Get L{LatLon3Tuple} once.
441 '''
442 T = _MODS.namedTuples.LatLon3Tuple
443 _GeoidBase._LL3T = T # overwrite poperty_RO
444 return T
446 def _llh3(self, lat, lon):
447 return self._LL3T(lat, lon, self._hGeoid(lat, lon), name=self.name)
449 def _llh3LL(self, llh, LatLon):
450 return llh if LatLon is None else self._xnamed(LatLon(*llh))
452 def _llh3minmax(self, highest=True, *unused):
453 hs, np = self._hs_y_x, self.numpy
454 # <https://docs.SciPy.org/doc/numpy/reference/generated/
455 # numpy.argmin.html#numpy.argmin>
456 arg = np.argmax if highest else np.argmin
457 y, x = np.unravel_index(arg(hs, axis=None), hs.shape)
458 return self._g2ll2(*self._gyx2g2(y, x)) + (float(hs[y, x]),)
460 def _load(self, g, dtype, n, offset=0):
461 # numpy.fromfile, like .frombuffer
462 g.seek(offset, _SEEK_SET)
463 return self.numpy.fromfile(g, dtype, n)
465 @Property_RO
466 def _lowerleft(self):
467 '''(INTERNAL) Cache for L{lowerleft}.
468 '''
469 return self._llh3(self._lat_lo, self._lon_lo)
471 def lowerleft(self, LatLon=None):
472 '''Return the lower-left location and height of this geoid.
474 @kwarg LatLon: Optional class to return the location
475 (C{LatLon}) and height or C{None}.
477 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
478 lon, height)} otherwise a B{C{LatLon}} instance
479 with the lat-, longitude and height of the lower-left,
480 SW grid corner.
481 '''
482 return self._llh3LL(self._lowerleft, LatLon)
484 @Property_RO
485 def _loweright(self):
486 '''(INTERNAL) Cache for L{loweright}.
487 '''
488 return self._llh3(self._lat_lo, self._lon_hi)
490 def loweright(self, LatLon=None):
491 '''Return the lower-right location and height of this geoid.
493 @kwarg LatLon: Optional class to return the location and height
494 (C{LatLon}) or C{None}.
496 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
497 lon, height)} otherwise a B{C{LatLon}} instance
498 with the lat-, longitude and height of the lower-right,
499 SE grid corner.
500 '''
502 return self._llh3LL(self._loweright, LatLon)
504 lowerright = loweright # synonymous
506 @Property_RO
507 def _lowest(self):
508 '''(INTERNAL) Cache for L{lowest}.
509 '''
510 return self._llh3minmax(False)
512 def lowest(self, LatLon=None, **unused):
513 '''Return the location and lowest height of this geoid.
515 @kwarg LatLon: Optional class to return the location and height
516 (C{LatLon}) or C{None}.
518 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
519 lon, height)} otherwise a B{C{LatLon}} instance
520 with the lat-, longitude and height of the lowest
521 grid location.
522 '''
523 return self._llh3LL(self._lowest, LatLon)
525 @Property_RO
526 def mean(self):
527 '''Get the mean of this geoid's heights (C{float}).
528 '''
529 if self._mean is None: # see GeoidKarney
530 self._mean = float(self.numpy.mean(self._hs_y_x))
531 return self._mean
533 @property_RO
534 def name(self):
535 '''Get the name of this geoid (C{str}).
536 '''
537 return _HeightBase.name.fget(self) or self._geoid # recursion
539 @Property_RO
540 def nBytes(self):
541 '''Get the grid in-memory size in bytes (C{int}).
542 '''
543 return self._nBytes
545 def _open(self, geoid, datum, kind, name, smooth):
546 # open the geoid file
547 try:
548 self._geoid = _os_path.basename(geoid)
549 self._sizeB = _os_path.getsize(geoid)
550 g = open(geoid, _rb_)
551 except (IOError, OSError) as x:
552 raise GeoidError(geoid=geoid, cause=x)
554 if datum not in (None, self._datum):
555 self._datum = _ellipsoidal_datum(datum, name=name)
556 self._kind = int(kind)
557 if name:
558 _HeightBase.name.fset(self, name) # rename
559 if smooth:
560 self._smooth = Int_(smooth=smooth, Error=GeoidError, low=0)
562 return g
564 def outside(self, lat, lon):
565 '''Check whether a location is outside this geoid's
566 lat-/longitude or crop range.
568 @arg lat: The latitude (C{degrees}).
569 @arg lon: The longitude (C{degrees}).
571 @return: A 1- or 2-character C{str} if outside or an
572 empty C{str} if inside.
573 '''
574 return (_S_ if lat < self._lat_lo else
575 (_N_ if lat > self._lat_hi else NN)) + \
576 (_W_ if lon < self._lon_lo else
577 (_E_ if lon > self._lon_hi else NN))
579 @Property_RO
580 def pgm(self):
581 '''Get the PGM attributes (C{_PGM} or C{None} if not available/applicable).
582 '''
583 return self._pgm
585 @Property_RO
586 def sizeB(self):
587 '''Get the geoid grid file size in bytes (C{int}).
588 '''
589 return self._sizeB
591 @Property_RO
592 def smooth(self):
593 '''Get the C{RectBivariateSpline} smoothing (C{int}).
594 '''
595 return self._smooth
597 @Property_RO
598 def stdev(self):
599 '''Get the standard deviation of this geoid's heights (C{float}) or C{None}.
600 '''
601 if self._stdev is None: # see GeoidKarney
602 self._stdev = float(self.numpy.std(self._hs_y_x))
603 return self._stdev
605 def _swne(self, crop):
606 # crop box to 4-tuple (s, w, n, e)
607 try:
608 if len(crop) == 2:
609 try: # sw, ne LatLons
610 swne = (crop[0].lat, crop[0].lon,
611 crop[1].lat, crop[1].lon)
612 except AttributeError: # (s, w), (n, e)
613 swne = tuple(crop[0]) + tuple(crop[1])
614 else: # (s, w, n, e)
615 swne = crop
616 if len(swne) == 4:
617 s, w, n, e = map(float, swne)
618 if -90 <= s <= (n - _1_0) <= 89 and \
619 -180 <= w <= (e - _1_0) <= 179:
620 return s, w, n, e
621 except (IndexError, TypeError, ValueError):
622 pass
623 raise GeoidError(crop=crop)
625 def toStr(self, prec=3, sep=_COMMASPACE_): # PYCHOK signature
626 '''This geoid and all geoid attributes as a string.
628 @kwarg prec: Number of decimal digits (0..9 or C{None} for
629 default). Trailing zero decimals are stripped
630 for B{C{prec}} values of 1 and above, but kept
631 for negative B{C{prec}} values.
632 @kwarg sep: Separator to join (C{str}).
634 @return: Geoid name and attributes (C{str}).
635 '''
636 s = 1 if self.kind < 0 else 2
637 t = tuple(Fmt.PAREN(m.__name__, fstr(m(), prec=prec)) for m in
638 (self.lowerleft, self.upperright,
639 self.center,
640 self.highest, self.lowest)) + \
641 attrs( _mean_, _stdev_, prec=prec, Nones=False) + \
642 attrs((_kind_, 'smooth')[:s], prec=prec, Nones=False) + \
643 attrs( 'cropped', 'dtype', _endian_, 'hits', _knots_, 'nBytes',
644 'sizeB', _scipy_, _numpy_, prec=prec, Nones=False)
645 return _COLONSPACE_(self, sep.join(t))
647 @Property_RO
648 def u2B(self):
649 '''Get the PGM itemsize in bytes (C{int}).
650 '''
651 return self._u2B
653 @Property_RO
654 def _upperleft(self):
655 '''(INTERNAL) Cache for method L{upperleft}.
656 '''
657 return self._llh3(self._lat_hi, self._lon_lo)
659 def upperleft(self, LatLon=None):
660 '''Return the upper-left location and height of this geoid.
662 @kwarg LatLon: Optional class to return the location and height
663 (C{LatLon}) or C{None}.
665 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
666 lon, height)} otherwise a B{C{LatLon}} instance
667 with the lat-, longitude and height of the upper-left,
668 NW grid corner.
669 '''
670 return self._llh3LL(self._upperleft, LatLon)
672 @Property_RO
673 def _upperright(self):
674 '''(INTERNAL) Cache for method L{upperright}.
675 '''
676 return self._llh3(self._lat_hi, self._lon_hi)
678 def upperright(self, LatLon=None):
679 '''Return the upper-right location and height of this geoid.
681 @kwarg LatLon: Optional class to return the location and height
682 (C{LatLon}) or C{None}.
684 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
685 lon, height)} otherwise a B{C{LatLon}} instance
686 with the lat-, longitude and height of the upper-right,
687 NE grid corner.
688 '''
689 return self._llh3LL(self._upperright, LatLon)
692class GeoidError(HeightError):
693 '''Geoid interpolator C{Geoid...} or interpolation issue.
694 '''
695 pass
698class GeoidG2012B(_GeoidBase):
699 '''Geoid height interpolator for U{GEOID12B Model
700 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/>} grids U{CONUS
701 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_CONUS.shtml>},
702 U{Alaska<https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_AK.shtml>},
703 U{Hawaii<https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_HI.shtml>},
704 U{Guam and Northern Mariana Islands
705 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_GMNI.shtml>},
706 U{Puerto Rico and U.S. Virgin Islands
707 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_PRVI.shtml>} and
708 U{American Samoa<https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_AS.shtml>}
709 based on C{SciPy} U{RectBivariateSpline<https://docs.SciPy.org/doc/
710 scipy/reference/generated/scipy.interpolate.RectBivariateSpline.html>}
711 or U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/
712 scipy.interpolate.interp2d.html>} interpolation.
714 Use any of the binary C{le} (little endian) or C{be} (big endian)
715 C{g2012b*.bin} grid files.
716 '''
717 def __init__(self, g2012b_bin, crop=None, datum=None, # NAD 83 Ellipsoid
718 kind=3, name=NN, smooth=0):
719 '''New L{GeoidG2012B} interpolator.
721 @arg g2012b_bin: A C{GEOID12B} grid file name (C{.bin}).
722 @kwarg crop: Optional crop box, not supported (C{None}).
723 @kwarg datum: Optional grid datum (L{Datum}, L{Ellipsoid}, L{Ellipsoid2}
724 or L{a_f2Tuple}), default C{WGS84}.
725 @kwarg kind: C{scipy.interpolate} order (C{int}), use 1..5 for
726 U{RectBivariateSpline<https://docs.SciPy.org/doc/scipy/
727 reference/generated/scipy.interpolate.RectBivariateSpline.html>},
728 -2 for U{interp2d linear<https://docs.SciPy.org/doc/scipy/
729 reference/generated/scipy.interpolate.interp2d.html>}, -3
730 for C{interp2d cubic} or -5 for C{interp2d quintic}.
731 @kwarg name: Optional geoid name (C{str}).
732 @kwarg smooth: Smoothing factor for U{RectBivariateSpline
733 <https://docs.SciPy.org/doc/scipy/reference/generated/
734 scipy.interpolate.RectBivariateSpline.html>}
735 only (C{int}).
737 @raise GeoidError: G2012B grid file B{C{g2012b_bin}} issue or invalid
738 B{C{crop}}, B{C{kind}} or B{C{smooth}}.
740 @raise ImportError: Package C{numpy} or C{scipy} not found or
741 not installed.
743 @raise LenError: Grid file B{C{g2012b_bin}} axis mismatch.
745 @raise SciPyError: A C{RectBivariateSpline} or C{inter2d} issue.
747 @raise SciPyWarning: A C{RectBivariateSpline} or C{inter2d}
748 warning as exception.
750 @raise TypeError: Invalid B{C{datum}}.
751 '''
752 if crop is not None:
753 raise GeoidError(crop=crop, txt=_not_(_supported_))
755 g = self._open(g2012b_bin, datum, kind, name, smooth)
756 _ = self.numpy # import numpy for ._load and
758 try:
759 p = _Gpars()
760 n = (self.sizeB // 4) - 11 # number of f4 heights
761 # U{numpy dtype formats are different from Python struct formats
762 # <https://docs.SciPy.org/doc/numpy-1.15.0/reference/arrays.dtypes.html>}
763 for en_ in ('<', '>'):
764 # skip 4xf8, get 3xi4
765 p.nlat, p.nlon, ien = map(int, self._load(g, en_+'i4', 3, 32))
766 if ien == 1: # correct endian
767 p.knots = p.nlat * p.nlon
768 if p.knots == n and 1 < p.nlat < n \
769 and 1 < p.nlon < n:
770 self._endian = en_+'f4'
771 break
772 else: # couldn't validate endian
773 raise GeoidError(_endian_)
775 # get the first 4xf8
776 p.slat, p.wlon, p.dlat, p.dlon = map(float, self._load(g, en_+'f8', 4))
777 # read all f4 heights, ignoring the first 4xf8 and 3xi4
778 hs = self._load(g, self._endian, n, 44).reshape(p.nlat, p.nlon)
779 p.wlon -= _360_0 # western-most East longitude to earth (..., lon)
780 _GeoidBase.__init__(self, hs, p)
782 except Exception as x:
783 raise _SciPyIssue(x, _in_, repr(g2012b_bin))
784 finally:
785 g.close()
787 def _g2ll2(self, lat, lon):
788 # convert grid (lat, lon) to earth (lat, lon)
789 return lat, lon
791 def _ll2g2(self, lat, lon):
792 # convert earth (lat, lon) to grid (lat, lon)
793 return lat, lon
795 if _FOR_DOCS:
796 __call__ = _GeoidBase.__call__
797 height = _GeoidBase.height
800class GeoidHeight5Tuple(_NamedTuple): # .geoids.py
801 '''5-Tuple C{(lat, lon, egm84, egm96, egm2008)} for U{GeoidHeights.dat
802 <https://SourceForge.net/projects/geographiclib/files/testdata/>}
803 tests with the heights for 3 different EGM grids at C{degrees90}
804 and C{degrees180} degrees (after converting C{lon} from original
805 C{0 <= EasterLon <= 360}).
806 '''
807 _Names_ = (_lat_, _lon_, 'egm84', 'egm96', 'egm2008')
808 _Units_ = ( Lat, Lon, Height, Height, Height)
811def _I(i):
812 '''(INTERNAL) Cache a single C{int} constant.
813 '''
814 return _intCs.setdefault(i, i) # PYCHOK undefined due to del _intCs
817def _T(*cs):
818 '''(INTERNAL) Cache a tuple of single C{int} constants.
819 '''
820 return map1(_I, *cs)
823class GeoidKarney(_GeoidBase):
824 '''Geoid height interpolator for I{Karney}'s U{GeographicLib Earth
825 Gravitational Model (EGM)<https://GeographicLib.SourceForge.io/C++/doc/
826 geoid.html>} geoid U{egm*.pgm<https://GeographicLib.SourceForge.io/
827 html/geoid.html#geoidinst>} datasets using bilinear or U{cubic
828 <https://dl.ACM.org/citation.cfm?id=368443>} interpolation and U{caching
829 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidcache>}
830 in pure Python, transcoded from I{Karney}'s U{C++ class Geoid
831 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinterp>}.
833 Use any of the geoid U{egm84-, egm96- or egm2008-*.pgm
834 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinst>}
835 datasets.
836 '''
837 _C0 = _F(372), _F(240), _F(372) # n, _ and s common denominators
838 # matrices c3n_, c3, c3s_, transposed from GeographicLib/Geoid.cpp
839 _C3 = ((_T(0, 0, 62, 124, 124, 62, 0, 0, 0, 0, 0, 0),
840 _T(0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0),
841 _T(-131, 7, -31, -62, -62, -31, 45, 216, 156, -45, -55, -7),
842 _T(0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0),
843 _T(138, -138, 0, 0, 0, 0, -183, 33, 153, -3, 48, -48), # PYCHOK indent
844 _T(144, 42, -62, -124, -124, -62, -9, 87, 99, 9, 42, -42),
845 _T(0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0),
846 _T(0, 0, 0, 0, 0, 0, 93, -93, -93, 93, 0, 0),
847 _T(-102, 102, 0, 0, 0, 0, 18, 12, -12, -18, -84, 84),
848 _T(-31, -31, 31, 62, 62, 31, 0, -93, -93, 0, 31, 31)), # PYCHOK indent
850 (_T(9, -9, 9, 186, 54, -9, -9, 54, -54, 9, -9, 9),
851 _T(-18, 18, -88, -42, 162, -32, 8, -78, 78, -8, 18, -18),
852 _T(-88, 8, -18, -42, -78, 18, 18, 162, 78, -18, -32, -8),
853 _T(0, 0, 90, -150, 30, 30, 30, -90, 90, -30, 0, 0),
854 _T(96, -96, 96, -96, -24, 24, -96, -24, 144, -24, 24, -24), # PYCHOK indent
855 _T(90, 30, 0, -150, -90, 0, 0, 30, 90, 0, 30, -30),
856 _T(0, 0, -20, 60, -60, 20, -20, 60, -60, 20, 0, 0),
857 _T(0, 0, -60, 60, 60, -60, 60, -60, -60, 60, 0, 0),
858 _T(-60, 60, 0, 60, -60, 0, 0, 60, -60, 0, -60, 60),
859 _T(-20, -20, 0, 60, 60, 0, 0, -60, -60, 0, 20, 20)),
861 (_T(18, -18, 36, 210, 162, -36, 0, 0, 0, 0, -18, 18), # PYCHOK indent
862 _T(-36, 36, -165, 45, 141, -21, 0, 0, 0, 0, 36, -36),
863 _T(-122, -2, -27, -111, -75, 27, 62, 124, 124, 62, -64, 2),
864 _T(0, 0, 93, -93, -93, 93, 0, 0, 0, 0, 0, 0),
865 _T(120, -120, 147, -57, -129, 39, 0, 0, 0, 0, 66, -66), # PYCHOK indent
866 _T(135, 51, -9, -192, -180, 9, 31, 62, 62, 31, 51, -51),
867 _T(0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0),
868 _T(0, 0, -93, 93, 93, -93, 0, 0, 0, 0, 0, 0),
869 _T(-84, 84, 18, 12, -12, -18, 0, 0, 0, 0, -102, 102),
870 _T(-31, -31, 0, 93, 93, 0, -31, -62, -62, -31, 31, 31)))
872 _BT = (_T(0, 0), # bilinear 4-tuple [i, j] indices
873 _T(1, 0),
874 _T(0, 1),
875 _T(1, 1))
877 _CM = (_T( 0, -1), # 10x12 cubic matrix [i, j] indices
878 _T( 1, -1),
879 _T(-1, 0),
880 _T( 0, 0),
881 _T( 1, 0),
882 _T( 2, 0),
883 _T(-1, 1),
884 _T( 0, 1),
885 _T( 1, 1),
886 _T( 2, 1),
887 _T( 0, 2),
888 _T( 1, 2))
890 _endian = '>H' # struct.unpack 1 ushort (big endian, unsigned short)
891 _4endian = '>4H' # struct.unpack 4 ushorts
892 _Rendian = NN # struct.unpack a row of ushorts
893# _highest = (-8.4, 147.367, 85.839) if egm2008-1.pgm else (
894# (-8.167, 147.25, 85.422) if egm96-5.pgm else
895# (-4.5, 148.75, 81.33)) # egm84-15.pgm
896# _lowest = (4.7, 78.767, -106.911) if egm2008-1.pgm else (
897# (4.667, 78.833, -107.043) if egm96-5.pgm else
898# (4.75, 79.25, -107.34)) # egm84-15.pgm
899 _mean = _F(-1.317) # from egm2008-1, -1.438 egm96-5, -0.855 egm84-15
900 _nBytes = None # not applicable
901 _nterms = len(_C3[0]) # columns length, number of row
902 _smooth = None # not applicable
903 _stdev = _F(29.244) # from egm2008-1, 29.227 egm96-5, 29.183 egm84-15
904 _u2B = _calcsize(_endian) # pixelsize_ in bytes
905 _4u2B = _calcsize(_4endian) # 4 pixelsize_s in bytes
906 _Ru2B = 0 # row of pixelsize_s in bytes
907 _yxH = () # cache (y, x) indices
908 _yxHt = () # cached 4- or 10-tuple for _ev2H resp. _ev3H
909 _yx_hits = 0 # cache hits
911 def __init__(self, egm_pgm, crop=None, datum=None, # WGS84
912 kind=3, name=NN, smooth=None):
913 '''New L{GeoidKarney} interpolator.
915 @arg egm_pgm: An U{EGM geoid dataset<https://GeographicLib.SourceForge.io/
916 html/geoid.html#geoidinst>} file name (C{egm*.pgm}), see
917 note below.
918 @kwarg crop: Optional box to limit geoid locations, a 4-tuple (C{south,
919 west, north, east}), 2-tuple (C{(south, west), (north,
920 east)}) or 2, in C{degrees90} lat- and C{degrees180}
921 longitudes or a 2-tuple (C{LatLonSW, LatLonNE}) of
922 C{LatLon} instances.
923 @kwarg datum: Optional grid datum (C{Datum}, L{Ellipsoid}, L{Ellipsoid2}
924 or L{a_f2Tuple}), default C{WGS84}.
925 @kwarg kind: Interpolation order (C{int}), 2 for C{bilinear} or 3
926 for C{cubic}.
927 @kwarg name: Optional geoid name (C{str}).
928 @kwarg smooth: Smoothing factor, unsupported (C{None}).
930 @raise GeoidError: EGM dataset B{C{egm_pgm}} issue or invalid
931 B{C{crop}}, B{C{kind}} or B{C{smooth}}.
933 @raise TypeError: Invalid B{C{datum}}.
935 @see: Class L{GeoidPGM} and function L{egmGeoidHeights}.
937 @note: Geoid file B{C{egm_pgm}} remains open and must be closed
938 by calling the C{close} method or by using this instance
939 in a C{with B{GeoidKarney}(...) as ...} context.
940 '''
941 if smooth is not None:
942 raise GeoidError(smooth=smooth, txt=_not_(_supported_))
944 if kind in (2,):
945 self._evH = self._ev2H
946 elif kind not in (3,):
947 raise GeoidError(kind=kind)
949 self._egm = g = self._open(egm_pgm, datum, kind, name, smooth)
950 self._pgm = p = _PGM(g, pgm=egm_pgm, itemsize=self.u2B, sizeB=self.sizeB)
952 self._Rendian = self._4endian.replace(_4_, str(p.nlon))
953 self._Ru2B = _calcsize(self._Rendian)
955 self._knots = p.knots # grid knots
956 self._lon_of = float(p.flon) # forward offset
957 self._lon_og = float(p.glon) # reverse offset
958 # set earth (lat, lon) limits (s, w, n, e)
959 self._lat_lo, \
960 self._lon_lo, \
961 self._lat_hi, \
962 self._lon_hi = self._swne(crop if crop else p.crop4)
963 self._cropped = True if crop else False
965 def __call__(self, *llis):
966 '''Interpolate the geoid height for one or several locations.
968 @arg llis: The location or locations (C{LatLon}, ... or
969 C{LatLon}s).
971 @return: A single interpolated geoid height (C{float}) or
972 a list or tuple of interpolated geoid heights
973 (C{float}s).
975 @raise GeoidError: Insufficient number of B{C{llis}}, an
976 invalid B{C{lli}} or the C{egm*.pgm}
977 geoid file is closed.
979 @raise RangeError: An B{C{lli}} is outside this geoid's lat-
980 or longitude range.
981 '''
982 return self._called(llis, False)
984 def _c0c3v(self, y, x):
985 # get the common denominator, the 10x12 cubic matrix and
986 # the 12 cubic v-coefficients around geoid index (y, x)
987 p = self._pgm
988 if 0 < x < (p.nlon - 2) and 0 < y < (p.nlat - 2):
989 # read 4x4 ushorts, drop the 4 corners
990 g = self._egm
991 e = self._4endian
992 n = self._4u2B
993 R = self._Ru2B
995 b = self._seek(y - 1, x - 1)
996 v = _unpack(e, g.read(n))[1:3]
997 b += R
998 g.seek(b, _SEEK_SET)
999 v += _unpack(e, g.read(n))
1000 b += R
1001 g.seek(b, _SEEK_SET)
1002 v += _unpack(e, g.read(n))
1003 b += R
1004 g.seek(b, _SEEK_SET)
1005 v += _unpack(e, g.read(n))[1:3]
1006 j = 1
1008 else: # likely some wrapped y and/or x's
1009 v = self._raws(y, x, GeoidKarney._CM)
1010 j = 0 if y < 1 else (1 if y < (p.nlat - 2) else 2)
1012 return GeoidKarney._C0[j], GeoidKarney._C3[j], v
1014 @Property_RO
1015 def dtype(self):
1016 '''Get the geoid's grid data type (C{str}).
1017 '''
1018 return 'ushort'
1020 def _ev(self, lat, lon): # PYCHOK expected
1021 # interpolate the geoid height at grid (lat, lon)
1022 fy, fx = self._g2yx2(lat, lon)
1023 y, x = int(floor(fy)), int(floor(fx))
1024 fy -= y
1025 fx -= x
1026 H = self._evH(fy, fx, y, x) # ._ev3H or ._ev2H
1027 H *= self._pgm.Scale # H.fmul(self._pgm.Scale)
1028 H += self._pgm.Offset # H.fadd(self._pgm.Offset)
1029 return H.fsum()
1031 def _ev2H(self, fy, fx, *yx):
1032 # compute the bilinear 4-tuple and interpolate raw H
1033 if self._yxH == yx:
1034 t = self._yxHt
1035 self._yx_hits += 1
1036 else:
1037 y, x = self._yxH = yx
1038 self._yxHt = t = self._raws(y, x, GeoidKarney._BT)
1039 v = _1_0, -fx, fx
1040 H = Fdot(v, t[0], t[0], t[1]).fmul(_1_0 - fy) # c = a * (1 - fy)
1041 H += Fdot(v, t[2], t[2], t[3]).fmul(fy) # c += b * fy
1042 return H
1044 def _ev3H(self, fy, fx, *yx):
1045 # compute the cubic 10-tuple and interpolate raw H
1046 if self._yxH == yx:
1047 t = self._yxHt
1048 self._yx_hits += 1
1049 else:
1050 self._yxH = yx
1051 c0, c3, v = self._c0c3v(*yx)
1052 t = [fdot(v, *c3[i]) / c0 for i in range(self._nterms)]
1053 self._yxHt = t = tuple(t)
1054 # GeographicLib/Geoid.cpp Geoid::height(lat, lon) ...
1055 # real h = t[0] + fx * (t[1] + fx * (t[3] + fx * t[6])) +
1056 # fy * (t[2] + fx * (t[4] + fx * t[7]) +
1057 # fy * (t[5] + fx * t[8] + fy * t[9]));
1058 v = _1_0, fx, fy
1059 H = Fdot(v, t[5], t[8], t[9])
1060 H *= fy
1061 H += Fhorner(fx, t[2], t[4], t[7])
1062 H *= fy
1063 H += Fhorner(fx, t[0], t[1], t[3], t[6])
1064 return H
1066 _evH = _ev3H # overriden for kind == 2
1068 def _g2ll2(self, lat, lon):
1069 # convert grid (lat, lon) to earth (lat, lon), uncropped
1070 while lon > _180_0:
1071 lon -= _360_0
1072 return lat, lon
1074 def _g2yx2(self, lat, lon):
1075 # convert grid (lat, lon) to grid (y, x) indices
1076 p = self._pgm
1077 # note, slat = +90, rlat < 0 makes y >=0
1078 return ((lat - p.slat) * p.rlat), ((lon - p.wlon) * p.rlon)
1080 def _gyx2g2(self, y, x):
1081 # convert grid (y, x) indices to grid (lat, lon)
1082 p = self._pgm
1083 return (p.slat + p.dlat * y), (p.wlon + p.dlon * x)
1085 def height(self, lats, lons):
1086 '''Interpolate the geoid height for one or several lat-/longitudes.
1088 @arg lats: Latitude or latitudes (C{degrees} or C{degrees}s).
1089 @arg lons: Longitude or longitudes (C{degrees} or C{degrees}s).
1091 @return: A single interpolated geoid height (C{float}) or a
1092 list of interpolated geoid heights (C{float}s).
1094 @raise GeoidError: Insufficient or non-matching number of
1095 B{C{lats}} and B{C{lons}} or the C{egm*.pgm}
1096 geoid file is closed.
1098 @raise RangeError: A B{C{lat}} or B{C{lon}} is outside this
1099 geoid's lat- or longitude range.
1100 '''
1101 return _HeightBase._height(self, lats, lons, Error=GeoidError)
1103 @Property_RO
1104 def _highest_ltd(self):
1105 '''(INTERNAL) Cache for L{highest} mesthod.
1106 '''
1107 return self._llh3minmax(True, -12, -4)
1109 def highest(self, LatLon=None, full=False): # PYCHOK full
1110 '''Return the location and largest height of this geoid.
1112 @kwarg LatLon: Optional class to return the location and height
1113 (C{LatLon}) or C{None}.
1114 @kwarg full: Search the full or limited latitude range (C{bool}).
1116 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
1117 lon, height)} otherwise a B{C{LatLon}} instance
1118 with the lat-, longitude and height of the highest
1119 grid location.
1120 '''
1121 llh = self._highest if full or self.cropped else self._highest_ltd
1122 return self._llh3LL(llh, LatLon)
1124 def _lat2y2(self, lat2):
1125 # convert earth lat(s) to min and max grid y indices
1126 ys, m = [], self._pgm.nlat - 1
1127 for lat in lat2:
1128 y, _ = self._g2yx2(*self._ll2g2(lat, 0))
1129 ys.append(max(min(int(y), m), 0))
1130 return min(ys), max(ys) + 1
1132 def _ll2g2(self, lat, lon):
1133 # convert earth (lat, lon) to grid (lat, lon), uncropped
1134 while lon < 0:
1135 lon += _360_0
1136 return lat, lon
1138 def _llh3minmax(self, highest=True, *lat2):
1139 # find highest or lowest, takes 10+ secs for egm2008-1.pgm geoid
1140 # (Python 2.7.16, macOS 10.13.6 High Sierra, iMac 3 GHz Core i3)
1141 y = x = 0
1142 h = self._raw(y, x)
1143 if highest:
1144 for j, r in self._raw2(*lat2):
1145 m = max(r)
1146 if m > h:
1147 h, y, x = m, j, r.index(m)
1148 else: # lowest
1149 for j, r in self._raw2(*lat2):
1150 m = min(r)
1151 if m < h:
1152 h, y, x = m, j, r.index(m)
1153 h *= self._pgm.Scale
1154 h += self._pgm.Offset
1155 return self._g2ll2(*self._gyx2g2(y, x)) + (h,)
1157 @Property_RO
1158 def _lowest_ltd(self):
1159 '''(INTERNAL) Cache for L{lowest}.
1160 '''
1161 return self._llh3minmax(False, 0, 8)
1163 def lowest(self, LatLon=None, full=False): # PYCHOK full
1164 '''Return the location and lowest height of this geoid.
1166 @kwarg LatLon: Optional class to return the location and height
1167 (C{LatLon}) or C{None}.
1168 @kwarg full: Search the full or limited latitude range (C{bool}).
1170 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
1171 lon, height)} otherwise a B{C{LatLon}} instance
1172 with the lat-, longitude and height of the lowest
1173 grid location.
1174 '''
1175 llh = self._lowest if full or self.cropped else self._lowest_ltd
1176 return self._llh3LL(llh, LatLon)
1178 def _raw(self, y, x):
1179 # get the ushort geoid height at geoid index (y, x),
1180 # like GeographicLib/Geoid.hpp real rawval(is, iy)
1181 p = self._pgm
1182 if x < 0:
1183 x += p.nlon
1184 elif x >= p.nlon:
1185 x -= p.nlon
1186 h = p.nlon // 2
1187 if y < 0:
1188 y = -y
1189 elif y >= p.nlat:
1190 y = (p.nlat - 1) * 2 - y
1191 else:
1192 h = 0
1193 x += h if x < h else -h
1194 self._seek(y, x)
1195 h = _unpack(self._endian, self._egm.read(self._u2B))
1196 return h[0]
1198 def _raws(self, y, x, ijs):
1199 # get bilinear 4-tuple or 10x12 cubic matrix
1200 return tuple(self._raw(y + j, x + i) for i, j in ijs)
1202 def _raw2(self, *lat2):
1203 # yield a 2-tuple (y, ushorts) for each row or for
1204 # the rows between two (or more) earth lat values
1205 p = self._pgm
1206 g = self._egm
1207 e = self._Rendian
1208 n = self._Ru2B
1209 # min(lat2) <= lat <= max(lat2) or 0 <= y < p.nlat
1210 s, t = self._lat2y2(lat2) if lat2 else (0, p.nlat)
1211 self._seek(s, 0) # to start of row s
1212 for y in range(s, t):
1213 yield y, _unpack(e, g.read(n))
1215 def _seek(self, y, x):
1216 # position geoid to grid index (y, x)
1217 p, g = self._pgm, self._egm
1218 if g:
1219 b = p.skip + (y * p.nlon + x) * self._u2B
1220 g.seek(b, _SEEK_SET)
1221 return b # position
1222 raise GeoidError('closed file: %r' % (p.egm,)) # IOError
1225class GeoidPGM(_GeoidBase):
1226 '''Geoid height interpolator for I{Karney}'s U{GeographicLib Earth
1227 Gravitational Model (EGM)<https://GeographicLib.SourceForge.io/C++/doc/
1228 geoid.html>} geoid U{egm*.pgm<https://GeographicLib.SourceForge.io/
1229 html/geoid.html#geoidinst>} datasets but based on C{SciPy}
1230 U{RectBivariateSpline<https://docs.SciPy.org/doc/scipy/reference/
1231 generated/scipy.interpolate.RectBivariateSpline.html>} or
1232 U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/
1233 scipy.interpolate.interp2d.html>} interpolation.
1235 Use any of the U{egm84-, egm96- or egm2008-*.pgm
1236 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinst>}
1237 datasets. However, unless cropped, an entire C{egm*.pgm} dataset
1238 is loaded into the C{SciPy} U{RectBivariateSpline<https://docs.SciPy.org/
1239 doc/scipy/reference/generated/scipy.interpolate.RectBivariateSpline.html>}
1240 or U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/
1241 scipy.interpolate.interp2d.html>} interpolator and converted from
1242 2-byte C{int} to 8-byte C{dtype float64}. Therefore, internal memory
1243 usage is 4x the U{egm*.pgm<https://GeographicLib.SourceForge.io/C++/doc/
1244 geoid.html#geoidinst>} file size and may exceed the available memory,
1245 especially with 32-bit Python, see properties C{.nBytes} and C{.sizeB}.
1246 '''
1247 _endian = '>u2'
1249 def __init__(self, egm_pgm, crop=None, datum=None, # WGS84
1250 kind=3, name=NN, smooth=0):
1251 '''New L{GeoidPGM} interpolator.
1253 @arg egm_pgm: An U{EGM geoid dataset<https://GeographicLib.SourceForge.io/
1254 html/geoid.html#geoidinst>} file name (C{egm*.pgm}).
1255 @kwarg crop: Optional box to crop B{C{egm_pgm}}, a 4-tuple (C{south, west,
1256 north, east}) or 2-tuple (C{(south, west), (north, east)}),
1257 in C{degrees90} lat- and C{degrees180} longitudes or a
1258 2-tuple (C{LatLonSW, LatLonNE}) of C{LatLon} instances.
1259 @kwarg datum: Optional grid datum (L{Datum}, L{Ellipsoid}, L{Ellipsoid2}
1260 or L{a_f2Tuple}), default C{WGS84}.
1261 @kwarg kind: C{scipy.interpolate} order (C{int}), use 1..5 for
1262 U{RectBivariateSpline<https://docs.SciPy.org/doc/scipy/
1263 reference/generated/scipy.interpolate.RectBivariateSpline.html>},
1264 -2 for U{interp2d linear<https://docs.SciPy.org/doc/scipy/
1265 reference/generated/scipy.interpolate.interp2d.html>}, -3
1266 for C{interp2d cubic} or -5 for C{interp2d quintic}.
1267 @kwarg name: Optional geoid name (C{str}).
1268 @kwarg smooth: Smoothing factor for U{RectBivariateSpline
1269 <https://docs.SciPy.org/doc/scipy/reference/generated/
1270 scipy.interpolate.RectBivariateSpline.html>}
1271 only (C{int}).
1273 @raise GeoidError: EGM dataset B{C{egm_pgm}} issue or invalid B{C{crop}},
1274 B{C{kind}} or B{C{smooth}}.
1276 @raise ImportError: Package C{numpy} or C{scipy} not found or not installed.
1278 @raise LenError: EGM dataset B{C{egm_pgm}} axis mismatch.
1280 @raise SciPyError: A C{RectBivariateSpline} or C{inter2d} issue.
1282 @raise SciPyWarning: A C{RectBivariateSpline} or C{inter2d}
1283 warning as exception.
1285 @raise TypeError: Invalid B{C{datum}}.
1287 @note: The U{GeographicLib egm*.pgm<https://GeographicLib.SourceForge.io/
1288 html/geoid.html#geoidinst>} file sizes are based on a 2-byte
1289 C{int} height converted to 8-byte C{dtype float64} for C{scipy}
1290 interpolators. Therefore, internal memory usage is 4 times the
1291 C{egm*.pgm} file size and may exceed the available memory,
1292 especially with 32-bit Python. To reduce memory usage, set
1293 keyword argument B{C{crop}} to the region of interest. For example
1294 C{B{crop}=(20, -125, 50, -65)} covers the U{conterminous US<https://
1295 www.NGS.NOAA.gov/GEOID/GEOID12B/maps/GEOID12B_CONUS_grids.png>}
1296 (CONUS), less than 3% of the entire C{egm2008-1.pgm} dataset.
1298 @see: Class L{GeoidKarney} and function L{egmGeoidHeights}.
1299 '''
1300 np = self.numpy
1301 self._u2B = np.dtype(self.endian).itemsize
1303 g = self._open(egm_pgm, datum, kind, name, smooth)
1304 self._pgm = p = _PGM(g, pgm=egm_pgm, itemsize=self.u2B, sizeB=self.sizeB)
1305 if crop:
1306 g = p._cropped(g, abs(kind) + 1, *self._swne(crop))
1307 self._g2ll2 = self._g2ll2_cropped
1308 self._ll2g2 = self._ll2g2_cropped
1309 if map2(int, np.__version__.split(_DOT_)[:2]) < (1, 9):
1310 g = open(g.name, _rb_) # reopen tempfile for numpy 1.8.0-
1311 self._cropped = True
1312 else:
1313 self._cropped = False
1315 try:
1316 # U{numpy dtype formats are different from Python struct formats
1317 # <https://docs.SciPy.org/doc/numpy-1.15.0/reference/arrays.dtypes.html>}
1318 # read all heights, skipping the PGM header lines, converted to float
1319 hs = self._load(g, self.endian, p.knots, p.skip).reshape(p.nlat, p.nlon) * p.Scale
1320 if p.Offset: # offset
1321 hs = p.Offset + hs
1322 if p.dlat < 0: # flip the rows
1323 hs = np.flipud(hs)
1324 _GeoidBase.__init__(self, hs, p)
1325 except Exception as x:
1326 raise _SciPyIssue(x, _in_, repr(egm_pgm))
1327 finally:
1328 g.close()
1330 def _g2ll2(self, lat, lon):
1331 # convert grid (lat, lon) to earth (lat, lon), uncropped
1332 while lon > _180_0:
1333 lon -= _360_0
1334 return lat, lon
1336 def _g2ll2_cropped(self, lat, lon):
1337 # convert cropped grid (lat, lon) to earth (lat, lon)
1338 return lat, lon - self._lon_og
1340 def _ll2g2(self, lat, lon):
1341 # convert earth (lat, lon) to grid (lat, lon), uncropped
1342 while lon < 0:
1343 lon += _360_0
1344 return lat, lon
1346 def _ll2g2_cropped(self, lat, lon):
1347 # convert earth (lat, lon) to cropped grid (lat, lon)
1348 return lat, lon + self._lon_of
1350 if _FOR_DOCS:
1351 __call__ = _GeoidBase.__call__
1352 height = _GeoidBase.height
1355class _Gpars(_Named):
1356 '''(INTERNAL) Basic geoid parameters.
1357 '''
1358 # interpolator parameters
1359 dlat = 0 # +/- latitude resolution in C{degrees}
1360 dlon = 0 # longitude resolution in C{degrees}
1361 nlat = 1 # number of latitude knots (C{int})
1362 nlon = 0 # number of longitude knots (C{int})
1363 rlat = 0 # +/- latitude resolution in C{float}, 1 / .dlat
1364 rlon = 0 # longitude resolution in C{float}, 1 / .dlon
1365 slat = 0 # nothern- or southern most latitude (C{degrees90})
1366 wlon = 0 # western-most longitude in Eastern lon (C{degrees360})
1368 flon = 0 # forward, earth to grid longitude offset
1369 glon = 0 # reverse, grid to earth longitude offset
1371 knots = 0 # number of knots, nlat * nlon (C{int})
1372 skip = 0 # header bytes to skip (C{int})
1374 def __repr__(self):
1375 t = _COMMASPACE_.join(pairs((a, getattr(self, a)) for
1376 a in dir(self.__class__)
1377 if a[:1].isupper()))
1378 return _COLONSPACE_(self, t)
1380 def __str__(self):
1381 return Fmt.PAREN(self.classname, repr(self.name))
1384class _PGM(_Gpars):
1385 '''(INTERNAL) Parse an C{egm*.pgm} geoid dataset file.
1387 # Geoid file in PGM format for the GeographicLib::Geoid class
1388 # Description WGS84 EGM96, 5-minute grid
1389 # URL https://Earth-Info.NGA.mil/GandG/wgs84/gravitymod/egm96/egm96.html
1390 # DateTime 2009-08-29 18:45:03
1391 # MaxBilinearError 0.140
1392 # RMSBilinearError 0.005
1393 # MaxCubicError 0.003
1394 # RMSCubicError 0.001
1395 # Offset -108
1396 # Scale 0.003
1397 # Origin 90N 0E
1398 # AREA_OR_POINT Point
1399 # Vertical_Datum WGS84
1400 <width> <height>
1401 <pixel>
1402 ...
1403 '''
1404 crop4 = () # 4-tuple (C{south, west, north, east}).
1405 egm = None
1406 glon = 180 # reverse offset, uncropped
1407# pgm = NN # name
1408 sizeB = 0
1409 u2B = 2 # item size of grid height (C{int}).
1411 @staticmethod
1412 def _llstr2floats(latlon):
1413 # llstr to (lat, lon) floats
1414 lat, lon = latlon.split()
1415 return parseDMS2(lat, lon)
1417 # PGM file attributes, CamelCase but not .istitle()
1418 AREA_OR_POINT = str
1419 DateTime = str
1420 Description = str # 'WGS84 EGM96, 5-minute grid'
1421 Geoid = str # 'file in PGM format for the GeographicLib::Geoid class'
1422 MaxBilinearError = float
1423 MaxCubicError = float
1424 Offset = float
1425 Origin = _llstr2floats
1426 Pixel = 0
1427 RMSBilinearError = float
1428 RMSCubicError = float
1429 Scale = float
1430 URL = str # 'https://Earth-Info.NGA.mil/GandG/wgs84/...'
1431 Vertical_Datum = str
1433 def __init__(self, g, pgm=NN, itemsize=0, sizeB=0): # MCCABE 22
1434 '''(INTERNAL) New C{_PGM} parsed C{egm*.pgm} geoid dataset.
1435 '''
1436 self.name = pgm # geoid file name
1437 if itemsize:
1438 self._u2B = itemsize
1439 if sizeB:
1440 self.sizeB = sizeB
1442 t = g.readline() # make sure newline == '\n'
1443 if t != b'P5\n' and t.strip() != b'P5':
1444 raise self._Errorf(_format_, _header_, t)
1446 while True: # read all # Attr ... lines,
1447 try: # ignore empty ones or comments
1448 t = g.readline().strip()
1449 if t.startswith(_bHASH_):
1450 t = t.lstrip(_bHASH_).lstrip()
1451 a, v = map(_ub2str, t.split(None, 1))
1452 f = getattr(_PGM, a, None)
1453 if callable(f) and a[:1].isupper():
1454 setattr(self, a, f(v))
1455 elif t:
1456 break
1457 except (TypeError, ValueError):
1458 raise self._Errorf(_format_, 'Attr', t)
1459 else: # should never get here
1460 raise self._Errorf(_format_, _header_, g.tell())
1462 try: # must be (even) width and (odd) height
1463 nlon, nlat = map(int, t.split())
1464 if nlon < 2 or nlon > (360 * 60) or isodd(nlon) or \
1465 nlat < 2 or nlat > (181 * 60) or not isodd(nlat):
1466 raise ValueError
1467 except (TypeError, ValueError):
1468 raise self._Errorf(_format_, _SPACE_(_width_, _height_), t)
1470 try: # must be 16 bit pixel height
1471 t = g.readline().strip()
1472 self.Pixel = int(t)
1473 if not 255 < self.Pixel < 65536: # >u2 or >H only
1474 raise ValueError
1475 except (TypeError, ValueError):
1476 raise self._Errorf(_format_, 'pixel', t)
1478 for a in dir(_PGM): # set undefined # Attr ... to None
1479 if a[:1].isupper() and callable(getattr(self, a)):
1480 setattr(self, a, None)
1482 if self.Origin is None:
1483 raise self._Errorf(_format_, 'Origin', self.Origin)
1484 if self.Offset is None or self.Offset > 0:
1485 raise self._Errorf(_format_, 'Offset', self.Offset)
1486 if self.Scale is None or self.Scale < EPS:
1487 raise self._Errorf(_format_, 'Scale', self.Scale)
1489 self.skip = g.tell()
1490 self.knots = nlat * nlon
1492 self.nlat, self.nlon = nlat, nlon
1493 self.slat, self.wlon = self.Origin
1494 # note, negative .dlat and .rlat since rows
1495 # are from .slat 90N down in decreasing lat
1496 self.dlat, self.dlon = _180_0 / (1 - nlat), _360_0 / nlon
1497 self.rlat, self.rlon = (1 - nlat) / _180_0, nlon / _360_0
1499 # grid corners in earth (lat, lon), .slat = 90, .dlat < 0
1500 n = float(self.slat)
1501 s = n + self.dlat * (nlat - 1)
1502 w = self.wlon - self.glon
1503 e = w + self.dlon * nlon
1504 self.crop4 = s, w, n, e
1506 n = self.sizeB - self.skip
1507 if n > 0 and n != (self.knots * self.u2B):
1508 raise self._Errorf('%s(%s x %s != %s)', _assert_, nlat, nlon, n)
1510 def _cropped(self, g, k1, south, west, north, east): # MCCABE 15
1511 '''Crop the geoid to (south, west, north, east) box.
1512 '''
1513 # flon offset for both west and east
1514 f = 360 if west < 0 else 0
1515 # earth (lat, lon) to grid indices (y, x),
1516 # note y is decreasing, i.e. n < s
1517 s, w = self._lle2yx2(south, west, f)
1518 n, e = self._lle2yx2(north, east, f)
1519 s += 1 # s > n
1520 e += 1 # e > w
1522 hi, wi = self.nlat, self.nlon
1523 # handle special cases
1524 if (s - n) > hi:
1525 n, s = 0, hi # entire lat range
1526 if (e - w) > wi:
1527 w, e, f = 0, wi, 180 # entire lon range
1528 if s == hi and w == n == 0 and e == wi:
1529 return g # use entire geoid as-is
1531 if (e - w) < k1 or (s - n) < (k1 + 1):
1532 raise self._Errorf(_format_, 'swne', (north - south, east - west))
1534 if e > wi > w: # wrap around
1535 # read w..wi and 0..e
1536 r, p = (wi - w), (e - wi)
1537 elif e > w:
1538 r, p = (e - w), 0
1539 else:
1540 raise self._Errorf('%s(%s < %s)', _assert_, w, e)
1542 # convert to bytes
1543 r *= self.u2B
1544 p *= self.u2B
1545 q = wi * self.u2B # stride
1546 # number of rows and cols to skip from
1547 # the original (.slat, .wlon) origin
1548 z = self.skip + (n * wi + w) * self.u2B
1549 # sanity check
1550 if r < 2 or p < 0 or q < 2 or z < self.skip \
1551 or z > self.sizeB:
1552 raise self._Errorf(_format_, _assert_, (r, p, q, z))
1554 # can't use _BytesIO since numpy
1555 # needs .fileno attr in .fromfile
1556 t, c = 0, self._tmpfile()
1557 # reading (s - n) rows, forward
1558 for y in range(n, s): # PYCHOK y unused
1559 g.seek(z, _SEEK_SET)
1560 # Python 2 tmpfile.write returns None
1561 t += c.write(g.read(r)) or r
1562 if p: # wrap around to start of row
1563 g.seek(-q, _SEEK_CUR)
1564 # assert(g.tell() == (z - w * self.u2B))
1565 # Python 2 tmpfile.write returns None
1566 t += c.write(g.read(p)) or p
1567 z += q
1568 c.flush()
1569 g.close()
1571 s -= n # nlat
1572 e -= w # nlon
1573 k = s * e # knots
1574 z = k * self.u2B
1575 if t != z:
1576 raise self._Errorf('%s(%s != %s) %s', _assert_, t, z, self)
1578 # update the _Gpars accordingly, note attributes
1579 # .dlat, .dlon, .rlat and .rlon remain unchanged
1580 self.slat += n * self.dlat
1581 self.wlon += w * self.dlon
1582 self.nlat = s
1583 self.nlon = e
1584 self.flon = self.glon = f
1586 self.crop4 = south, west, north, east
1587 self.knots = k
1588 self.skip = 0 # no header lines in c
1590 c.seek(0, _SEEK_SET)
1591 # c = open(c.name, _rb_) # reopen for numpy 1.8.0-
1592 return c
1594 def _Errorf(self, fmt, *args): # PYCHOK no cover
1595 t = fmt % args
1596 e = self.pgm or NN
1597 if e:
1598 t = _SPACE_(t, _in_, repr(e))
1599 return PGMError(t)
1601 def _lle2yx2(self, lat, lon, flon):
1602 # earth (lat, lon) to grid indices (y, x)
1603 # with .dlat decreasing from 90N .slat
1604 lat -= self.slat
1605 lon += flon - self.wlon
1606 return (min(self.nlat - 1, max(0, int(lat * self.rlat))),
1607 max(0, int(lon * self.rlon)))
1609 def _tmpfile(self):
1610 # create a tmpfile to hold the cropped geoid grid
1611 try:
1612 from tempfile import NamedTemporaryFile as tmpfile
1613 except ImportError: # Python 2.7.16-
1614 from os import tmpfile
1615 t = _os_path.splitext(_os_path.basename(self.pgm))[0]
1616 f = tmpfile(mode='w+b', prefix=t or 'egm')
1617 f.seek(0, _SEEK_SET) # force overwrite
1618 return f
1620 @Property_RO
1621 def pgm(self):
1622 '''Get the geoid file name (C{str}).
1623 '''
1624 return self.name
1627class PGMError(GeoidError):
1628 '''Issue parsing or cropping an C{egm*.pgm} geoid dataset.
1629 '''
1630 pass
1633def egmGeoidHeights(GeoidHeights_dat):
1634 '''Generate geoid U{egm*.pgm<https://GeographicLib.SourceForge.io/
1635 html/geoid.html#geoidinst>} height tests from U{GeoidHeights.dat
1636 <https://SourceForge.net/projects/geographiclib/files/testdata/>}
1637 U{Test data for Geoids<https://GeographicLib.SourceForge.io/C++/doc/
1638 geoid.html#testgeoid>}.
1640 @arg GeoidHeights_dat: The un-gz-ed C{GeoidHeights.dat} file
1641 (C{str} or C{file} handle).
1643 @return: For each test, yield a L{GeoidHeight5Tuple}C{(lat, lon,
1644 egm84, egm96, egm2008)}.
1646 @raise GeoidError: Invalid B{C{GeoidHeights_dat}}.
1648 @note: Function L{egmGeoidHeights} is used to test the geoids
1649 L{GeoidKarney} and L{GeoidPGM}, see PyGeodesy module
1650 C{test/testGeoids.py}.
1651 '''
1652 dat = GeoidHeights_dat
1653 if isinstance(dat, bytes):
1654 dat = _BytesIO(dat)
1656 try:
1657 dat.seek(0, _SEEK_SET) # reset
1658 except AttributeError as x:
1659 raise GeoidError(GeoidHeights_dat=type(dat), cause=x)
1661 for t in dat.readlines():
1662 t = t.strip()
1663 if t and not t.startswith(_bHASH_):
1664 lat, lon, egm84, egm96, egm2008 = map(float, t.split())
1665 while lon > _180_0: # EasternLon to earth lon
1666 lon -= _360_0
1667 yield GeoidHeight5Tuple(lat, lon, egm84, egm96, egm2008)
1670__all__ += _ALL_DOCS(_GeoidBase)
1672if __name__ == '__main__':
1674 from pygeodesy.lazily import printf, _sys
1676 _crop = ()
1677 _GeoidEGM = GeoidKarney
1678 _kind = 3
1680 geoids = _sys.argv[1:]
1681 while geoids:
1682 geoid = geoids.pop(0)
1684 if '-crop'.startswith(geoid.lower()):
1685 _crop = 20, -125, 50, -65 # CONUS
1687 elif '-karney'.startswith(geoid.lower()):
1688 _GeoidEGM = GeoidKarney
1690 elif '-kind'.startswith(geoid.lower()):
1691 _kind = int(geoids.pop(0))
1693 elif '-pgm'.startswith(geoid.lower()):
1694 _GeoidEGM = GeoidPGM
1696 elif geoid[-4:].lower() in ('.pgm',):
1697 g = _GeoidEGM(geoid, crop=_crop, kind=_kind)
1698 printf(g.toStr(), nt=1, nl=1)
1699 printf(repr(g.pgm), nt=1)
1700 # <https://GeographicLib.SourceForge.io/cgi-bin/GeoidEval>:
1701 # The height of the EGM96 geoid at Timbuktu
1702 # echo 16:46:33N 3:00:34W | GeoidEval
1703 # => 28.7068 -0.02e-6 -1.73e-6
1704 # The 1st number is the height of the geoid, the 2nd and
1705 # 3rd are its slopes in northerly and easterly direction
1706 t = 'Timbuktu %s' % (g,)
1707 k = {'egm84-15.pgm': '31.2979',
1708 'egm96-5.pgm': '28.7067',
1709 'egm2008-1.pgm': '28.7880'}.get(g.name.lower(), '28.7880')
1710 ll = parseDMS2('16:46:33N', '3:00:34W', sep=':')
1711 for ll in (ll, (16.776, -3.009),):
1712 try:
1713 h, ll = g.height(*ll), fstr(ll, prec=6)
1714 printf('%s.height(%s): %.4F vs %s', t, ll, h, k)
1715 except (GeoidError, RangeError) as x:
1716 printf(_COLONSPACE_(t, str(x)))
1718 elif geoid[-4:].lower() in ('.bin',):
1719 g = GeoidG2012B(geoid, kind=_kind)
1720 print(g.toStr())
1722 else:
1723 raise GeoidError(grid=repr(geoid))
1725_I = int # PYCHOK unused _I
1726del _intCs # trash ints cache
1728# **) MIT License
1729#
1730# Copyright (C) 2016-2023 -- mrJean1 at Gmail -- All Rights Reserved.
1731#
1732# Permission is hereby granted, free of charge, to any person obtaining a
1733# copy of this software and associated documentation files (the "Software"),
1734# to deal in the Software without restriction, including without limitation
1735# the rights to use, copy, modify, merge, publish, distribute, sublicense,
1736# and/or sell copies of the Software, and to permit persons to whom the
1737# Software is furnished to do so, subject to the following conditions:
1738#
1739# The above copyright notice and this permission notice shall be included
1740# in all copies or substantial portions of the Software.
1741#
1742# THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
1743# OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
1744# FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
1745# THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR
1746# OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE,
1747# ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR
1748# OTHER DEALINGS IN THE SOFTWARE.
1750# <https://GeographicLib.SourceForge.io/cgi-bin/GeoidEval>
1751# _lowerleft = -90, -179, -30.1500 # egm2008-1.pgm
1752# _lowerleft = -90, -179, -29.5350 # egm96-5.pgm
1753# _lowerleft = -90, -179, -29.7120 # egm84-15.pgm
1755# _center = 0, 0, 17.2260 # egm2008-1.pgm
1756# _center = 0, 0, 17.1630 # egm96-5.pgm
1757# _center = 0, 0, 18.3296 # egm84-15.pgm
1759# _upperright = 90, 180, 14.8980 # egm2008-1.pgm
1760# _upperright = 90, 180, 13.6050 # egm96-5.pgm
1761# _upperright = 90, 180, 13.0980 # egm84-15.pgm
1764# % python3 -m pygeodesy.geoids [-Karney] ../testGeoids/egm*.pgm
1765#
1766# 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)
1767#
1768# _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'
1769#
1770# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1771# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1772#
1773# 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)
1774#
1775# _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'
1776#
1777# Timbuktu GeoidKarney('egm84-15.pgm').height(16.775833, -3.009444): 31.2983 vs 31.2979
1778# Timbuktu GeoidKarney('egm84-15.pgm').height(16.776, -3.009): 31.2979 vs 31.2979
1779#
1780# 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)
1781#
1782# _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'
1783#
1784# Timbuktu GeoidKarney('egm96-5.pgm').height(16.775833, -3.009444): 28.7068 vs 28.7067
1785# Timbuktu GeoidKarney('egm96-5.pgm').height(16.776, -3.009): 28.7067 vs 28.7067
1788# % python3 -m pygeodesy.geoids -Karney ../testGeoids/egm*.pgm
1789#
1790# 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)
1791#
1792# _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'
1793#
1794# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1795# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1796#
1797# 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)
1798#
1799# _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'
1800#
1801# Timbuktu GeoidKarney('egm84-15.pgm').height(16.775833, -3.009444): 31.2983 vs 31.2979
1802# Timbuktu GeoidKarney('egm84-15.pgm').height(16.776, -3.009): 31.2979 vs 31.2979
1803#
1804# 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)
1805#
1806# _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'
1807#
1808# Timbuktu GeoidKarney('egm96-5.pgm').height(16.775833, -3.009444): 28.7068 vs 28.7067
1809# Timbuktu GeoidKarney('egm96-5.pgm').height(16.776, -3.009): 28.7067 vs 28.7067
1812# % python2 -m pygeodesy.geoids -PGM ../testGeoids/egm*.pgm
1813#
1814# 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)
1815#
1816# _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'
1817#
1818# Timbuktu GeoidPGM('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1819# Timbuktu GeoidPGM('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1820#
1821# 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)
1822#
1823# _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'
1824#
1825# Timbuktu GeoidPGM('egm84-15.pgm').height(16.775833, -3.009444): 31.2979 vs 31.2979
1826# Timbuktu GeoidPGM('egm84-15.pgm').height(16.776, -3.009): 31.2975 vs 31.2979
1827#
1828# 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)
1829#
1830# _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'
1831#
1832# Timbuktu GeoidPGM('egm96-5.pgm').height(16.775833, -3.009444): 28.7065 vs 28.7067
1833# Timbuktu GeoidPGM('egm96-5.pgm').height(16.776, -3.009): 28.7064 vs 28.7067