GeographicLib  1.37
Geoid.cpp
Go to the documentation of this file.
1 /**
2  * \file Geoid.cpp
3  * \brief Implementation for GeographicLib::Geoid class
4  *
5  * Copyright (c) Charles Karney (2009-2012) <charles@karney.com> and licensed
6  * under the MIT/X11 License. For more information, see
7  * http://geographiclib.sourceforge.net/
8  **********************************************************************/
9 
10 #include <GeographicLib/Geoid.hpp>
11 #include <cstdlib>
13 
14 #if !defined(GEOGRAPHICLIB_DATA)
15 # if defined(_WIN32)
16 # define GEOGRAPHICLIB_DATA "C:/ProgramData/GeographicLib"
17 # else
18 # define GEOGRAPHICLIB_DATA "/usr/local/share/GeographicLib"
19 # endif
20 #endif
21 
22 #if !defined(GEOGRAPHICLIB_GEOID_DEFAULT_NAME)
23 # define GEOGRAPHICLIB_GEOID_DEFAULT_NAME "egm96-5"
24 #endif
25 
26 #if defined(_MSC_VER)
27 // Squelch warnings about unsafe use of getenv
28 # pragma warning (disable: 4996)
29 #endif
30 
31 namespace GeographicLib {
32 
33  using namespace std;
34 
35  // This is the transfer matrix for a 3rd order fit with a 12-point stencil
36  // with weights
37  //
38  // \x -1 0 1 2
39  // y
40  // -1 . 1 1 .
41  // 0 1 2 2 1
42  // 1 1 2 2 1
43  // 2 . 1 1 .
44  //
45  // A algorithm for n-dimensional polynomial fits is described in
46  // F. H. Lesh,
47  // Multi-dimensional least-squares polynomial curve fitting,
48  // CACM 2, 29-30 (1959).
49  //
50  // Here's the Maxima code to generate this matrix:
51  //
52  // /* The stencil and the weights */
53  // xarr:[
54  // 0, 1,
55  // -1, 0, 1, 2,
56  // -1, 0, 1, 2,
57  // 0, 1]$
58  // yarr:[
59  // -1,-1,
60  // 0, 0, 0, 0,
61  // 1, 1, 1, 1,
62  // 2, 2]$
63  // warr:[
64  // 1, 1,
65  // 1, 2, 2, 1,
66  // 1, 2, 2, 1,
67  // 1, 1]$
68  //
69  // /* [x exponent, y exponent] for cubic fit */
70  // pows:[
71  // [0,0],
72  // [1,0],[0,1],
73  // [2,0],[1,1],[0,2],
74  // [3,0],[2,1],[1,2],[0,3]]$
75  //
76  // basisvec(x,y,pows):=map(lambda([ex],(if ex[1]=0 then 1 else x^ex[1])*
77  // (if ex[2]=0 then 1 else y^ex[2])),pows)$
78  // addterm(x,y,f,w,pows):=block([a,b,bb:basisvec(x,y,pows)],
79  // a:w*(transpose(bb).bb),
80  // b:(w*f) * bb,
81  // [a,b])$
82  //
83  // c3row(k):=block([a,b,c,pows:pows,n],
84  // n:length(pows),
85  // a:zeromatrix(n,n),
86  // b:copylist(part(a,1)),
87  // c:[a,b],
88  // for i:1 thru length(xarr) do
89  // c:c+addterm(xarr[i],yarr[i],if i=k then 1 else 0,warr[i],pows),
90  // a:c[1],b:c[2],
91  // part(transpose( a^^-1 . transpose(b)),1))$
92  // c3:[]$
93  // for k:1 thru length(warr) do c3:endcons(c3row(k),c3)$
94  // c3:apply(matrix,c3)$
95  // c0:part(ratsimp(
96  // genmatrix(yc,1,length(warr)).abs(c3).genmatrix(yd,length(pows),1)),2)$
97  // c3:c0*c3$
98 
99  const int Geoid::c0_ = 240; // Common denominator
100  const int Geoid::c3_[stencilsize_ * nterms_] = {
101  9, -18, -88, 0, 96, 90, 0, 0, -60, -20,
102  -9, 18, 8, 0, -96, 30, 0, 0, 60, -20,
103  9, -88, -18, 90, 96, 0, -20, -60, 0, 0,
104  186, -42, -42, -150, -96, -150, 60, 60, 60, 60,
105  54, 162, -78, 30, -24, -90, -60, 60, -60, 60,
106  -9, -32, 18, 30, 24, 0, 20, -60, 0, 0,
107  -9, 8, 18, 30, -96, 0, -20, 60, 0, 0,
108  54, -78, 162, -90, -24, 30, 60, -60, 60, -60,
109  -54, 78, 78, 90, 144, 90, -60, -60, -60, -60,
110  9, -8, -18, -30, -24, 0, 20, 60, 0, 0,
111  -9, 18, -32, 0, 24, 30, 0, 0, -60, 20,
112  9, -18, -8, 0, -24, -30, 0, 0, 60, 20,
113  };
114 
115  // Like c3, but with the coeffs of x, x^2, and x^3 constrained to be zero.
116  // Use this at the N pole so that the height in independent of the longitude
117  // there.
118  //
119  // Here's the Maxima code to generate this matrix (continued from above).
120  //
121  // /* figure which terms to exclude so that fit is indep of x at y=0 */
122  // mask:part(zeromatrix(1,length(pows)),1)+1$
123  // for i:1 thru length(pows) do
124  // if pows[i][1]>0 and pows[i][2]=0 then mask[i]:0$
125  //
126  // /* Same as c3row but with masked pows. */
127  // c3nrow(k):=block([a,b,c,powsa:[],n,d,e],
128  // for i:1 thru length(mask) do if mask[i]>0 then
129  // powsa:endcons(pows[i],powsa),
130  // n:length(powsa),
131  // a:zeromatrix(n,n),
132  // b:copylist(part(a,1)),
133  // c:[a,b],
134  // for i:1 thru length(xarr) do
135  // c:c+addterm(xarr[i],yarr[i],if i=k then 1 else 0,warr[i],powsa),
136  // a:c[1],b:c[2],
137  // d:part(transpose( a^^-1 . transpose(b)),1),
138  // e:[],
139  // for i:1 thru length(mask) do
140  // if mask[i]>0 then (e:endcons(first(d),e),d:rest(d)) else e:endcons(0,e),
141  // e)$
142  // c3n:[]$
143  // for k:1 thru length(warr) do c3n:endcons(c3nrow(k),c3n)$
144  // c3n:apply(matrix,c3n)$
145  // c0n:part(ratsimp(
146  // genmatrix(yc,1,length(warr)).abs(c3n).genmatrix(yd,length(pows),1)),2)$
147  // c3n:c0n*c3n$
148 
149  const int Geoid::c0n_ = 372; // Common denominator
150  const int Geoid::c3n_[stencilsize_ * nterms_] = {
151  0, 0, -131, 0, 138, 144, 0, 0, -102, -31,
152  0, 0, 7, 0, -138, 42, 0, 0, 102, -31,
153  62, 0, -31, 0, 0, -62, 0, 0, 0, 31,
154  124, 0, -62, 0, 0, -124, 0, 0, 0, 62,
155  124, 0, -62, 0, 0, -124, 0, 0, 0, 62,
156  62, 0, -31, 0, 0, -62, 0, 0, 0, 31,
157  0, 0, 45, 0, -183, -9, 0, 93, 18, 0,
158  0, 0, 216, 0, 33, 87, 0, -93, 12, -93,
159  0, 0, 156, 0, 153, 99, 0, -93, -12, -93,
160  0, 0, -45, 0, -3, 9, 0, 93, -18, 0,
161  0, 0, -55, 0, 48, 42, 0, 0, -84, 31,
162  0, 0, -7, 0, -48, -42, 0, 0, 84, 31,
163  };
164 
165  // Like c3n, but y -> 1-y so that h is independent of x at y = 1. Use this
166  // at the S pole so that the height in independent of the longitude there.
167  //
168  // Here's the Maxima code to generate this matrix (continued from above).
169  //
170  // /* Transform c3n to c3s by transforming y -> 1-y */
171  // vv:[
172  // v[11],v[12],
173  // v[7],v[8],v[9],v[10],
174  // v[3],v[4],v[5],v[6],
175  // v[1],v[2]]$
176  // poly:expand(vv.(c3n/c0n).transpose(basisvec(x,1-y,pows)))$
177  // c3sf[i,j]:=coeff(coeff(coeff(poly,v[i]),x,pows[j][1]),y,pows[j][2])$
178  // c3s:genmatrix(c3sf,length(vv),length(pows))$
179  // c0s:part(ratsimp(
180  // genmatrix(yc,1,length(warr)).abs(c3s).genmatrix(yd,length(pows),1)),2)$
181  // c3s:c0s*c3s$
182 
183  const int Geoid::c0s_ = 372; // Common denominator
184  const int Geoid::c3s_[stencilsize_ * nterms_] = {
185  18, -36, -122, 0, 120, 135, 0, 0, -84, -31,
186  -18, 36, -2, 0, -120, 51, 0, 0, 84, -31,
187  36, -165, -27, 93, 147, -9, 0, -93, 18, 0,
188  210, 45, -111, -93, -57, -192, 0, 93, 12, 93,
189  162, 141, -75, -93, -129, -180, 0, 93, -12, 93,
190  -36, -21, 27, 93, 39, 9, 0, -93, -18, 0,
191  0, 0, 62, 0, 0, 31, 0, 0, 0, -31,
192  0, 0, 124, 0, 0, 62, 0, 0, 0, -62,
193  0, 0, 124, 0, 0, 62, 0, 0, 0, -62,
194  0, 0, 62, 0, 0, 31, 0, 0, 0, -31,
195  -18, 36, -64, 0, 66, 51, 0, 0, -102, 31,
196  18, -36, 2, 0, -66, -51, 0, 0, 102, 31,
197  };
198 
199  Geoid::Geoid(const std::string& name, const std::string& path, bool cubic,
200  bool threadsafe)
201  : _name(name)
202  , _dir(path)
203  , _cubic(cubic)
204  , _a( Constants::WGS84_a() )
205  , _e2( (2 - Constants::WGS84_f()) * Constants::WGS84_f() )
206  , _degree( Math::degree() )
207  , _eps( sqrt(numeric_limits<real>::epsilon()) )
208  , _threadsafe(false) // Set after cache is read
209  {
210  GEOGRAPHICLIB_STATIC_ASSERT(sizeof(pixel_t) == pixel_size_,
211  "pixel_t has the wrong size");
212  if (_dir.empty())
213  _dir = DefaultGeoidPath();
214  _filename = _dir + "/" + _name + (pixel_size_ != 4 ? ".pgm" : ".pgm4");
215  _file.open(_filename.c_str(), ios::binary);
216  if (!(_file.good()))
217  throw GeographicErr("File not readable " + _filename);
218  string s;
219  if (!(getline(_file, s) && s == "P5"))
220  throw GeographicErr("File not in PGM format " + _filename);
221  _offset = numeric_limits<real>::max();
222  _scale = 0;
223  _maxerror = _rmserror = -1;
224  _description = "NONE";
225  _datetime = "UNKNOWN";
226  while (getline(_file, s)) {
227  if (s.empty())
228  continue;
229  if (s[0] == '#') {
230  istringstream is(s);
231  string commentid, key;
232  if (!(is >> commentid >> key) || commentid != "#")
233  continue;
234  if (key == "Description" || key =="DateTime") {
235  string::size_type p =
236  s.find_first_not_of(" \t", unsigned(is.tellg()));
237  if (p != string::npos)
238  (key == "Description" ? _description : _datetime) = s.substr(p);
239  } else if (key == "Offset") {
240  if (!(is >> _offset))
241  throw GeographicErr("Error reading offset " + _filename);
242  } else if (key == "Scale") {
243  if (!(is >> _scale))
244  throw GeographicErr("Error reading scale " + _filename);
245  } else if (key == (_cubic ? "MaxCubicError" : "MaxBilinearError")) {
246  // It's not an error if the error can't be read
247  is >> _maxerror;
248  } else if (key == (_cubic ? "RMSCubicError" : "RMSBilinearError")) {
249  // It's not an error if the error can't be read
250  is >> _rmserror;
251  }
252  } else {
253  istringstream is(s);
254  if (!(is >> _width >> _height))
255  throw GeographicErr("Error reading raster size " + _filename);
256  break;
257  }
258  }
259  {
260  unsigned maxval;
261  if (!(_file >> maxval))
262  throw GeographicErr("Error reading maxval " + _filename);
263  if (maxval != pixel_max_)
264  throw GeographicErr("Incorrect value of maxval " + _filename);
265  // Add 1 for whitespace after maxval
266  _datastart = (unsigned long long)(_file.tellg()) + 1ULL;
267  _swidth = (unsigned long long)(_width);
268  }
269  if (_offset == numeric_limits<real>::max())
270  throw GeographicErr("Offset not set " + _filename);
271  if (_scale == 0)
272  throw GeographicErr("Scale not set " + _filename);
273  if (_scale < 0)
274  throw GeographicErr("Scale must be positive " + _filename);
275  if (_height < 2 || _width < 2)
276  // Coarsest grid spacing is 180deg.
277  throw GeographicErr("Raster size too small " + _filename);
278  if (_width & 1)
279  // This is so that longitude grids can be extended thru the poles.
280  throw GeographicErr("Raster width is odd " + _filename);
281  if (!(_height & 1))
282  // This is so that latitude grid includes the equator.
283  throw GeographicErr("Raster height is even " + _filename);
284  _file.seekg(0, ios::end);
285  if (!_file.good() ||
286  _datastart + pixel_size_ * _swidth * (unsigned long long)(_height) !=
287  (unsigned long long)(_file.tellg()))
288  // Possibly this test should be "<" because the file contains, e.g., a
289  // second image. However, for now we are more strict.
290  throw GeographicErr("File has the wrong length " + _filename);
291  _rlonres = _width / real(360);
292  _rlatres = (_height - 1) / real(180);
293  _cache = false;
294  _ix = _width;
295  _iy = _height;
296  // Ensure that file errors throw exceptions
297  _file.exceptions(ifstream::eofbit | ifstream::failbit | ifstream::badbit);
298  if (threadsafe) {
299  CacheAll();
300  _file.close();
301  _threadsafe = true;
302  }
303  }
304 
305  Math::real Geoid::height(real lat, real lon, bool gradp,
306  real& gradn, real& grade) const {
307  if (Math::isnan(lat) || Math::isnan(lon)) {
308  if (gradp) gradn = grade = Math::NaN();
309  return Math::NaN();
310  }
311  lon = Math::AngNormalize(lon);
312  real
313  fx = lon * _rlonres,
314  fy = -lat * _rlatres;
315  int
316  ix = int(floor(fx)),
317  iy = min((_height - 1)/2 - 1, int(floor(fy)));
318  fx -= ix;
319  fy -= iy;
320  iy += (_height - 1)/2;
321  ix += ix < 0 ? _width : (ix >= _width ? -_width : 0);
322  real v00 = 0, v01 = 0, v10 = 0, v11 = 0;
323  real t[nterms_];
324 
325  if (_threadsafe || !(ix == _ix && iy == _iy)) {
326  if (!_cubic) {
327  v00 = rawval(ix , iy );
328  v01 = rawval(ix + 1, iy );
329  v10 = rawval(ix , iy + 1);
330  v11 = rawval(ix + 1, iy + 1);
331  } else {
332  real v[stencilsize_];
333  int k = 0;
334  v[k++] = rawval(ix , iy - 1);
335  v[k++] = rawval(ix + 1, iy - 1);
336  v[k++] = rawval(ix - 1, iy );
337  v[k++] = rawval(ix , iy );
338  v[k++] = rawval(ix + 1, iy );
339  v[k++] = rawval(ix + 2, iy );
340  v[k++] = rawval(ix - 1, iy + 1);
341  v[k++] = rawval(ix , iy + 1);
342  v[k++] = rawval(ix + 1, iy + 1);
343  v[k++] = rawval(ix + 2, iy + 1);
344  v[k++] = rawval(ix , iy + 2);
345  v[k++] = rawval(ix + 1, iy + 2);
346 
347  const int* c3x = iy == 0 ? c3n_ : (iy == _height - 2 ? c3s_ : c3_);
348  int c0x = iy == 0 ? c0n_ : (iy == _height - 2 ? c0s_ : c0_);
349  for (unsigned i = 0; i < nterms_; ++i) {
350  t[i] = 0;
351  for (unsigned j = 0; j < stencilsize_; ++j)
352  t[i] += v[j] * c3x[nterms_ * j + i];
353  t[i] /= c0x;
354  }
355  }
356  } else { // same cell; used cached coefficients
357  if (!_cubic) {
358  v00 = _v00;
359  v01 = _v01;
360  v10 = _v10;
361  v11 = _v11;
362  } else
363  copy(_t, _t + nterms_, t);
364  }
365  if (!_cubic) {
366  real
367  a = (1 - fx) * v00 + fx * v01,
368  b = (1 - fx) * v10 + fx * v11,
369  c = (1 - fy) * a + fy * b,
370  h = _offset + _scale * c;
371  if (gradp) {
372  real
373  phi = lat * _degree,
374  cosphi = cos(phi),
375  sinphi = sin(phi),
376  n = 1/sqrt(1 - _e2 * sinphi * sinphi);
377  gradn = ((1 - fx) * (v00 - v10) + fx * (v01 - v11)) *
378  _rlatres / (_degree * _a * (1 - _e2) * n * n * n);
379  grade = (cosphi > _eps ?
380  ((1 - fy) * (v01 - v00) + fy * (v11 - v10)) / cosphi :
381  (sinphi > 0 ? v11 - v10 : v01 - v00) *
382  _rlatres / _degree ) *
383  _rlonres / (_degree * _a * n);
384  gradn *= _scale;
385  grade *= _scale;
386  }
387  if (!_threadsafe) {
388  _ix = ix;
389  _iy = iy;
390  _v00 = v00;
391  _v01 = v01;
392  _v10 = v10;
393  _v11 = v11;
394  }
395  return h;
396  } else {
397  real h = t[0] + fx * (t[1] + fx * (t[3] + fx * t[6])) +
398  fy * (t[2] + fx * (t[4] + fx * t[7]) +
399  fy * (t[5] + fx * t[8] + fy * t[9]));
400  h = _offset + _scale * h;
401  if (gradp) {
402  // Avoid 0/0 at the poles by backing off 1/100 of a cell size
403  lat = min(lat, 90 - 1/(100 * _rlatres));
404  lat = max(lat, -90 + 1/(100 * _rlatres));
405  fy = (90 - lat) * _rlatres;
406  fy -= int(fy);
407  real
408  phi = lat * _degree,
409  cosphi = cos(phi),
410  sinphi = sin(phi),
411  n = 1/sqrt(1 - _e2 * sinphi * sinphi);
412  gradn = t[2] + fx * (t[4] + fx * t[7]) +
413  fy * (2 * t[5] + fx * 2 * t[8] + 3 * fy * t[9]);
414  grade = t[1] + fx * (2 * t[3] + fx * 3 * t[6]) +
415  fy * (t[4] + fx * 2 * t[7] + fy * t[8]);
416  gradn *= - _rlatres / (_degree * _a * (1 - _e2) * n * n * n) * _scale;
417  grade *= _rlonres / (_degree * _a * n * cosphi) * _scale;
418  }
419  if (!_threadsafe) {
420  _ix = ix;
421  _iy = iy;
422  copy(t, t + nterms_, _t);
423  }
424  return h;
425  }
426  }
427 
428  void Geoid::CacheClear() const {
429  if (!_threadsafe) {
430  _cache = false;
431  try {
432  _data.clear();
433  // Use swap to release memory back to system
434  vector< vector<pixel_t> >().swap(_data);
435  }
436  catch (const exception&) {
437  }
438  }
439  }
440 
441  void Geoid::CacheArea(real south, real west, real north, real east) const {
442  if (_threadsafe)
443  throw GeographicErr("Attempt to change cache of threadsafe Geoid");
444  if (south > north) {
445  CacheClear();
446  return;
447  }
448  west = Math::AngNormalize(west); // west in [-180, 180)
449  east = Math::AngNormalize(east);
450  if (east <= west)
451  east += 360; // east - west in (0, 360]
452  int
453  iw = int(floor(west * _rlonres)),
454  ie = int(floor(east * _rlonres)),
455  in = int(floor(-north * _rlatres)) + (_height - 1)/2,
456  is = int(floor(-south * _rlatres)) + (_height - 1)/2;
457  in = max(0, min(_height - 2, in));
458  is = max(0, min(_height - 2, is));
459  is += 1;
460  ie += 1;
461  if (_cubic) {
462  in -= 1;
463  is += 1;
464  iw -= 1;
465  ie += 1;
466  }
467  if (ie - iw >= _width - 1) {
468  // Include entire longitude range
469  iw = 0;
470  ie = _width - 1;
471  } else {
472  ie += iw < 0 ? _width : (iw >= _width ? -_width : 0);
473  iw += iw < 0 ? _width : (iw >= _width ? -_width : 0);
474  }
475  int oysize = int(_data.size());
476  _xsize = ie - iw + 1;
477  _ysize = is - in + 1;
478  _xoffset = iw;
479  _yoffset = in;
480 
481  try {
482  _data.resize(_ysize, vector<pixel_t>(_xsize));
483  for (int iy = min(oysize, _ysize); iy--;)
484  _data[iy].resize(_xsize);
485  }
486  catch (const bad_alloc&) {
487  CacheClear();
488  throw GeographicErr("Insufficient memory for caching " + _filename);
489  }
490 
491  try {
492  for (int iy = in; iy <= is; ++iy) {
493  int iy1 = iy, iw1 = iw;
494  if (iy < 0 || iy >= _height) {
495  // Allow points "beyond" the poles to support interpolation
496  iy1 = iy1 < 0 ? -iy1 : 2 * (_height - 1) - iy1;
497  iw1 += _width/2;
498  if (iw1 >= _width)
499  iw1 -= _width;
500  }
501  int xs1 = min(_width - iw1, _xsize);
502  filepos(iw1, iy1);
503  Utility::readarray<pixel_t, pixel_t, true>
504  (_file, &(_data[iy - in][0]), xs1);
505  if (xs1 < _xsize) {
506  // Wrap around longitude = 0
507  filepos(0, iy1);
508  Utility::readarray<pixel_t, pixel_t, true>
509  (_file, &(_data[iy - in][xs1]), _xsize - xs1);
510  }
511  }
512  _cache = true;
513  }
514  catch (const exception& e) {
515  CacheClear();
516  throw GeographicErr(string("Error filling cache ") + e.what());
517  }
518  }
519 
520  std::string Geoid::DefaultGeoidPath() {
521  string path;
522  char* geoidpath = getenv("GEOGRAPHICLIB_GEOID_PATH");
523  if (geoidpath)
524  path = string(geoidpath);
525  if (!path.empty())
526  return path;
527  char* datapath = getenv("GEOGRAPHICLIB_DATA");
528  if (datapath)
529  path = string(datapath);
530  return (!path.empty() ? path : string(GEOGRAPHICLIB_DATA)) + "/geoids";
531  }
532 
533  std::string Geoid::DefaultGeoidName() {
534  string name;
535  char* geoidname = getenv("GEOGRAPHICLIB_GEOID_NAME");
536  if (geoidname)
537  name = string(geoidname);
538  return !name.empty() ? name : string(GEOGRAPHICLIB_GEOID_DEFAULT_NAME);
539  }
540 
541 } // namespace GeographicLib
static T AngNormalize(T x)
Definition: Math.hpp:399
static T NaN()
Definition: Math.hpp:460
GeographicLib::Math::real real
Definition: GeodSolve.cpp:40
Header for GeographicLib::Utility class.
static bool isnan(T x)
Definition: Math.hpp:477
Mathematical functions needed by GeographicLib.
Definition: Math.hpp:101
void CacheClear() const
Definition: Geoid.cpp:428
static std::string DefaultGeoidPath()
Definition: Geoid.cpp:520
#define GEOGRAPHICLIB_GEOID_DEFAULT_NAME
Definition: Geoid.cpp:23
void CacheArea(real south, real west, real north, real east) const
Definition: Geoid.cpp:441
#define GEOGRAPHICLIB_DATA
Definition: Geoid.cpp:18
Namespace for GeographicLib.
Definition: Accumulator.cpp:12
Constants needed by GeographicLib
Definition: Constants.hpp:95
Exception handling for GeographicLib.
Definition: Constants.hpp:362
static std::string DefaultGeoidName()
Definition: Geoid.cpp:533
Header for GeographicLib::Geoid class.
void CacheAll() const
Definition: Geoid.hpp:262