Coverage for pygeodesy/geoids.py: 95%
673 statements
« prev ^ index » next coverage.py v7.2.2, created at 2024-02-07 13:12 -0500
« prev ^ index » next coverage.py v7.2.2, created at 2024-02-07 13:12 -0500
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
86# from pygeodesy.datums import _ellipsoidal_datum # from .heights
87# from pygeodesy.dms import parseDMS2 # _MODS
88from pygeodesy.errors import _incompatible, LenError, RangeError, SciPyError, \
89 _SciPyIssue
90from pygeodesy.fmath import favg, Fdot, fdot, Fhorner, frange
91# from pygoedesy.formy import heightOrthometric # _MODS
92from pygeodesy.heights import _as_llis2, _ascalar, _height_called, HeightError, \
93 _HeightsBase, _ellipsoidal_datum, _Wrap
94from pygeodesy.interns import MISSING, NN, _4_, _COLONSPACE_, _COMMASPACE_, \
95 _cubic_, _DOT_, _E_, _height_, _in_, _kind_, \
96 _knots_, _lat_, _linear_, _lon_, _mean_, _N_, \
97 _n_a_, _not_, _numpy_, _on_, _outside_, _S_, \
98 _s_, _scipy_, _SPACE_, _stdev_, _supported_, \
99 _tbd_, _W_, _width_
100from pygeodesy.lazily import _ALL_DOCS, _ALL_LAZY, _ALL_MODS as _MODS, _FOR_DOCS
101from pygeodesy.named import _Named, _NamedTuple, notOverloaded
102# from pygeodesy.namedTuples import LatLon3Tuple # _MODS
103from pygeodesy.props import deprecated_method, Property_RO, property_RO
104from pygeodesy.streprs import attrs, Fmt, fstr, pairs
105from pygeodesy.units import Height, Int_, Lat, Lon
106# from pygeodesy.utily import _Wrap # from .heights
108from math import floor
109import os.path as _os_path
110from os import SEEK_CUR as _SEEK_CUR, SEEK_SET as _SEEK_SET
111from struct import calcsize as _calcsize, unpack as _unpack
113try:
114 from StringIO import StringIO as _BytesIO # reads bytes
115 _ub2str = str # PYCHOK convert bytes to str for egm*.pgm text
117except ImportError: # Python 3+
118 from io import BytesIO as _BytesIO # PYCHOK expected
120__all__ = _ALL_LAZY.geoids
121__version__ = '23.12.14'
123_assert_ = 'assert'
124_bHASH_ = b'#'
125_endian_ = 'endian'
126_format_ = '%s %r'
127_header_ = 'header'
128# temporarily hold a single instance for each int value
129_intCs = {}
130_interp2d_ks = {-2: _linear_,
131 -3: _cubic_,
132 -5: 'quintic'}
133_lli_ = 'lli'
134_non_increasing_ = 'non-increasing'
135_rb_ = 'rb'
138class _GeoidBase(_HeightsBase):
139 '''(INTERNAL) Base class for C{Geoid...}s.
140 '''
141 _cropped = None
142# _datum = _WGS84 # from _HeightsBase
143 _egm = None # open C{egm*.pgm} geoid file
144 _endian = _tbd_
145 _geoid = _n_a_
146 _hs_y_x = None # numpy 2darray, row-major order
147 _interp2d = None # interp2d interpolation
148 _kind = 3 # order for interp2d, RectBivariateSpline
149# _kmin = 2 # min number of knots
150 _knots = 0 # nlat * nlon
151 _mean = None # fixed in GeoidKarney
152# _name = NN # _Named
153 _nBytes = 0 # numpy size in bytes, float64
154 _pgm = None # PGM attributes, C{_PGM} or C{None}
155 _sizeB = 0 # geoid file size in bytes
156 _smooth = 0 # used only for RectBivariateSpline
157 _stdev = None # fixed in GeoidKarney
158 _u2B = 0 # np.itemsize or undefined
160 _lat_d = _0_0 # increment, +tive
161 _lat_lo = _0_0 # lower lat, south
162 _lat_hi = _0_0 # upper lat, noth
163 _lon_d = _0_0 # increment, +tive
164 _lon_lo = _0_0 # left lon, west
165 _lon_hi = _0_0 # right lon, east
166 _lon_of = _0_0 # forward lon offset
167 _lon_og = _0_0 # reverse lon offset
169 _center = None # (lat, lon, height)
170 _yx_hits = None # cache hits, ala Karney
172 def __init__(self, hs, p):
173 '''(INTERNAL) Set up the grid axes, the C{SciPy} interpolator
174 and several internal geoid attributes.
176 @arg hs: Grid knots with known height (C{numpy 2darray}).
177 @arg p: The C{slat, wlon, nlat, nlon, dlat, dlon} and
178 other geoid parameters (C{INTERNAL}).
180 @raise GeoidError: Incompatible grid B{C{hs}} shape or
181 invalid B{C{kind}}.
183 @raise LenError: Mismatch grid B{C{hs}} axis.
185 @raise SciPyError: A C{scipy.interpolate.inter2d} or
186 C{-.RectBivariateSpline} issue.
188 @raise SciPyWarning: A C{scipy.interpolate.inter2d} or
189 C{-.RectBivariateSpline} warning as
190 exception.
192 @note: C{scipy.interpolate.interp2d} has been C{DEPRECATED},
193 specify keyword argument C{B{kind}=1..5} to use
194 C{scipy.interpolate.RectBivariateSpline}.
195 '''
196 spi = self.scipy_interpolate
197 # for 2d scipy.interpolate.interp2d(xs, ys, hs, ...) and
198 # scipy.interpolate.RectBivariateSpline(ys, xs, hs, ...)
199 # require the shape of hs to be (len(ys), len(xs)), note
200 # the different (xs, ys, ...) and (ys, xs, ...) orders
201 if (p.nlat, p.nlon) != hs.shape:
202 raise GeoidError(shape=hs.shape, txt=_incompatible((p.nlat, p.nlon)))
204 # both axes and bounding box
205 ys, self._lat_d = self._gaxis2(p.slat, p.dlat, p.nlat, _lat_ + _s_)
206 xs, self._lon_d = self._gaxis2(p.wlon, p.dlon, p.nlon, _lon_ + _s_)
208 bb = ys[0], ys[-1], xs[0], xs[-1] + p.dlon # fudge lon_hi
209 # geoid grids are typically stored in row-major order, some
210 # with rows (90..-90) reversed and columns (0..360) wrapped
211 # to Easten longitude, 0 <= east < 180 and 180 <= west < 360
212 k = self.kind
213 if k in _interp2d_ks: # .interp2d DEPRECATED since scipy 1.10
214 if self._scipy_version() < (1, 10):
215 self._interp2d = spi.interp2d(xs, ys, hs, kind=_interp2d_ks[k])
216 else: # call and overwrite the deprecated .interp2d
217 self._interp2d = self._interp2d(xs, ys, hs, k)
218 elif 1 <= k <= 5:
219 self._ev = spi.RectBivariateSpline(ys, xs, hs, bbox=bb, ky=k, kx=k,
220 s=self._smooth).ev
221 else:
222 raise GeoidError(kind=k)
224 self._hs_y_x = hs # numpy 2darray, row-major
225 self._nBytes = hs.nbytes # numpy size in bytes
226 self._knots = p.knots # grid knots
227 self._lon_of = float(p.flon) # forward offset
228 self._lon_og = float(p.glon) # reverse offset
229 # shrink the box by 1 unit on every side
230 # bb += self._lat_d, -self._lat_d, self._lon_d, -self._lon_d
231 self._lat_lo = float(bb[0])
232 self._lat_hi = float(bb[1])
233 self._lon_lo = float(bb[2] - p.glon)
234 self._lon_hi = float(bb[3] - p.glon)
236 def __call__(self, *llis, **wrap_H):
237 '''Interpolate the geoid height for one or several locations.
239 @arg llis: One or more locations (C{LatLon}s), all positional.
240 @kwarg wrap_H: Keyword arguments C{B{wrap}=False, B{H}=False}.
241 If C{B{wrap} is True}, wrap or I{normalize} all
242 B{C{llis}} locations (C{bool}). If C{B{H} is True},
243 return the I{orthometric} height instead of the
244 I{geoid} height at each location (C{bool}).
246 @return: A single interpolated geoid (or orthometric) height
247 (C{float}) or a list or tuple of interpolated geoid
248 (or orthometric) heights (C{float}s).
250 @raise GeoidError: Insufficient number of B{C{llis}}, an
251 invalid B{C{lli}} or the C{egm*.pgm}
252 geoid file is closed.
254 @raise RangeError: An B{C{lli}} is outside this geoid's lat-
255 or longitude range.
257 @raise SciPyError: A C{scipy.interpolate.inter2d} or
258 C{-.RectBivariateSpline} issue.
260 @raise SciPyWarning: A C{scipy.interpolate.inter2d} or
261 C{-.RectBivariateSpline} warning as
262 exception.
264 @note: To obtain I{orthometric} heights, each B{C{llis}}
265 location must have an ellipsoid C{height} or C{h}
266 attribute, otherwise C{height=0} is used.
268 @see: Function L{pygeodesy.heightOrthometric}.
269 '''
270 return self._called(llis, True, **wrap_H)
272 def __enter__(self):
273 '''Open context.
274 '''
275 return self
277 def __exit__(self, *unused): # PYCHOK exc_type, exc_value, exc_traceback)
278 '''Close context.
279 '''
280 self.close()
281 # return None # XXX False
283 def __repr__(self):
284 return self.toStr()
286 def __str__(self):
287 return Fmt.PAREN(self.classname, repr(self.name))
289 def _called(self, llis, scipy, wrap=False, H=False):
290 # handle __call__
291 _H = _MODS.formy.heightOrthometric if H else None
292 _as, llis = _as_llis2(llis, Error=GeoidError)
293 try:
294 hs, _w = [], _Wrap._latlonop(wrap)
295 for i, lli in enumerate(llis):
296 N = self._hGeoid(*_w(lli.lat, lli.lon))
297 if _H: # convert to orthometric
298 N = _H(lli, N)
299 hs.append(N)
300 return _as(hs)
302 except (GeoidError, RangeError) as x:
303 # XXX avoid str(LatLon()) degree symbols
304 t = _lli_ if _as is _ascalar else Fmt.SQUARE(llis=i)
305 lli = fstr((lli.lat, lli.lon), strepr=repr)
306 raise type(x)(t, lli, wrap=wrap, H=H, cause=x)
307 except Exception as x:
308 if scipy and self.scipy:
309 raise _SciPyIssue(x)
310 else:
311 raise
313 @Property_RO
314 def _center(self):
315 ''' Cache for method L{center}.
316 '''
317 return self._llh3(favg(self._lat_lo, self._lat_hi),
318 favg(self._lon_lo, self._lon_hi))
320 def center(self, LatLon=None):
321 '''Return the center location and height of this geoid.
323 @kwarg LatLon: Optional class to return the location and height
324 (C{LatLon}) or C{None}.
326 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
327 lon, height)} otherwise a B{C{LatLon}} instance
328 with the lat-, longitude and geoid height of the
329 center grid location.
330 '''
331 return self._llh3LL(self._center, LatLon)
333 def close(self):
334 '''Close the C{egm*.pgm} geoid file if open (and applicable).
335 '''
336 if not self.closed:
337 self._egm.close()
338 self._egm = None
340 @property_RO
341 def closed(self):
342 '''Get the C{egm*.pgm} geoid file status.
343 '''
344 return self._egm is None
346 @Property_RO
347 def cropped(self):
348 '''Is geoid cropped (C{bool} or C{None} if crop not supported).
349 '''
350 return self._cropped
352 @Property_RO
353 def dtype(self):
354 '''Get the grid C{scipy} U{dtype<https://docs.SciPy.org/doc/numpy/
355 reference/generated/numpy.ndarray.dtype.html>} (C{numpy.dtype}).
356 '''
357 return self._hs_y_x.dtype
359 @Property_RO
360 def endian(self):
361 '''Get the geoid endianess and U{dtype<https://docs.SciPy.org/
362 doc/numpy/reference/generated/numpy.dtype.html>} (C{str}).
363 '''
364 return self._endian
366 def _ev(self, y, x): # PYCHOK expected
367 # only used for .interpolate.interp2d, but
368 # overwritten for .RectBivariateSpline,
369 # note (y, x) must be flipped!
370 return self._interp2d(x, y)
372 def _gaxis2(self, lo, d, n, name):
373 # build grid axis, hi = lo + (n - 1) * d
374 m, a = len2(frange(lo, n, d))
375 if m != n:
376 raise LenError(self.__class__, grid=m, **{name: n})
377 if d < 0:
378 d, a = -d, list(reversed(a))
379 for i in range(1, m):
380 e = a[i] - a[i-1]
381 if e < EPS: # non-increasing axis
382 i = Fmt.SQUARE(name, i)
383 raise GeoidError(i, e, txt=_non_increasing_)
384 return self.numpy.array(a), d
386 def _g2ll2(self, lat, lon): # PYCHOK no cover
387 '''(INTERNAL) I{Must be overloaded}.'''
388 notOverloaded(self, lat, lon)
390 def _gyx2g2(self, y, x):
391 # convert grid (y, x) indices to grid (lat, lon)
392 return ((self._lat_lo + self._lat_d * y),
393 (self._lon_lo + self._lon_of + self._lon_d * x))
395 def height(self, lats, lons, **wrap):
396 '''Interpolate the geoid height for one or several lat-/longitudes.
398 @arg lats: Latitude or latitudes (C{degrees} or C{degrees}s).
399 @arg lons: Longitude or longitudes (C{degrees} or C{degrees}s).
400 @kwarg wrap: If C{True}, wrap or I{normalize} all B{C{lats}}
401 and B{C{lons}} locations (C{bool}).
403 @return: A single interpolated geoid height (C{float}) or a
404 list of interpolated geoid heights (C{float}s).
406 @raise GeoidError: Insufficient or non-matching number of
407 B{C{lats}} and B{C{lons}}.
409 @raise RangeError: A B{C{lat}} or B{C{lon}} is outside this
410 geoid's lat- or longitude range.
412 @raise SciPyError: A C{scipy.interpolate.inter2d} or
413 C{-.RectBivariateSpline} issue.
415 @raise SciPyWarning: A C{scipy.interpolate.inter2d} or
416 C{-.RectBivariateSpline} warning as
417 exception.
418 '''
419 return _height_called(self, lats, lons, Error=GeoidError, **wrap)
421 def _hGeoid(self, lat, lon):
422 out = self.outside(lat, lon)
423 if out:
424 lli = fstr((lat, lon), strepr=repr)
425 raise RangeError(lli=lli, txt=_SPACE_(_outside_, _on_, out))
426 return float(self._ev(*self._ll2g2(lat, lon)))
428 @Property_RO
429 def _highest(self):
430 '''(INTERNAL) Cache for L{highest} method.
431 '''
432 return self._llh3minmax(True)
434 def highest(self, LatLon=None, **unused):
435 '''Return the location and largest height of this geoid.
437 @kwarg LatLon: Optional class to return the location and height
438 (C{LatLon}) or C{None}.
440 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
441 lon, height)} otherwise a B{C{LatLon}} instance
442 with the lat-, longitude and geoid height of the
443 highest grid location.
444 '''
445 return self._llh3LL(self._highest, LatLon)
447 @Property_RO
448 def hits(self):
449 '''Get the number of cache hits (C{int} or C{None}).
450 '''
451 return self._yx_hits
453 @deprecated_method
454 def _interp2d(self, xs, ys, hs=(), k=0): # overwritten in .__init__ above
455 '''DEPRECATED on 23.01.06, use keyword argument C{B{kind}=1..5}.'''
456 # assert k in _interp2d_ks # and len(hs) == len(xs) == len(ys)
457 try:
458 return self.scipy_interpolate.interp2d(xs, ys, hs, kind=_interp2d_ks[k])
459 except AttributeError as x:
460 raise SciPyError(interp2d=MISSING, kind=k, cause=x)
462 @Property_RO
463 def kind(self):
464 '''Get the interpolator kind and order (C{int}).
465 '''
466 return self._kind
468 @Property_RO
469 def knots(self):
470 '''Get the number of grid knots (C{int}).
471 '''
472 return self._knots
474 def _ll2g2(self, lat, lon): # PYCHOK no cover
475 '''(INTERNAL) I{Must be overloaded}.'''
476 notOverloaded(self, lat, lon)
478 @property_RO
479 def _LL3T(self):
480 '''(INTERNAL) Get L{LatLon3Tuple} once.
481 '''
482 T = _MODS.namedTuples.LatLon3Tuple
483 _GeoidBase._LL3T = T # overwrite poperty_RO
484 return T
486 def _llh3(self, lat, lon):
487 return self._LL3T(lat, lon, self._hGeoid(lat, lon), name=self.name)
489 def _llh3LL(self, llh, LatLon):
490 return llh if LatLon is None else self._xnamed(LatLon(*llh))
492 def _llh3minmax(self, highest=True, *unused):
493 hs, np = self._hs_y_x, self.numpy
494 # <https://docs.SciPy.org/doc/numpy/reference/generated/
495 # numpy.argmin.html#numpy.argmin>
496 arg = np.argmax if highest else np.argmin
497 y, x = np.unravel_index(arg(hs, axis=None), hs.shape)
498 return self._g2ll2(*self._gyx2g2(y, x)) + (float(hs[y, x]),)
500 def _load(self, g, dtype, n, offset=0):
501 # numpy.fromfile, like .frombuffer
502 g.seek(offset, _SEEK_SET)
503 return self.numpy.fromfile(g, dtype, n)
505 @Property_RO
506 def _lowerleft(self):
507 '''(INTERNAL) Cache for L{lowerleft}.
508 '''
509 return self._llh3(self._lat_lo, self._lon_lo)
511 def lowerleft(self, LatLon=None):
512 '''Return the lower-left location and height of this geoid.
514 @kwarg LatLon: Optional class to return the location
515 (C{LatLon}) and height or C{None}.
517 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
518 lon, height)} otherwise a B{C{LatLon}} instance
519 with the lat-, longitude and geoid height of the
520 lower-left, SW grid corner.
521 '''
522 return self._llh3LL(self._lowerleft, LatLon)
524 @Property_RO
525 def _loweright(self):
526 '''(INTERNAL) Cache for L{loweright}.
527 '''
528 return self._llh3(self._lat_lo, self._lon_hi)
530 def loweright(self, LatLon=None):
531 '''Return the lower-right location and height of this geoid.
533 @kwarg LatLon: Optional class to return the location and height
534 (C{LatLon}) or C{None}.
536 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
537 lon, height)} otherwise a B{C{LatLon}} instance
538 with the lat-, longitude and geoid height of the
539 lower-right, SE grid corner.
540 '''
542 return self._llh3LL(self._loweright, LatLon)
544 lowerright = loweright # synonymous
546 @Property_RO
547 def _lowest(self):
548 '''(INTERNAL) Cache for L{lowest}.
549 '''
550 return self._llh3minmax(False)
552 def lowest(self, LatLon=None, **unused):
553 '''Return the location and lowest height of this geoid.
555 @kwarg LatLon: Optional class to return the location and height
556 (C{LatLon}) or C{None}.
558 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
559 lon, height)} otherwise a B{C{LatLon}} instance
560 with the lat-, longitude and geoid height of the
561 lowest grid location.
562 '''
563 return self._llh3LL(self._lowest, LatLon)
565 @Property_RO
566 def mean(self):
567 '''Get the mean of this geoid's heights (C{float}).
568 '''
569 if self._mean is None: # see GeoidKarney
570 self._mean = float(self.numpy.mean(self._hs_y_x))
571 return self._mean
573 @property_RO
574 def name(self):
575 '''Get the name of this geoid (C{str}).
576 '''
577 return _HeightsBase.name.fget(self) or self._geoid # recursion
579 @Property_RO
580 def nBytes(self):
581 '''Get the grid in-memory size in bytes (C{int}).
582 '''
583 return self._nBytes
585 def _open(self, geoid, datum, kind, name, smooth):
586 # open the geoid file
587 try:
588 self._geoid = _os_path.basename(geoid)
589 self._sizeB = _os_path.getsize(geoid)
590 g = open(geoid, _rb_)
591 except (IOError, OSError) as x:
592 raise GeoidError(geoid=geoid, cause=x)
594 if datum not in (None, self._datum):
595 self._datum = _ellipsoidal_datum(datum, name=name)
596 self._kind = int(kind)
597 if name:
598 _HeightsBase.name.fset(self, name) # rename
599 if smooth:
600 self._smooth = Int_(smooth=smooth, Error=GeoidError, low=0)
602 return g
604 def outside(self, lat, lon):
605 '''Check whether a location is outside this geoid's
606 lat-/longitude or crop range.
608 @arg lat: The latitude (C{degrees}).
609 @arg lon: The longitude (C{degrees}).
611 @return: A 1- or 2-character C{str} if outside or an
612 empty C{str} if inside.
613 '''
614 return (_S_ if lat < self._lat_lo else
615 (_N_ if lat > self._lat_hi else NN)) + \
616 (_W_ if lon < self._lon_lo else
617 (_E_ if lon > self._lon_hi else NN))
619 @Property_RO
620 def pgm(self):
621 '''Get the PGM attributes (C{_PGM} or C{None} if not available/applicable).
622 '''
623 return self._pgm
625 @Property_RO
626 def sizeB(self):
627 '''Get the geoid grid file size in bytes (C{int}).
628 '''
629 return self._sizeB
631 @Property_RO
632 def smooth(self):
633 '''Get the C{RectBivariateSpline} smoothing (C{int}).
634 '''
635 return self._smooth
637 @Property_RO
638 def stdev(self):
639 '''Get the standard deviation of this geoid's heights (C{float}) or C{None}.
640 '''
641 if self._stdev is None: # see GeoidKarney
642 self._stdev = float(self.numpy.std(self._hs_y_x))
643 return self._stdev
645 def _swne(self, crop):
646 # crop box to 4-tuple (s, w, n, e)
647 try:
648 if len(crop) == 2:
649 try: # sw, ne LatLons
650 swne = (crop[0].lat, crop[0].lon,
651 crop[1].lat, crop[1].lon)
652 except AttributeError: # (s, w), (n, e)
653 swne = tuple(crop[0]) + tuple(crop[1])
654 else: # (s, w, n, e)
655 swne = crop
656 if len(swne) == 4:
657 s, w, n, e = map(float, swne)
658 if -90 <= s <= (n - _1_0) <= 89 and \
659 -180 <= w <= (e - _1_0) <= 179:
660 return s, w, n, e
661 except (IndexError, TypeError, ValueError):
662 pass
663 raise GeoidError(crop=crop)
665 def toStr(self, prec=3, sep=_COMMASPACE_): # PYCHOK signature
666 '''This geoid and all geoid attributes as a string.
668 @kwarg prec: Number of decimal digits (0..9 or C{None} for
669 default). Trailing zero decimals are stripped
670 for B{C{prec}} values of 1 and above, but kept
671 for negative B{C{prec}} values.
672 @kwarg sep: Separator to join (C{str}).
674 @return: Geoid name and attributes (C{str}).
675 '''
676 s = 1 if self.kind < 0 else 2
677 t = tuple(Fmt.PAREN(m.__name__, fstr(m(), prec=prec)) for m in
678 (self.lowerleft, self.upperright,
679 self.center,
680 self.highest, self.lowest)) + \
681 attrs( _mean_, _stdev_, prec=prec, Nones=False) + \
682 attrs((_kind_, 'smooth')[:s], prec=prec, Nones=False) + \
683 attrs( 'cropped', 'dtype', _endian_, 'hits', _knots_, 'nBytes',
684 'sizeB', _scipy_, _numpy_, prec=prec, Nones=False)
685 return _COLONSPACE_(self, sep.join(t))
687 @Property_RO
688 def u2B(self):
689 '''Get the PGM itemsize in bytes (C{int}).
690 '''
691 return self._u2B
693 @Property_RO
694 def _upperleft(self):
695 '''(INTERNAL) Cache for method L{upperleft}.
696 '''
697 return self._llh3(self._lat_hi, self._lon_lo)
699 def upperleft(self, LatLon=None):
700 '''Return the upper-left location and height of this geoid.
702 @kwarg LatLon: Optional class to return the location and height
703 (C{LatLon}) or C{None}.
705 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
706 lon, height)} otherwise a B{C{LatLon}} instance
707 with the lat-, longitude and geoid height of the
708 upper-left, NW grid corner.
709 '''
710 return self._llh3LL(self._upperleft, LatLon)
712 @Property_RO
713 def _upperright(self):
714 '''(INTERNAL) Cache for method L{upperright}.
715 '''
716 return self._llh3(self._lat_hi, self._lon_hi)
718 def upperright(self, LatLon=None):
719 '''Return the upper-right location and height of this geoid.
721 @kwarg LatLon: Optional class to return the location and height
722 (C{LatLon}) or C{None}.
724 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
725 lon, height)} otherwise a B{C{LatLon}} instance
726 with the lat-, longitude and geoid height of the
727 upper-right, NE grid corner.
728 '''
729 return self._llh3LL(self._upperright, LatLon)
732class GeoidError(HeightError):
733 '''Geoid interpolator C{Geoid...} or interpolation issue.
734 '''
735 pass
738class GeoidG2012B(_GeoidBase):
739 '''Geoid height interpolator for U{GEOID12B Model
740 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/>} grids U{CONUS
741 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_CONUS.shtml>},
742 U{Alaska<https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_AK.shtml>},
743 U{Hawaii<https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_HI.shtml>},
744 U{Guam and Northern Mariana Islands
745 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_GMNI.shtml>},
746 U{Puerto Rico and U.S. Virgin Islands
747 <https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_PRVI.shtml>} and
748 U{American Samoa<https://www.NGS.NOAA.gov/GEOID/GEOID12B/GEOID12B_AS.shtml>}
749 based on C{SciPy} U{RectBivariateSpline<https://docs.SciPy.org/doc/
750 scipy/reference/generated/scipy.interpolate.RectBivariateSpline.html>}
751 or U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/
752 scipy.interpolate.interp2d.html>} interpolation.
754 Use any of the binary C{le} (little endian) or C{be} (big endian)
755 C{g2012b*.bin} grid files.
756 '''
757 def __init__(self, g2012b_bin, crop=None, datum=None, # NAD 83 Ellipsoid
758 kind=3, name=NN, smooth=0):
759 '''New L{GeoidG2012B} interpolator.
761 @arg g2012b_bin: A C{GEOID12B} grid file name (C{.bin}).
762 @kwarg crop: Optional crop box, not supported (C{None}).
763 @kwarg datum: Optional grid datum (L{Datum}, L{Ellipsoid}, L{Ellipsoid2}
764 or L{a_f2Tuple}), default C{WGS84}.
765 @kwarg kind: C{scipy.interpolate} order (C{int}), use 1..5 for
766 U{RectBivariateSpline<https://docs.SciPy.org/doc/scipy/
767 reference/generated/scipy.interpolate.RectBivariateSpline.html>},
768 -2 for U{interp2d linear<https://docs.SciPy.org/doc/scipy/
769 reference/generated/scipy.interpolate.interp2d.html>}, -3
770 for C{interp2d cubic} or -5 for C{interp2d quintic}.
771 @kwarg name: Optional geoid name (C{str}).
772 @kwarg smooth: Smoothing factor for U{RectBivariateSpline
773 <https://docs.SciPy.org/doc/scipy/reference/generated/
774 scipy.interpolate.RectBivariateSpline.html>}
775 only (C{int}).
777 @raise GeoidError: G2012B grid file B{C{g2012b_bin}} issue or invalid
778 B{C{crop}}, B{C{kind}} or B{C{smooth}}.
780 @raise ImportError: Package C{numpy} or C{scipy} not found or
781 not installed.
783 @raise LenError: Grid file B{C{g2012b_bin}} axis mismatch.
785 @raise SciPyError: A C{RectBivariateSpline} or C{inter2d} issue.
787 @raise SciPyWarning: A C{RectBivariateSpline} or C{inter2d}
788 warning as exception.
790 @raise TypeError: Invalid B{C{datum}}.
792 @note: C{scipy.interpolate.interp2d} has been C{DEPRECATED}, specify
793 C{B{kind}=1..5} for C{scipy.interpolate.RectBivariateSpline}.
794 '''
795 if crop is not None:
796 raise GeoidError(crop=crop, txt=_not_(_supported_))
798 g = self._open(g2012b_bin, datum, kind, name, smooth)
799 _ = self.numpy # import numpy for ._load and
801 try:
802 p = _Gpars()
803 n = (self.sizeB // 4) - 11 # number of f4 heights
804 # U{numpy dtype formats are different from Python struct formats
805 # <https://docs.SciPy.org/doc/numpy-1.15.0/reference/arrays.dtypes.html>}
806 for en_ in ('<', '>'):
807 # skip 4xf8, get 3xi4
808 p.nlat, p.nlon, ien = map(int, self._load(g, en_+'i4', 3, 32))
809 if ien == 1: # correct endian
810 p.knots = p.nlat * p.nlon
811 if p.knots == n and 1 < p.nlat < n \
812 and 1 < p.nlon < n:
813 self._endian = en_+'f4'
814 break
815 else: # couldn't validate endian
816 raise GeoidError(_endian_)
818 # get the first 4xf8
819 p.slat, p.wlon, p.dlat, p.dlon = map(float, self._load(g, en_+'f8', 4))
820 # read all f4 heights, ignoring the first 4xf8 and 3xi4
821 hs = self._load(g, self._endian, n, 44).reshape(p.nlat, p.nlon)
822 p.wlon -= _360_0 # western-most East longitude to earth (..., lon)
823 _GeoidBase.__init__(self, hs, p)
825 except Exception as x:
826 raise _SciPyIssue(x, _in_, repr(g2012b_bin))
827 finally:
828 g.close()
830 def _g2ll2(self, lat, lon):
831 # convert grid (lat, lon) to earth (lat, lon)
832 return lat, lon
834 def _ll2g2(self, lat, lon):
835 # convert earth (lat, lon) to grid (lat, lon)
836 return lat, lon
838 if _FOR_DOCS:
839 __call__ = _GeoidBase.__call__
840 height = _GeoidBase.height
843class GeoidHeight5Tuple(_NamedTuple): # .geoids.py
844 '''5-Tuple C{(lat, lon, egm84, egm96, egm2008)} for U{GeoidHeights.dat
845 <https://SourceForge.net/projects/geographiclib/files/testdata/>}
846 tests with the heights for 3 different EGM grids at C{degrees90}
847 and C{degrees180} degrees (after converting C{lon} from original
848 C{0 <= EasterLon <= 360}).
849 '''
850 _Names_ = (_lat_, _lon_, 'egm84', 'egm96', 'egm2008')
851 _Units_ = ( Lat, Lon, Height, Height, Height)
854def _I(i):
855 '''(INTERNAL) Cache a single C{int} constant.
856 '''
857 return _intCs.setdefault(i, i) # PYCHOK undefined due to del _intCs
860def _T(*cs):
861 '''(INTERNAL) Cache a tuple of single C{int} constants.
862 '''
863 return map1(_I, *cs)
866class GeoidKarney(_GeoidBase):
867 '''Geoid height interpolator for I{Karney}'s U{GeographicLib Earth
868 Gravitational Model (EGM)<https://GeographicLib.SourceForge.io/C++/doc/
869 geoid.html>} geoid U{egm*.pgm<https://GeographicLib.SourceForge.io/
870 html/geoid.html#geoidinst>} datasets using bilinear or U{cubic
871 <https://dl.ACM.org/citation.cfm?id=368443>} interpolation and U{caching
872 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidcache>}
873 in pure Python, transcoded from I{Karney}'s U{C++ class Geoid
874 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinterp>}.
876 Use any of the geoid U{egm84-, egm96- or egm2008-*.pgm
877 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinst>}
878 datasets.
879 '''
880 _C0 = _F(372), _F(240), _F(372) # n, _ and s common denominators
881 # matrices c3n_, c3, c3s_, transposed from GeographicLib/Geoid.cpp
882 _C3 = ((_T(0, 0, 62, 124, 124, 62, 0, 0, 0, 0, 0, 0),
883 _T(0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0),
884 _T(-131, 7, -31, -62, -62, -31, 45, 216, 156, -45, -55, -7),
885 _T(0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0),
886 _T(138, -138, 0, 0, 0, 0, -183, 33, 153, -3, 48, -48), # PYCHOK indent
887 _T(144, 42, -62, -124, -124, -62, -9, 87, 99, 9, 42, -42),
888 _T(0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0),
889 _T(0, 0, 0, 0, 0, 0, 93, -93, -93, 93, 0, 0),
890 _T(-102, 102, 0, 0, 0, 0, 18, 12, -12, -18, -84, 84),
891 _T(-31, -31, 31, 62, 62, 31, 0, -93, -93, 0, 31, 31)), # PYCHOK indent
893 (_T(9, -9, 9, 186, 54, -9, -9, 54, -54, 9, -9, 9),
894 _T(-18, 18, -88, -42, 162, -32, 8, -78, 78, -8, 18, -18),
895 _T(-88, 8, -18, -42, -78, 18, 18, 162, 78, -18, -32, -8),
896 _T(0, 0, 90, -150, 30, 30, 30, -90, 90, -30, 0, 0),
897 _T(96, -96, 96, -96, -24, 24, -96, -24, 144, -24, 24, -24), # PYCHOK indent
898 _T(90, 30, 0, -150, -90, 0, 0, 30, 90, 0, 30, -30),
899 _T(0, 0, -20, 60, -60, 20, -20, 60, -60, 20, 0, 0),
900 _T(0, 0, -60, 60, 60, -60, 60, -60, -60, 60, 0, 0),
901 _T(-60, 60, 0, 60, -60, 0, 0, 60, -60, 0, -60, 60),
902 _T(-20, -20, 0, 60, 60, 0, 0, -60, -60, 0, 20, 20)),
904 (_T(18, -18, 36, 210, 162, -36, 0, 0, 0, 0, -18, 18), # PYCHOK indent
905 _T(-36, 36, -165, 45, 141, -21, 0, 0, 0, 0, 36, -36),
906 _T(-122, -2, -27, -111, -75, 27, 62, 124, 124, 62, -64, 2),
907 _T(0, 0, 93, -93, -93, 93, 0, 0, 0, 0, 0, 0),
908 _T(120, -120, 147, -57, -129, 39, 0, 0, 0, 0, 66, -66), # PYCHOK indent
909 _T(135, 51, -9, -192, -180, 9, 31, 62, 62, 31, 51, -51),
910 _T(0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0),
911 _T(0, 0, -93, 93, 93, -93, 0, 0, 0, 0, 0, 0),
912 _T(-84, 84, 18, 12, -12, -18, 0, 0, 0, 0, -102, 102),
913 _T(-31, -31, 0, 93, 93, 0, -31, -62, -62, -31, 31, 31)))
915 _BT = (_T(0, 0), # bilinear 4-tuple [i, j] indices
916 _T(1, 0),
917 _T(0, 1),
918 _T(1, 1))
920 _CM = (_T( 0, -1), # 10x12 cubic matrix [i, j] indices
921 _T( 1, -1),
922 _T(-1, 0),
923 _T( 0, 0),
924 _T( 1, 0),
925 _T( 2, 0),
926 _T(-1, 1),
927 _T( 0, 1),
928 _T( 1, 1),
929 _T( 2, 1),
930 _T( 0, 2),
931 _T( 1, 2))
933 _endian = '>H' # struct.unpack 1 ushort (big endian, unsigned short)
934 _4endian = '>4H' # struct.unpack 4 ushorts
935 _Rendian = NN # struct.unpack a row of ushorts
936# _highest = (-8.4, 147.367, 85.839) if egm2008-1.pgm else (
937# (-8.167, 147.25, 85.422) if egm96-5.pgm else
938# (-4.5, 148.75, 81.33)) # egm84-15.pgm
939# _lowest = (4.7, 78.767, -106.911) if egm2008-1.pgm else (
940# (4.667, 78.833, -107.043) if egm96-5.pgm else
941# (4.75, 79.25, -107.34)) # egm84-15.pgm
942 _mean = _F(-1.317) # from egm2008-1, -1.438 egm96-5, -0.855 egm84-15
943 _nBytes = None # not applicable
944 _nterms = len(_C3[0]) # columns length, number of row
945 _smooth = None # not applicable
946 _stdev = _F(29.244) # from egm2008-1, 29.227 egm96-5, 29.183 egm84-15
947 _u2B = _calcsize(_endian) # pixelsize_ in bytes
948 _4u2B = _calcsize(_4endian) # 4 pixelsize_s in bytes
949 _Ru2B = 0 # row of pixelsize_s in bytes
950 _yxH = () # cache (y, x) indices
951 _yxHt = () # cached 4- or 10-tuple for _ev2H resp. _ev3H
952 _yx_hits = 0 # cache hits
954 def __init__(self, egm_pgm, crop=None, datum=None, # WGS84
955 kind=3, name=NN, smooth=None):
956 '''New L{GeoidKarney} interpolator.
958 @arg egm_pgm: An U{EGM geoid dataset<https://GeographicLib.SourceForge.io/
959 html/geoid.html#geoidinst>} file name (C{egm*.pgm}), see
960 note below.
961 @kwarg crop: Optional box to limit geoid locations, a 4-tuple (C{south,
962 west, north, east}), 2-tuple (C{(south, west), (north,
963 east)}) or 2, in C{degrees90} lat- and C{degrees180}
964 longitudes or a 2-tuple (C{LatLonSW, LatLonNE}) of
965 C{LatLon} instances.
966 @kwarg datum: Optional grid datum (C{Datum}, L{Ellipsoid}, L{Ellipsoid2}
967 or L{a_f2Tuple}), default C{WGS84}.
968 @kwarg kind: Interpolation order (C{int}), 2 for C{bilinear} or 3
969 for C{cubic}.
970 @kwarg name: Optional geoid name (C{str}).
971 @kwarg smooth: Smoothing factor, unsupported (C{None}).
973 @raise GeoidError: EGM dataset B{C{egm_pgm}} issue or invalid
974 B{C{crop}}, B{C{kind}} or B{C{smooth}}.
976 @raise TypeError: Invalid B{C{datum}}.
978 @see: Class L{GeoidPGM} and function L{egmGeoidHeights}.
980 @note: Geoid file B{C{egm_pgm}} remains open and must be closed
981 by calling the C{close} method or by using this instance
982 in a C{with B{GeoidKarney}(...) as ...} context.
983 '''
984 if smooth is not None:
985 raise GeoidError(smooth=smooth, txt=_not_(_supported_))
987 if kind in (2,):
988 self._evH = self._ev2H
989 elif kind not in (3,):
990 raise GeoidError(kind=kind)
992 self._egm = g = self._open(egm_pgm, datum, kind, name, smooth)
993 self._pgm = p = _PGM(g, pgm=egm_pgm, itemsize=self.u2B, sizeB=self.sizeB)
995 self._Rendian = self._4endian.replace(_4_, str(p.nlon))
996 self._Ru2B = _calcsize(self._Rendian)
998 self._knots = p.knots # grid knots
999 self._lon_of = float(p.flon) # forward offset
1000 self._lon_og = float(p.glon) # reverse offset
1001 # set earth (lat, lon) limits (s, w, n, e)
1002 self._lat_lo, self._lon_lo, \
1003 self._lat_hi, self._lon_hi = self._swne(crop if crop else p.crop4)
1004 self._cropped = True if crop else False
1006 def __call__(self, *llis, **wrap_H):
1007 '''Interpolate the geoid height for one or several locations.
1009 @arg llis: One or more locations (C{LatLon}s), all positional.
1010 @kwarg wrap_H: Keyword arguments C{B{wrap}=False, B{H}=False}.
1011 If C{B{wrap} is True}, wrap or I{normalize} all
1012 B{C{llis}} locations (C{bool}). If C{B{H} is True},
1013 return the I{orthometric} height instead of the
1014 I{geoid} height at each location (C{bool}).
1016 @return: A single interpolated geoid (or orthometric) height
1017 (C{float}) or a list or tuple of interpolated geoid
1018 (or orthometric) heights (C{float}s).
1020 @raise GeoidError: Insufficient number of B{C{llis}}, an
1021 invalid B{C{lli}} or the C{egm*.pgm}
1022 geoid file is closed.
1024 @raise RangeError: An B{C{lli}} is outside this geoid's lat-
1025 or longitude range.
1027 @note: To obtain I{orthometric} heights, each B{C{llis}}
1028 location must have an ellipsoid C{height} or C{h}
1029 attribute, otherwise C{height=0} is used.
1031 @see: Function L{pygeodesy.heightOrthometric}.
1032 '''
1033 return self._called(llis, False, **wrap_H)
1035 def _c0c3v(self, y, x):
1036 # get the common denominator, the 10x12 cubic matrix and
1037 # the 12 cubic v-coefficients around geoid index (y, x)
1038 p = self._pgm
1039 if 0 < x < (p.nlon - 2) and 0 < y < (p.nlat - 2):
1040 # read 4x4 ushorts, drop the 4 corners
1041 g = self._egm
1042 e = self._4endian
1043 n = self._4u2B
1044 R = self._Ru2B
1046 b = self._seek(y - 1, x - 1)
1047 v = _unpack(e, g.read(n))[1:3]
1048 b += R
1049 g.seek(b, _SEEK_SET)
1050 v += _unpack(e, g.read(n))
1051 b += R
1052 g.seek(b, _SEEK_SET)
1053 v += _unpack(e, g.read(n))
1054 b += R
1055 g.seek(b, _SEEK_SET)
1056 v += _unpack(e, g.read(n))[1:3]
1057 j = 1
1059 else: # likely some wrapped y and/or x's
1060 v = self._raws(y, x, GeoidKarney._CM)
1061 j = 0 if y < 1 else (1 if y < (p.nlat - 2) else 2)
1063 return GeoidKarney._C0[j], GeoidKarney._C3[j], v
1065 @Property_RO
1066 def dtype(self):
1067 '''Get the geoid's grid data type (C{str}).
1068 '''
1069 return 'ushort'
1071 def _ev(self, lat, lon): # PYCHOK expected
1072 # interpolate the geoid height at grid (lat, lon)
1073 fy, fx = self._g2yx2(lat, lon)
1074 y, x = int(floor(fy)), int(floor(fx))
1075 fy -= y
1076 fx -= x
1077 H = self._evH(fy, fx, y, x) # ._ev3H or ._ev2H
1078 H *= self._pgm.Scale # H.fmul(self._pgm.Scale)
1079 H += self._pgm.Offset # H.fadd(self._pgm.Offset)
1080 return H.fsum()
1082 def _ev2H(self, fy, fx, *yx):
1083 # compute the bilinear 4-tuple and interpolate raw H
1084 if self._yxH == yx:
1085 t = self._yxHt
1086 self._yx_hits += 1
1087 else:
1088 y, x = self._yxH = yx
1089 self._yxHt = t = self._raws(y, x, GeoidKarney._BT)
1090 v = _1_0, -fx, fx
1091 H = Fdot(v, t[0], t[0], t[1]).fmul(_1_0 - fy) # c = a * (1 - fy)
1092 H += Fdot(v, t[2], t[2], t[3]).fmul(fy) # c += b * fy
1093 return H
1095 def _ev3H(self, fy, fx, *yx):
1096 # compute the cubic 10-tuple and interpolate raw H
1097 if self._yxH == yx:
1098 t = self._yxHt
1099 self._yx_hits += 1
1100 else:
1101 self._yxH = yx
1102 c0, c3, v = self._c0c3v(*yx)
1103 t = [fdot(v, *c3[i]) / c0 for i in range(self._nterms)]
1104 self._yxHt = t = tuple(t)
1105 # GeographicLib/Geoid.cpp Geoid::height(lat, lon) ...
1106 # real h = t[0] + fx * (t[1] + fx * (t[3] + fx * t[6])) +
1107 # fy * (t[2] + fx * (t[4] + fx * t[7]) +
1108 # fy * (t[5] + fx * t[8] + fy * t[9]));
1109 v = _1_0, fx, fy
1110 H = Fdot(v, t[5], t[8], t[9])
1111 H *= fy
1112 H += Fhorner(fx, t[2], t[4], t[7])
1113 H *= fy
1114 H += Fhorner(fx, t[0], t[1], t[3], t[6])
1115 return H
1117 _evH = _ev3H # overriden for kind == 2
1119 def _g2ll2(self, lat, lon):
1120 # convert grid (lat, lon) to earth (lat, lon), uncropped
1121 while lon > _180_0:
1122 lon -= _360_0
1123 return lat, lon
1125 def _g2yx2(self, lat, lon):
1126 # convert grid (lat, lon) to grid (y, x) indices
1127 p = self._pgm
1128 # note, slat = +90, rlat < 0 makes y >=0
1129 return ((lat - p.slat) * p.rlat), ((lon - p.wlon) * p.rlon)
1131 def _gyx2g2(self, y, x):
1132 # convert grid (y, x) indices to grid (lat, lon)
1133 p = self._pgm
1134 return (p.slat + p.dlat * y), (p.wlon + p.dlon * x)
1136 def height(self, lats, lons, **wrap):
1137 '''Interpolate the geoid height for one or several lat-/longitudes.
1139 @arg lats: Latitude or latitudes (C{degrees} or C{degrees}s).
1140 @arg lons: Longitude or longitudes (C{degrees} or C{degrees}s).
1141 @kwarg wrap: If C{True}, wrap or I{normalize} all B{C{lats}}
1142 and B{C{lons}} locations (C{bool}).
1144 @return: A single interpolated geoid height (C{float}) or a
1145 list of interpolated geoid heights (C{float}s).
1147 @raise GeoidError: Insufficient or non-matching number of
1148 B{C{lats}} and B{C{lons}} or the C{egm*.pgm}
1149 geoid file is closed.
1151 @raise RangeError: A B{C{lat}} or B{C{lon}} is outside this
1152 geoid's lat- or longitude range.
1153 '''
1154 return _height_called(self, lats, lons, Error=GeoidError, **wrap)
1156 @Property_RO
1157 def _highest_ltd(self):
1158 '''(INTERNAL) Cache for L{highest} mesthod.
1159 '''
1160 return self._llh3minmax(True, -12, -4)
1162 def highest(self, LatLon=None, full=False): # PYCHOK full
1163 '''Return the location and largest height of this geoid.
1165 @kwarg LatLon: Optional class to return the location and height
1166 (C{LatLon}) or C{None}.
1167 @kwarg full: Search the full or limited latitude range (C{bool}).
1169 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
1170 lon, height)} otherwise a B{C{LatLon}} instance
1171 with the lat-, longitude and geoid height at the
1172 highest grid location.
1173 '''
1174 llh = self._highest if full or self.cropped else self._highest_ltd
1175 return self._llh3LL(llh, LatLon)
1177 def _lat2y2(self, lat2):
1178 # convert earth lat(s) to min and max grid y indices
1179 ys, m = [], self._pgm.nlat - 1
1180 for lat in lat2:
1181 y, _ = self._g2yx2(*self._ll2g2(lat, 0))
1182 ys.append(max(min(int(y), m), 0))
1183 return min(ys), max(ys) + 1
1185 def _ll2g2(self, lat, lon):
1186 # convert earth (lat, lon) to grid (lat, lon), uncropped
1187 while lon < 0:
1188 lon += _360_0
1189 return lat, lon
1191 def _llh3minmax(self, highest=True, *lat2):
1192 # find highest or lowest, takes 10+ secs for egm2008-1.pgm geoid
1193 # (Python 2.7.16, macOS 10.13.6 High Sierra, iMac 3 GHz Core i3)
1194 y = x = 0
1195 h = self._raw(y, x)
1196 if highest:
1197 for j, r in self._raw2(*lat2):
1198 m = max(r)
1199 if m > h:
1200 h, y, x = m, j, r.index(m)
1201 else: # lowest
1202 for j, r in self._raw2(*lat2):
1203 m = min(r)
1204 if m < h:
1205 h, y, x = m, j, r.index(m)
1206 h *= self._pgm.Scale
1207 h += self._pgm.Offset
1208 return self._g2ll2(*self._gyx2g2(y, x)) + (h,)
1210 @Property_RO
1211 def _lowest_ltd(self):
1212 '''(INTERNAL) Cache for L{lowest}.
1213 '''
1214 return self._llh3minmax(False, 0, 8)
1216 def lowest(self, LatLon=None, full=False): # PYCHOK full
1217 '''Return the location and lowest height of this geoid.
1219 @kwarg LatLon: Optional class to return the location and height
1220 (C{LatLon}) or C{None}.
1221 @kwarg full: Search the full or limited latitude range (C{bool}).
1223 @return: If B{C{LatLon}} is C{None}, a L{LatLon3Tuple}C{(lat,
1224 lon, height)} otherwise a B{C{LatLon}} instance
1225 with the lat-, longitude and geoid height of the
1226 lowest grid location.
1227 '''
1228 llh = self._lowest if full or self.cropped else self._lowest_ltd
1229 return self._llh3LL(llh, LatLon)
1231 def _raw(self, y, x):
1232 # get the ushort geoid height at geoid index (y, x),
1233 # like GeographicLib/Geoid.hpp real rawval(is, iy)
1234 p = self._pgm
1235 if x < 0:
1236 x += p.nlon
1237 elif x >= p.nlon:
1238 x -= p.nlon
1239 h = p.nlon // 2
1240 if y < 0:
1241 y = -y
1242 elif y >= p.nlat:
1243 y = (p.nlat - 1) * 2 - y
1244 else:
1245 h = 0
1246 x += h if x < h else -h
1247 self._seek(y, x)
1248 h = _unpack(self._endian, self._egm.read(self._u2B))
1249 return h[0]
1251 def _raws(self, y, x, ijs):
1252 # get bilinear 4-tuple or 10x12 cubic matrix
1253 return tuple(self._raw(y + j, x + i) for i, j in ijs)
1255 def _raw2(self, *lat2):
1256 # yield a 2-tuple (y, ushorts) for each row or for
1257 # the rows between two (or more) earth lat values
1258 p = self._pgm
1259 g = self._egm
1260 e = self._Rendian
1261 n = self._Ru2B
1262 # min(lat2) <= lat <= max(lat2) or 0 <= y < p.nlat
1263 s, t = self._lat2y2(lat2) if lat2 else (0, p.nlat)
1264 self._seek(s, 0) # to start of row s
1265 for y in range(s, t):
1266 yield y, _unpack(e, g.read(n))
1268 def _seek(self, y, x):
1269 # position geoid to grid index (y, x)
1270 p, g = self._pgm, self._egm
1271 if g:
1272 b = p.skip + (y * p.nlon + x) * self._u2B
1273 g.seek(b, _SEEK_SET)
1274 return b # position
1275 raise GeoidError('closed file: %r' % (p.egm,)) # IOError
1278class GeoidPGM(_GeoidBase):
1279 '''Geoid height interpolator for I{Karney}'s U{GeographicLib Earth
1280 Gravitational Model (EGM)<https://GeographicLib.SourceForge.io/C++/doc/
1281 geoid.html>} geoid U{egm*.pgm<https://GeographicLib.SourceForge.io/
1282 html/geoid.html#geoidinst>} datasets but based on C{SciPy}
1283 U{RectBivariateSpline<https://docs.SciPy.org/doc/scipy/reference/
1284 generated/scipy.interpolate.RectBivariateSpline.html>} or
1285 U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/
1286 scipy.interpolate.interp2d.html>} interpolation.
1288 Use any of the U{egm84-, egm96- or egm2008-*.pgm
1289 <https://GeographicLib.SourceForge.io/C++/doc/geoid.html#geoidinst>}
1290 datasets. However, unless cropped, an entire C{egm*.pgm} dataset
1291 is loaded into the C{SciPy} U{RectBivariateSpline<https://docs.SciPy.org/
1292 doc/scipy/reference/generated/scipy.interpolate.RectBivariateSpline.html>}
1293 or U{interp2d<https://docs.SciPy.org/doc/scipy/reference/generated/
1294 scipy.interpolate.interp2d.html>} interpolator and converted from
1295 2-byte C{int} to 8-byte C{dtype float64}. Therefore, internal memory
1296 usage is 4x the U{egm*.pgm<https://GeographicLib.SourceForge.io/C++/doc/
1297 geoid.html#geoidinst>} file size and may exceed the available memory,
1298 especially with 32-bit Python, see properties C{.nBytes} and C{.sizeB}.
1299 '''
1300 _endian = '>u2'
1302 def __init__(self, egm_pgm, crop=None, datum=None, # WGS84
1303 kind=3, name=NN, smooth=0):
1304 '''New L{GeoidPGM} interpolator.
1306 @arg egm_pgm: An U{EGM geoid dataset<https://GeographicLib.SourceForge.io/
1307 html/geoid.html#geoidinst>} file name (C{egm*.pgm}).
1308 @kwarg crop: Optional box to crop B{C{egm_pgm}}, a 4-tuple (C{south, west,
1309 north, east}) or 2-tuple (C{(south, west), (north, east)}),
1310 in C{degrees90} lat- and C{degrees180} longitudes or a
1311 2-tuple (C{LatLonSW, LatLonNE}) of C{LatLon} instances.
1312 @kwarg datum: Optional grid datum (L{Datum}, L{Ellipsoid}, L{Ellipsoid2}
1313 or L{a_f2Tuple}), default C{WGS84}.
1314 @kwarg kind: C{scipy.interpolate} order (C{int}), use 1..5 for
1315 U{RectBivariateSpline<https://docs.SciPy.org/doc/scipy/
1316 reference/generated/scipy.interpolate.RectBivariateSpline.html>},
1317 -2 for U{interp2d linear<https://docs.SciPy.org/doc/scipy/
1318 reference/generated/scipy.interpolate.interp2d.html>}, -3
1319 for C{interp2d cubic} or -5 for C{interp2d quintic}.
1320 @kwarg name: Optional geoid name (C{str}).
1321 @kwarg smooth: Smoothing factor for U{RectBivariateSpline
1322 <https://docs.SciPy.org/doc/scipy/reference/generated/
1323 scipy.interpolate.RectBivariateSpline.html>}
1324 only (C{int}).
1326 @raise GeoidError: EGM dataset B{C{egm_pgm}} issue or invalid B{C{crop}},
1327 B{C{kind}} or B{C{smooth}}.
1329 @raise ImportError: Package C{numpy} or C{scipy} not found or not installed.
1331 @raise LenError: EGM dataset B{C{egm_pgm}} axis mismatch.
1333 @raise SciPyError: A C{RectBivariateSpline} or C{inter2d} issue.
1335 @raise SciPyWarning: A C{RectBivariateSpline} or C{inter2d}
1336 warning as exception.
1338 @raise TypeError: Invalid B{C{datum}}.
1340 @note: C{scipy.interpolate.interp2d} has been C{DEPRECATED}, specify
1341 C{B{kind}=1..5} for C{scipy.interpolate.RectBivariateSpline}.
1343 @note: The U{GeographicLib egm*.pgm<https://GeographicLib.SourceForge.io/
1344 html/geoid.html#geoidinst>} file sizes are based on a 2-byte
1345 C{int} height converted to 8-byte C{dtype float64} for C{scipy}
1346 interpolators. Therefore, internal memory usage is 4 times the
1347 C{egm*.pgm} file size and may exceed the available memory,
1348 especially with 32-bit Python. To reduce memory usage, set
1349 keyword argument B{C{crop}} to the region of interest. For example
1350 C{B{crop}=(20, -125, 50, -65)} covers the U{conterminous US<https://
1351 www.NGS.NOAA.gov/GEOID/GEOID12B/maps/GEOID12B_CONUS_grids.png>}
1352 (CONUS), less than 3% of the entire C{egm2008-1.pgm} dataset.
1354 @see: Class L{GeoidKarney} and function L{egmGeoidHeights}.
1355 '''
1356 np = self.numpy
1357 self._u2B = np.dtype(self.endian).itemsize
1359 g = self._open(egm_pgm, datum, kind, name, smooth)
1360 self._pgm = p = _PGM(g, pgm=egm_pgm, itemsize=self.u2B, sizeB=self.sizeB)
1361 if crop:
1362 g = p._cropped(g, abs(kind) + 1, *self._swne(crop))
1363 self._g2ll2 = self._g2ll2_cropped
1364 self._ll2g2 = self._ll2g2_cropped
1365 if map2(int, np.__version__.split(_DOT_)[:2]) < (1, 9):
1366 g = open(g.name, _rb_) # reopen tempfile for numpy 1.8.0-
1367 self._cropped = True
1368 else:
1369 self._cropped = False
1371 try:
1372 # U{numpy dtype formats are different from Python struct formats
1373 # <https://docs.SciPy.org/doc/numpy-1.15.0/reference/arrays.dtypes.html>}
1374 # read all heights, skipping the PGM header lines, converted to float
1375 hs = self._load(g, self.endian, p.knots, p.skip).reshape(p.nlat, p.nlon) * p.Scale
1376 if p.Offset: # offset
1377 hs = p.Offset + hs
1378 if p.dlat < 0: # flip the rows
1379 hs = np.flipud(hs)
1380 _GeoidBase.__init__(self, hs, p)
1381 except Exception as x:
1382 raise _SciPyIssue(x, _in_, repr(egm_pgm))
1383 finally:
1384 g.close()
1386 def _g2ll2(self, lat, lon):
1387 # convert grid (lat, lon) to earth (lat, lon), uncropped
1388 while lon > _180_0:
1389 lon -= _360_0
1390 return lat, lon
1392 def _g2ll2_cropped(self, lat, lon):
1393 # convert cropped grid (lat, lon) to earth (lat, lon)
1394 return lat, lon - self._lon_og
1396 def _ll2g2(self, lat, lon):
1397 # convert earth (lat, lon) to grid (lat, lon), uncropped
1398 while lon < 0:
1399 lon += _360_0
1400 return lat, lon
1402 def _ll2g2_cropped(self, lat, lon):
1403 # convert earth (lat, lon) to cropped grid (lat, lon)
1404 return lat, lon + self._lon_of
1406 if _FOR_DOCS:
1407 __call__ = _GeoidBase.__call__
1408 height = _GeoidBase.height
1411class _Gpars(_Named):
1412 '''(INTERNAL) Basic geoid parameters.
1413 '''
1414 # interpolator parameters
1415 dlat = 0 # +/- latitude resolution in C{degrees}
1416 dlon = 0 # longitude resolution in C{degrees}
1417 nlat = 1 # number of latitude knots (C{int})
1418 nlon = 0 # number of longitude knots (C{int})
1419 rlat = 0 # +/- latitude resolution in C{float}, 1 / .dlat
1420 rlon = 0 # longitude resolution in C{float}, 1 / .dlon
1421 slat = 0 # nothern- or southern most latitude (C{degrees90})
1422 wlon = 0 # western-most longitude in Eastern lon (C{degrees360})
1424 flon = 0 # forward, earth to grid longitude offset
1425 glon = 0 # reverse, grid to earth longitude offset
1427 knots = 0 # number of knots, nlat * nlon (C{int})
1428 skip = 0 # header bytes to skip (C{int})
1430 def __repr__(self):
1431 t = _COMMASPACE_.join(pairs((a, getattr(self, a)) for
1432 a in dir(self.__class__)
1433 if a[:1].isupper()))
1434 return _COLONSPACE_(self, t)
1436 def __str__(self):
1437 return Fmt.PAREN(self.classname, repr(self.name))
1440class _PGM(_Gpars):
1441 '''(INTERNAL) Parse an C{egm*.pgm} geoid dataset file.
1443 # Geoid file in PGM format for the GeographicLib::Geoid class
1444 # Description WGS84 EGM96, 5-minute grid
1445 # URL https://Earth-Info.NGA.mil/GandG/wgs84/gravitymod/egm96/egm96.html
1446 # DateTime 2009-08-29 18:45:03
1447 # MaxBilinearError 0.140
1448 # RMSBilinearError 0.005
1449 # MaxCubicError 0.003
1450 # RMSCubicError 0.001
1451 # Offset -108
1452 # Scale 0.003
1453 # Origin 90N 0E
1454 # AREA_OR_POINT Point
1455 # Vertical_Datum WGS84
1456 <width> <height>
1457 <pixel>
1458 ...
1459 '''
1460 crop4 = () # 4-tuple (C{south, west, north, east}).
1461 egm = None
1462 glon = 180 # reverse offset, uncropped
1463# pgm = NN # name
1464 sizeB = 0
1465 u2B = 2 # item size of grid height (C{int}).
1467 @staticmethod
1468 def _llstr2floats(latlon):
1469 # llstr to (lat, lon) floats
1470 lat, lon = latlon.split()
1471 return _MODS.dms.parseDMS2(lat, lon)
1473 # PGM file attributes, CamelCase but not .istitle()
1474 AREA_OR_POINT = str
1475 DateTime = str
1476 Description = str # 'WGS84 EGM96, 5-minute grid'
1477 Geoid = str # 'file in PGM format for the GeographicLib::Geoid class'
1478 MaxBilinearError = float
1479 MaxCubicError = float
1480 Offset = float
1481 Origin = _llstr2floats
1482 Pixel = 0
1483 RMSBilinearError = float
1484 RMSCubicError = float
1485 Scale = float
1486 URL = str # 'https://Earth-Info.NGA.mil/GandG/wgs84/...'
1487 Vertical_Datum = str
1489 def __init__(self, g, pgm=NN, itemsize=0, sizeB=0): # MCCABE 22
1490 '''(INTERNAL) New C{_PGM} parsed C{egm*.pgm} geoid dataset.
1491 '''
1492 self.name = pgm # geoid file name
1493 if itemsize:
1494 self._u2B = itemsize
1495 if sizeB:
1496 self.sizeB = sizeB
1498 t = g.readline() # make sure newline == '\n'
1499 if t != b'P5\n' and t.strip() != b'P5':
1500 raise self._Errorf(_format_, _header_, t)
1502 while True: # read all # Attr ... lines,
1503 try: # ignore empty ones or comments
1504 t = g.readline().strip()
1505 if t.startswith(_bHASH_):
1506 t = t.lstrip(_bHASH_).lstrip()
1507 a, v = map(_ub2str, t.split(None, 1))
1508 f = getattr(_PGM, a, None)
1509 if callable(f) and a[:1].isupper():
1510 setattr(self, a, f(v))
1511 elif t:
1512 break
1513 except (TypeError, ValueError):
1514 raise self._Errorf(_format_, 'Attr', t)
1515 else: # should never get here
1516 raise self._Errorf(_format_, _header_, g.tell())
1518 try: # must be (even) width and (odd) height
1519 nlon, nlat = map(int, t.split())
1520 if nlon < 2 or nlon > (360 * 60) or isodd(nlon) or \
1521 nlat < 2 or nlat > (181 * 60) or not isodd(nlat):
1522 raise ValueError
1523 except (TypeError, ValueError):
1524 raise self._Errorf(_format_, _SPACE_(_width_, _height_), t)
1526 try: # must be 16 bit pixel height
1527 t = g.readline().strip()
1528 self.Pixel = int(t)
1529 if not 255 < self.Pixel < 65536: # >u2 or >H only
1530 raise ValueError
1531 except (TypeError, ValueError):
1532 raise self._Errorf(_format_, 'pixel', t)
1534 for a in dir(_PGM): # set undefined # Attr ... to None
1535 if a[:1].isupper() and callable(getattr(self, a)):
1536 setattr(self, a, None)
1538 if self.Origin is None:
1539 raise self._Errorf(_format_, 'Origin', self.Origin)
1540 if self.Offset is None or self.Offset > 0:
1541 raise self._Errorf(_format_, 'Offset', self.Offset)
1542 if self.Scale is None or self.Scale < EPS:
1543 raise self._Errorf(_format_, 'Scale', self.Scale)
1545 self.skip = g.tell()
1546 self.knots = nlat * nlon
1548 self.nlat, self.nlon = nlat, nlon
1549 self.slat, self.wlon = self.Origin
1550 # note, negative .dlat and .rlat since rows
1551 # are from .slat 90N down in decreasing lat
1552 self.dlat, self.dlon = _180_0 / (1 - nlat), _360_0 / nlon
1553 self.rlat, self.rlon = (1 - nlat) / _180_0, nlon / _360_0
1555 # grid corners in earth (lat, lon), .slat = 90, .dlat < 0
1556 n = float(self.slat)
1557 s = n + self.dlat * (nlat - 1)
1558 w = self.wlon - self.glon
1559 e = w + self.dlon * nlon
1560 self.crop4 = s, w, n, e
1562 n = self.sizeB - self.skip
1563 if n > 0 and n != (self.knots * self.u2B):
1564 raise self._Errorf('%s(%s x %s != %s)', _assert_, nlat, nlon, n)
1566 def _cropped(self, g, k1, south, west, north, east): # MCCABE 15
1567 '''Crop the geoid to (south, west, north, east) box.
1568 '''
1569 # flon offset for both west and east
1570 f = 360 if west < 0 else 0
1571 # earth (lat, lon) to grid indices (y, x),
1572 # note y is decreasing, i.e. n < s
1573 s, w = self._lle2yx2(south, west, f)
1574 n, e = self._lle2yx2(north, east, f)
1575 s += 1 # s > n
1576 e += 1 # e > w
1578 hi, wi = self.nlat, self.nlon
1579 # handle special cases
1580 if (s - n) > hi:
1581 n, s = 0, hi # entire lat range
1582 if (e - w) > wi:
1583 w, e, f = 0, wi, 180 # entire lon range
1584 if s == hi and w == n == 0 and e == wi:
1585 return g # use entire geoid as-is
1587 if (e - w) < k1 or (s - n) < (k1 + 1):
1588 raise self._Errorf(_format_, 'swne', (north - south, east - west))
1590 if e > wi > w: # wrap around
1591 # read w..wi and 0..e
1592 r, p = (wi - w), (e - wi)
1593 elif e > w:
1594 r, p = (e - w), 0
1595 else:
1596 raise self._Errorf('%s(%s < %s)', _assert_, w, e)
1598 # convert to bytes
1599 r *= self.u2B
1600 p *= self.u2B
1601 q = wi * self.u2B # stride
1602 # number of rows and cols to skip from
1603 # the original (.slat, .wlon) origin
1604 z = self.skip + (n * wi + w) * self.u2B
1605 # sanity check
1606 if r < 2 or p < 0 or q < 2 or z < self.skip \
1607 or z > self.sizeB:
1608 raise self._Errorf(_format_, _assert_, (r, p, q, z))
1610 # can't use _BytesIO since numpy
1611 # needs .fileno attr in .fromfile
1612 t, c = 0, self._tmpfile()
1613 # reading (s - n) rows, forward
1614 for y in range(n, s): # PYCHOK y unused
1615 g.seek(z, _SEEK_SET)
1616 # Python 2 tmpfile.write returns None
1617 t += c.write(g.read(r)) or r
1618 if p: # wrap around to start of row
1619 g.seek(-q, _SEEK_CUR)
1620 # assert(g.tell() == (z - w * self.u2B))
1621 # Python 2 tmpfile.write returns None
1622 t += c.write(g.read(p)) or p
1623 z += q
1624 c.flush()
1625 g.close()
1627 s -= n # nlat
1628 e -= w # nlon
1629 k = s * e # knots
1630 z = k * self.u2B
1631 if t != z:
1632 raise self._Errorf('%s(%s != %s) %s', _assert_, t, z, self)
1634 # update the _Gpars accordingly, note attributes
1635 # .dlat, .dlon, .rlat and .rlon remain unchanged
1636 self.slat += n * self.dlat
1637 self.wlon += w * self.dlon
1638 self.nlat = s
1639 self.nlon = e
1640 self.flon = self.glon = f
1642 self.crop4 = south, west, north, east
1643 self.knots = k
1644 self.skip = 0 # no header lines in c
1646 c.seek(0, _SEEK_SET)
1647 # c = open(c.name, _rb_) # reopen for numpy 1.8.0-
1648 return c
1650 def _Errorf(self, fmt, *args): # PYCHOK no cover
1651 t = fmt % args
1652 e = self.pgm or NN
1653 if e:
1654 t = _SPACE_(t, _in_, repr(e))
1655 return PGMError(t)
1657 def _lle2yx2(self, lat, lon, flon):
1658 # earth (lat, lon) to grid indices (y, x)
1659 # with .dlat decreasing from 90N .slat
1660 lat -= self.slat
1661 lon += flon - self.wlon
1662 return (min(self.nlat - 1, max(0, int(lat * self.rlat))),
1663 max(0, int(lon * self.rlon)))
1665 def _tmpfile(self):
1666 # create a tmpfile to hold the cropped geoid grid
1667 try:
1668 from tempfile import NamedTemporaryFile as tmpfile
1669 except ImportError: # Python 2.7.16-
1670 from os import tmpfile
1671 t = _os_path.splitext(_os_path.basename(self.pgm))[0]
1672 f = tmpfile(mode='w+b', prefix=t or 'egm')
1673 f.seek(0, _SEEK_SET) # force overwrite
1674 return f
1676 @Property_RO
1677 def pgm(self):
1678 '''Get the geoid file name (C{str}).
1679 '''
1680 return self.name
1683class PGMError(GeoidError):
1684 '''Issue parsing or cropping an C{egm*.pgm} geoid dataset.
1685 '''
1686 pass
1689def egmGeoidHeights(GeoidHeights_dat):
1690 '''Generate geoid U{egm*.pgm<https://GeographicLib.SourceForge.io/
1691 html/geoid.html#geoidinst>} height tests from U{GeoidHeights.dat
1692 <https://SourceForge.net/projects/geographiclib/files/testdata/>}
1693 U{Test data for Geoids<https://GeographicLib.SourceForge.io/C++/doc/
1694 geoid.html#testgeoid>}.
1696 @arg GeoidHeights_dat: The un-gz-ed C{GeoidHeights.dat} file
1697 (C{str} or C{file} handle).
1699 @return: For each test, yield a L{GeoidHeight5Tuple}C{(lat, lon,
1700 egm84, egm96, egm2008)}.
1702 @raise GeoidError: Invalid B{C{GeoidHeights_dat}}.
1704 @note: Function L{egmGeoidHeights} is used to test the geoids
1705 L{GeoidKarney} and L{GeoidPGM}, see PyGeodesy module
1706 C{test/testGeoids.py}.
1707 '''
1708 dat = GeoidHeights_dat
1709 if isinstance(dat, bytes):
1710 dat = _BytesIO(dat)
1712 try:
1713 dat.seek(0, _SEEK_SET) # reset
1714 except AttributeError as x:
1715 raise GeoidError(GeoidHeights_dat=type(dat), cause=x)
1717 for t in dat.readlines():
1718 t = t.strip()
1719 if t and not t.startswith(_bHASH_):
1720 lat, lon, egm84, egm96, egm2008 = map(float, t.split())
1721 while lon > _180_0: # EasternLon to earth lon
1722 lon -= _360_0
1723 yield GeoidHeight5Tuple(lat, lon, egm84, egm96, egm2008)
1726__all__ += _ALL_DOCS(_GeoidBase)
1728if __name__ == '__main__':
1730 from pygeodesy.lazily import printf, _sys
1732 _crop = ()
1733 _GeoidEGM = GeoidKarney
1734 _kind = 3
1736 geoids = _sys.argv[1:]
1737 while geoids:
1738 geoid = geoids.pop(0)
1740 if '-crop'.startswith(geoid.lower()):
1741 _crop = 20, -125, 50, -65 # CONUS
1743 elif '-karney'.startswith(geoid.lower()):
1744 _GeoidEGM = GeoidKarney
1746 elif '-kind'.startswith(geoid.lower()):
1747 _kind = int(geoids.pop(0))
1749 elif '-pgm'.startswith(geoid.lower()):
1750 _GeoidEGM = GeoidPGM
1752 elif geoid[-4:].lower() in ('.pgm',):
1753 g = _GeoidEGM(geoid, crop=_crop, kind=_kind)
1754 printf(g.toStr(), nt=1, nl=1)
1755 printf(repr(g.pgm), nt=1)
1756 # <https://GeographicLib.SourceForge.io/cgi-bin/GeoidEval>:
1757 # The height of the EGM96 geoid at Timbuktu
1758 # echo 16:46:33N 3:00:34W | GeoidEval
1759 # => 28.7068 -0.02e-6 -1.73e-6
1760 # The 1st number is the height of the geoid, the 2nd and
1761 # 3rd are its slopes in northerly and easterly direction
1762 t = 'Timbuktu %s' % (g,)
1763 k = {'egm84-15.pgm': '31.2979',
1764 'egm96-5.pgm': '28.7067',
1765 'egm2008-1.pgm': '28.7880'}.get(g.name.lower(), '28.7880')
1766 ll = _MODS.dms.parseDMS2('16:46:33N', '3:00:34W', sep=':')
1767 for ll in (ll, (16.776, -3.009),):
1768 try:
1769 h, ll = g.height(*ll), fstr(ll, prec=6)
1770 printf('%s.height(%s): %.4F vs %s', t, ll, h, k)
1771 except (GeoidError, RangeError) as x:
1772 printf(_COLONSPACE_(t, str(x)))
1774 elif geoid[-4:].lower() in ('.bin',):
1775 g = GeoidG2012B(geoid, kind=_kind)
1776 printf(g.toStr())
1778 else:
1779 raise GeoidError(grid=repr(geoid))
1781_I = int # PYCHOK unused _I
1782del _intCs # trash ints cache
1784# **) MIT License
1785#
1786# Copyright (C) 2016-2024 -- mrJean1 at Gmail -- All Rights Reserved.
1787#
1788# Permission is hereby granted, free of charge, to any person obtaining a
1789# copy of this software and associated documentation files (the "Software"),
1790# to deal in the Software without restriction, including without limitation
1791# the rights to use, copy, modify, merge, publish, distribute, sublicense,
1792# and/or sell copies of the Software, and to permit persons to whom the
1793# Software is furnished to do so, subject to the following conditions:
1794#
1795# The above copyright notice and this permission notice shall be included
1796# in all copies or substantial portions of the Software.
1797#
1798# THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
1799# OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
1800# FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
1801# THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR
1802# OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE,
1803# ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR
1804# OTHER DEALINGS IN THE SOFTWARE.
1806# <https://GeographicLib.SourceForge.io/cgi-bin/GeoidEval>
1807# _lowerleft = -90, -179, -30.1500 # egm2008-1.pgm
1808# _lowerleft = -90, -179, -29.5350 # egm96-5.pgm
1809# _lowerleft = -90, -179, -29.7120 # egm84-15.pgm
1811# _center = 0, 0, 17.2260 # egm2008-1.pgm
1812# _center = 0, 0, 17.1630 # egm96-5.pgm
1813# _center = 0, 0, 18.3296 # egm84-15.pgm
1815# _upperright = 90, 180, 14.8980 # egm2008-1.pgm
1816# _upperright = 90, 180, 13.6050 # egm96-5.pgm
1817# _upperright = 90, 180, 13.0980 # egm84-15.pgm
1820# % python3 -m pygeodesy.geoids [-Karney] ../testGeoids/egm*.pgm
1821#
1822# 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)
1823#
1824# _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'
1825#
1826# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1827# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1828#
1829# 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)
1830#
1831# _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'
1832#
1833# Timbuktu GeoidKarney('egm84-15.pgm').height(16.775833, -3.009444): 31.2983 vs 31.2979
1834# Timbuktu GeoidKarney('egm84-15.pgm').height(16.776, -3.009): 31.2979 vs 31.2979
1835#
1836# 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)
1837#
1838# _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'
1839#
1840# Timbuktu GeoidKarney('egm96-5.pgm').height(16.775833, -3.009444): 28.7068 vs 28.7067
1841# Timbuktu GeoidKarney('egm96-5.pgm').height(16.776, -3.009): 28.7067 vs 28.7067
1844# % python3 -m pygeodesy.geoids -Karney ../testGeoids/egm*.pgm
1845#
1846# 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)
1847#
1848# _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'
1849#
1850# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1851# Timbuktu GeoidKarney('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1852#
1853# 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)
1854#
1855# _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'
1856#
1857# Timbuktu GeoidKarney('egm84-15.pgm').height(16.775833, -3.009444): 31.2983 vs 31.2979
1858# Timbuktu GeoidKarney('egm84-15.pgm').height(16.776, -3.009): 31.2979 vs 31.2979
1859#
1860# 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)
1861#
1862# _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'
1863#
1864# Timbuktu GeoidKarney('egm96-5.pgm').height(16.775833, -3.009444): 28.7068 vs 28.7067
1865# Timbuktu GeoidKarney('egm96-5.pgm').height(16.776, -3.009): 28.7067 vs 28.7067
1868# % python2 -m pygeodesy.geoids -PGM ../testGeoids/egm*.pgm
1869#
1870# 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)
1871#
1872# _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'
1873#
1874# Timbuktu GeoidPGM('egm2008-1.pgm').height(16.775833, -3.009444): 28.7881 vs 28.7880
1875# Timbuktu GeoidPGM('egm2008-1.pgm').height(16.776, -3.009): 28.7880 vs 28.7880
1876#
1877# 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)
1878#
1879# _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'
1880#
1881# Timbuktu GeoidPGM('egm84-15.pgm').height(16.775833, -3.009444): 31.2979 vs 31.2979
1882# Timbuktu GeoidPGM('egm84-15.pgm').height(16.776, -3.009): 31.2975 vs 31.2979
1883#
1884# 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)
1885#
1886# _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'
1887#
1888# Timbuktu GeoidPGM('egm96-5.pgm').height(16.775833, -3.009444): 28.7065 vs 28.7067
1889# Timbuktu GeoidPGM('egm96-5.pgm').height(16.776, -3.009): 28.7064 vs 28.7067