/usr/share/doc/geographiclib/html/TransverseMercator_8cpp_source.html is in geographiclib-tools 1.8-2.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 | <!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml">
<head>
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
<title>GeographicLib: TransverseMercator.cpp Source File</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<link href="doxygen.css" rel="stylesheet" type="text/css"/>
</head>
<body>
<!-- Generated by Doxygen 1.7.1 -->
<div class="navigation" id="top">
<div class="tabs">
<ul class="tablist">
<li><a href="index.html"><span>Main Page</span></a></li>
<li><a href="pages.html"><span>Related Pages</span></a></li>
<li><a href="namespaces.html"><span>Namespaces</span></a></li>
<li><a href="annotated.html"><span>Classes</span></a></li>
<li class="current"><a href="files.html"><span>Files</span></a></li>
<li><a href="dirs.html"><span>Directories</span></a></li>
</ul>
</div>
<div class="tabs2">
<ul class="tablist">
<li><a href="files.html"><span>File List</span></a></li>
<li><a href="globals.html"><span>File Members</span></a></li>
</ul>
</div>
<div class="navpath">
<ul>
<li><a class="el" href="dir_b6dff8fdfd5d70bd61ba562840823df9.html">src</a> </li>
</ul>
</div>
</div>
<div class="header">
<div class="headertitle">
<h1>TransverseMercator.cpp</h1> </div>
</div>
<div class="contents">
<a href="TransverseMercator_8cpp.html">Go to the documentation of this file.</a><div class="fragment"><pre class="fragment"><a name="l00001"></a>00001 <span class="comment">/**</span>
<a name="l00002"></a>00002 <span class="comment"> * \file TransverseMercator.cpp</span>
<a name="l00003"></a>00003 <span class="comment"> * \brief Implementation for GeographicLib::TransverseMercator class</span>
<a name="l00004"></a>00004 <span class="comment"> *</span>
<a name="l00005"></a>00005 <span class="comment"> * Copyright (c) Charles Karney (2008, 2009, 2010, 2011) <charles@karney.com></span>
<a name="l00006"></a>00006 <span class="comment"> * and licensed under the LGPL. For more information, see</span>
<a name="l00007"></a>00007 <span class="comment"> * http://geographiclib.sourceforge.net/</span>
<a name="l00008"></a>00008 <span class="comment"> *</span>
<a name="l00009"></a>00009 <span class="comment"> * This implementation follows closely</span>
<a name="l00010"></a>00010 <span class="comment"> * <a href="http://www.jhs-suositukset.fi/suomi/jhs154"> JHS 154, ETRS89 -</span>
<a name="l00011"></a>00011 <span class="comment"> * j&auml;rjestelm&auml;&auml;n liittyv&auml;t karttaprojektiot,</span>
<a name="l00012"></a>00012 <span class="comment"> * tasokoordinaatistot ja karttalehtijako</a> (Map projections, plane</span>
<a name="l00013"></a>00013 <span class="comment"> * coordinates, and map sheet index for ETRS89), published by JUHTA, Finnish</span>
<a name="l00014"></a>00014 <span class="comment"> * Geodetic Institute, and the National Land Survey of Finland (2006).</span>
<a name="l00015"></a>00015 <span class="comment"> *</span>
<a name="l00016"></a>00016 <span class="comment"> * The relevant section is available as the 2008 PDF file</span>
<a name="l00017"></a>00017 <span class="comment"> * http://docs.jhs-suositukset.fi/jhs-suositukset/JHS154/JHS154_liite1.pdf</span>
<a name="l00018"></a>00018 <span class="comment"> *</span>
<a name="l00019"></a>00019 <span class="comment"> * This is a straight transcription of the formulas in this paper with the</span>
<a name="l00020"></a>00020 <span class="comment"> * following exceptions:</span>
<a name="l00021"></a>00021 <span class="comment"> * - use of 6th order series instead of 4th order series. This reduces the</span>
<a name="l00022"></a>00022 <span class="comment"> * error to about 5nm for the UTM range of coordinates (instead of 200nm),</span>
<a name="l00023"></a>00023 <span class="comment"> * with a speed penalty of only 1%;</span>
<a name="l00024"></a>00024 <span class="comment"> * - use Newton's method instead of plain iteration to solve for latitude in</span>
<a name="l00025"></a>00025 <span class="comment"> * terms of isometric latitude in the Reverse method;</span>
<a name="l00026"></a>00026 <span class="comment"> * - use of Horner's representation for evaluating polynomials and Clenshaw's</span>
<a name="l00027"></a>00027 <span class="comment"> * method for summing trigonometric series;</span>
<a name="l00028"></a>00028 <span class="comment"> * - several modifications of the formulas to improve the numerical accuracy;</span>
<a name="l00029"></a>00029 <span class="comment"> * - evaluating the convergence and scale using the expression for the</span>
<a name="l00030"></a>00030 <span class="comment"> * projection or its inverse.</span>
<a name="l00031"></a>00031 <span class="comment"> *</span>
<a name="l00032"></a>00032 <span class="comment"> * If the preprocessor variable TM_TX_MAXPOW is set to an integer between 4 and</span>
<a name="l00033"></a>00033 <span class="comment"> * 8, then this specifies the order of the series used for the forward and</span>
<a name="l00034"></a>00034 <span class="comment"> * reverse transformations. The default value is 6. (The series accurate to</span>
<a name="l00035"></a>00035 <span class="comment"> * 12th order is given in \ref tmseries.)</span>
<a name="l00036"></a>00036 <span class="comment"> *</span>
<a name="l00037"></a>00037 <span class="comment"> * Other equivalent implementations are given in</span>
<a name="l00038"></a>00038 <span class="comment"> * - http://www.ign.fr/DISPLAY/000/526/702/5267021/NTG_76.pdf</span>
<a name="l00039"></a>00039 <span class="comment"> * - http://www.lantmateriet.se/upload/filer/kartor/geodesi_gps_och_detaljmatning/geodesi/Formelsamling/Gauss_Conformal_Projection.pdf</span>
<a name="l00040"></a>00040 <span class="comment"> **********************************************************************/</span>
<a name="l00041"></a>00041
<a name="l00042"></a>00042 <span class="preprocessor">#include "<a class="code" href="TransverseMercator_8hpp.html" title="Header for GeographicLib::TransverseMercator class.">GeographicLib/TransverseMercator.hpp</a>"</span>
<a name="l00043"></a>00043
<a name="l00044"></a><a class="code" href="TransverseMercator_8cpp.html#ab09ba83d5531a4ca91dc8ce82ae98084">00044</a> <span class="preprocessor">#define GEOGRAPHICLIB_TRANSVERSEMERCATOR_CPP "$Id: TransverseMercator.cpp 6937 2011-02-01 20:17:13Z karney $"</span>
<a name="l00045"></a>00045 <span class="preprocessor"></span>
<a name="l00046"></a>00046 <a class="code" href="Constants_8hpp.html#af90fa899707a2ac513d5e4c76853bbf5">RCSID_DECL</a>(<a class="code" href="TransverseMercator_8cpp.html#ab09ba83d5531a4ca91dc8ce82ae98084">GEOGRAPHICLIB_TRANSVERSEMERCATOR_CPP</a>)
<a name="l00047"></a>00047 <a class="code" href="Constants_8hpp.html#af90fa899707a2ac513d5e4c76853bbf5">RCSID_DECL</a>(<a class="code" href="TransverseMercator_8hpp.html#abbbe2faaff284db3183ddf3afcbf5623">GEOGRAPHICLIB_TRANSVERSEMERCATOR_HPP</a>)
<a name="l00048"></a>00048
<a name="l00049"></a>00049 namespace GeographicLib {
<a name="l00050"></a>00050
<a name="l00051"></a>00051 <span class="keyword">using namespace </span>std;
<a name="l00052"></a>00052
<a name="l00053"></a>00053 <span class="keyword">const</span> <a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">Math::real</a> TransverseMercator::tol =
<a name="l00054"></a>00054 <a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(0.1)*sqrt(numeric_limits<real>::epsilon());
<a name="l00055"></a>00055 <span class="comment">// Overflow value s.t. atan(overflow) = pi/2</span>
<a name="l00056"></a>00056 <span class="keyword">const</span> <a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">Math::real</a> TransverseMercator::overflow =
<a name="l00057"></a>00057 1 / sq(numeric_limits<real>::epsilon());
<a name="l00058"></a>00058
<a name="l00059"></a><a class="code" href="classGeographicLib_1_1TransverseMercator.html#a1322ba37e8f01d724f9da1dfd94d91d7">00059</a> TransverseMercator::TransverseMercator(real a, real r, real k0)
<a name="l00060"></a>00060 : _a(a)
<a name="l00061"></a>00061 , _r(r)
<a name="l00062"></a>00062 , _f(_r != 0 ? 1 / _r : 0)
<a name="l00063"></a>00063 , _k0(k0)
<a name="l00064"></a>00064 , _e2(_f * (2 - _f))
<a name="l00065"></a>00065 , _e(sqrt(abs(_e2)))
<a name="l00066"></a>00066 , _e2m(1 - _e2)
<a name="l00067"></a>00067 <span class="comment">// _c = sqrt( pow(1 + _e, 1 + _e) * pow(1 - _e, 1 - _e) ) )</span>
<a name="l00068"></a>00068 <span class="comment">// See, for example, Lee (1976), p 100.</span>
<a name="l00069"></a>00069 , _c( sqrt(_e2m) * exp(eatanhe(real(1))) )
<a name="l00070"></a>00070 , _n(_f / (2 - _f))
<a name="l00071"></a>00071 {
<a name="l00072"></a>00072 <span class="keywordflow">if</span> (!(_a > 0))
<a name="l00073"></a>00073 <span class="keywordflow">throw</span> <a class="code" href="classGeographicLib_1_1GeographicErr.html" title="Exception handling for GeographicLib">GeographicErr</a>(<span class="stringliteral">"Major radius is not positive"</span>);
<a name="l00074"></a>00074 <span class="keywordflow">if</span> (!(_f < 1))
<a name="l00075"></a>00075 <span class="keywordflow">throw</span> <a class="code" href="classGeographicLib_1_1GeographicErr.html" title="Exception handling for GeographicLib">GeographicErr</a>(<span class="stringliteral">"Minor radius is not positive"</span>);
<a name="l00076"></a>00076 <span class="keywordflow">if</span> (!(_k0 > 0))
<a name="l00077"></a>00077 <span class="keywordflow">throw</span> <a class="code" href="classGeographicLib_1_1GeographicErr.html" title="Exception handling for GeographicLib">GeographicErr</a>(<span class="stringliteral">"Scale is not positive"</span>);
<a name="l00078"></a>00078 <span class="comment">// If coefficents might overflow an int, convert them to double (and they</span>
<a name="l00079"></a>00079 <span class="comment">// are all exactly representable as doubles).</span>
<a name="l00080"></a>00080 real nx = sq(_n);
<a name="l00081"></a>00081 <span class="keywordflow">switch</span> (maxpow) {
<a name="l00082"></a>00082 <span class="keywordflow">case</span> 4:
<a name="l00083"></a>00083 _b1 = 1/(1+_n)*(nx*(nx+16)+64)/64;
<a name="l00084"></a>00084 _alp[1] = _n*(_n*(_n*(164*_n+225)-480)+360)/720;
<a name="l00085"></a>00085 _bet[1] = _n*(_n*((555-4*_n)*_n-960)+720)/1440;
<a name="l00086"></a>00086 _alp[2] = nx*(_n*(557*_n-864)+390)/1440;
<a name="l00087"></a>00087 _bet[2] = nx*((96-437*_n)*_n+30)/1440;
<a name="l00088"></a>00088 nx *= _n;
<a name="l00089"></a>00089 _alp[3] = (427-1236*_n)*nx/1680;
<a name="l00090"></a>00090 _bet[3] = (119-148*_n)*nx/3360;
<a name="l00091"></a>00091 nx *= _n;
<a name="l00092"></a>00092 _alp[4] = 49561*nx/161280;
<a name="l00093"></a>00093 _bet[4] = 4397*nx/161280;
<a name="l00094"></a>00094 <span class="keywordflow">break</span>;
<a name="l00095"></a>00095 <span class="keywordflow">case</span> 5:
<a name="l00096"></a>00096 _b1 = 1/(1+_n)*(nx*(nx+16)+64)/64;
<a name="l00097"></a>00097 _alp[1] = _n*(_n*(_n*((328-635*_n)*_n+450)-960)+720)/1440;
<a name="l00098"></a>00098 _bet[1] = _n*(_n*(_n*((-3645*_n-64)*_n+8880)-15360)+11520)/23040;
<a name="l00099"></a>00099 _alp[2] = nx*(_n*(_n*(4496*_n+3899)-6048)+2730)/10080;
<a name="l00100"></a>00100 _bet[2] = nx*(_n*(_n*(4416*_n-3059)+672)+210)/10080;
<a name="l00101"></a>00101 nx *= _n;
<a name="l00102"></a>00102 _alp[3] = nx*(_n*(15061*_n-19776)+6832)/26880;
<a name="l00103"></a>00103 _bet[3] = nx*((-627*_n-592)*_n+476)/13440;
<a name="l00104"></a>00104 nx *= _n;
<a name="l00105"></a>00105 _alp[4] = (49561-171840*_n)*nx/161280;
<a name="l00106"></a>00106 _bet[4] = (4397-3520*_n)*nx/161280;
<a name="l00107"></a>00107 nx *= _n;
<a name="l00108"></a>00108 _alp[5] = 34729*nx/80640;
<a name="l00109"></a>00109 _bet[5] = 4583*nx/161280;
<a name="l00110"></a>00110 <span class="keywordflow">break</span>;
<a name="l00111"></a>00111 <span class="keywordflow">case</span> 6:
<a name="l00112"></a>00112 _b1 = 1/(1+_n)*(nx*(nx*(nx+4)+64)+256)/256;
<a name="l00113"></a>00113 _alp[1] = _n*(_n*(_n*(_n*(_n*(31564*_n-66675)+34440)+47250)-100800)+
<a name="l00114"></a>00114 75600)/151200;
<a name="l00115"></a>00115 _bet[1] = _n*(_n*(_n*(_n*(_n*(384796*_n-382725)-6720)+932400)-1612800)+
<a name="l00116"></a>00116 1209600)/2419200;
<a name="l00117"></a>00117 _alp[2] = nx*(_n*(_n*((863232-1983433*_n)*_n+748608)-1161216)+524160)/
<a name="l00118"></a>00118 1935360;
<a name="l00119"></a>00119 _bet[2] = nx*(_n*(_n*((1695744-1118711*_n)*_n-1174656)+258048)+80640)/
<a name="l00120"></a>00120 3870720;
<a name="l00121"></a>00121 nx *= _n;
<a name="l00122"></a>00122 _alp[3] = nx*(_n*(_n*(670412*_n+406647)-533952)+184464)/725760;
<a name="l00123"></a>00123 _bet[3] = nx*(_n*(_n*(22276*_n-16929)-15984)+12852)/362880;
<a name="l00124"></a>00124 nx *= _n;
<a name="l00125"></a>00125 _alp[4] = nx*(_n*(6601661*_n-7732800)+2230245)/7257600;
<a name="l00126"></a>00126 _bet[4] = nx*((-830251*_n-158400)*_n+197865)/7257600;
<a name="l00127"></a>00127 nx *= _n;
<a name="l00128"></a>00128 _alp[5] = (3438171-13675556*_n)*nx/7983360;
<a name="l00129"></a>00129 _bet[5] = (453717-435388*_n)*nx/15966720;
<a name="l00130"></a>00130 nx *= _n;
<a name="l00131"></a>00131 _alp[6] = 212378941*nx/319334400;
<a name="l00132"></a>00132 _bet[6] = 20648693*nx/638668800;
<a name="l00133"></a>00133 <span class="keywordflow">break</span>;
<a name="l00134"></a>00134 <span class="keywordflow">case</span> 7:
<a name="l00135"></a>00135 _b1 = 1/(1+_n)*(nx*(nx*(nx+4)+64)+256)/256;
<a name="l00136"></a>00136 _alp[1] = _n*(_n*(_n*(_n*(_n*(_n*(1804025*_n+2020096)-4267200)+2204160)+
<a name="l00137"></a>00137 3024000)-6451200)+4838400)/9676800;
<a name="l00138"></a>00138 _bet[1] = _n*(_n*(_n*(_n*(_n*((6156736-5406467*_n)*_n-6123600)-107520)+
<a name="l00139"></a>00139 14918400)-25804800)+19353600)/38707200;
<a name="l00140"></a>00140 _alp[2] = nx*(_n*(_n*(_n*(_n*(4626384*_n-9917165)+4316160)+3743040)-
<a name="l00141"></a>00141 5806080)+2620800)/9676800;
<a name="l00142"></a>00142 _bet[2] = nx*(_n*(_n*(_n*(_n*(829456*_n-5593555)+8478720)-5873280)+
<a name="l00143"></a>00143 1290240)+403200)/19353600;
<a name="l00144"></a>00144 nx *= _n;
<a name="l00145"></a>00145 _alp[3] = nx*(_n*(_n*((26816480-67102379*_n)*_n+16265880)-21358080)+
<a name="l00146"></a>00146 7378560)/29030400;
<a name="l00147"></a>00147 _bet[3] = nx*(_n*(_n*(_n*(9261899*_n+3564160)-2708640)-2557440)+
<a name="l00148"></a>00148 2056320)/58060800;
<a name="l00149"></a>00149 nx *= _n;
<a name="l00150"></a>00150 _alp[4] = nx*(_n*(_n*(155912000*_n+72618271)-85060800)+24532695)/
<a name="l00151"></a>00151 79833600;
<a name="l00152"></a>00152 _bet[4] = nx*(_n*(_n*(14928352*_n-9132761)-1742400)+2176515)/79833600;
<a name="l00153"></a>00153 nx *= _n;
<a name="l00154"></a>00154 _alp[5] = nx*(_n*(102508609*_n-109404448)+27505368)/63866880;
<a name="l00155"></a>00155 _bet[5] = nx*((-8005831*_n-1741552)*_n+1814868)/63866880;
<a name="l00156"></a>00156 nx *= _n;
<a name="l00157"></a>00157 _alp[6] = (2760926233.0-12282192400.0*_n)*nx/4151347200.0;
<a name="l00158"></a>00158 _bet[6] = (268433009-261810608*_n)*nx/8302694400.0;
<a name="l00159"></a>00159 nx *= _n;
<a name="l00160"></a>00160 _alp[7] = 1522256789.0*nx/1383782400.0;
<a name="l00161"></a>00161 _bet[7] = 219941297*nx/5535129600.0;
<a name="l00162"></a>00162 <span class="keywordflow">break</span>;
<a name="l00163"></a>00163 <span class="keywordflow">case</span> 8:
<a name="l00164"></a>00164 _b1 = 1/(1+_n)*(nx*(nx*(nx*(25*nx+64)+256)+4096)+16384)/16384;
<a name="l00165"></a>00165 _alp[1] = _n*(_n*(_n*(_n*(_n*(_n*((37884525-75900428*_n)*_n+42422016)-
<a name="l00166"></a>00166 89611200)+46287360)+63504000)-135475200)+
<a name="l00167"></a>00167 101606400)/203212800;
<a name="l00168"></a>00168 _bet[1] = _n*(_n*(_n*(_n*(_n*(_n*(_n*(31777436*_n-37845269)+43097152)-
<a name="l00169"></a>00169 42865200)-752640)+104428800)-180633600)+
<a name="l00170"></a>00170 135475200)/270950400;
<a name="l00171"></a>00171 _alp[2] = nx*(_n*(_n*(_n*(_n*(_n*(148003883*_n+83274912)-178508970)+
<a name="l00172"></a>00172 77690880)+67374720)-104509440)+47174400)/
<a name="l00173"></a>00173 174182400;
<a name="l00174"></a>00174 _bet[2] = nx*(_n*(_n*(_n*(_n*(_n*(24749483*_n+14930208)-100683990)+
<a name="l00175"></a>00175 152616960)-105719040)+23224320)+7257600)/
<a name="l00176"></a>00176 348364800;
<a name="l00177"></a>00177 nx *= _n;
<a name="l00178"></a>00178 _alp[3] = nx*(_n*(_n*(_n*(_n*(318729724*_n-738126169)+294981280)+
<a name="l00179"></a>00179 178924680)-234938880)+81164160)/319334400;
<a name="l00180"></a>00180 _bet[3] = nx*(_n*(_n*(_n*((101880889-232468668*_n)*_n+39205760)-
<a name="l00181"></a>00181 29795040)-28131840)+22619520)/638668800;
<a name="l00182"></a>00182 nx *= _n;
<a name="l00183"></a>00183 _alp[4] = nx*(_n*(_n*((14967552000.0-40176129013.0*_n)*_n+6971354016.0)-
<a name="l00184"></a>00184 8165836800.0)+2355138720.0)/7664025600.0;
<a name="l00185"></a>00185 _bet[4] = nx*(_n*(_n*(_n*(324154477*_n+1433121792.0)-876745056)-
<a name="l00186"></a>00186 167270400)+208945440)/7664025600.0;
<a name="l00187"></a>00187 nx *= _n;
<a name="l00188"></a>00188 _alp[5] = nx*(_n*(_n*(10421654396.0*_n+3997835751.0)-4266773472.0)+
<a name="l00189"></a>00189 1072709352.0)/2490808320.0;
<a name="l00190"></a>00190 _bet[5] = nx*(_n*(_n*(457888660*_n-312227409)-67920528)+70779852)/
<a name="l00191"></a>00191 2490808320.0;
<a name="l00192"></a>00192 nx *= _n;
<a name="l00193"></a>00193 _alp[6] = nx*(_n*(175214326799.0*_n-171950693600.0)+38652967262.0)/
<a name="l00194"></a>00194 58118860800.0;
<a name="l00195"></a>00195 _bet[6] = nx*((-19841813847.0*_n-3665348512.0)*_n+3758062126.0)/
<a name="l00196"></a>00196 116237721600.0;
<a name="l00197"></a>00197 nx *= _n;
<a name="l00198"></a>00198 _alp[7] = (13700311101.0-67039739596.0*_n)*nx/12454041600.0;
<a name="l00199"></a>00199 _bet[7] = (1979471673.0-1989295244.0*_n)*nx/49816166400.0;
<a name="l00200"></a>00200 nx *= _n;
<a name="l00201"></a>00201 _alp[8] = 1424729850961.0*nx/743921418240.0;
<a name="l00202"></a>00202 _bet[8] = 191773887257.0*nx/3719607091200.0;
<a name="l00203"></a>00203 <span class="keywordflow">break</span>;
<a name="l00204"></a>00204 <span class="keywordflow">default</span>:
<a name="l00205"></a>00205 <a class="code" href="Constants_8hpp.html#a8f24445c1bccd69b63e365aa5d5bb129">STATIC_ASSERT</a>(maxpow >= 4 && maxpow <= 8, <span class="stringliteral">"Bad value of maxpow"</span>);
<a name="l00206"></a>00206 }
<a name="l00207"></a>00207 <span class="comment">// _a1 is the equivalent radius for computing the circumference of</span>
<a name="l00208"></a>00208 <span class="comment">// ellipse.</span>
<a name="l00209"></a>00209 _a1 = _b1 * _a;
<a name="l00210"></a>00210 }
<a name="l00211"></a>00211
<a name="l00212"></a>00212 <span class="keyword">const</span> <a class="code" href="classGeographicLib_1_1TransverseMercator.html" title="Transverse Mercator Projection.">TransverseMercator</a>
<a name="l00213"></a>00213 <a class="code" href="classGeographicLib_1_1TransverseMercator.html#aa25b52e35bc54c368e6c8b17e02b0542">TransverseMercator::UTM</a>(Constants::WGS84_a<real>(),
<a name="l00214"></a>00214 Constants::WGS84_r<real>(),
<a name="l00215"></a>00215 Constants::UTM_k0<real>());
<a name="l00216"></a>00216
<a name="l00217"></a><a class="code" href="classGeographicLib_1_1TransverseMercator.html#a07d73a6d94e2434cbb937978d61a5ae7">00217</a> <span class="keywordtype">void</span> <a class="code" href="classGeographicLib_1_1TransverseMercator.html#a07d73a6d94e2434cbb937978d61a5ae7">TransverseMercator::Forward</a>(real lon0, real lat, real lon,
<a name="l00218"></a>00218 real& x, real& y, real& gamma, real& k)
<a name="l00219"></a>00219 <span class="keyword">const</span> <span class="keywordflow">throw</span>() {
<a name="l00220"></a>00220 <span class="comment">// Avoid losing a bit of accuracy in lon (assuming lon0 is an integer)</span>
<a name="l00221"></a>00221 <span class="keywordflow">if</span> (lon - lon0 > 180)
<a name="l00222"></a>00222 lon -= lon0 + 360;
<a name="l00223"></a>00223 <span class="keywordflow">else</span> <span class="keywordflow">if</span> (lon - lon0 <= -180)
<a name="l00224"></a>00224 lon -= lon0 - 360;
<a name="l00225"></a>00225 <span class="keywordflow">else</span>
<a name="l00226"></a>00226 lon -= lon0;
<a name="l00227"></a>00227 <span class="comment">// Now lon in (-180, 180]</span>
<a name="l00228"></a>00228 <span class="comment">// Explicitly enforce the parity</span>
<a name="l00229"></a>00229 <span class="keywordtype">int</span>
<a name="l00230"></a>00230 latsign = lat < 0 ? -1 : 1,
<a name="l00231"></a>00231 lonsign = lon < 0 ? -1 : 1;
<a name="l00232"></a>00232 lon *= lonsign;
<a name="l00233"></a>00233 lat *= latsign;
<a name="l00234"></a>00234 <span class="keywordtype">bool</span> backside = lon > 90;
<a name="l00235"></a>00235 <span class="keywordflow">if</span> (backside) {
<a name="l00236"></a>00236 <span class="keywordflow">if</span> (lat == 0)
<a name="l00237"></a>00237 latsign = -1;
<a name="l00238"></a>00238 lon = 180 - lon;
<a name="l00239"></a>00239 }
<a name="l00240"></a>00240 real
<a name="l00241"></a>00241 phi = lat * Math::degree<real>(),
<a name="l00242"></a>00242 lam = lon * Math::degree<real>();
<a name="l00243"></a>00243 <span class="comment">// phi = latitude</span>
<a name="l00244"></a>00244 <span class="comment">// phi' = conformal latitude</span>
<a name="l00245"></a>00245 <span class="comment">// psi = isometric latitude</span>
<a name="l00246"></a>00246 <span class="comment">// tau = tan(phi)</span>
<a name="l00247"></a>00247 <span class="comment">// tau' = tan(phi')</span>
<a name="l00248"></a>00248 <span class="comment">// [xi', eta'] = Gauss-Schreiber TM coordinates</span>
<a name="l00249"></a>00249 <span class="comment">// [xi, eta] = Gauss-Krueger TM coordinates</span>
<a name="l00250"></a>00250 <span class="comment">//</span>
<a name="l00251"></a>00251 <span class="comment">// We use</span>
<a name="l00252"></a>00252 <span class="comment">// tan(phi') = sinh(psi)</span>
<a name="l00253"></a>00253 <span class="comment">// sin(phi') = tanh(psi)</span>
<a name="l00254"></a>00254 <span class="comment">// cos(phi') = sech(psi)</span>
<a name="l00255"></a>00255 <span class="comment">// denom^2 = 1-cos(phi')^2*sin(lam)^2 = 1-sech(psi)^2*sin(lam)^2</span>
<a name="l00256"></a>00256 <span class="comment">// sin(xip) = sin(phi')/denom = tanh(psi)/denom</span>
<a name="l00257"></a>00257 <span class="comment">// cos(xip) = cos(phi')*cos(lam)/denom = sech(psi)*cos(lam)/denom</span>
<a name="l00258"></a>00258 <span class="comment">// cosh(etap) = 1/denom = 1/denom</span>
<a name="l00259"></a>00259 <span class="comment">// sinh(etap) = cos(phi')*sin(lam)/denom = sech(psi)*sin(lam)/denom</span>
<a name="l00260"></a>00260 real etap, xip;
<a name="l00261"></a>00261 <span class="keywordflow">if</span> (lat != 90) {
<a name="l00262"></a>00262 real
<a name="l00263"></a>00263 c = max(<a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(0), cos(lam)), <span class="comment">// cos(pi/2) might be negative</span>
<a name="l00264"></a>00264 tau = tan(phi),
<a name="l00265"></a>00265 secphi = <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(<a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(1), tau),
<a name="l00266"></a>00266 sig = sinh( eatanhe(tau / secphi) ),
<a name="l00267"></a>00267 taup = <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(<a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(1), sig) * tau - sig * secphi;
<a name="l00268"></a>00268 xip = atan2(taup, c);
<a name="l00269"></a>00269 <span class="comment">// Used to be</span>
<a name="l00270"></a>00270 <span class="comment">// etap = Math::atanh(sin(lam) / cosh(psi));</span>
<a name="l00271"></a>00271 etap = <a class="code" href="classGeographicLib_1_1Math.html#ab0998a80c8946d1c016c1bc4810a0698">Math::asinh</a>(sin(lam) / <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(taup, c));
<a name="l00272"></a>00272 <span class="comment">// convergence and scale for Gauss-Schreiber TM (xip, etap) -- gamma0 =</span>
<a name="l00273"></a>00273 <span class="comment">// atan(tan(xip) * tanh(etap)) = atan(tan(lam) * sin(phi'));</span>
<a name="l00274"></a>00274 <span class="comment">// sin(phi') = tau'/sqrt(1 + tau'^2)</span>
<a name="l00275"></a>00275 gamma = atan(tanx(lam) *
<a name="l00276"></a>00276 taup / <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(<a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(1), taup)); <span class="comment">// Krueger p 22 (44)</span>
<a name="l00277"></a>00277 <span class="comment">// k0 = sqrt(1 - _e2 * sin(phi)^2) * (cos(phi') / cos(phi)) * cosh(etap)</span>
<a name="l00278"></a>00278 <span class="comment">// Note 1/cos(phi) = cosh(psip);</span>
<a name="l00279"></a>00279 <span class="comment">// and cos(phi') * cosh(etap) = 1/hypot(sinh(psi), cos(lam))</span>
<a name="l00280"></a>00280 <span class="comment">//</span>
<a name="l00281"></a>00281 <span class="comment">// This form has cancelling errors. This property is lost if cosh(psip)</span>
<a name="l00282"></a>00282 <span class="comment">// is replaced by 1/cos(phi), even though it's using "primary" data (phi</span>
<a name="l00283"></a>00283 <span class="comment">// instead of psip).</span>
<a name="l00284"></a>00284 k = sqrt(_e2m + _e2 * sq(cos(phi))) * secphi / <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(taup, c);
<a name="l00285"></a>00285 } <span class="keywordflow">else</span> {
<a name="l00286"></a>00286 xip = Math::pi<real>()/2;
<a name="l00287"></a>00287 etap = 0;
<a name="l00288"></a>00288 gamma = lam;
<a name="l00289"></a>00289 k = _c;
<a name="l00290"></a>00290 }
<a name="l00291"></a>00291 <span class="comment">// {xi',eta'} is {northing,easting} for Gauss-Schreiber transverse Mercator</span>
<a name="l00292"></a>00292 <span class="comment">// (for eta' = 0, xi' = bet). {xi,eta} is {northing,easting} for transverse</span>
<a name="l00293"></a>00293 <span class="comment">// Mercator with constant scale on the central meridian (for eta = 0, xip =</span>
<a name="l00294"></a>00294 <span class="comment">// rectifying latitude). Define</span>
<a name="l00295"></a>00295 <span class="comment">//</span>
<a name="l00296"></a>00296 <span class="comment">// zeta = xi + i*eta</span>
<a name="l00297"></a>00297 <span class="comment">// zeta' = xi' + i*eta'</span>
<a name="l00298"></a>00298 <span class="comment">//</span>
<a name="l00299"></a>00299 <span class="comment">// The conversion from conformal to rectifying latitude can be expresses as</span>
<a name="l00300"></a>00300 <span class="comment">// a series in _n:</span>
<a name="l00301"></a>00301 <span class="comment">//</span>
<a name="l00302"></a>00302 <span class="comment">// zeta = zeta' + sum(h[j-1]' * sin(2 * j * zeta'), j = 1..maxpow)</span>
<a name="l00303"></a>00303 <span class="comment">//</span>
<a name="l00304"></a>00304 <span class="comment">// where h[j]' = O(_n^j). The reversion of this series gives</span>
<a name="l00305"></a>00305 <span class="comment">//</span>
<a name="l00306"></a>00306 <span class="comment">// zeta' = zeta - sum(h[j-1] * sin(2 * j * zeta), j = 1..maxpow)</span>
<a name="l00307"></a>00307 <span class="comment">//</span>
<a name="l00308"></a>00308 <span class="comment">// which is used in Reverse.</span>
<a name="l00309"></a>00309 <span class="comment">//</span>
<a name="l00310"></a>00310 <span class="comment">// Evaluate sums via Clenshaw method. See</span>
<a name="l00311"></a>00311 <span class="comment">// http://mathworld.wolfram.com/ClenshawRecurrenceFormula.html</span>
<a name="l00312"></a>00312 <span class="comment">//</span>
<a name="l00313"></a>00313 <span class="comment">// Let</span>
<a name="l00314"></a>00314 <span class="comment">//</span>
<a name="l00315"></a>00315 <span class="comment">// S = sum(c[k] * F[k](x), k = 0..N)</span>
<a name="l00316"></a>00316 <span class="comment">// F[n+1](x) = alpha(n,x) * F[n](x) + beta(n,x) * F[n-1](x)</span>
<a name="l00317"></a>00317 <span class="comment">//</span>
<a name="l00318"></a>00318 <span class="comment">// Evaluate S with</span>
<a name="l00319"></a>00319 <span class="comment">//</span>
<a name="l00320"></a>00320 <span class="comment">// y[N+2] = y[N+1] = 0</span>
<a name="l00321"></a>00321 <span class="comment">// y[k] = alpha(k,x) * y[k+1] + beta(k+1,x) * y[k+2] + c[k]</span>
<a name="l00322"></a>00322 <span class="comment">// S = c[0] * F[0](x) + y[1] * F[1](x) + beta(1,x) * F[0](x) * y[2]</span>
<a name="l00323"></a>00323 <span class="comment">//</span>
<a name="l00324"></a>00324 <span class="comment">// Here we have</span>
<a name="l00325"></a>00325 <span class="comment">//</span>
<a name="l00326"></a>00326 <span class="comment">// x = 2 * zeta'</span>
<a name="l00327"></a>00327 <span class="comment">// F[n](x) = sin(n * x)</span>
<a name="l00328"></a>00328 <span class="comment">// a(n, x) = 2 * cos(x)</span>
<a name="l00329"></a>00329 <span class="comment">// b(n, x) = -1</span>
<a name="l00330"></a>00330 <span class="comment">// [ sin(A+B) - 2*cos(B)*sin(A) + sin(A-B) = 0, A = n*x, B = x ]</span>
<a name="l00331"></a>00331 <span class="comment">// N = maxpow</span>
<a name="l00332"></a>00332 <span class="comment">// c[k] = _alp[k]</span>
<a name="l00333"></a>00333 <span class="comment">// S = y[1] * sin(x)</span>
<a name="l00334"></a>00334 <span class="comment">//</span>
<a name="l00335"></a>00335 <span class="comment">// For the derivative we have</span>
<a name="l00336"></a>00336 <span class="comment">//</span>
<a name="l00337"></a>00337 <span class="comment">// x = 2 * zeta'</span>
<a name="l00338"></a>00338 <span class="comment">// F[n](x) = cos(n * x)</span>
<a name="l00339"></a>00339 <span class="comment">// a(n, x) = 2 * cos(x)</span>
<a name="l00340"></a>00340 <span class="comment">// b(n, x) = -1</span>
<a name="l00341"></a>00341 <span class="comment">// [ cos(A+B) - 2*cos(B)*cos(A) + cos(A-B) = 0, A = n*x, B = x ]</span>
<a name="l00342"></a>00342 <span class="comment">// c[0] = 1; c[k] = 2*k*_alp[k]</span>
<a name="l00343"></a>00343 <span class="comment">// S = (c[0] - y[2]) + y[1] * cos(x)</span>
<a name="l00344"></a>00344 real
<a name="l00345"></a>00345 c0 = cos(2 * xip), ch0 = cosh(2 * etap),
<a name="l00346"></a>00346 s0 = sin(2 * xip), sh0 = sinh(2 * etap),
<a name="l00347"></a>00347 ar = 2 * c0 * ch0, ai = -2 * s0 * sh0; <span class="comment">// 2 * cos(2*zeta')</span>
<a name="l00348"></a>00348 <span class="keywordtype">int</span> n = maxpow;
<a name="l00349"></a>00349 real
<a name="l00350"></a>00350 xi0 = (n & 1 ? _alp[n] : 0), eta0 = 0,
<a name="l00351"></a>00351 xi1 = 0, eta1 = 0;
<a name="l00352"></a>00352 real <span class="comment">// Accumulators for dzeta/dzeta'</span>
<a name="l00353"></a>00353 yr0 = (n & 1 ? 2 * maxpow * _alp[n--] : 0), yi0 = 0,
<a name="l00354"></a>00354 yr1 = 0, yi1 = 0;
<a name="l00355"></a>00355 <span class="keywordflow">while</span> (n) {
<a name="l00356"></a>00356 xi1 = ar * xi0 - ai * eta0 - xi1 + _alp[n];
<a name="l00357"></a>00357 eta1 = ai * xi0 + ar * eta0 - eta1;
<a name="l00358"></a>00358 yr1 = ar * yr0 - ai * yi0 - yr1 + 2 * n * _alp[n];
<a name="l00359"></a>00359 yi1 = ai * yr0 + ar * yi0 - yi1;
<a name="l00360"></a>00360 --n;
<a name="l00361"></a>00361 xi0 = ar * xi1 - ai * eta1 - xi0 + _alp[n];
<a name="l00362"></a>00362 eta0 = ai * xi1 + ar * eta1 - eta0;
<a name="l00363"></a>00363 yr0 = ar * yr1 - ai * yi1 - yr0 + 2 * n * _alp[n];
<a name="l00364"></a>00364 yi0 = ai * yr1 + ar * yi1 - yi0;
<a name="l00365"></a>00365 --n;
<a name="l00366"></a>00366 }
<a name="l00367"></a>00367 ar /= 2; ai /= 2; <span class="comment">// cos(2*zeta')</span>
<a name="l00368"></a>00368 yr1 = 1 - yr1 + ar * yr0 - ai * yi0;
<a name="l00369"></a>00369 yi1 = - yi1 + ai * yr0 + ar * yi0;
<a name="l00370"></a>00370 ar = s0 * ch0; ai = c0 * sh0; <span class="comment">// sin(2*zeta')</span>
<a name="l00371"></a>00371 real
<a name="l00372"></a>00372 xi = xip + ar * xi0 - ai * eta0,
<a name="l00373"></a>00373 eta = etap + ai * xi0 + ar * eta0;
<a name="l00374"></a>00374 <span class="comment">// Fold in change in convergence and scale for Gauss-Schreiber TM to</span>
<a name="l00375"></a>00375 <span class="comment">// Gauss-Krueger TM.</span>
<a name="l00376"></a>00376 gamma -= atan2(yi1, yr1);
<a name="l00377"></a>00377 k *= _b1 * <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(yr1, yi1);
<a name="l00378"></a>00378 gamma /= Math::degree<real>();
<a name="l00379"></a>00379 y = _a1 * _k0 * (backside ? Math::pi<real>() - xi : xi) * latsign;
<a name="l00380"></a>00380 x = _a1 * _k0 * eta * lonsign;
<a name="l00381"></a>00381 <span class="keywordflow">if</span> (backside)
<a name="l00382"></a>00382 gamma = 180 - gamma;
<a name="l00383"></a>00383 gamma *= latsign * lonsign;
<a name="l00384"></a>00384 k *= _k0;
<a name="l00385"></a>00385 }
<a name="l00386"></a>00386
<a name="l00387"></a><a class="code" href="classGeographicLib_1_1TransverseMercator.html#a15611aee4e3707e155278bab02403a07">00387</a> <span class="keywordtype">void</span> <a class="code" href="classGeographicLib_1_1TransverseMercator.html#a15611aee4e3707e155278bab02403a07">TransverseMercator::Reverse</a>(real lon0, real x, real y,
<a name="l00388"></a>00388 real& lat, real& lon, real& gamma, real& k)
<a name="l00389"></a>00389 <span class="keyword">const</span> <span class="keywordflow">throw</span>() {
<a name="l00390"></a>00390 <span class="comment">// This undoes the steps in Forward. The wrinkles are: (1) Use of the</span>
<a name="l00391"></a>00391 <span class="comment">// reverted series to express zeta' in terms of zeta. (2) Newton's method</span>
<a name="l00392"></a>00392 <span class="comment">// to solve for phi in terms of tan(phi).</span>
<a name="l00393"></a>00393 real
<a name="l00394"></a>00394 xi = y / (_a1 * _k0),
<a name="l00395"></a>00395 eta = x / (_a1 * _k0);
<a name="l00396"></a>00396 <span class="comment">// Explicitly enforce the parity</span>
<a name="l00397"></a>00397 <span class="keywordtype">int</span>
<a name="l00398"></a>00398 xisign = xi < 0 ? -1 : 1,
<a name="l00399"></a>00399 etasign = eta < 0 ? -1 : 1;
<a name="l00400"></a>00400 xi *= xisign;
<a name="l00401"></a>00401 eta *= etasign;
<a name="l00402"></a>00402 <span class="keywordtype">bool</span> backside = xi > Math::pi<real>()/2;
<a name="l00403"></a>00403 <span class="keywordflow">if</span> (backside)
<a name="l00404"></a>00404 xi = Math::pi<real>() - xi;
<a name="l00405"></a>00405 real
<a name="l00406"></a>00406 c0 = cos(2 * xi), ch0 = cosh(2 * eta),
<a name="l00407"></a>00407 s0 = sin(2 * xi), sh0 = sinh(2 * eta),
<a name="l00408"></a>00408 ar = 2 * c0 * ch0, ai = -2 * s0 * sh0; <span class="comment">// 2 * cos(2*zeta)</span>
<a name="l00409"></a>00409 <span class="keywordtype">int</span> n = maxpow;
<a name="l00410"></a>00410 real <span class="comment">// Accumulators for zeta'</span>
<a name="l00411"></a>00411 xip0 = (n & 1 ? -_bet[n] : 0), etap0 = 0,
<a name="l00412"></a>00412 xip1 = 0, etap1 = 0;
<a name="l00413"></a>00413 real <span class="comment">// Accumulators for dzeta'/dzeta</span>
<a name="l00414"></a>00414 yr0 = (n & 1 ? - 2 * maxpow * _bet[n--] : 0), yi0 = 0,
<a name="l00415"></a>00415 yr1 = 0, yi1 = 0;
<a name="l00416"></a>00416 <span class="keywordflow">while</span> (n) {
<a name="l00417"></a>00417 xip1 = ar * xip0 - ai * etap0 - xip1 - _bet[n];
<a name="l00418"></a>00418 etap1 = ai * xip0 + ar * etap0 - etap1;
<a name="l00419"></a>00419 yr1 = ar * yr0 - ai * yi0 - yr1 - 2 * n * _bet[n];
<a name="l00420"></a>00420 yi1 = ai * yr0 + ar * yi0 - yi1;
<a name="l00421"></a>00421 --n;
<a name="l00422"></a>00422 xip0 = ar * xip1 - ai * etap1 - xip0 - _bet[n];
<a name="l00423"></a>00423 etap0 = ai * xip1 + ar * etap1 - etap0;
<a name="l00424"></a>00424 yr0 = ar * yr1 - ai * yi1 - yr0 - 2 * n * _bet[n];
<a name="l00425"></a>00425 yi0 = ai * yr1 + ar * yi1 - yi0;
<a name="l00426"></a>00426 --n;
<a name="l00427"></a>00427 }
<a name="l00428"></a>00428 ar /= 2; ai /= 2; <span class="comment">// cos(2*zeta')</span>
<a name="l00429"></a>00429 yr1 = 1 - yr1 + ar * yr0 - ai * yi0;
<a name="l00430"></a>00430 yi1 = - yi1 + ai * yr0 + ar * yi0;
<a name="l00431"></a>00431 ar = s0 * ch0; ai = c0 * sh0; <span class="comment">// sin(2*zeta)</span>
<a name="l00432"></a>00432 real
<a name="l00433"></a>00433 xip = xi + ar * xip0 - ai * etap0,
<a name="l00434"></a>00434 etap = eta + ai * xip0 + ar * etap0;
<a name="l00435"></a>00435 <span class="comment">// Convergence and scale for Gauss-Schreiber TM to Gauss-Krueger TM.</span>
<a name="l00436"></a>00436 gamma = atan2(yi1, yr1);
<a name="l00437"></a>00437 k = _b1 / <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(yr1, yi1);
<a name="l00438"></a>00438 <span class="comment">// JHS 154 has</span>
<a name="l00439"></a>00439 <span class="comment">//</span>
<a name="l00440"></a>00440 <span class="comment">// phi' = asin(sin(xi') / cosh(eta')) (Krueger p 17 (25))</span>
<a name="l00441"></a>00441 <span class="comment">// lam = asin(tanh(eta') / cos(phi')</span>
<a name="l00442"></a>00442 <span class="comment">// psi = asinh(tan(phi'))</span>
<a name="l00443"></a>00443 real lam, phi;
<a name="l00444"></a>00444 real
<a name="l00445"></a>00445 s = sinh(etap),
<a name="l00446"></a>00446 c = max(<a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(0), cos(xip)), <span class="comment">// cos(pi/2) might be negative</span>
<a name="l00447"></a>00447 r = <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(s, c);
<a name="l00448"></a>00448 <span class="keywordflow">if</span> (r != 0) {
<a name="l00449"></a>00449 lam = atan2(s, c); <span class="comment">// Krueger p 17 (25)</span>
<a name="l00450"></a>00450 <span class="comment">// Use Newton's method to solve for tau</span>
<a name="l00451"></a>00451 real
<a name="l00452"></a>00452 taup = sin(xip)/r,
<a name="l00453"></a>00453 tau = taup,
<a name="l00454"></a>00454 stol = tol * max(<a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(1), abs(taup));
<a name="l00455"></a>00455 <span class="comment">// min iterations = 1, max iterations = 2; mean = 1.99</span>
<a name="l00456"></a>00456 <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i = 0; i < numit; ++i) {
<a name="l00457"></a>00457 real
<a name="l00458"></a>00458 tau1 = <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(<a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(1), tau),
<a name="l00459"></a>00459 sig = sinh( eatanhe( tau / tau1 ) ),
<a name="l00460"></a>00460 taupa = <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(<a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(1), sig) * tau - sig * tau1,
<a name="l00461"></a>00461 dtau = (taup - taupa) * (1 + _e2m * sq(tau)) /
<a name="l00462"></a>00462 ( _e2m * tau1 * <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(<a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(1), taupa) );
<a name="l00463"></a>00463 tau += dtau;
<a name="l00464"></a>00464 <span class="keywordflow">if</span> (!(abs(dtau) >= stol))
<a name="l00465"></a>00465 <span class="keywordflow">break</span>;
<a name="l00466"></a>00466 }
<a name="l00467"></a>00467 phi = atan(tau);
<a name="l00468"></a>00468 gamma += atan(tanx(xip) * tanh(etap)); <span class="comment">// Krueger p 19 (31)</span>
<a name="l00469"></a>00469 <span class="comment">// Note cos(phi') * cosh(eta') = r</span>
<a name="l00470"></a>00470 k *= sqrt(_e2m + _e2 * sq(cos(phi))) * <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(<a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(1), tau) * r;
<a name="l00471"></a>00471 } <span class="keywordflow">else</span> {
<a name="l00472"></a>00472 phi = Math::pi<real>()/2;
<a name="l00473"></a>00473 lam = 0;
<a name="l00474"></a>00474 k *= _c;
<a name="l00475"></a>00475 }
<a name="l00476"></a>00476 lat = phi / Math::degree<real>() * xisign;
<a name="l00477"></a>00477 lon = lam / Math::degree<real>();
<a name="l00478"></a>00478 <span class="keywordflow">if</span> (backside)
<a name="l00479"></a>00479 lon = 180 - lon;
<a name="l00480"></a>00480 lon *= etasign;
<a name="l00481"></a>00481 <span class="comment">// Avoid losing a bit of accuracy in lon (assuming lon0 is an integer)</span>
<a name="l00482"></a>00482 <span class="keywordflow">if</span> (lon + lon0 >= 180)
<a name="l00483"></a>00483 lon += lon0 - 360;
<a name="l00484"></a>00484 <span class="keywordflow">else</span> <span class="keywordflow">if</span> (lon + lon0 < -180)
<a name="l00485"></a>00485 lon += lon0 + 360;
<a name="l00486"></a>00486 <span class="keywordflow">else</span>
<a name="l00487"></a>00487 lon += lon0;
<a name="l00488"></a>00488 gamma /= Math::degree<real>();
<a name="l00489"></a>00489 <span class="keywordflow">if</span> (backside)
<a name="l00490"></a>00490 gamma = 180 - gamma;
<a name="l00491"></a>00491 gamma *= xisign * etasign;
<a name="l00492"></a>00492 k *= _k0;
<a name="l00493"></a>00493 }
<a name="l00494"></a>00494
<a name="l00495"></a>00495 } <span class="comment">// namespace GeographicLib</span>
</pre></div></div>
</div>
<hr class="footer"/><address class="footer"><small>Generated on Tue Feb 22 2011 for GeographicLib by
<a href="http://www.doxygen.org/index.html">
<img class="footer" src="doxygen.png" alt="doxygen"/></a> 1.7.1 </small></address>
</body>
</html>
|