src/TransverseMercator.cpp
Go to the documentation of this file.
1 
37 #include <iostream>
38 #include <complex>
40 
41 namespace GeographicLib {
42 
43  using namespace std;
44 
46  : _a(a)
47  , _f(f)
48  , _k0(k0)
49  , _e2(_f * (2 - _f))
50  , _es((f < 0 ? -1 : 1) * sqrt(abs(_e2)))
51  , _e2m(1 - _e2)
52  // _c = sqrt( pow(1 + _e, 1 + _e) * pow(1 - _e, 1 - _e) ) )
53  // See, for example, Lee (1976), p 100.
54  , _c( sqrt(_e2m) * exp(Math::eatanhe(real(1), _es)) )
55  , _n(_f / (2 - _f))
56  {
57  if (!(Math::isfinite(_a) && _a > 0))
58  throw GeographicErr("Equatorial radius is not positive");
59  if (!(Math::isfinite(_f) && _f < 1))
60  throw GeographicErr("Polar semi-axis is not positive");
61  if (!(Math::isfinite(_k0) && _k0 > 0))
62  throw GeographicErr("Scale is not positive");
63 
64  // Generated by Maxima on 2015-05-14 22:55:13-04:00
65 #if GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER/2 == 2
66  static const real b1coeff[] = {
67  // b1*(n+1), polynomial in n2 of order 2
68  1, 16, 64, 64,
69  }; // count = 4
70 #elif GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER/2 == 3
71  static const real b1coeff[] = {
72  // b1*(n+1), polynomial in n2 of order 3
73  1, 4, 64, 256, 256,
74  }; // count = 5
75 #elif GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER/2 == 4
76  static const real b1coeff[] = {
77  // b1*(n+1), polynomial in n2 of order 4
78  25, 64, 256, 4096, 16384, 16384,
79  }; // count = 6
80 #else
81 #error "Bad value for GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER"
82 #endif
83 
84 #if GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER == 4
85  static const real alpcoeff[] = {
86  // alp[1]/n^1, polynomial in n of order 3
87  164, 225, -480, 360, 720,
88  // alp[2]/n^2, polynomial in n of order 2
89  557, -864, 390, 1440,
90  // alp[3]/n^3, polynomial in n of order 1
91  -1236, 427, 1680,
92  // alp[4]/n^4, polynomial in n of order 0
93  49561, 161280,
94  }; // count = 14
95 #elif GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER == 5
96  static const real alpcoeff[] = {
97  // alp[1]/n^1, polynomial in n of order 4
98  -635, 328, 450, -960, 720, 1440,
99  // alp[2]/n^2, polynomial in n of order 3
100  4496, 3899, -6048, 2730, 10080,
101  // alp[3]/n^3, polynomial in n of order 2
102  15061, -19776, 6832, 26880,
103  // alp[4]/n^4, polynomial in n of order 1
104  -171840, 49561, 161280,
105  // alp[5]/n^5, polynomial in n of order 0
106  34729, 80640,
107  }; // count = 20
108 #elif GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER == 6
109  static const real alpcoeff[] = {
110  // alp[1]/n^1, polynomial in n of order 5
111  31564, -66675, 34440, 47250, -100800, 75600, 151200,
112  // alp[2]/n^2, polynomial in n of order 4
113  -1983433, 863232, 748608, -1161216, 524160, 1935360,
114  // alp[3]/n^3, polynomial in n of order 3
115  670412, 406647, -533952, 184464, 725760,
116  // alp[4]/n^4, polynomial in n of order 2
117  6601661, -7732800, 2230245, 7257600,
118  // alp[5]/n^5, polynomial in n of order 1
119  -13675556, 3438171, 7983360,
120  // alp[6]/n^6, polynomial in n of order 0
121  212378941, 319334400,
122  }; // count = 27
123 #elif GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER == 7
124  static const real alpcoeff[] = {
125  // alp[1]/n^1, polynomial in n of order 6
126  1804025, 2020096, -4267200, 2204160, 3024000, -6451200, 4838400, 9676800,
127  // alp[2]/n^2, polynomial in n of order 5
128  4626384, -9917165, 4316160, 3743040, -5806080, 2620800, 9676800,
129  // alp[3]/n^3, polynomial in n of order 4
130  -67102379, 26816480, 16265880, -21358080, 7378560, 29030400,
131  // alp[4]/n^4, polynomial in n of order 3
132  155912000, 72618271, -85060800, 24532695, 79833600,
133  // alp[5]/n^5, polynomial in n of order 2
134  102508609, -109404448, 27505368, 63866880,
135  // alp[6]/n^6, polynomial in n of order 1
136  -12282192400LL, 2760926233LL, 4151347200LL,
137  // alp[7]/n^7, polynomial in n of order 0
138  1522256789, 1383782400,
139  }; // count = 35
140 #elif GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER == 8
141  static const real alpcoeff[] = {
142  // alp[1]/n^1, polynomial in n of order 7
143  -75900428, 37884525, 42422016, -89611200, 46287360, 63504000, -135475200,
144  101606400, 203212800,
145  // alp[2]/n^2, polynomial in n of order 6
146  148003883, 83274912, -178508970, 77690880, 67374720, -104509440,
147  47174400, 174182400,
148  // alp[3]/n^3, polynomial in n of order 5
149  318729724, -738126169, 294981280, 178924680, -234938880, 81164160,
150  319334400,
151  // alp[4]/n^4, polynomial in n of order 4
152  -40176129013LL, 14967552000LL, 6971354016LL, -8165836800LL, 2355138720LL,
153  7664025600LL,
154  // alp[5]/n^5, polynomial in n of order 3
155  10421654396LL, 3997835751LL, -4266773472LL, 1072709352, 2490808320LL,
156  // alp[6]/n^6, polynomial in n of order 2
157  175214326799LL, -171950693600LL, 38652967262LL, 58118860800LL,
158  // alp[7]/n^7, polynomial in n of order 1
159  -67039739596LL, 13700311101LL, 12454041600LL,
160  // alp[8]/n^8, polynomial in n of order 0
161  1424729850961LL, 743921418240LL,
162  }; // count = 44
163 #else
164 #error "Bad value for GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER"
165 #endif
166 
167 #if GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER == 4
168  static const real betcoeff[] = {
169  // bet[1]/n^1, polynomial in n of order 3
170  -4, 555, -960, 720, 1440,
171  // bet[2]/n^2, polynomial in n of order 2
172  -437, 96, 30, 1440,
173  // bet[3]/n^3, polynomial in n of order 1
174  -148, 119, 3360,
175  // bet[4]/n^4, polynomial in n of order 0
176  4397, 161280,
177  }; // count = 14
178 #elif GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER == 5
179  static const real betcoeff[] = {
180  // bet[1]/n^1, polynomial in n of order 4
181  -3645, -64, 8880, -15360, 11520, 23040,
182  // bet[2]/n^2, polynomial in n of order 3
183  4416, -3059, 672, 210, 10080,
184  // bet[3]/n^3, polynomial in n of order 2
185  -627, -592, 476, 13440,
186  // bet[4]/n^4, polynomial in n of order 1
187  -3520, 4397, 161280,
188  // bet[5]/n^5, polynomial in n of order 0
189  4583, 161280,
190  }; // count = 20
191 #elif GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER == 6
192  static const real betcoeff[] = {
193  // bet[1]/n^1, polynomial in n of order 5
194  384796, -382725, -6720, 932400, -1612800, 1209600, 2419200,
195  // bet[2]/n^2, polynomial in n of order 4
196  -1118711, 1695744, -1174656, 258048, 80640, 3870720,
197  // bet[3]/n^3, polynomial in n of order 3
198  22276, -16929, -15984, 12852, 362880,
199  // bet[4]/n^4, polynomial in n of order 2
200  -830251, -158400, 197865, 7257600,
201  // bet[5]/n^5, polynomial in n of order 1
202  -435388, 453717, 15966720,
203  // bet[6]/n^6, polynomial in n of order 0
204  20648693, 638668800,
205  }; // count = 27
206 #elif GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER == 7
207  static const real betcoeff[] = {
208  // bet[1]/n^1, polynomial in n of order 6
209  -5406467, 6156736, -6123600, -107520, 14918400, -25804800, 19353600,
210  38707200,
211  // bet[2]/n^2, polynomial in n of order 5
212  829456, -5593555, 8478720, -5873280, 1290240, 403200, 19353600,
213  // bet[3]/n^3, polynomial in n of order 4
214  9261899, 3564160, -2708640, -2557440, 2056320, 58060800,
215  // bet[4]/n^4, polynomial in n of order 3
216  14928352, -9132761, -1742400, 2176515, 79833600,
217  // bet[5]/n^5, polynomial in n of order 2
218  -8005831, -1741552, 1814868, 63866880,
219  // bet[6]/n^6, polynomial in n of order 1
220  -261810608, 268433009, 8302694400LL,
221  // bet[7]/n^7, polynomial in n of order 0
222  219941297, 5535129600LL,
223  }; // count = 35
224 #elif GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER == 8
225  static const real betcoeff[] = {
226  // bet[1]/n^1, polynomial in n of order 7
227  31777436, -37845269, 43097152, -42865200, -752640, 104428800, -180633600,
228  135475200, 270950400,
229  // bet[2]/n^2, polynomial in n of order 6
230  24749483, 14930208, -100683990, 152616960, -105719040, 23224320, 7257600,
231  348364800,
232  // bet[3]/n^3, polynomial in n of order 5
233  -232468668, 101880889, 39205760, -29795040, -28131840, 22619520,
234  638668800,
235  // bet[4]/n^4, polynomial in n of order 4
236  324154477, 1433121792, -876745056, -167270400, 208945440, 7664025600LL,
237  // bet[5]/n^5, polynomial in n of order 3
238  457888660, -312227409, -67920528, 70779852, 2490808320LL,
239  // bet[6]/n^6, polynomial in n of order 2
240  -19841813847LL, -3665348512LL, 3758062126LL, 116237721600LL,
241  // bet[7]/n^7, polynomial in n of order 1
242  -1989295244, 1979471673, 49816166400LL,
243  // bet[8]/n^8, polynomial in n of order 0
244  191773887257LL, 3719607091200LL,
245  }; // count = 44
246 #else
247 #error "Bad value for GEOGRAPHICLIB_TRANSVERSEMERCATOR_ORDER"
248 #endif
249 
250  GEOGRAPHICLIB_STATIC_ASSERT(sizeof(b1coeff) / sizeof(real) ==
251  maxpow_/2 + 2,
252  "Coefficient array size mismatch for b1");
253  GEOGRAPHICLIB_STATIC_ASSERT(sizeof(alpcoeff) / sizeof(real) ==
254  (maxpow_ * (maxpow_ + 3))/2,
255  "Coefficient array size mismatch for alp");
256  GEOGRAPHICLIB_STATIC_ASSERT(sizeof(betcoeff) / sizeof(real) ==
257  (maxpow_ * (maxpow_ + 3))/2,
258  "Coefficient array size mismatch for bet");
259  int m = maxpow_/2;
260  _b1 = Math::polyval(m, b1coeff, Math::sq(_n)) / (b1coeff[m + 1] * (1+_n));
261  // _a1 is the equivalent radius for computing the circumference of
262  // ellipse.
263  _a1 = _b1 * _a;
264  int o = 0;
265  real d = _n;
266  for (int l = 1; l <= maxpow_; ++l) {
267  m = maxpow_ - l;
268  _alp[l] = d * Math::polyval(m, alpcoeff + o, _n) / alpcoeff[o + m + 1];
269  _bet[l] = d * Math::polyval(m, betcoeff + o, _n) / betcoeff[o + m + 1];
270  o += m + 2;
271  d *= _n;
272  }
273  // Post condition: o == sizeof(alpcoeff) / sizeof(real) &&
274  // o == sizeof(betcoeff) / sizeof(real)
275  }
276 
278  static const TransverseMercator utm(Constants::WGS84_a(),
281  return utm;
282  }
283 
284  // Engsager and Poder (2007) use trigonometric series to convert between phi
285  // and phip. Here are the series...
286  //
287  // Conversion from phi to phip:
288  //
289  // phip = phi + sum(c[j] * sin(2*j*phi), j, 1, 6)
290  //
291  // c[1] = - 2 * n
292  // + 2/3 * n^2
293  // + 4/3 * n^3
294  // - 82/45 * n^4
295  // + 32/45 * n^5
296  // + 4642/4725 * n^6;
297  // c[2] = 5/3 * n^2
298  // - 16/15 * n^3
299  // - 13/9 * n^4
300  // + 904/315 * n^5
301  // - 1522/945 * n^6;
302  // c[3] = - 26/15 * n^3
303  // + 34/21 * n^4
304  // + 8/5 * n^5
305  // - 12686/2835 * n^6;
306  // c[4] = 1237/630 * n^4
307  // - 12/5 * n^5
308  // - 24832/14175 * n^6;
309  // c[5] = - 734/315 * n^5
310  // + 109598/31185 * n^6;
311  // c[6] = 444337/155925 * n^6;
312  //
313  // Conversion from phip to phi:
314  //
315  // phi = phip + sum(d[j] * sin(2*j*phip), j, 1, 6)
316  //
317  // d[1] = 2 * n
318  // - 2/3 * n^2
319  // - 2 * n^3
320  // + 116/45 * n^4
321  // + 26/45 * n^5
322  // - 2854/675 * n^6;
323  // d[2] = 7/3 * n^2
324  // - 8/5 * n^3
325  // - 227/45 * n^4
326  // + 2704/315 * n^5
327  // + 2323/945 * n^6;
328  // d[3] = 56/15 * n^3
329  // - 136/35 * n^4
330  // - 1262/105 * n^5
331  // + 73814/2835 * n^6;
332  // d[4] = 4279/630 * n^4
333  // - 332/35 * n^5
334  // - 399572/14175 * n^6;
335  // d[5] = 4174/315 * n^5
336  // - 144838/6237 * n^6;
337  // d[6] = 601676/22275 * n^6;
338  //
339  // In order to maintain sufficient relative accuracy close to the pole use
340  //
341  // S = sum(c[i]*sin(2*i*phi),i,1,6)
342  // taup = (tau + tan(S)) / (1 - tau * tan(S))
343 
344  // In Math::taupf and Math::tauf we evaluate the forward transform explicitly
345  // and solve the reverse one by Newton's method.
346  //
347  // There are adapted from TransverseMercatorExact (taup and taupinv). tau =
348  // tan(phi), taup = sinh(psi)
349 
351  real& x, real& y,
352  real& gamma, real& k) const {
353  lat = Math::LatFix(lat);
354  lon = Math::AngDiff(lon0, lon);
355  // Explicitly enforce the parity
356  int
357  latsign = (lat < 0) ? -1 : 1,
358  lonsign = (lon < 0) ? -1 : 1;
359  lon *= lonsign;
360  lat *= latsign;
361  bool backside = lon > 90;
362  if (backside) {
363  if (lat == 0)
364  latsign = -1;
365  lon = 180 - lon;
366  }
367  real sphi, cphi, slam, clam;
368  Math::sincosd(lat, sphi, cphi);
369  Math::sincosd(lon, slam, clam);
370  // phi = latitude
371  // phi' = conformal latitude
372  // psi = isometric latitude
373  // tau = tan(phi)
374  // tau' = tan(phi')
375  // [xi', eta'] = Gauss-Schreiber TM coordinates
376  // [xi, eta] = Gauss-Krueger TM coordinates
377  //
378  // We use
379  // tan(phi') = sinh(psi)
380  // sin(phi') = tanh(psi)
381  // cos(phi') = sech(psi)
382  // denom^2 = 1-cos(phi')^2*sin(lam)^2 = 1-sech(psi)^2*sin(lam)^2
383  // sin(xip) = sin(phi')/denom = tanh(psi)/denom
384  // cos(xip) = cos(phi')*cos(lam)/denom = sech(psi)*cos(lam)/denom
385  // cosh(etap) = 1/denom = 1/denom
386  // sinh(etap) = cos(phi')*sin(lam)/denom = sech(psi)*sin(lam)/denom
387  real etap, xip;
388  if (lat != 90) {
389  real
390  tau = sphi / cphi,
391  taup = Math::taupf(tau, _es);
392  xip = atan2(taup, clam);
393  // Used to be
394  // etap = Math::atanh(sin(lam) / cosh(psi));
395  etap = Math::asinh(slam / Math::hypot(taup, clam));
396  // convergence and scale for Gauss-Schreiber TM (xip, etap) -- gamma0 =
397  // atan(tan(xip) * tanh(etap)) = atan(tan(lam) * sin(phi'));
398  // sin(phi') = tau'/sqrt(1 + tau'^2)
399  // Krueger p 22 (44)
400  gamma = Math::atan2d(slam * taup, clam * Math::hypot(real(1), taup));
401  // k0 = sqrt(1 - _e2 * sin(phi)^2) * (cos(phi') / cos(phi)) * cosh(etap)
402  // Note 1/cos(phi) = cosh(psip);
403  // and cos(phi') * cosh(etap) = 1/hypot(sinh(psi), cos(lam))
404  //
405  // This form has cancelling errors. This property is lost if cosh(psip)
406  // is replaced by 1/cos(phi), even though it's using "primary" data (phi
407  // instead of psip).
408  k = sqrt(_e2m + _e2 * Math::sq(cphi)) * Math::hypot(real(1), tau)
409  / Math::hypot(taup, clam);
410  } else {
411  xip = Math::pi()/2;
412  etap = 0;
413  gamma = lon;
414  k = _c;
415  }
416  // {xi',eta'} is {northing,easting} for Gauss-Schreiber transverse Mercator
417  // (for eta' = 0, xi' = bet). {xi,eta} is {northing,easting} for transverse
418  // Mercator with constant scale on the central meridian (for eta = 0, xip =
419  // rectifying latitude). Define
420  //
421  // zeta = xi + i*eta
422  // zeta' = xi' + i*eta'
423  //
424  // The conversion from conformal to rectifying latitude can be expressed as
425  // a series in _n:
426  //
427  // zeta = zeta' + sum(h[j-1]' * sin(2 * j * zeta'), j = 1..maxpow_)
428  //
429  // where h[j]' = O(_n^j). The reversion of this series gives
430  //
431  // zeta' = zeta - sum(h[j-1] * sin(2 * j * zeta), j = 1..maxpow_)
432  //
433  // which is used in Reverse.
434  //
435  // Evaluate sums via Clenshaw method. See
436  // https://en.wikipedia.org/wiki/Clenshaw_algorithm
437  //
438  // Let
439  //
440  // S = sum(a[k] * phi[k](x), k = 0..n)
441  // phi[k+1](x) = alpha[k](x) * phi[k](x) + beta[k](x) * phi[k-1](x)
442  //
443  // Evaluate S with
444  //
445  // b[n+2] = b[n+1] = 0
446  // b[k] = alpha[k](x) * b[k+1] + beta[k+1](x) * b[k+2] + a[k]
447  // S = (a[0] + beta[1](x) * b[2]) * phi[0](x) + b[1] * phi[1](x)
448  //
449  // Here we have
450  //
451  // x = 2 * zeta'
452  // phi[k](x) = sin(k * x)
453  // alpha[k](x) = 2 * cos(x)
454  // beta[k](x) = -1
455  // [ sin(A+B) - 2*cos(B)*sin(A) + sin(A-B) = 0, A = k*x, B = x ]
456  // n = maxpow_
457  // a[k] = _alp[k]
458  // S = b[1] * sin(x)
459  //
460  // For the derivative we have
461  //
462  // x = 2 * zeta'
463  // phi[k](x) = cos(k * x)
464  // alpha[k](x) = 2 * cos(x)
465  // beta[k](x) = -1
466  // [ cos(A+B) - 2*cos(B)*cos(A) + cos(A-B) = 0, A = k*x, B = x ]
467  // a[0] = 1; a[k] = 2*k*_alp[k]
468  // S = (a[0] - b[2]) + b[1] * cos(x)
469  //
470  // Matrix formulation (not used here):
471  // phi[k](x) = [sin(k * x); k * cos(k * x)]
472  // alpha[k](x) = 2 * [cos(x), 0; -sin(x), cos(x)]
473  // beta[k](x) = -1 * [1, 0; 0, 1]
474  // a[k] = _alp[k] * [1, 0; 0, 1]
475  // b[n+2] = b[n+1] = [0, 0; 0, 0]
476  // b[k] = alpha[k](x) * b[k+1] + beta[k+1](x) * b[k+2] + a[k]
477  // N.B., for all k: b[k](1,2) = 0; b[k](1,1) = b[k](2,2)
478  // S = (a[0] + beta[1](x) * b[2]) * phi[0](x) + b[1] * phi[1](x)
479  // phi[0](x) = [0; 0]
480  // phi[1](x) = [sin(x); cos(x)]
481  real
482  c0 = cos(2 * xip), ch0 = cosh(2 * etap),
483  s0 = sin(2 * xip), sh0 = sinh(2 * etap);
484  complex<real> a(2 * c0 * ch0, -2 * s0 * sh0); // 2 * cos(2*zeta')
485  int n = maxpow_;
487  y0(n & 1 ? _alp[n] : 0), y1, // default initializer is 0+i0
488  z0(n & 1 ? 2*n * _alp[n] : 0), z1;
489  if (n & 1) --n;
490  while (n) {
491  y1 = a * y0 - y1 + _alp[n];
492  z1 = a * z0 - z1 + 2*n * _alp[n];
493  --n;
494  y0 = a * y1 - y0 + _alp[n];
495  z0 = a * z1 - z0 + 2*n * _alp[n];
496  --n;
497  }
498  a /= real(2); // cos(2*zeta')
499  z1 = real(1) - z1 + a * z0;
500  a = complex<real>(s0 * ch0, c0 * sh0); // sin(2*zeta')
501  y1 = complex<real>(xip, etap) + a * y0;
502  // Fold in change in convergence and scale for Gauss-Schreiber TM to
503  // Gauss-Krueger TM.
504  gamma -= Math::atan2d(z1.imag(), z1.real());
505  k *= _b1 * abs(z1);
506  real xi = y1.real(), eta = y1.imag();
507  y = _a1 * _k0 * (backside ? Math::pi() - xi : xi) * latsign;
508  x = _a1 * _k0 * eta * lonsign;
509  if (backside)
510  gamma = 180 - gamma;
511  gamma *= latsign * lonsign;
512  gamma = Math::AngNormalize(gamma);
513  k *= _k0;
514  }
515 
517  real& lat, real& lon,
518  real& gamma, real& k) const {
519  // This undoes the steps in Forward. The wrinkles are: (1) Use of the
520  // reverted series to express zeta' in terms of zeta. (2) Newton's method
521  // to solve for phi in terms of tan(phi).
522  real
523  xi = y / (_a1 * _k0),
524  eta = x / (_a1 * _k0);
525  // Explicitly enforce the parity
526  int
527  xisign = (xi < 0) ? -1 : 1,
528  etasign = (eta < 0) ? -1 : 1;
529  xi *= xisign;
530  eta *= etasign;
531  bool backside = xi > Math::pi()/2;
532  if (backside)
533  xi = Math::pi() - xi;
534  real
535  c0 = cos(2 * xi), ch0 = cosh(2 * eta),
536  s0 = sin(2 * xi), sh0 = sinh(2 * eta);
537  complex<real> a(2 * c0 * ch0, -2 * s0 * sh0); // 2 * cos(2*zeta)
538  int n = maxpow_;
540  y0(n & 1 ? -_bet[n] : 0), y1, // default initializer is 0+i0
541  z0(n & 1 ? -2*n * _bet[n] : 0), z1;
542  if (n & 1) --n;
543  while (n) {
544  y1 = a * y0 - y1 - _bet[n];
545  z1 = a * z0 - z1 - 2*n * _bet[n];
546  --n;
547  y0 = a * y1 - y0 - _bet[n];
548  z0 = a * z1 - z0 - 2*n * _bet[n];
549  --n;
550  }
551  a /= real(2); // cos(2*zeta)
552  z1 = real(1) - z1 + a * z0;
553  a = complex<real>(s0 * ch0, c0 * sh0); // sin(2*zeta)
554  y1 = complex<real>(xi, eta) + a * y0;
555  // Convergence and scale for Gauss-Schreiber TM to Gauss-Krueger TM.
556  gamma = Math::atan2d(z1.imag(), z1.real());
557  k = _b1 / abs(z1);
558  // JHS 154 has
559  //
560  // phi' = asin(sin(xi') / cosh(eta')) (Krueger p 17 (25))
561  // lam = asin(tanh(eta') / cos(phi')
562  // psi = asinh(tan(phi'))
563  real
564  xip = y1.real(), etap = y1.imag(),
565  s = sinh(etap),
566  c = max(real(0), cos(xip)), // cos(pi/2) might be negative
567  r = Math::hypot(s, c);
568  if (r != 0) {
569  lon = Math::atan2d(s, c); // Krueger p 17 (25)
570  // Use Newton's method to solve for tau
571  real
572  sxip = sin(xip),
573  tau = Math::tauf(sxip/r, _es);
574  gamma += Math::atan2d(sxip * tanh(etap), c); // Krueger p 19 (31)
575  lat = Math::atand(tau);
576  // Note cos(phi') * cosh(eta') = r
577  k *= sqrt(_e2m + _e2 / (1 + Math::sq(tau))) *
578  Math::hypot(real(1), tau) * r;
579  } else {
580  lat = 90;
581  lon = 0;
582  k *= _c;
583  }
584  lat *= xisign;
585  if (backside)
586  lon = 180 - lon;
587  lon *= etasign;
588  lon = Math::AngNormalize(lon + lon0);
589  if (backside)
590  gamma = 180 - gamma;
591  gamma *= xisign * etasign;
592  gamma = Math::AngNormalize(gamma);
593  k *= _k0;
594  }
595 
596 } // namespace GeographicLib
static T AngNormalize(T x)
Definition: Math.hpp:440
Matrix3f m
#define max(a, b)
Definition: datatypes.h:20
static T pi()
Definition: Math.hpp:202
Scalar * y
static T atand(T x)
Definition: Math.hpp:723
static const double lat
EIGEN_DEVICE_FUNC const ExpReturnType exp() const
static bool isfinite(T x)
Definition: Math.hpp:806
int n
static T LatFix(T x)
Definition: Math.hpp:467
Scalar Scalar * c
Definition: benchVecAdd.cpp:17
EIGEN_DEVICE_FUNC const TanhReturnType tanh() const
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
EIGEN_DEVICE_FUNC const CoshReturnType cosh() const
Mathematical functions needed by GeographicLib.
Definition: Math.hpp:102
Definition: Half.h:150
static void sincosd(T x, T &sinx, T &cosx)
Definition: Math.hpp:558
void Forward(real lon0, real lat, real lon, real &x, real &y, real &gamma, real &k) const
static T AngDiff(T x, T y, T &e)
Definition: Math.hpp:486
Transverse Mercator projection.
static const TransverseMercator & UTM()
Header for GeographicLib::TransverseMercator class.
Array33i a
static T asinh(T x)
Definition: Math.hpp:311
EIGEN_DEVICE_FUNC const CosReturnType cos() const
TransverseMercator(real a, real f, real k0)
static const Line3 l(Rot3(), 1, 1)
static T hypot(T x, T y)
Definition: Math.hpp:243
const mpreal gamma(const mpreal &x, mp_rnd_t r=mpreal::get_default_rnd())
Definition: mpreal.h:2262
EIGEN_DEVICE_FUNC const SinhReturnType sinh() const
static T sq(T x)
Definition: Math.hpp:232
static T atan2d(T y, T x)
Definition: Math.hpp:691
static T polyval(int N, const T p[], T x)
Definition: Math.hpp:425
Point2(* f)(const Point3 &, OptionalJacobian< 2, 3 >)
Namespace for GeographicLib.
RealScalar s
const double lon0
void Reverse(real lon0, real x, real y, real &lat, real &lon, real &gamma, real &k) const
Jet< T, N > atan2(const Jet< T, N > &g, const Jet< T, N > &f)
Definition: jet.h:556
Vector xi
Definition: testPose2.cpp:150
Point2 z1
Exception handling for GeographicLib.
Definition: Constants.hpp:389
static const double lon
EIGEN_DEVICE_FUNC const SinReturnType sin() const
static T tauf(T taup, T es)
Definition: Math.cpp:31
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
#define abs(x)
Definition: datatypes.h:17
static T taupf(T tau, T es)
Definition: Math.cpp:25


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:51:10