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