This file is indexed.

/usr/share/doc/geographiclib/html/PolarStereographic_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
<!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: PolarStereographic.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&nbsp;Page</span></a></li>
      <li><a href="pages.html"><span>Related&nbsp;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&nbsp;List</span></a></li>
      <li><a href="globals.html"><span>File&nbsp;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>PolarStereographic.cpp</h1>  </div>
</div>
<div class="contents">
<a href="PolarStereographic_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 PolarStereographic.cpp</span>
<a name="l00003"></a>00003 <span class="comment"> * \brief Implementation for GeographicLib::PolarStereographic 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) &lt;charles@karney.com&gt;</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 
<a name="l00010"></a>00010 <span class="preprocessor">#include &quot;<a class="code" href="PolarStereographic_8hpp.html" title="Header for GeographicLib::PolarStereographic class.">GeographicLib/PolarStereographic.hpp</a>&quot;</span>
<a name="l00011"></a>00011 
<a name="l00012"></a><a class="code" href="PolarStereographic_8cpp.html#a2887f1d8fd2633373698a2d1468befc2">00012</a> <span class="preprocessor">#define GEOGRAPHICLIB_POLARSTEREOGRAPHIC_CPP &quot;$Id: PolarStereographic.cpp 6937 2011-02-01 20:17:13Z karney $&quot;</span>
<a name="l00013"></a>00013 <span class="preprocessor"></span>
<a name="l00014"></a>00014 <a class="code" href="Constants_8hpp.html#af90fa899707a2ac513d5e4c76853bbf5">RCSID_DECL</a>(<a class="code" href="PolarStereographic_8cpp.html#a2887f1d8fd2633373698a2d1468befc2">GEOGRAPHICLIB_POLARSTEREOGRAPHIC_CPP</a>)
<a name="l00015"></a>00015 <a class="code" href="Constants_8hpp.html#af90fa899707a2ac513d5e4c76853bbf5">RCSID_DECL</a>(<a class="code" href="PolarStereographic_8hpp.html#a7450566a6ed1475dc10f7c0dc85f5dd9">GEOGRAPHICLIB_POLARSTEREOGRAPHIC_HPP</a>)
<a name="l00016"></a>00016 
<a name="l00017"></a>00017 namespace GeographicLib {
<a name="l00018"></a>00018 
<a name="l00019"></a>00019   <span class="keyword">using namespace </span>std;
<a name="l00020"></a>00020 
<a name="l00021"></a>00021   <span class="keyword">const</span> <a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">Math::real</a> PolarStereographic::tol =
<a name="l00022"></a>00022     <a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(0.1)*sqrt(numeric_limits&lt;real&gt;::epsilon());
<a name="l00023"></a>00023   <span class="comment">// Overflow value s.t. atan(overflow) = pi/2</span>
<a name="l00024"></a>00024   <span class="keyword">const</span> <a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">Math::real</a> PolarStereographic::overflow =
<a name="l00025"></a>00025     1 / sq(numeric_limits&lt;real&gt;::epsilon());
<a name="l00026"></a>00026 
<a name="l00027"></a><a class="code" href="classGeographicLib_1_1PolarStereographic.html#a56558071847ee7d71674ac29a23ec653">00027</a>   PolarStereographic::PolarStereographic(real a, real r, real k0)
<a name="l00028"></a>00028     : _a(a)
<a name="l00029"></a>00029     , _r(r)
<a name="l00030"></a>00030     , _f(_r != 0 ? 1 / _r : 0)
<a name="l00031"></a>00031     , _e2(_f * (2 - _f))
<a name="l00032"></a>00032     , _e(sqrt(abs(_e2)))
<a name="l00033"></a>00033     , _e2m(1 - _e2)
<a name="l00034"></a>00034     , _Cx(exp(eatanhe(real(1))))
<a name="l00035"></a>00035     , _c( (1 - _f) * _Cx )
<a name="l00036"></a>00036     , _k0(k0)
<a name="l00037"></a>00037   {
<a name="l00038"></a>00038     <span class="keywordflow">if</span> (!(_a &gt; 0))
<a name="l00039"></a>00039       <span class="keywordflow">throw</span> <a class="code" href="classGeographicLib_1_1GeographicErr.html" title="Exception handling for GeographicLib">GeographicErr</a>(<span class="stringliteral">&quot;Major radius is not positive&quot;</span>);
<a name="l00040"></a>00040     <span class="keywordflow">if</span> (!(_f &lt; 1))
<a name="l00041"></a>00041       <span class="keywordflow">throw</span> <a class="code" href="classGeographicLib_1_1GeographicErr.html" title="Exception handling for GeographicLib">GeographicErr</a>(<span class="stringliteral">&quot;Minor radius is not positive&quot;</span>);
<a name="l00042"></a>00042     <span class="keywordflow">if</span> (!(_k0 &gt; 0))
<a name="l00043"></a>00043       <span class="keywordflow">throw</span> <a class="code" href="classGeographicLib_1_1GeographicErr.html" title="Exception handling for GeographicLib">GeographicErr</a>(<span class="stringliteral">&quot;Scale is not positive&quot;</span>);
<a name="l00044"></a>00044   }
<a name="l00045"></a>00045 
<a name="l00046"></a>00046   <span class="keyword">const</span> <a class="code" href="classGeographicLib_1_1PolarStereographic.html" title="Polar Stereographic Projection.">PolarStereographic</a>
<a name="l00047"></a>00047   <a class="code" href="classGeographicLib_1_1PolarStereographic.html#a2db6bcb1b59a6ddc6087ee04c64c9825">PolarStereographic::UPS</a>(Constants::WGS84_a&lt;real&gt;(),
<a name="l00048"></a>00048                           Constants::WGS84_r&lt;real&gt;(),
<a name="l00049"></a>00049                           Constants::UPS_k0&lt;real&gt;());
<a name="l00050"></a>00050 
<a name="l00051"></a>00051   <span class="comment">// This formulation converts to conformal coordinates by tau = tan(phi) and</span>
<a name="l00052"></a>00052   <span class="comment">// tau&#39; = tan(phi&#39;) where phi&#39; is the conformal latitude.  The formulas are:</span>
<a name="l00053"></a>00053   <span class="comment">//    tau = tan(phi)</span>
<a name="l00054"></a>00054   <span class="comment">//    secphi = hypot(1, tau)</span>
<a name="l00055"></a>00055   <span class="comment">//    sig = sinh(e * atanh(e * tau / secphi))</span>
<a name="l00056"></a>00056   <span class="comment">//    taup = tan(phip) = tau * hypot(1, sig) - sig * hypot(1, tau)</span>
<a name="l00057"></a>00057   <span class="comment">//    c = (1 - f) * exp(e * atanh(e))</span>
<a name="l00058"></a>00058   <span class="comment">//</span>
<a name="l00059"></a>00059   <span class="comment">// Forward:</span>
<a name="l00060"></a>00060   <span class="comment">//   rho = (2*k0*a/c) / (hypot(1, taup) + taup)  (taup &gt;= 0)</span>
<a name="l00061"></a>00061   <span class="comment">//       = (2*k0*a/c) * (hypot(1, taup) - taup)  (taup &lt;  0)</span>
<a name="l00062"></a>00062   <span class="comment">//</span>
<a name="l00063"></a>00063   <span class="comment">// Reverse:</span>
<a name="l00064"></a>00064   <span class="comment">//   taup = ((2*k0*a/c) / rho - rho / (2*k0*a/c))/2</span>
<a name="l00065"></a>00065   <span class="comment">//</span>
<a name="l00066"></a>00066   <span class="comment">// Scale:</span>
<a name="l00067"></a>00067   <span class="comment">//   k = (rho/a) * secphi * sqrt((1-e2) + e2 / secphi^2)</span>
<a name="l00068"></a>00068   <span class="comment">//</span>
<a name="l00069"></a>00069   <span class="comment">// In limit rho -&gt; 0, tau -&gt; inf, taup -&gt; inf, secphi -&gt; inf, secphip -&gt; inf</span>
<a name="l00070"></a>00070   <span class="comment">//   secphip = taup = exp(-e * atanh(e)) * tau = exp(-e * atanh(e)) * secphi</span>
<a name="l00071"></a>00071 
<a name="l00072"></a><a class="code" href="classGeographicLib_1_1PolarStereographic.html#a30ef1a1f906ee389e2a5f7e5fd7d8fa4">00072</a>   <span class="keywordtype">void</span> <a class="code" href="classGeographicLib_1_1PolarStereographic.html#a30ef1a1f906ee389e2a5f7e5fd7d8fa4">PolarStereographic::Forward</a>(<span class="keywordtype">bool</span> northp, real lat, real lon,
<a name="l00073"></a>00073                                    real&amp; x, real&amp; y, real&amp; gamma, real&amp; k)
<a name="l00074"></a>00074     <span class="keyword">const</span> <span class="keywordflow">throw</span>() {
<a name="l00075"></a>00075     lat *= northp ? 1 : -1;
<a name="l00076"></a>00076     real
<a name="l00077"></a>00077       phi = lat * Math::degree&lt;real&gt;(),
<a name="l00078"></a>00078       tau = lat != -90 ? tan(phi) : -overflow,
<a name="l00079"></a>00079       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="l00080"></a>00080       sig = sinh( eatanhe(tau / secphi) ),
<a name="l00081"></a>00081       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="l00082"></a>00082       rho =  <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(<a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(1), taup) + abs(taup);
<a name="l00083"></a>00083     rho = taup &gt;= 0 ? (lat != 90 ? 1/rho : 0) : rho;
<a name="l00084"></a>00084     rho *= 2 * _k0 * _a / _c;
<a name="l00085"></a>00085     k = lat != 90 ? (rho / _a) * secphi * sqrt(_e2m + _e2 / sq(secphi)) :
<a name="l00086"></a>00086       _k0;
<a name="l00087"></a>00087     lon = lon &gt;= 180 ? lon - 360 : lon &lt; -180 ? lon + 360 : lon;
<a name="l00088"></a>00088     real
<a name="l00089"></a>00089       lam = lon * Math::degree&lt;real&gt;();
<a name="l00090"></a>00090     x = rho * (lon == -180 ? 0 : sin(lam));
<a name="l00091"></a>00091     y = (northp ? -rho : rho) * (abs(lon) == 90 ? 0 : cos(lam));
<a name="l00092"></a>00092     gamma = northp ? lon : -lon;
<a name="l00093"></a>00093   }
<a name="l00094"></a>00094 
<a name="l00095"></a><a class="code" href="classGeographicLib_1_1PolarStereographic.html#a01302b8dba43c57e3c3849f94123a157">00095</a>   <span class="keywordtype">void</span> <a class="code" href="classGeographicLib_1_1PolarStereographic.html#a01302b8dba43c57e3c3849f94123a157">PolarStereographic::Reverse</a>(<span class="keywordtype">bool</span> northp, real x, real y,
<a name="l00096"></a>00096                                    real&amp; lat, real&amp; lon, real&amp; gamma, real&amp; k)
<a name="l00097"></a>00097     <span class="keyword">const</span> <span class="keywordflow">throw</span>() {
<a name="l00098"></a>00098     real
<a name="l00099"></a>00099       rho = <a class="code" href="classGeographicLib_1_1Math.html#a0d422863198d4bec2aae6b187a60760c">Math::hypot</a>(x, y),
<a name="l00100"></a>00100       t = rho / (2 * _k0 * _a / _c),
<a name="l00101"></a>00101       taup = (1 / t - t) / 2,
<a name="l00102"></a>00102       tau = taup * _Cx,
<a name="l00103"></a>00103       stol = tol * max(<a class="code" href="Geod_8cpp.html#a5caf95d46b184d9ca1d3764b3781b3c9">real</a>(1), abs(taup));
<a name="l00104"></a>00104     <span class="comment">// min iterations = 1, max iterations = 2; mean = 1.99</span>
<a name="l00105"></a>00105     <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i = 0; i &lt; numit; ++i) {
<a name="l00106"></a>00106       real
<a name="l00107"></a>00107         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="l00108"></a>00108         sig = sinh( eatanhe( tau / tau1 ) ),
<a name="l00109"></a>00109         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="l00110"></a>00110         dtau = (taup - taupa) * (1 + _e2m * sq(tau)) /
<a name="l00111"></a>00111         ( _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="l00112"></a>00112       tau += dtau;
<a name="l00113"></a>00113       <span class="keywordflow">if</span> (!(abs(dtau) &gt;= stol))
<a name="l00114"></a>00114         <span class="keywordflow">break</span>;
<a name="l00115"></a>00115     }
<a name="l00116"></a>00116     real
<a name="l00117"></a>00117       phi = atan(tau),
<a name="l00118"></a>00118       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="l00119"></a>00119     k = rho != 0 ? (rho / _a) * secphi * sqrt(_e2m + _e2 / sq(secphi)) : _k0;
<a name="l00120"></a>00120     lat = (northp ? 1 : -1) * (rho != 0 ? phi / Math::degree&lt;real&gt;() : 90);
<a name="l00121"></a>00121     lon = -atan2( -x, northp ? -y : y ) / Math::degree&lt;real&gt;();
<a name="l00122"></a>00122     gamma = northp ? lon : -lon;
<a name="l00123"></a>00123   }
<a name="l00124"></a>00124 
<a name="l00125"></a><a class="code" href="classGeographicLib_1_1PolarStereographic.html#a3f957214eb1d1248277a680e4c4ceed5">00125</a>   <span class="keywordtype">void</span> <a class="code" href="classGeographicLib_1_1PolarStereographic.html#a3f957214eb1d1248277a680e4c4ceed5">PolarStereographic::SetScale</a>(real lat, real k) {
<a name="l00126"></a>00126     <span class="keywordflow">if</span> (!(k &gt; 0))
<a name="l00127"></a>00127       <span class="keywordflow">throw</span> <a class="code" href="classGeographicLib_1_1GeographicErr.html" title="Exception handling for GeographicLib">GeographicErr</a>(<span class="stringliteral">&quot;Scale is not positive&quot;</span>);
<a name="l00128"></a>00128     <span class="keywordflow">if</span> (!(-90 &lt; lat &amp;&amp; lat &lt;= 90))
<a name="l00129"></a>00129       <span class="keywordflow">throw</span> <a class="code" href="classGeographicLib_1_1GeographicErr.html" title="Exception handling for GeographicLib">GeographicErr</a>(<span class="stringliteral">&quot;Latitude must be in (-90d, 90d]&quot;</span>);
<a name="l00130"></a>00130     real x, y, gamma, kold;
<a name="l00131"></a>00131     _k0 = 1;
<a name="l00132"></a>00132     <a class="code" href="classGeographicLib_1_1PolarStereographic.html#a30ef1a1f906ee389e2a5f7e5fd7d8fa4">Forward</a>(<span class="keyword">true</span>, lat, 0, x, y, gamma, kold);
<a name="l00133"></a>00133     _k0 *= k/kold;
<a name="l00134"></a>00134   }
<a name="l00135"></a>00135 
<a name="l00136"></a>00136 } <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&nbsp;
<a href="http://www.doxygen.org/index.html">
<img class="footer" src="doxygen.png" alt="doxygen"/></a> 1.7.1 </small></address>
</body>
</html>