Fortran library for Geodesics  1.52
geodesic.for
Go to the documentation of this file.
1 * The subroutines in this files are documented at
2 * https://geographiclib.sourceforge.io/html/Fortran/
3 *
4 *> @file geodesic.for
5 *! @brief Implementation of geodesic routines in Fortran
6 *!
7 *! This is a Fortran implementation of the geodesic algorithms described
8 *! in
9 *! - C. F. F. Karney,
10 *! <a href="https://doi.org/10.1007/s00190-012-0578-z">
11 *! Algorithms for geodesics</a>,
12 *! J. Geodesy <b>87</b>, 43--55 (2013);
13 *! DOI: <a href="https://doi.org/10.1007/s00190-012-0578-z">
14 *! 10.1007/s00190-012-0578-z</a>;
15 *! addenda: <a href="https://geographiclib.sourceforge.io/geod-addenda.html">
16 *! geod-addenda.html</a>.
17 *! .
18 *! The principal advantages of these algorithms over previous ones
19 *! (e.g., Vincenty, 1975) are
20 *! - accurate to round off for |<i>f</i>| &lt; 1/50;
21 *! - the solution of the inverse problem is always found;
22 *! - differential and integral properties of geodesics are computed.
23 *!
24 *! The shortest path between two points on the ellipsoid at (\e lat1, \e
25 *! lon1) and (\e lat2, \e lon2) is called the geodesic. Its length is
26 *! \e s12 and the geodesic from point 1 to point 2 has forward azimuths
27 *! \e azi1 and \e azi2 at the two end points.
28 *!
29 *! Traditionally two geodesic problems are considered:
30 *! - the direct problem -- given \e lat1, \e lon1, \e s12, and \e azi1,
31 *! determine \e lat2, \e lon2, and \e azi2. This is solved by the
32 *! subroutine direct().
33 *! - the inverse problem -- given \e lat1, \e lon1, \e lat2, \e lon2,
34 *! determine \e s12, \e azi1, and \e azi2. This is solved by the
35 *! subroutine invers().
36 *!
37 *! The ellipsoid is specified by its equatorial radius \e a (typically
38 *! in meters) and flattening \e f. The routines are accurate to round
39 *! off with double precision arithmetic provided that |<i>f</i>| &lt;
40 *! 1/50; for the WGS84 ellipsoid, the errors are less than 15
41 *! nanometers. (Reasonably accurate results are obtained for |<i>f</i>|
42 *! &lt; 1/5.) For a prolate ellipsoid, specify \e f &lt; 0.
43 *!
44 *! The routines also calculate several other quantities of interest
45 *! - \e SS12 is the area between the geodesic from point 1 to point 2
46 *! and the equator; i.e., it is the area, measured counter-clockwise,
47 *! of the geodesic quadrilateral with corners (\e lat1,\e lon1), (0,\e
48 *! lon1), (0,\e lon2), and (\e lat2,\e lon2).
49 *! - \e m12, the reduced length of the geodesic is defined such that if
50 *! the initial azimuth is perturbed by \e dazi1 (radians) then the
51 *! second point is displaced by \e m12 \e dazi1 in the direction
52 *! perpendicular to the geodesic. On a curved surface the reduced
53 *! length obeys a symmetry relation, \e m12 + \e m21 = 0. On a flat
54 *! surface, we have \e m12 = \e s12.
55 *! - \e MM12 and \e MM21 are geodesic scales. If two geodesics are
56 *! parallel at point 1 and separated by a small distance \e dt, then
57 *! they are separated by a distance \e MM12 \e dt at point 2. \e MM21
58 *! is defined similarly (with the geodesics being parallel to one
59 *! another at point 2). On a flat surface, we have \e MM12 = \e MM21
60 *! = 1.
61 *! - \e a12 is the arc length on the auxiliary sphere. This is a
62 *! construct for converting the problem to one in spherical
63 *! trigonometry. \e a12 is measured in degrees. The spherical arc
64 *! length from one equator crossing to the next is always 180&deg;.
65 *!
66 *! If points 1, 2, and 3 lie on a single geodesic, then the following
67 *! addition rules hold:
68 *! - \e s13 = \e s12 + \e s23
69 *! - \e a13 = \e a12 + \e a23
70 *! - \e SS13 = \e SS12 + \e SS23
71 *! - \e m13 = \e m12 \e MM23 + \e m23 \e MM21
72 *! - \e MM13 = \e MM12 \e MM23 &minus; (1 &minus; \e MM12 \e MM21) \e
73 *! m23 / \e m12
74 *! - \e MM31 = \e MM32 \e MM21 &minus; (1 &minus; \e MM23 \e MM32) \e
75 *! m12 / \e m23
76 *!
77 *! The shortest distance returned by the solution of the inverse problem
78 *! is (obviously) uniquely defined. However, in a few special cases
79 *! there are multiple azimuths which yield the same shortest distance.
80 *! Here is a catalog of those cases:
81 *! - \e lat1 = &minus;\e lat2 (with neither point at a pole). If \e
82 *! azi1 = \e azi2, the geodesic is unique. Otherwise there are two
83 *! geodesics and the second one is obtained by setting [\e azi1, \e
84 *! azi2] &rarr; [\e azi2, \e azi1], [\e MM12, \e MM21] &rarr; [\e
85 *! MM21, \e MM12], \e SS12 &rarr; &minus;\e SS12. (This occurs when
86 *! the longitude difference is near &plusmn;180&deg; for oblate
87 *! ellipsoids.)
88 *! - \e lon2 = \e lon1 &plusmn; 180&deg; (with neither point at a pole).
89 *! If \e azi1 = 0&deg; or &plusmn;180&deg;, the geodesic is unique.
90 *! Otherwise there are two geodesics and the second one is obtained by
91 *! setting [\e azi1, \e azi2] &rarr; [&minus;\e azi1, &minus;\e azi2],
92 *! \e SS12 &rarr; &minus;\e SS12. (This occurs when \e lat2 is near
93 *! &minus;\e lat1 for prolate ellipsoids.)
94 *! - Points 1 and 2 at opposite poles. There are infinitely many
95 *! geodesics which can be generated by setting [\e azi1, \e azi2]
96 *! &rarr; [\e azi1, \e azi2] + [\e d, &minus;\e d], for arbitrary \e
97 *! d. (For spheres, this prescription applies when points 1 and 2 are
98 *! antipodal.)
99 *! - \e s12 = 0 (coincident points). There are infinitely many
100 *! geodesics which can be generated by setting [\e azi1, \e azi2]
101 *! &rarr; [\e azi1, \e azi2] + [\e d, \e d], for arbitrary \e d.
102 *!
103 *! These routines are a simple transcription of the corresponding C++
104 *! classes in <a href="https://geographiclib.sourceforge.io">
105 *! GeographicLib</a>. Because of the limitations of Fortran 77, the
106 *! classes have been replaced by simple subroutines with no attempt to
107 *! save "state" across subroutine calls. Most of the internal comments
108 *! have been retained. However, in the process of transcription some
109 *! documentation has been lost and the documentation for the C++
110 *! classes, GeographicLib::Geodesic, GeographicLib::GeodesicLine, and
111 *! GeographicLib::PolygonAreaT, should be consulted. The C++ code
112 *! remains the "reference implementation". Think twice about
113 *! restructuring the internals of the Fortran code since this may make
114 *! porting fixes from the C++ code more difficult.
115 *!
116 *! Copyright (c) Charles Karney (2012-2021) <charles@karney.com> and
117 *! licensed under the MIT/X11 License. For more information, see
118 *! https://geographiclib.sourceforge.io/
119 *!
120 *! This library was distributed with
121 *! <a href="../index.html">GeographicLib</a> 1.52.
122 
123 *> Solve the direct geodesic problem
124 *!
125 *! @param[in] a the equatorial radius (meters).
126 *! @param[in] f the flattening of the ellipsoid. Setting \e f = 0 gives
127 *! a sphere. Negative \e f gives a prolate ellipsoid.
128 *! @param[in] lat1 latitude of point 1 (degrees).
129 *! @param[in] lon1 longitude of point 1 (degrees).
130 *! @param[in] azi1 azimuth at point 1 (degrees).
131 *! @param[in] s12a12 if \e arcmode is not set, this is the distance
132 *! from point 1 to point 2 (meters); otherwise it is the arc
133 *! length from point 1 to point 2 (degrees); it can be negative.
134 *! @param[in] flags a bitor'ed combination of the \e arcmode and \e
135 *! unroll flags.
136 *! @param[out] lat2 latitude of point 2 (degrees).
137 *! @param[out] lon2 longitude of point 2 (degrees).
138 *! @param[out] azi2 (forward) azimuth at point 2 (degrees).
139 *! @param[in] omask a bitor'ed combination of mask values
140 *! specifying which of the following parameters should be set.
141 *! @param[out] a12s12 if \e arcmode is not set, this is the arc length
142 *! from point 1 to point 2 (degrees); otherwise it is the distance
143 *! from point 1 to point 2 (meters).
144 *! @param[out] m12 reduced length of geodesic (meters).
145 *! @param[out] MM12 geodesic scale of point 2 relative to point 1
146 *! (dimensionless).
147 *! @param[out] MM21 geodesic scale of point 1 relative to point 2
148 *! (dimensionless).
149 *! @param[out] SS12 area under the geodesic (meters<sup>2</sup>).
150 *!
151 *! \e flags is an integer in [0, 4) whose binary bits are interpreted
152 *! as follows
153 *! - 1 the \e arcmode flag
154 *! - 2 the \e unroll flag
155 *! .
156 *! If \e arcmode is not set, \e s12a12 is \e s12 and \e a12s12 is \e
157 *! a12; otherwise, \e s12a12 is \e a12 and \e a12s12 is \e s12. It \e
158 *! unroll is not set, the value \e lon2 returned is in the range
159 *! [&minus;180&deg;, 180&deg;]; if unroll is set, the longitude variable
160 *! is "unrolled" so that \e lon2 &minus; \e lon1 indicates how many
161 *! times and in what sense the geodesic encircles the ellipsoid.
162 *!
163 *! \e omask is an integer in [0, 16) whose binary bits are interpreted
164 *! as follows
165 *! - 1 return \e a12
166 *! - 2 return \e m12
167 *! - 4 return \e MM12 and \e MM21
168 *! - 8 return \e SS12
169 *!
170 *! \e lat1 should be in the range [&minus;90&deg;, 90&deg;]. The value
171 *! \e azi2 returned is in the range [&minus;180&deg;, 180&deg;].
172 *!
173 *! If either point is at a pole, the azimuth is defined by keeping the
174 *! longitude fixed, writing \e lat = \e lat = &plusmn;(90&deg; &minus;
175 *! &epsilon;), and taking the limit &epsilon; &rarr; 0+. An arc length
176 *! greater that 180&deg; signifies a geodesic which is not a shortest
177 *! path. (For a prolate ellipsoid, an additional condition is necessary
178 *! for a shortest path: the longitudinal extent must not exceed of
179 *! 180&deg;.)
180 *!
181 *! Example of use:
182 *! \include geoddirect.for
183 
184  subroutine direct(a, f, lat1, lon1, azi1, s12a12, flags,
185  + lat2, lon2, azi2, omask, a12s12, m12, MM12, MM21, SS12)
186 * input
187  double precision a, f, lat1, lon1, azi1, s12a12
188  integer flags, omask
189 * output
190  double precision lat2, lon2, azi2
191 * optional output
192  double precision a12s12, m12, MM12, MM21, SS12
193 
194  integer ord, nC1, nC1p, nC2, nA3, nA3x, nC3, nC3x, nC4, nC4x
195  parameter(ord = 6, nc1 = ord, nc1p = ord,
196  + nc2 = ord, na3 = ord, na3x = na3,
197  + nc3 = ord, nc3x = (nc3 * (nc3 - 1)) / 2,
198  + nc4 = ord, nc4x = (nc4 * (nc4 + 1)) / 2)
199  double precision A3x(0:nA3x-1), C3x(0:nC3x-1), C4x(0:nC4x-1),
200  + c1a(nc1), c1pa(nc1p), c2a(nc2), c3a(nc3-1), c4a(0:nc4-1)
201 
202  double precision atanhx, hypotx,
203  + angnm, angrnd, trgsum, a1m1f, a2m1f, a3f, atn2dx, latfix
204  logical arcmod, unroll, arcp, redlp, scalp, areap
205  double precision e2, f1, ep2, n, b, c2,
206  + salp0, calp0, k2, eps,
207  + salp1, calp1, ssig1, csig1, cbet1, sbet1, dn1, somg1, comg1,
208  + salp2, calp2, ssig2, csig2, sbet2, cbet2, dn2, somg2, comg2,
209  + ssig12, csig12, salp12, calp12, omg12, lam12, lon12,
210  + sig12, stau1, ctau1, tau12, t, s, c, serr, e,
211  + a1m1, a2m1, a3c, a4, ab1, ab2,
212  + b11, b12, b21, b22, b31, b41, b42, j12
213 
214  double precision dblmin, dbleps, pi, degree, tiny,
215  + tol0, tol1, tol2, tolb, xthrsh
216  integer digits, maxit1, maxit2
217  logical init
218  common /geocom/ dblmin, dbleps, pi, degree, tiny,
219  + tol0, tol1, tol2, tolb, xthrsh, digits, maxit1, maxit2, init
220 
221  if (.not.init) call geoini
222 
223  e2 = f * (2 - f)
224  ep2 = e2 / (1 - e2)
225  f1 = 1 - f
226  n = f / (2 - f)
227  b = a * f1
228  c2 = 0
229 
230  arcmod = mod(flags/1, 2) .eq. 1
231  unroll = mod(flags/2, 2) .eq. 1
232 
233  arcp = mod(omask/1, 2) .eq. 1
234  redlp = mod(omask/2, 2) .eq. 1
235  scalp = mod(omask/4, 2) .eq. 1
236  areap = mod(omask/8, 2) .eq. 1
237 
238  if (areap) then
239  if (e2 .eq. 0) then
240  c2 = a**2
241  else if (e2 .gt. 0) then
242  c2 = (a**2 + b**2 * atanhx(sqrt(e2)) / sqrt(e2)) / 2
243  else
244  c2 = (a**2 + b**2 * atan(sqrt(abs(e2))) / sqrt(abs(e2))) / 2
245  end if
246  end if
247 
248  call a3cof(n, a3x)
249  call c3cof(n, c3x)
250  if (areap) call c4cof(n, c4x)
251 
252 * Guard against underflow in salp0
253  call sncsdx(angrnd(azi1), salp1, calp1)
254 
255  call sncsdx(angrnd(latfix(lat1)), sbet1, cbet1)
256  sbet1 = f1 * sbet1
257  call norm2x(sbet1, cbet1)
258 * Ensure cbet1 = +dbleps at poles
259  cbet1 = max(tiny, cbet1)
260  dn1 = sqrt(1 + ep2 * sbet1**2)
261 
262 * Evaluate alp0 from sin(alp1) * cos(bet1) = sin(alp0),
263 * alp0 in [0, pi/2 - |bet1|]
264  salp0 = salp1 * cbet1
265 * Alt: calp0 = hypot(sbet1, calp1 * cbet1). The following
266 * is slightly better (consider the case salp1 = 0).
267  calp0 = hypotx(calp1, salp1 * sbet1)
268 * Evaluate sig with tan(bet1) = tan(sig1) * cos(alp1).
269 * sig = 0 is nearest northward crossing of equator.
270 * With bet1 = 0, alp1 = pi/2, we have sig1 = 0 (equatorial line).
271 * With bet1 = pi/2, alp1 = -pi, sig1 = pi/2
272 * With bet1 = -pi/2, alp1 = 0 , sig1 = -pi/2
273 * Evaluate omg1 with tan(omg1) = sin(alp0) * tan(sig1).
274 * With alp0 in (0, pi/2], quadrants for sig and omg coincide.
275 * No atan2(0,0) ambiguity at poles since cbet1 = +dbleps.
276 * With alp0 = 0, omg1 = 0 for alp1 = 0, omg1 = pi for alp1 = pi.
277  ssig1 = sbet1
278  somg1 = salp0 * sbet1
279  if (sbet1 .ne. 0 .or. calp1 .ne. 0) then
280  csig1 = cbet1 * calp1
281  else
282  csig1 = 1
283  end if
284  comg1 = csig1
285 * sig1 in (-pi, pi]
286  call norm2x(ssig1, csig1)
287 * norm2x(somg1, comg1); -- don't need to normalize!
288 
289  k2 = calp0**2 * ep2
290  eps = k2 / (2 * (1 + sqrt(1 + k2)) + k2)
291 
292  a1m1 = a1m1f(eps)
293  call c1f(eps, c1a)
294  b11 = trgsum(.true., ssig1, csig1, c1a, nc1)
295  s = sin(b11)
296  c = cos(b11)
297 * tau1 = sig1 + B11
298  stau1 = ssig1 * c + csig1 * s
299  ctau1 = csig1 * c - ssig1 * s
300 * Not necessary because C1pa reverts C1a
301 * B11 = -TrgSum(true, stau1, ctau1, C1pa, nC1p)
302 
303  if (.not. arcmod) call c1pf(eps, c1pa)
304 
305  if (redlp .or. scalp) then
306  a2m1 = a2m1f(eps)
307  call c2f(eps, c2a)
308  b21 = trgsum(.true., ssig1, csig1, c2a, nc2)
309  else
310 * Suppress bogus warnings about unitialized variables
311  a2m1 = 0
312  b21 = 0
313  end if
314 
315  call c3f(eps, c3x, c3a)
316  a3c = -f * salp0 * a3f(eps, a3x)
317  b31 = trgsum(.true., ssig1, csig1, c3a, nc3-1)
318 
319  if (areap) then
320  call c4f(eps, c4x, c4a)
321 * Multiplier = a^2 * e^2 * cos(alpha0) * sin(alpha0)
322  a4 = a**2 * calp0 * salp0 * e2
323  b41 = trgsum(.false., ssig1, csig1, c4a, nc4)
324  else
325 * Suppress bogus warnings about unitialized variables
326  a4 = 0
327  b41 = 0
328  end if
329 
330  if (arcmod) then
331 * Interpret s12a12 as spherical arc length
332  sig12 = s12a12 * degree
333  call sncsdx(s12a12, ssig12, csig12)
334 * Suppress bogus warnings about unitialized variables
335  b12 = 0
336  else
337 * Interpret s12a12 as distance
338  tau12 = s12a12 / (b * (1 + a1m1))
339  s = sin(tau12)
340  c = cos(tau12)
341 * tau2 = tau1 + tau12
342  b12 = - trgsum(.true.,
343  + stau1 * c + ctau1 * s, ctau1 * c - stau1 * s, c1pa, nc1p)
344  sig12 = tau12 - (b12 - b11)
345  ssig12 = sin(sig12)
346  csig12 = cos(sig12)
347  if (abs(f) .gt. 0.01d0) then
348 * Reverted distance series is inaccurate for |f| > 1/100, so correct
349 * sig12 with 1 Newton iteration. The following table shows the
350 * approximate maximum error for a = WGS_a() and various f relative to
351 * GeodesicExact.
352 * erri = the error in the inverse solution (nm)
353 * errd = the error in the direct solution (series only) (nm)
354 * errda = the error in the direct solution (series + 1 Newton) (nm)
355 *
356 * f erri errd errda
357 * -1/5 12e6 1.2e9 69e6
358 * -1/10 123e3 12e6 765e3
359 * -1/20 1110 108e3 7155
360 * -1/50 18.63 200.9 27.12
361 * -1/100 18.63 23.78 23.37
362 * -1/150 18.63 21.05 20.26
363 * 1/150 22.35 24.73 25.83
364 * 1/100 22.35 25.03 25.31
365 * 1/50 29.80 231.9 30.44
366 * 1/20 5376 146e3 10e3
367 * 1/10 829e3 22e6 1.5e6
368 * 1/5 157e6 3.8e9 280e6
369  ssig2 = ssig1 * csig12 + csig1 * ssig12
370  csig2 = csig1 * csig12 - ssig1 * ssig12
371  b12 = trgsum(.true., ssig2, csig2, c1a, nc1)
372  serr = (1 + a1m1) * (sig12 + (b12 - b11)) - s12a12 / b
373  sig12 = sig12 - serr / sqrt(1 + k2 * ssig2**2)
374  ssig12 = sin(sig12)
375  csig12 = cos(sig12)
376 * Update B12 below
377  end if
378  end if
379 
380 * sig2 = sig1 + sig12
381  ssig2 = ssig1 * csig12 + csig1 * ssig12
382  csig2 = csig1 * csig12 - ssig1 * ssig12
383  dn2 = sqrt(1 + k2 * ssig2**2)
384  if (arcmod .or. abs(f) .gt. 0.01d0)
385  + b12 = trgsum(.true., ssig2, csig2, c1a, nc1)
386  ab1 = (1 + a1m1) * (b12 - b11)
387 
388 * sin(bet2) = cos(alp0) * sin(sig2)
389  sbet2 = calp0 * ssig2
390 * Alt: cbet2 = hypot(csig2, salp0 * ssig2)
391  cbet2 = hypotx(salp0, calp0 * csig2)
392  if (cbet2 .eq. 0) then
393 * I.e., salp0 = 0, csig2 = 0. Break the degeneracy in this case
394  cbet2 = tiny
395  csig2 = cbet2
396  end if
397 * tan(omg2) = sin(alp0) * tan(sig2)
398 * No need to normalize
399  somg2 = salp0 * ssig2
400  comg2 = csig2
401 * tan(alp0) = cos(sig2)*tan(alp2)
402 * No need to normalize
403  salp2 = salp0
404  calp2 = calp0 * csig2
405 * East or west going?
406  e = sign(1d0, salp0)
407 * omg12 = omg2 - omg1
408  if (unroll) then
409  omg12 = e * (sig12
410  + - (atan2( ssig2, csig2) - atan2( ssig1, csig1))
411  + + (atan2(e * somg2, comg2) - atan2(e * somg1, comg1)))
412  else
413  omg12 = atan2(somg2 * comg1 - comg2 * somg1,
414  + comg2 * comg1 + somg2 * somg1)
415  end if
416 
417  lam12 = omg12 + a3c *
418  + ( sig12 + (trgsum(.true., ssig2, csig2, c3a, nc3-1)
419  + - b31))
420  lon12 = lam12 / degree
421  if (unroll) then
422  lon2 = lon1 + lon12
423  else
424  lon2 = angnm(angnm(lon1) + angnm(lon12))
425  end if
426  lat2 = atn2dx(sbet2, f1 * cbet2)
427  azi2 = atn2dx(salp2, calp2)
428 
429  if (redlp .or. scalp) then
430  b22 = trgsum(.true., ssig2, csig2, c2a, nc2)
431  ab2 = (1 + a2m1) * (b22 - b21)
432  j12 = (a1m1 - a2m1) * sig12 + (ab1 - ab2)
433  end if
434 * Add parens around (csig1 * ssig2) and (ssig1 * csig2) to ensure
435 * accurate cancellation in the case of coincident points.
436  if (redlp) m12 = b * ((dn2 * (csig1 * ssig2) -
437  + dn1 * (ssig1 * csig2)) - csig1 * csig2 * j12)
438  if (scalp) then
439  t = k2 * (ssig2 - ssig1) * (ssig2 + ssig1) / (dn1 + dn2)
440  mm12 = csig12 + (t * ssig2 - csig2 * j12) * ssig1 / dn1
441  mm21 = csig12 - (t * ssig1 - csig1 * j12) * ssig2 / dn2
442  end if
443 
444  if (areap) then
445  b42 = trgsum(.false., ssig2, csig2, c4a, nc4)
446  if (calp0 .eq. 0 .or. salp0 .eq. 0) then
447 * alp12 = alp2 - alp1, used in atan2 so no need to normalize
448  salp12 = salp2 * calp1 - calp2 * salp1
449  calp12 = calp2 * calp1 + salp2 * salp1
450  else
451 * tan(alp) = tan(alp0) * sec(sig)
452 * tan(alp2-alp1) = (tan(alp2) -tan(alp1)) / (tan(alp2)*tan(alp1)+1)
453 * = calp0 * salp0 * (csig1-csig2) / (salp0^2 + calp0^2 * csig1*csig2)
454 * If csig12 > 0, write
455 * csig1 - csig2 = ssig12 * (csig1 * ssig12 / (1 + csig12) + ssig1)
456 * else
457 * csig1 - csig2 = csig1 * (1 - csig12) + ssig12 * ssig1
458 * No need to normalize
459  if (csig12 .le. 0) then
460  salp12 = csig1 * (1 - csig12) + ssig12 * ssig1
461  else
462  salp12 = ssig12 * (csig1 * ssig12 / (1 + csig12) + ssig1)
463  end if
464  salp12 = calp0 * salp0 * salp12
465  calp12 = salp0**2 + calp0**2 * csig1 * csig2
466  end if
467  ss12 = c2 * atan2(salp12, calp12) + a4 * (b42 - b41)
468  end if
469 
470  if (arcp) then
471  if (arcmod) then
472  a12s12 = b * ((1 + a1m1) * sig12 + ab1)
473  else
474  a12s12 = sig12 / degree
475  end if
476  end if
477 
478  return
479  end
480 
481 *> Solve the inverse geodesic problem.
482 *!
483 *! @param[in] a the equatorial radius (meters).
484 *! @param[in] f the flattening of the ellipsoid. Setting \e f = 0 gives
485 *! a sphere. Negative \e f gives a prolate ellipsoid.
486 *! @param[in] lat1 latitude of point 1 (degrees).
487 *! @param[in] lon1 longitude of point 1 (degrees).
488 *! @param[in] lat2 latitude of point 2 (degrees).
489 *! @param[in] lon2 longitude of point 2 (degrees).
490 *! @param[out] s12 distance from point 1 to point 2 (meters).
491 *! @param[out] azi1 azimuth at point 1 (degrees).
492 *! @param[out] azi2 (forward) azimuth at point 2 (degrees).
493 *! @param[in] omask a bitor'ed combination of mask values
494 *! specifying which of the following parameters should be set.
495 *! @param[out] a12 arc length from point 1 to point 2 (degrees).
496 *! @param[out] m12 reduced length of geodesic (meters).
497 *! @param[out] MM12 geodesic scale of point 2 relative to point 1
498 *! (dimensionless).
499 *! @param[out] MM21 geodesic scale of point 1 relative to point 2
500 *! (dimensionless).
501 *! @param[out] SS12 area under the geodesic (meters<sup>2</sup>).
502 *!
503 *! \e omask is an integer in [0, 16) whose binary bits are interpreted
504 *! as follows
505 *! - 1 return \e a12
506 *! - 2 return \e m12
507 *! - 4 return \e MM12 and \e MM21
508 *! - 8 return \e SS12
509 *!
510 *! \e lat1 and \e lat2 should be in the range [&minus;90&deg;, 90&deg;].
511 *! The values of \e azi1 and \e azi2 returned are in the range
512 *! [&minus;180&deg;, 180&deg;].
513 *!
514 *! If either point is at a pole, the azimuth is defined by keeping the
515 *! longitude fixed, writing \e lat = &plusmn;(90&deg; &minus;
516 *! &epsilon;), and taking the limit &epsilon; &rarr; 0+.
517 *!
518 *! The solution to the inverse problem is found using Newton's method.
519 *! If this fails to converge (this is very unlikely in geodetic
520 *! applications but does occur for very eccentric ellipsoids), then the
521 *! bisection method is used to refine the solution.
522 *!
523 *! Example of use:
524 *! \include geodinverse.for
525 
526  subroutine invers(a, f, lat1, lon1, lat2, lon2,
527  + s12, azi1, azi2, omask, a12, m12, MM12, MM21, SS12)
528 * input
529  double precision a, f, lat1, lon1, lat2, lon2
530  integer omask
531 * output
532  double precision s12, azi1, azi2
533 * optional output
534  double precision a12, m12, MM12, MM21, SS12
535 
536  integer ord, nA3, nA3x, nC3, nC3x, nC4, nC4x, nC
537  parameter (ord = 6, na3 = ord, na3x = na3,
538  + nc3 = ord, nc3x = (nc3 * (nc3 - 1)) / 2,
539  + nc4 = ord, nc4x = (nc4 * (nc4 + 1)) / 2,
540  + nc = ord)
541  double precision A3x(0:nA3x-1), C3x(0:nC3x-1), C4x(0:nC4x-1),
542  + Ca(nC)
543 
544  double precision atanhx, hypotx,
545  + AngDif, AngRnd, TrgSum, Lam12f, InvSta, atn2dx, LatFix
546  integer latsgn, lonsgn, swapp, numit
547  logical arcp, redlp, scalp, areap, merid, tripn, tripb
548 
549  double precision e2, f1, ep2, n, b, c2,
550  + lat1x, lat2x, salp0, calp0, k2, eps,
551  + salp1, calp1, ssig1, csig1, cbet1, sbet1, dbet1, dn1,
552  + salp2, calp2, ssig2, csig2, sbet2, cbet2, dbet2, dn2,
553  + slam12, clam12, salp12, calp12, omg12, lam12, lon12, lon12s,
554  + salp1a, calp1a, salp1b, calp1b,
555  + dalp1, sdalp1, cdalp1, nsalp1, alp12, somg12, comg12, domg12,
556  + sig12, v, dv, dnm, dummy,
557  + a4, b41, b42, s12x, m12x, a12x, sdomg12, cdomg12
558 
559  double precision dblmin, dbleps, pi, degree, tiny,
560  + tol0, tol1, tol2, tolb, xthrsh
561  integer digits, maxit1, maxit2, lmask
562  logical init
563  common /geocom/ dblmin, dbleps, pi, degree, tiny,
564  + tol0, tol1, tol2, tolb, xthrsh, digits, maxit1, maxit2, init
565 
566  if (.not.init) call geoini
567 
568  f1 = 1 - f
569  e2 = f * (2 - f)
570  ep2 = e2 / f1**2
571  n = f / ( 2 - f)
572  b = a * f1
573  c2 = 0
574 
575  arcp = mod(omask/1, 2) .eq. 1
576  redlp = mod(omask/2, 2) .eq. 1
577  scalp = mod(omask/4, 2) .eq. 1
578  areap = mod(omask/8, 2) .eq. 1
579  if (scalp) then
580  lmask = 16 + 2 + 4
581  else
582  lmask = 16 + 2
583  end if
584 
585  if (areap) then
586  if (e2 .eq. 0) then
587  c2 = a**2
588  else if (e2 .gt. 0) then
589  c2 = (a**2 + b**2 * atanhx(sqrt(e2)) / sqrt(e2)) / 2
590  else
591  c2 = (a**2 + b**2 * atan(sqrt(abs(e2))) / sqrt(abs(e2))) / 2
592  end if
593  end if
594 
595  call a3cof(n, a3x)
596  call c3cof(n, c3x)
597  if (areap) call c4cof(n, c4x)
598 
599 * Compute longitude difference (AngDiff does this carefully). Result is
600 * in [-180, 180] but -180 is only for west-going geodesics. 180 is for
601 * east-going and meridional geodesics.
602 * If very close to being on the same half-meridian, then make it so.
603  lon12 = angdif(lon1, lon2, lon12s)
604 * Make longitude difference positive.
605  if (lon12 .ge. 0) then
606  lonsgn = 1
607  else
608  lonsgn = -1
609  end if
610  lon12 = lonsgn * angrnd(lon12)
611  lon12s = angrnd((180 - lon12) - lonsgn * lon12s)
612  lam12 = lon12 * degree
613  if (lon12 .gt. 90) then
614  call sncsdx(lon12s, slam12, clam12)
615  clam12 = -clam12
616  else
617  call sncsdx(lon12, slam12, clam12)
618  end if
619 
620 * If really close to the equator, treat as on equator.
621  lat1x = angrnd(latfix(lat1))
622  lat2x = angrnd(latfix(lat2))
623 * Swap points so that point with higher (abs) latitude is point 1
624 * If one latitude is a nan, then it becomes lat1.
625  if (abs(lat1x) .lt. abs(lat2x)) then
626  swapp = -1
627  else
628  swapp = 1
629  end if
630  if (swapp .lt. 0) then
631  lonsgn = -lonsgn
632  call swap(lat1x, lat2x)
633  end if
634 * Make lat1 <= 0
635  if (lat1x .lt. 0) then
636  latsgn = 1
637  else
638  latsgn = -1
639  end if
640  lat1x = lat1x * latsgn
641  lat2x = lat2x * latsgn
642 * Now we have
643 *
644 * 0 <= lon12 <= 180
645 * -90 <= lat1 <= 0
646 * lat1 <= lat2 <= -lat1
647 *
648 * longsign, swapp, latsgn register the transformation to bring the
649 * coordinates to this canonical form. In all cases, 1 means no change
650 * was made. We make these transformations so that there are few cases
651 * to check, e.g., on verifying quadrants in atan2. In addition, this
652 * enforces some symmetries in the results returned.
653 
654  call sncsdx(lat1x, sbet1, cbet1)
655  sbet1 = f1 * sbet1
656  call norm2x(sbet1, cbet1)
657 * Ensure cbet1 = +dbleps at poles
658  cbet1 = max(tiny, cbet1)
659 
660  call sncsdx(lat2x, sbet2, cbet2)
661  sbet2 = f1 * sbet2
662  call norm2x(sbet2, cbet2)
663 * Ensure cbet2 = +dbleps at poles
664  cbet2 = max(tiny, cbet2)
665 
666 * If cbet1 < -sbet1, then cbet2 - cbet1 is a sensitive measure of the
667 * |bet1| - |bet2|. Alternatively (cbet1 >= -sbet1), abs(sbet2) + sbet1
668 * is a better measure. This logic is used in assigning calp2 in
669 * Lambda12. Sometimes these quantities vanish and in that case we force
670 * bet2 = +/- bet1 exactly. An example where is is necessary is the
671 * inverse problem 48.522876735459 0 -48.52287673545898293
672 * 179.599720456223079643 which failed with Visual Studio 10 (Release and
673 * Debug)
674 
675  if (cbet1 .lt. -sbet1) then
676  if (cbet2 .eq. cbet1) sbet2 = sign(sbet1, sbet2)
677  else
678  if (abs(sbet2) .eq. -sbet1) cbet2 = cbet1
679  end if
680 
681  dn1 = sqrt(1 + ep2 * sbet1**2)
682  dn2 = sqrt(1 + ep2 * sbet2**2)
683 
684 * Suppress bogus warnings about unitialized variables
685  a12x = 0
686  merid = lat1x .eq. -90 .or. slam12 .eq. 0
687 
688  if (merid) then
689 
690 * Endpoints are on a single full meridian, so the geodesic might lie on
691 * a meridian.
692 
693 * Head to the target longitude
694  calp1 = clam12
695  salp1 = slam12
696 * At the target we're heading north
697  calp2 = 1
698  salp2 = 0
699 
700 * tan(bet) = tan(sig) * cos(alp)
701  ssig1 = sbet1
702  csig1 = calp1 * cbet1
703  ssig2 = sbet2
704  csig2 = calp2 * cbet2
705 
706 * sig12 = sig2 - sig1
707  sig12 = atan2(0d0 + max(0d0, csig1 * ssig2 - ssig1 * csig2),
708  + csig1 * csig2 + ssig1 * ssig2)
709  call lengs(n, sig12, ssig1, csig1, dn1, ssig2, csig2, dn2,
710  + cbet1, cbet2, lmask,
711  + s12x, m12x, dummy, mm12, mm21, ep2, ca)
712 
713 * Add the check for sig12 since zero length geodesics might yield m12 <
714 * 0. Test case was
715 *
716 * echo 20.001 0 20.001 0 | GeodSolve -i
717 *
718 * In fact, we will have sig12 > pi/2 for meridional geodesic which is
719 * not a shortest path.
720  if (sig12 .lt. 1 .or. m12x .ge. 0) then
721  if (sig12 .lt. 3 * tiny .or.
722  + (sig12 .lt. tol0 .and.
723  + (s12x .lt. 0 .or. m12x .lt. 0))) then
724 * Prevent negative s12 or m12 for short lines
725  sig12 = 0
726  m12x = 0
727  s12x = 0
728  end if
729  m12x = m12x * b
730  s12x = s12x * b
731  a12x = sig12 / degree
732  else
733 * m12 < 0, i.e., prolate and too close to anti-podal
734  merid = .false.
735  end if
736  end if
737 
738  omg12 = 0
739 * somg12 > 1 marks that it needs to be calculated
740  somg12 = 2
741  comg12 = 0
742  if (.not. merid .and. sbet1 .eq. 0 .and.
743  + (f .le. 0 .or. lon12s .ge. f * 180)) then
744 
745 * Geodesic runs along equator
746  calp1 = 0
747  calp2 = 0
748  salp1 = 1
749  salp2 = 1
750  s12x = a * lam12
751  sig12 = lam12 / f1
752  omg12 = sig12
753  m12x = b * sin(sig12)
754  if (scalp) then
755  mm12 = cos(sig12)
756  mm21 = mm12
757  end if
758  a12x = lon12 / f1
759  else if (.not. merid) then
760 * Now point1 and point2 belong within a hemisphere bounded by a
761 * meridian and geodesic is neither meridional or equatorial.
762 
763 * Figure a starting point for Newton's method
764  sig12 = invsta(sbet1, cbet1, dn1, sbet2, cbet2, dn2, lam12,
765  + slam12, clam12, f, a3x, salp1, calp1, salp2, calp2, dnm, ca)
766 
767  if (sig12 .ge. 0) then
768 * Short lines (InvSta sets salp2, calp2, dnm)
769  s12x = sig12 * b * dnm
770  m12x = dnm**2 * b * sin(sig12 / dnm)
771  if (scalp) then
772  mm12 = cos(sig12 / dnm)
773  mm21 = mm12
774  end if
775  a12x = sig12 / degree
776  omg12 = lam12 / (f1 * dnm)
777  else
778 
779 * Newton's method. This is a straightforward solution of f(alp1) =
780 * lambda12(alp1) - lam12 = 0 with one wrinkle. f(alp) has exactly one
781 * root in the interval (0, pi) and its derivative is positive at the
782 * root. Thus f(alp) is positive for alp > alp1 and negative for alp <
783 * alp1. During the course of the iteration, a range (alp1a, alp1b) is
784 * maintained which brackets the root and with each evaluation of
785 * f(alp) the range is shrunk, if possible. Newton's method is
786 * restarted whenever the derivative of f is negative (because the new
787 * value of alp1 is then further from the solution) or if the new
788 * estimate of alp1 lies outside (0,pi); in this case, the new starting
789 * guess is taken to be (alp1a + alp1b) / 2.
790 
791 * Bracketing range
792  salp1a = tiny
793  calp1a = 1
794  salp1b = tiny
795  calp1b = -1
796  tripn = .false.
797  tripb = .false.
798  do 10 numit = 0, maxit2-1
799 * the WGS84 test set: mean = 1.47, sd = 1.25, max = 16
800 * WGS84 and random input: mean = 2.85, sd = 0.60
801  v = lam12f(sbet1, cbet1, dn1, sbet2, cbet2, dn2,
802  + salp1, calp1, slam12, clam12, f, a3x, c3x, salp2, calp2,
803  + sig12, ssig1, csig1, ssig2, csig2,
804  + eps, domg12, numit .lt. maxit1, dv, ca)
805 * Reversed test to allow escape with NaNs
806  if (tripn) then
807  dummy = 8
808  else
809  dummy = 1
810  end if
811  if (tripb .or. .not. (abs(v) .ge. dummy * tol0))
812  + go to 20
813 * Update bracketing values
814  if (v .gt. 0 .and. (numit .gt. maxit1 .or.
815  + calp1/salp1 .gt. calp1b/salp1b)) then
816  salp1b = salp1
817  calp1b = calp1
818  else if (v .lt. 0 .and. (numit .gt. maxit1 .or.
819  + calp1/salp1 .lt. calp1a/salp1a)) then
820  salp1a = salp1
821  calp1a = calp1
822  end if
823  if (numit .lt. maxit1 .and. dv .gt. 0) then
824  dalp1 = -v/dv
825  sdalp1 = sin(dalp1)
826  cdalp1 = cos(dalp1)
827  nsalp1 = salp1 * cdalp1 + calp1 * sdalp1
828  if (nsalp1 .gt. 0 .and. abs(dalp1) .lt. pi) then
829  calp1 = calp1 * cdalp1 - salp1 * sdalp1
830  salp1 = nsalp1
831  call norm2x(salp1, calp1)
832 * In some regimes we don't get quadratic convergence because
833 * slope -> 0. So use convergence conditions based on dbleps
834 * instead of sqrt(dbleps).
835  tripn = abs(v) .le. 16 * tol0
836  go to 10
837  end if
838  end if
839 * Either dv was not positive or updated value was outside legal
840 * range. Use the midpoint of the bracket as the next estimate.
841 * This mechanism is not needed for the WGS84 ellipsoid, but it does
842 * catch problems with more eccentric ellipsoids. Its efficacy is
843 * such for the WGS84 test set with the starting guess set to alp1 =
844 * 90deg:
845 * the WGS84 test set: mean = 5.21, sd = 3.93, max = 24
846 * WGS84 and random input: mean = 4.74, sd = 0.99
847  salp1 = (salp1a + salp1b)/2
848  calp1 = (calp1a + calp1b)/2
849  call norm2x(salp1, calp1)
850  tripn = .false.
851  tripb = abs(salp1a - salp1) + (calp1a - calp1) .lt. tolb
852  + .or. abs(salp1 - salp1b) + (calp1 - calp1b) .lt. tolb
853  10 continue
854  20 continue
855  call lengs(eps, sig12, ssig1, csig1, dn1, ssig2, csig2, dn2,
856  + cbet1, cbet2, lmask,
857  + s12x, m12x, dummy, mm12, mm21, ep2, ca)
858  m12x = m12x * b
859  s12x = s12x * b
860  a12x = sig12 / degree
861  if (areap) then
862  sdomg12 = sin(domg12)
863  cdomg12 = cos(domg12)
864  somg12 = slam12 * cdomg12 - clam12 * sdomg12
865  comg12 = clam12 * cdomg12 + slam12 * sdomg12
866  end if
867  end if
868  end if
869 
870 * Convert -0 to 0
871  s12 = 0 + s12x
872  if (redlp) m12 = 0 + m12x
873 
874  if (areap) then
875 * From Lambda12: sin(alp1) * cos(bet1) = sin(alp0)
876  salp0 = salp1 * cbet1
877  calp0 = hypotx(calp1, salp1 * sbet1)
878  if (calp0 .ne. 0 .and. salp0 .ne. 0) then
879 * From Lambda12: tan(bet) = tan(sig) * cos(alp)
880  ssig1 = sbet1
881  csig1 = calp1 * cbet1
882  ssig2 = sbet2
883  csig2 = calp2 * cbet2
884  k2 = calp0**2 * ep2
885  eps = k2 / (2 * (1 + sqrt(1 + k2)) + k2)
886 * Multiplier = a^2 * e^2 * cos(alpha0) * sin(alpha0).
887  a4 = a**2 * calp0 * salp0 * e2
888  call norm2x(ssig1, csig1)
889  call norm2x(ssig2, csig2)
890  call c4f(eps, c4x, ca)
891  b41 = trgsum(.false., ssig1, csig1, ca, nc4)
892  b42 = trgsum(.false., ssig2, csig2, ca, nc4)
893  ss12 = a4 * (b42 - b41)
894  else
895 * Avoid problems with indeterminate sig1, sig2 on equator
896  ss12 = 0
897  end if
898 
899  if (.not. merid .and. somg12 .gt. 1) then
900  somg12 = sin(omg12)
901  comg12 = cos(omg12)
902  end if
903 
904  if (.not. merid .and. comg12 .ge. 0.7071d0
905  + .and. sbet2 - sbet1 .lt. 1.75d0) then
906 * Use tan(Gamma/2) = tan(omg12/2)
907 * * (tan(bet1/2)+tan(bet2/2))/(1+tan(bet1/2)*tan(bet2/2))
908 * with tan(x/2) = sin(x)/(1+cos(x))
909  domg12 = 1 + comg12
910  dbet1 = 1 + cbet1
911  dbet2 = 1 + cbet2
912  alp12 = 2 * atan2(somg12 * (sbet1 * dbet2 + sbet2 * dbet1),
913  + domg12 * ( sbet1 * sbet2 + dbet1 * dbet2 ) )
914  else
915 * alp12 = alp2 - alp1, used in atan2 so no need to normalize
916  salp12 = salp2 * calp1 - calp2 * salp1
917  calp12 = calp2 * calp1 + salp2 * salp1
918 * The right thing appears to happen if alp1 = +/-180 and alp2 = 0, viz
919 * salp12 = -0 and alp12 = -180. However this depends on the sign
920 * being attached to 0 correctly. The following ensures the correct
921 * behavior.
922  if (salp12 .eq. 0 .and. calp12 .lt. 0) then
923  salp12 = tiny * calp1
924  calp12 = -1
925  end if
926  alp12 = atan2(salp12, calp12)
927  end if
928  ss12 = ss12 + c2 * alp12
929  ss12 = ss12 * swapp * lonsgn * latsgn
930 * Convert -0 to 0
931  ss12 = 0 + ss12
932  end if
933 
934 * Convert calp, salp to azimuth accounting for lonsgn, swapp, latsgn.
935  if (swapp .lt. 0) then
936  call swap(salp1, salp2)
937  call swap(calp1, calp2)
938  if (scalp) call swap(mm12, mm21)
939  end if
940 
941  salp1 = salp1 * swapp * lonsgn
942  calp1 = calp1 * swapp * latsgn
943  salp2 = salp2 * swapp * lonsgn
944  calp2 = calp2 * swapp * latsgn
945 
946  azi1 = atn2dx(salp1, calp1)
947  azi2 = atn2dx(salp2, calp2)
948 
949  if (arcp) a12 = a12x
950 
951  return
952  end
953 
954 *> Determine the area of a geodesic polygon
955 *!
956 *! @param[in] a the equatorial radius (meters).
957 *! @param[in] f the flattening of the ellipsoid. Setting \e f = 0 gives
958 *! a sphere. Negative \e f gives a prolate ellipsoid.
959 *! @param[in] lats an array of the latitudes of the vertices (degrees).
960 *! @param[in] lons an array of the longitudes of the vertices (degrees).
961 *! @param[in] n the number of vertices.
962 *! @param[out] AA the (signed) area of the polygon (meters<sup>2</sup>).
963 *! @param[out] PP the perimeter of the polygon.
964 *!
965 *! \e lats should be in the range [&minus;90&deg;, 90&deg;].
966 *!
967 *! Arbitrarily complex polygons are allowed. In the case of
968 *! self-intersecting polygons the area is accumulated "algebraically",
969 *! e.g., the areas of the 2 loops in a figure-8 polygon will partially
970 *! cancel. There's no need to "close" the polygon by repeating the
971 *! first vertex. The area returned is signed with counter-clockwise
972 *! traversal being treated as positive.
973 
974  subroutine area(a, f, lats, lons, n, AA, PP)
975 * input
976  integer n
977  double precision a, f, lats(n), lons(n)
978 * output
979  double precision AA, PP
980 
981  integer i, omask, cross, trnsit
982  double precision s12, azi1, azi2, dummy, SS12, b, e2, c2, area0,
983  + atanhx, aacc(2), pacc(2)
984 
985  double precision dblmin, dbleps, pi, degree, tiny,
986  + tol0, tol1, tol2, tolb, xthrsh
987  integer digits, maxit1, maxit2
988  logical init
989  common /geocom/ dblmin, dbleps, pi, degree, tiny,
990  + tol0, tol1, tol2, tolb, xthrsh, digits, maxit1, maxit2, init
991 
992  omask = 8
993  call accini(aacc)
994  call accini(pacc)
995  cross = 0
996  do 10 i = 0, n-1
997  call invers(a, f, lats(i+1), lons(i+1),
998  + lats(mod(i+1,n)+1), lons(mod(i+1,n)+1),
999  + s12, azi1, azi2, omask, dummy, dummy, dummy, dummy, ss12)
1000  call accadd(pacc, s12)
1001  call accadd(aacc, -ss12)
1002  cross = cross + trnsit(lons(i+1), lons(mod(i+1,n)+1))
1003  10 continue
1004  pp = pacc(1)
1005  b = a * (1 - f)
1006  e2 = f * (2 - f)
1007  if (e2 .eq. 0) then
1008  c2 = a**2
1009  else if (e2 .gt. 0) then
1010  c2 = (a**2 + b**2 * atanhx(sqrt(e2)) / sqrt(e2)) / 2
1011  else
1012  c2 = (a**2 + b**2 * atan(sqrt(abs(e2))) / sqrt(abs(e2))) / 2
1013  end if
1014  area0 = 4 * pi * c2
1015  call accrem(aacc, area0)
1016  if (mod(abs(cross), 2) .eq. 1) then
1017  if (aacc(1) .lt. 0) then
1018  call accadd(aacc, +area0/2)
1019  else
1020  call accadd(aacc, -area0/2)
1021  end if
1022  end if
1023  if (aacc(1) .gt. area0/2) then
1024  call accadd(aacc, -area0)
1025  else if (aacc(1) .le. -area0/2) then
1026  call accadd(aacc, +area0)
1027  end if
1028  aa = aacc(1)
1029 
1030  return
1031  end
1032 
1033 *> Return the version numbers for this package.
1034 *!
1035 *! @param[out] major the major version number.
1036 *! @param[out] minor the minor version number.
1037 *! @param[out] patch the patch number.
1038 *!
1039 *! This subroutine was added with version 1.44.
1040 
1041  subroutine geover(major, minor, patch)
1042 * output
1043  integer major, minor, patch
1044 
1045  major = 1
1046  minor = 52
1047  patch = 0
1048 
1049  return
1050  end
1051 
1052 *> @cond SKIP
1053 
1054  block data geodat
1055  double precision dblmin, dbleps, pi, degree, tiny,
1056  + tol0, tol1, tol2, tolb, xthrsh
1057  integer digits, maxit1, maxit2
1058  logical init
1059  data init /.false./
1060  common /geocom/ dblmin, dbleps, pi, degree, tiny,
1061  + tol0, tol1, tol2, tolb, xthrsh, digits, maxit1, maxit2, init
1062  end
1063 
1064  subroutine geoini
1065  double precision dblmin, dbleps, pi, degree, tiny,
1066  + tol0, tol1, tol2, tolb, xthrsh
1067  integer digits, maxit1, maxit2
1068  logical init
1069  common /geocom/ dblmin, dbleps, pi, degree, tiny,
1070  + tol0, tol1, tol2, tolb, xthrsh, digits, maxit1, maxit2, init
1071 
1072  digits = 53
1073  dblmin = 0.5d0**1022
1074  dbleps = 0.5d0**(digits-1)
1075 
1076  pi = atan2(0d0, -1d0)
1077  degree = pi/180
1078 * This is about cbrt(dblmin). With other implementations, sqrt(dblmin)
1079 * is used. The larger value is used here to avoid complaints about a
1080 * IEEE_UNDERFLOW_FLAG IEEE_DENORMAL signal. This is triggered when
1081 * invers is called with points at opposite poles.
1082  tiny = 0.5d0**((1022+1)/3)
1083  tol0 = dbleps
1084 * Increase multiplier in defn of tol1 from 100 to 200 to fix inverse
1085 * case 52.784459512564 0 -52.784459512563990912 179.634407464943777557
1086 * which otherwise failed for Visual Studio 10 (Release and Debug)
1087  tol1 = 200 * tol0
1088  tol2 = sqrt(tol0)
1089 * Check on bisection interval
1090  tolb = tol0 * tol2
1091  xthrsh = 1000 * tol2
1092  maxit1 = 20
1093  maxit2 = maxit1 + digits + 10
1094 
1095  init = .true.
1096 
1097  return
1098  end
1099 
1100  subroutine lengs(eps, sig12, ssig1, csig1, dn1, ssig2, csig2, dn2,
1101  + cbet1, cbet2, omask,
1102  + s12b, m12b, m0, MM12, MM21, ep2, Ca)
1103 * input
1104  double precision eps, sig12, ssig1, csig1, dn1, ssig2, csig2, dn2,
1105  + cbet1, cbet2, ep2
1106  integer omask
1107 * optional output
1108  double precision s12b, m12b, m0, MM12, MM21
1109 * temporary storage
1110  double precision Ca(*)
1111 
1112  integer ord, nC1, nC2
1113  parameter (ord = 6, nc1 = ord, nc2 = ord)
1114 
1115  double precision A1m1f, A2m1f, TrgSum
1116  double precision m0x, J12, A1, A2, B1, B2, csig12, t, Cb(nC2)
1117  logical distp, redlp, scalp
1118  integer l
1119 
1120 * Return m12b = (reduced length)/b; also calculate s12b = distance/b,
1121 * and m0 = coefficient of secular term in expression for reduced length.
1122 
1123  distp = (mod(omask/16, 2) .eq. 1)
1124  redlp = (mod(omask/2, 2) .eq. 1)
1125  scalp = (mod(omask/4, 2) .eq. 1)
1126 
1127 * Suppress compiler warnings
1128  m0x = 0
1129  j12 = 0
1130  a1 = 0
1131  a2 = 0
1132  if (distp .or. redlp .or. scalp) then
1133  a1 = a1m1f(eps)
1134  call c1f(eps, ca)
1135  if (redlp .or. scalp) then
1136  a2 = a2m1f(eps)
1137  call c2f(eps, cb)
1138  m0x = a1 - a2
1139  a2 = 1 + a2
1140  end if
1141  a1 = 1 + a1
1142  end if
1143  if (distp) then
1144  b1 = trgsum(.true., ssig2, csig2, ca, nc1) -
1145  + trgsum(.true., ssig1, csig1, ca, nc1)
1146 * Missing a factor of b
1147  s12b = a1 * (sig12 + b1)
1148  if (redlp .or. scalp) then
1149  b2 = trgsum(.true., ssig2, csig2, cb, nc2) -
1150  + trgsum(.true., ssig1, csig1, cb, nc2)
1151  j12 = m0x * sig12 + (a1 * b1 - a2 * b2)
1152  end if
1153  else if (redlp .or. scalp) then
1154 * Assume here that nC1 >= nC2
1155  do 10 l = 1, nc2
1156  cb(l) = a1 * ca(l) - a2 * cb(l)
1157  10 continue
1158  j12 = m0x * sig12 + (trgsum(.true., ssig2, csig2, cb, nc2) -
1159  + trgsum(.true., ssig1, csig1, cb, nc2))
1160  end if
1161  if (redlp) then
1162  m0 = m0x
1163 * Missing a factor of b.
1164 * Add parens around (csig1 * ssig2) and (ssig1 * csig2) to ensure
1165 * accurate cancellation in the case of coincident points.
1166  m12b = dn2 * (csig1 * ssig2) - dn1 * (ssig1 * csig2) -
1167  + csig1 * csig2 * j12
1168  end if
1169  if (scalp) then
1170  csig12 = csig1 * csig2 + ssig1 * ssig2
1171  t = ep2 * (cbet1 - cbet2) * (cbet1 + cbet2) / (dn1 + dn2)
1172  mm12 = csig12 + (t * ssig2 - csig2 * j12) * ssig1 / dn1
1173  mm21 = csig12 - (t * ssig1 - csig1 * j12) * ssig2 / dn2
1174  end if
1175 
1176  return
1177  end
1178 
1179  double precision function astrd(x, y)
1180 * Solve k^4+2*k^3-(x^2+y^2-1)*k^2-2*y^2*k-y^2 = 0 for positive root k.
1181 * This solution is adapted from Geocentric::Reverse.
1182 * input
1183  double precision x, y
1184 
1185  double precision cbrt
1186  double precision k, p, q, r, S, r2, r3, disc, u,
1187  + t3, t, ang, v, uv, w
1188 
1189  p = x**2
1190  q = y**2
1191  r = (p + q - 1) / 6
1192  if ( .not. (q .eq. 0 .and. r .lt. 0) ) then
1193 * Avoid possible division by zero when r = 0 by multiplying equations
1194 * for s and t by r^3 and r, resp.
1195 * S = r^3 * s
1196  s = p * q / 4
1197  r2 = r**2
1198  r3 = r * r2
1199 * The discriminant of the quadratic equation for T3. This is zero on
1200 * the evolute curve p^(1/3)+q^(1/3) = 1
1201  disc = s * (s + 2 * r3)
1202  u = r
1203  if (disc .ge. 0) then
1204  t3 = s + r3
1205 * Pick the sign on the sqrt to maximize abs(T3). This minimizes loss
1206 * of precision due to cancellation. The result is unchanged because
1207 * of the way the T is used in definition of u.
1208 * T3 = (r * t)^3
1209  if (t3 .lt. 0) then
1210  disc = -sqrt(disc)
1211  else
1212  disc = sqrt(disc)
1213  end if
1214  t3 = t3 + disc
1215 * N.B. cbrt always returns the real root. cbrt(-8) = -2.
1216 * T = r * t
1217  t = cbrt(t3)
1218 * T can be zero; but then r2 / T -> 0.
1219  if (t .ne. 0) u = u + t + r2 / t
1220  else
1221 * T is complex, but the way u is defined the result is real.
1222  ang = atan2(sqrt(-disc), -(s + r3))
1223 * There are three possible cube roots. We choose the root which
1224 * avoids cancellation. Note that disc < 0 implies that r < 0.
1225  u = u + 2 * r * cos(ang / 3)
1226  end if
1227 * guaranteed positive
1228  v = sqrt(u**2 + q)
1229 * Avoid loss of accuracy when u < 0.
1230 * u+v, guaranteed positive
1231  if (u .lt. 0) then
1232  uv = q / (v - u)
1233  else
1234  uv = u + v
1235  end if
1236 * positive?
1237  w = (uv - q) / (2 * v)
1238 * Rearrange expression for k to avoid loss of accuracy due to
1239 * subtraction. Division by 0 not possible because uv > 0, w >= 0.
1240 * guaranteed positive
1241  k = uv / (sqrt(uv + w**2) + w)
1242  else
1243 * q == 0 && r <= 0
1244 * y = 0 with |x| <= 1. Handle this case directly.
1245 * for y small, positive root is k = abs(y)/sqrt(1-x^2)
1246  k = 0
1247  end if
1248  astrd = k
1249 
1250  return
1251  end
1252 
1253  double precision function invsta(sbet1, cbet1, dn1,
1254  + sbet2, cbet2, dn2, lam12, slam12, clam12, f, A3x,
1255  + salp1, calp1, salp2, calp2, dnm,
1256  + Ca)
1257 * Return a starting point for Newton's method in salp1 and calp1
1258 * (function value is -1). If Newton's method doesn't need to be used,
1259 * return also salp2, calp2, and dnm and function value is sig12.
1260 * input
1261  double precision sbet1, cbet1, dn1, sbet2, cbet2, dn2,
1262  + lam12, slam12, clam12, f, A3x(*)
1263 * output
1264  double precision salp1, calp1, salp2, calp2, dnm
1265 * temporary
1266  double precision Ca(*)
1267 
1268  double precision hypotx, A3f, Astrd
1269  logical shortp
1270  double precision f1, e2, ep2, n, etol2, k2, eps, sig12,
1271  + sbet12, cbet12, sbt12a, omg12, somg12, comg12, ssig12, csig12,
1272  + x, y, lamscl, betscl, cbt12a, bt12a, m12b, m0, dummy,
1273  + k, omg12a, sbetm2, lam12x
1274 
1275  double precision dblmin, dbleps, pi, degree, tiny,
1276  + tol0, tol1, tol2, tolb, xthrsh
1277  integer digits, maxit1, maxit2
1278  logical init
1279  common /geocom/ dblmin, dbleps, pi, degree, tiny,
1280  + tol0, tol1, tol2, tolb, xthrsh, digits, maxit1, maxit2, init
1281 
1282  f1 = 1 - f
1283  e2 = f * (2 - f)
1284  ep2 = e2 / (1 - e2)
1285  n = f / (2 - f)
1286 * The sig12 threshold for "really short". Using the auxiliary sphere
1287 * solution with dnm computed at (bet1 + bet2) / 2, the relative error in
1288 * the azimuth consistency check is sig12^2 * abs(f) * min(1, 1-f/2) / 2.
1289 * (Error measured for 1/100 < b/a < 100 and abs(f) >= 1/1000. For a
1290 * given f and sig12, the max error occurs for lines near the pole. If
1291 * the old rule for computing dnm = (dn1 + dn2)/2 is used, then the error
1292 * increases by a factor of 2.) Setting this equal to epsilon gives
1293 * sig12 = etol2. Here 0.1 is a safety factor (error decreased by 100)
1294 * and max(0.001, abs(f)) stops etol2 getting too large in the nearly
1295 * spherical case.
1296  etol2 = 0.1d0 * tol2 /
1297  + sqrt( max(0.001d0, abs(f)) * min(1d0, 1 - f/2) / 2 )
1298 
1299 * Return value
1300  sig12 = -1
1301 * bet12 = bet2 - bet1 in [0, pi); bt12a = bet2 + bet1 in (-pi, 0]
1302  sbet12 = sbet2 * cbet1 - cbet2 * sbet1
1303  cbet12 = cbet2 * cbet1 + sbet2 * sbet1
1304  sbt12a = sbet2 * cbet1 + cbet2 * sbet1
1305 
1306  shortp = cbet12 .ge. 0 .and. sbet12 .lt. 0.5d0 .and.
1307  + cbet2 * lam12 .lt. 0.5d0
1308 
1309  if (shortp) then
1310  sbetm2 = (sbet1 + sbet2)**2
1311 * sin((bet1+bet2)/2)^2
1312 * = (sbet1 + sbet2)^2 / ((sbet1 + sbet2)^2 + (cbet1 + cbet2)^2)
1313  sbetm2 = sbetm2 / (sbetm2 + (cbet1 + cbet2)**2)
1314  dnm = sqrt(1 + ep2 * sbetm2)
1315  omg12 = lam12 / (f1 * dnm)
1316  somg12 = sin(omg12)
1317  comg12 = cos(omg12)
1318  else
1319  somg12 = slam12
1320  comg12 = clam12
1321  end if
1322 
1323  salp1 = cbet2 * somg12
1324  if (comg12 .ge. 0) then
1325  calp1 = sbet12 + cbet2 * sbet1 * somg12**2 / (1 + comg12)
1326  else
1327  calp1 = sbt12a - cbet2 * sbet1 * somg12**2 / (1 - comg12)
1328  end if
1329 
1330  ssig12 = hypotx(salp1, calp1)
1331  csig12 = sbet1 * sbet2 + cbet1 * cbet2 * comg12
1332 
1333  if (shortp .and. ssig12 .lt. etol2) then
1334 * really short lines
1335  salp2 = cbet1 * somg12
1336  if (comg12 .ge. 0) then
1337  calp2 = somg12**2 / (1 + comg12)
1338  else
1339  calp2 = 1 - comg12
1340  end if
1341  calp2 = sbet12 - cbet1 * sbet2 * calp2
1342  call norm2x(salp2, calp2)
1343 * Set return value
1344  sig12 = atan2(ssig12, csig12)
1345  else if (abs(n) .gt. 0.1d0 .or. csig12 .ge. 0 .or.
1346  + ssig12 .ge. 6 * abs(n) * pi * cbet1**2) then
1347 * Nothing to do, zeroth order spherical approximation is OK
1348  continue
1349  else
1350 * lam12 - pi
1351  lam12x = atan2(-slam12, -clam12)
1352 * Scale lam12 and bet2 to x, y coordinate system where antipodal point
1353 * is at origin and singular point is at y = 0, x = -1.
1354  if (f .ge. 0) then
1355 * x = dlong, y = dlat
1356  k2 = sbet1**2 * ep2
1357  eps = k2 / (2 * (1 + sqrt(1 + k2)) + k2)
1358  lamscl = f * cbet1 * a3f(eps, a3x) * pi
1359  betscl = lamscl * cbet1
1360  x = lam12x / lamscl
1361  y = sbt12a / betscl
1362  else
1363 * f < 0: x = dlat, y = dlong
1364  cbt12a = cbet2 * cbet1 - sbet2 * sbet1
1365  bt12a = atan2(sbt12a, cbt12a)
1366 * In the case of lon12 = 180, this repeats a calculation made in
1367 * Inverse.
1368  call lengs(n, pi + bt12a,
1369  + sbet1, -cbet1, dn1, sbet2, cbet2, dn2, cbet1, cbet2, 2,
1370  + dummy, m12b, m0, dummy, dummy, ep2, ca)
1371  x = -1 + m12b / (cbet1 * cbet2 * m0 * pi)
1372  if (x .lt. -0.01d0) then
1373  betscl = sbt12a / x
1374  else
1375  betscl = -f * cbet1**2 * pi
1376  end if
1377  lamscl = betscl / cbet1
1378  y = lam12x / lamscl
1379  end if
1380 
1381  if (y .gt. -tol1 .and. x .gt. -1 - xthrsh) then
1382 * strip near cut
1383  if (f .ge. 0) then
1384  salp1 = min(1d0, -x)
1385  calp1 = - sqrt(1 - salp1**2)
1386  else
1387  if (x .gt. -tol1) then
1388  calp1 = 0
1389  else
1390  calp1 = 1
1391  end if
1392  calp1 = max(calp1, x)
1393  salp1 = sqrt(1 - calp1**2)
1394  end if
1395  else
1396 * Estimate alp1, by solving the astroid problem.
1397 *
1398 * Could estimate alpha1 = theta + pi/2, directly, i.e.,
1399 * calp1 = y/k; salp1 = -x/(1+k); for f >= 0
1400 * calp1 = x/(1+k); salp1 = -y/k; for f < 0 (need to check)
1401 *
1402 * However, it's better to estimate omg12 from astroid and use
1403 * spherical formula to compute alp1. This reduces the mean number of
1404 * Newton iterations for astroid cases from 2.24 (min 0, max 6) to 2.12
1405 * (min 0 max 5). The changes in the number of iterations are as
1406 * follows:
1407 *
1408 * change percent
1409 * 1 5
1410 * 0 78
1411 * -1 16
1412 * -2 0.6
1413 * -3 0.04
1414 * -4 0.002
1415 *
1416 * The histogram of iterations is (m = number of iterations estimating
1417 * alp1 directly, n = number of iterations estimating via omg12, total
1418 * number of trials = 148605):
1419 *
1420 * iter m n
1421 * 0 148 186
1422 * 1 13046 13845
1423 * 2 93315 102225
1424 * 3 36189 32341
1425 * 4 5396 7
1426 * 5 455 1
1427 * 6 56 0
1428 *
1429 * Because omg12 is near pi, estimate work with omg12a = pi - omg12
1430  k = astrd(x, y)
1431  if (f .ge. 0) then
1432  omg12a = -x * k/(1 + k)
1433  else
1434  omg12a = -y * (1 + k)/k
1435  end if
1436  omg12a = lamscl * omg12a
1437  somg12 = sin(omg12a)
1438  comg12 = -cos(omg12a)
1439 * Update spherical estimate of alp1 using omg12 instead of lam12
1440  salp1 = cbet2 * somg12
1441  calp1 = sbt12a - cbet2 * sbet1 * somg12**2 / (1 - comg12)
1442  end if
1443  end if
1444 * Sanity check on starting guess. Backwards check allows NaN through.
1445  if (.not. (salp1 .le. 0)) then
1446  call norm2x(salp1, calp1)
1447  else
1448  salp1 = 1
1449  calp1 = 0
1450  end if
1451  invsta = sig12
1452 
1453  return
1454  end
1455 
1456  double precision function lam12f(sbet1, cbet1, dn1,
1457  + sbet2, cbet2, dn2, salp1, calp1, slm120, clm120, f, A3x, C3x,
1458  + salp2, calp2, sig12, ssig1, csig1, ssig2, csig2, eps,
1459  + domg12, diffp, dlam12, Ca)
1460 * input
1461  double precision sbet1, cbet1, dn1, sbet2, cbet2, dn2,
1462  + salp1, calp1, slm120, clm120, f, A3x(*), C3x(*)
1463  logical diffp
1464 * output
1465  double precision salp2, calp2, sig12, ssig1, csig1, ssig2, csig2,
1466  + eps, domg12
1467 * optional output
1468  double precision dlam12
1469 * temporary
1470  double precision Ca(*)
1471 
1472  integer ord, nC3
1473  parameter(ord = 6, nc3 = ord)
1474 
1475  double precision hypotx, A3f, TrgSum
1476 
1477  double precision f1, e2, ep2, salp0, calp0,
1478  + somg1, comg1, somg2, comg2, somg12, comg12,
1479  + lam12, eta, b312, k2, dummy
1480 
1481  double precision dblmin, dbleps, pi, degree, tiny,
1482  + tol0, tol1, tol2, tolb, xthrsh
1483  integer digits, maxit1, maxit2
1484  logical init
1485  common /geocom/ dblmin, dbleps, pi, degree, tiny,
1486  + tol0, tol1, tol2, tolb, xthrsh, digits, maxit1, maxit2, init
1487 
1488  f1 = 1 - f
1489  e2 = f * (2 - f)
1490  ep2 = e2 / (1 - e2)
1491 * Break degeneracy of equatorial line. This case has already been
1492 * handled.
1493  if (sbet1 .eq. 0 .and. calp1 .eq. 0) calp1 = -tiny
1494 
1495 * sin(alp1) * cos(bet1) = sin(alp0)
1496  salp0 = salp1 * cbet1
1497 * calp0 > 0
1498  calp0 = hypotx(calp1, salp1 * sbet1)
1499 
1500 * tan(bet1) = tan(sig1) * cos(alp1)
1501 * tan(omg1) = sin(alp0) * tan(sig1) = tan(omg1)=tan(alp1)*sin(bet1)
1502  ssig1 = sbet1
1503  somg1 = salp0 * sbet1
1504  csig1 = calp1 * cbet1
1505  comg1 = csig1
1506  call norm2x(ssig1, csig1)
1507 * norm2x(somg1, comg1); -- don't need to normalize!
1508 
1509 * Enforce symmetries in the case abs(bet2) = -bet1. Need to be careful
1510 * about this case, since this can yield singularities in the Newton
1511 * iteration.
1512 * sin(alp2) * cos(bet2) = sin(alp0)
1513  if (cbet2 .ne. cbet1) then
1514  salp2 = salp0 / cbet2
1515  else
1516  salp2 = salp1
1517  end if
1518 * calp2 = sqrt(1 - sq(salp2))
1519 * = sqrt(sq(calp0) - sq(sbet2)) / cbet2
1520 * and subst for calp0 and rearrange to give (choose positive sqrt
1521 * to give alp2 in [0, pi/2]).
1522  if (cbet2 .ne. cbet1 .or. abs(sbet2) .ne. -sbet1) then
1523  if (cbet1 .lt. -sbet1) then
1524  calp2 = (cbet2 - cbet1) * (cbet1 + cbet2)
1525  else
1526  calp2 = (sbet1 - sbet2) * (sbet1 + sbet2)
1527  end if
1528  calp2 = sqrt((calp1 * cbet1)**2 + calp2) / cbet2
1529  else
1530  calp2 = abs(calp1)
1531  end if
1532 * tan(bet2) = tan(sig2) * cos(alp2)
1533 * tan(omg2) = sin(alp0) * tan(sig2).
1534  ssig2 = sbet2
1535  somg2 = salp0 * sbet2
1536  csig2 = calp2 * cbet2
1537  comg2 = csig2
1538  call norm2x(ssig2, csig2)
1539 * norm2x(somg2, comg2); -- don't need to normalize!
1540 
1541 * sig12 = sig2 - sig1, limit to [0, pi]
1542  sig12 = atan2(0d0 + max(0d0, csig1 * ssig2 - ssig1 * csig2),
1543  + csig1 * csig2 + ssig1 * ssig2)
1544 
1545 * omg12 = omg2 - omg1, limit to [0, pi]
1546  somg12 = 0d0 + max(0d0, comg1 * somg2 - somg1 * comg2)
1547  comg12 = comg1 * comg2 + somg1 * somg2
1548 * eta = omg12 - lam120
1549  eta = atan2(somg12 * clm120 - comg12 * slm120,
1550  + comg12 * clm120 + somg12 * slm120)
1551  k2 = calp0**2 * ep2
1552  eps = k2 / (2 * (1 + sqrt(1 + k2)) + k2)
1553  call c3f(eps, c3x, ca)
1554  b312 = (trgsum(.true., ssig2, csig2, ca, nc3-1) -
1555  + trgsum(.true., ssig1, csig1, ca, nc3-1))
1556  domg12 = -f * a3f(eps, a3x) * salp0 * (sig12 + b312)
1557  lam12 = eta + domg12
1558 
1559  if (diffp) then
1560  if (calp2 .eq. 0) then
1561  dlam12 = - 2 * f1 * dn1 / sbet1
1562  else
1563  call lengs(eps, sig12, ssig1, csig1, dn1, ssig2, csig2, dn2,
1564  + cbet1, cbet2, 2,
1565  + dummy, dlam12, dummy, dummy, dummy, ep2, ca)
1566  dlam12 = dlam12 * f1 / (calp2 * cbet2)
1567  end if
1568  end if
1569  lam12f = lam12
1570 
1571  return
1572  end
1573 
1574  double precision function a3f(eps, A3x)
1575 * Evaluate A3
1576  integer ord, nA3, nA3x
1577  parameter(ord = 6, na3 = ord, na3x = na3)
1578 
1579 * input
1580  double precision eps
1581 * output
1582  double precision A3x(0: nA3x-1)
1583 
1584  double precision polval
1585  A3f = polval(na3 - 1, a3x, eps)
1586 
1587  return
1588  end
1589 
1590  subroutine c3f(eps, C3x, c)
1591 * Evaluate C3 coeffs
1592 * Elements c[1] thru c[nC3-1] are set
1593  integer ord, nC3, nC3x
1594  parameter(ord = 6, nc3 = ord, nc3x = (nc3 * (nc3 - 1)) / 2)
1595 
1596 * input
1597  double precision eps, C3x(0:nC3x-1)
1598 * output
1599  double precision c(nC3-1)
1600 
1601  integer o, m, l
1602  double precision mult, polval
1603 
1604  mult = 1
1605  o = 0
1606  do 10 l = 1, nc3 - 1
1607  m = nc3 - l - 1
1608  mult = mult * eps
1609  c(l) = mult * polval(m, c3x(o), eps)
1610  o = o + m + 1
1611  10 continue
1612 
1613  return
1614  end
1615 
1616  subroutine c4f(eps, C4x, c)
1617 * Evaluate C4
1618 * Elements c[0] thru c[nC4-1] are set
1619  integer ord, nC4, nC4x
1620  parameter(ord = 6, nc4 = ord, nc4x = (nc4 * (nc4 + 1)) / 2)
1621 
1622 * input
1623  double precision eps, C4x(0:nC4x-1)
1624 *output
1625  double precision c(0:nC4-1)
1626 
1627  integer o, m, l
1628  double precision mult, polval
1629 
1630  mult = 1
1631  o = 0
1632  do 10 l = 0, nc4 - 1
1633  m = nc4 - l - 1
1634  c(l) = mult * polval(m, c4x(o), eps)
1635  o = o + m + 1
1636  mult = mult * eps
1637  10 continue
1638 
1639  return
1640  end
1641 
1642  double precision function a1m1f(eps)
1643 * The scale factor A1-1 = mean value of (d/dsigma)I1 - 1
1644 * input
1645  double precision eps
1646 
1647  double precision t
1648  integer ord, nA1, o, m
1649  parameter(ord = 6, na1 = ord)
1650  double precision polval, coeff(nA1/2 + 2)
1651  data coeff /1, 4, 64, 0, 256/
1652 
1653  o = 1
1654  m = na1/2
1655  t = polval(m, coeff(o), eps**2) / coeff(o + m + 1)
1656  a1m1f = (t + eps) / (1 - eps)
1657 
1658  return
1659  end
1660 
1661  subroutine c1f(eps, c)
1662 * The coefficients C1[l] in the Fourier expansion of B1
1663  integer ord, nC1
1664  parameter(ord = 6, nc1 = ord)
1665 
1666 * input
1667  double precision eps
1668 * output
1669  double precision c(nC1)
1670 
1671  double precision eps2, d
1672  integer o, m, l
1673  double precision polval, coeff((nC1**2 + 7*nC1 - 2*(nC1/2))/4)
1674  data coeff /
1675  + -1, 6, -16, 32,
1676  + -9, 64, -128, 2048,
1677  + 9, -16, 768,
1678  + 3, -5, 512,
1679  + -7, 1280,
1680  + -7, 2048/
1681 
1682  eps2 = eps**2
1683  d = eps
1684  o = 1
1685  do 10 l = 1, nc1
1686  m = (nc1 - l) / 2
1687  c(l) = d * polval(m, coeff(o), eps2) / coeff(o + m + 1)
1688  o = o + m + 2
1689  d = d * eps
1690  10 continue
1691 
1692  return
1693  end
1694 
1695  subroutine c1pf(eps, c)
1696 * The coefficients C1p[l] in the Fourier expansion of B1p
1697  integer ord, nC1p
1698  parameter(ord = 6, nc1p = ord)
1699 
1700 * input
1701  double precision eps
1702 * output
1703  double precision c(nC1p)
1704 
1705  double precision eps2, d
1706  integer o, m, l
1707  double precision polval, coeff((nC1p**2 + 7*nC1p - 2*(nC1p/2))/4)
1708  data coeff /
1709  + 205, -432, 768, 1536,
1710  + 4005, -4736, 3840, 12288,
1711  + -225, 116, 384,
1712  + -7173, 2695, 7680,
1713  + 3467, 7680,
1714  + 38081, 61440/
1715 
1716  eps2 = eps**2
1717  d = eps
1718  o = 1
1719  do 10 l = 1, nc1p
1720  m = (nc1p - l) / 2
1721  c(l) = d * polval(m, coeff(o), eps2) / coeff(o + m + 1)
1722  o = o + m + 2
1723  d = d * eps
1724  10 continue
1725 
1726  return
1727  end
1728 
1729 * The scale factor A2-1 = mean value of (d/dsigma)I2 - 1
1730  double precision function a2m1f(eps)
1731 * input
1732  double precision eps
1733 
1734  double precision t
1735  integer ord, nA2, o, m
1736  parameter(ord = 6, na2 = ord)
1737  double precision polval, coeff(nA2/2 + 2)
1738  data coeff /-11, -28, -192, 0, 256/
1739 
1740  o = 1
1741  m = na2/2
1742  t = polval(m, coeff(o), eps**2) / coeff(o + m + 1)
1743  a2m1f = (t - eps) / (1 + eps)
1744 
1745  return
1746  end
1747 
1748  subroutine c2f(eps, c)
1749 * The coefficients C2[l] in the Fourier expansion of B2
1750  integer ord, nC2
1751  parameter(ord = 6, nc2 = ord)
1752 
1753 * input
1754  double precision eps
1755 * output
1756  double precision c(nC2)
1757 
1758  double precision eps2, d
1759  integer o, m, l
1760  double precision polval, coeff((nC2**2 + 7*nC2 - 2*(nC2/2))/4)
1761  data coeff /
1762  + 1, 2, 16, 32,
1763  + 35, 64, 384, 2048,
1764  + 15, 80, 768,
1765  + 7, 35, 512,
1766  + 63, 1280,
1767  + 77, 2048/
1768 
1769  eps2 = eps**2
1770  d = eps
1771  o = 1
1772  do 10 l = 1, nc2
1773  m = (nc2 - l) / 2
1774  c(l) = d * polval(m, coeff(o), eps2) / coeff(o + m + 1)
1775  o = o + m + 2
1776  d = d * eps
1777  10 continue
1778 
1779  return
1780  end
1781 
1782  subroutine a3cof(n, A3x)
1783 * The scale factor A3 = mean value of (d/dsigma)I3
1784  integer ord, nA3, nA3x
1785  parameter(ord = 6, na3 = ord, na3x = na3)
1786 
1787 * input
1788  double precision n
1789 * output
1790  double precision A3x(0:nA3x-1)
1791 
1792  integer o, m, k, j
1793  double precision polval, coeff((nA3**2 + 7*nA3 - 2*(nA3/2))/4)
1794  data coeff /
1795  + -3, 128,
1796  + -2, -3, 64,
1797  + -1, -3, -1, 16,
1798  + 3, -1, -2, 8,
1799  + 1, -1, 2,
1800  + 1, 1/
1801 
1802  o = 1
1803  k = 0
1804  do 10 j = na3 - 1, 0, -1
1805  m = min(na3 - j - 1, j)
1806  a3x(k) = polval(m, coeff(o), n) / coeff(o + m + 1)
1807  k = k + 1
1808  o = o + m + 2
1809  10 continue
1810 
1811  return
1812  end
1813 
1814  subroutine c3cof(n, C3x)
1815 * The coefficients C3[l] in the Fourier expansion of B3
1816  integer ord, nC3, nC3x
1817  parameter(ord = 6, nc3 = ord, nc3x = (nc3 * (nc3 - 1)) / 2)
1818 
1819 * input
1820  double precision n
1821 * output
1822  double precision C3x(0:nC3x-1)
1823 
1824  integer o, m, l, j, k
1825  double precision polval,
1826  + coeff(((nc3-1)*(nc3**2 + 7*nc3 - 2*(nc3/2)))/8)
1827  data coeff /
1828  + 3, 128,
1829  + 2, 5, 128,
1830  + -1, 3, 3, 64,
1831  + -1, 0, 1, 8,
1832  + -1, 1, 4,
1833  + 5, 256,
1834  + 1, 3, 128,
1835  + -3, -2, 3, 64,
1836  + 1, -3, 2, 32,
1837  + 7, 512,
1838  + -10, 9, 384,
1839  + 5, -9, 5, 192,
1840  + 7, 512,
1841  + -14, 7, 512,
1842  + 21, 2560/
1843 
1844  o = 1
1845  k = 0
1846  do 20 l = 1, nc3 - 1
1847  do 10 j = nc3 - 1, l, -1
1848  m = min(nc3 - j - 1, j)
1849  c3x(k) = polval(m, coeff(o), n) / coeff(o + m + 1)
1850  k = k + 1
1851  o = o + m + 2
1852  10 continue
1853  20 continue
1854 
1855  return
1856  end
1857 
1858  subroutine c4cof(n, C4x)
1859 * The coefficients C4[l] in the Fourier expansion of I4
1860  integer ord, nC4, nC4x
1861  parameter(ord = 6, nc4 = ord, nc4x = (nc4 * (nc4 + 1)) / 2)
1862 
1863 * input
1864  double precision n
1865 * output
1866  double precision C4x(0:nC4x-1)
1867 
1868  integer o, m, l, j, k
1869  double precision polval, coeff((nC4 * (nC4 + 1) * (nC4 + 5)) / 6)
1870  data coeff /
1871  + 97, 15015, 1088, 156, 45045, -224, -4784, 1573, 45045,
1872  + -10656, 14144, -4576, -858, 45045,
1873  + 64, 624, -4576, 6864, -3003, 15015,
1874  + 100, 208, 572, 3432, -12012, 30030, 45045,
1875  + 1, 9009, -2944, 468, 135135, 5792, 1040, -1287, 135135,
1876  + 5952, -11648, 9152, -2574, 135135,
1877  + -64, -624, 4576, -6864, 3003, 135135,
1878  + 8, 10725, 1856, -936, 225225, -8448, 4992, -1144, 225225,
1879  + -1440, 4160, -4576, 1716, 225225,
1880  + -136, 63063, 1024, -208, 105105,
1881  + 3584, -3328, 1144, 315315,
1882  + -128, 135135, -2560, 832, 405405, 128, 99099/
1883 
1884  o = 1
1885  k = 0
1886  do 20 l = 0, nc4 - 1
1887  do 10 j = nc4 - 1, l, -1
1888  m = nc4 - j - 1
1889  c4x(k) = polval(m, coeff(o), n) / coeff(o + m + 1)
1890  k = k + 1
1891  o = o + m + 2
1892  10 continue
1893  20 continue
1894 
1895  return
1896  end
1897 
1898  double precision function sumx(u, v, t)
1899 * input
1900  double precision u, v
1901 * output
1902  double precision t
1903 
1904  double precision up, vpp
1905  sumx = u + v
1906  up = sumx - v
1907  vpp = sumx - up
1908  up = up - u
1909  vpp = vpp - v
1910  t = -(up + vpp)
1911 
1912  return
1913  end
1914 
1915  double precision function remx(x, y)
1916 * the remainder function but not worrying how ties are handled
1917 * y must be positive
1918 * input
1919  double precision x, y
1920 
1921  remx = mod(x, y)
1922  if (remx .lt. -y/2) then
1923  remx = remx + y
1924  else if (remx .gt. +y/2) then
1925  remx = remx - y
1926  end if
1927 
1928  return
1929  end
1930 
1931  double precision function angnm(x)
1932 * input
1933  double precision x
1934 
1935  double precision remx
1936  angnm = remx(x, 360d0)
1937  if (angnm .eq. -180) then
1938  angnm = 180
1939  end if
1940 
1941  return
1942  end
1943 
1944  double precision function latfix(x)
1945 * input
1946  double precision x
1947 
1948  latfix = x
1949  if (.not. (abs(x) .gt. 90)) return
1950 * concoct a NaN
1951  latfix = sqrt(90 - abs(x))
1952 
1953  return
1954  end
1955 
1956  double precision function angdif(x, y, e)
1957 * Compute y - x. x and y must both lie in [-180, 180]. The result is
1958 * equivalent to computing the difference exactly, reducing it to (-180,
1959 * 180] and rounding the result. Note that this prescription allows -180
1960 * to be returned (e.g., if x is tiny and negative and y = 180). The
1961 * error in the difference is returned in e
1962 * input
1963  double precision x, y
1964 * output
1965  double precision e
1966 
1967  double precision d, t, sumx, AngNm
1968  d = angnm(sumx(angnm(-x), angnm(y), t))
1969  if (d .eq. 180 .and. t .gt. 0) then
1970  d = -180
1971  end if
1972  angdif = sumx(d, t, e)
1973 
1974  return
1975  end
1976 
1977  double precision function angrnd(x)
1978 * The makes the smallest gap in x = 1/16 - nextafter(1/16, 0) = 1/2^57
1979 * for reals = 0.7 pm on the earth if x is an angle in degrees. (This is
1980 * about 1000 times more resolution than we get with angles around 90
1981 * degrees.) We use this to avoid having to deal with near singular
1982 * cases when x is non-zero but tiny (e.g., 1.0e-200). This also
1983 * converts -0 to +0.
1984 * input
1985  double precision x
1986 
1987  double precision y, z
1988  z = 1/16d0
1989  y = abs(x)
1990 * The compiler mustn't "simplify" z - (z - y) to y
1991  if (y .lt. z) y = z - (z - y)
1992  angrnd = sign(y, x)
1993  if (x .eq. 0) angrnd = 0
1994 
1995  return
1996  end
1997 
1998  subroutine swap(x, y)
1999 * input/output
2000  double precision x, y
2001 
2002  double precision z
2003  z = x
2004  x = y
2005  y = z
2006 
2007  return
2008  end
2009 
2010  double precision function hypotx(x, y)
2011 * input
2012  double precision x, y
2013 
2014 * With Fortran 2008, this becomes: hypotx = hypot(x, y)
2015  hypotx = sqrt(x**2 + y**2)
2016 
2017  return
2018  end
2019 
2020  subroutine norm2x(x, y)
2021 * input/output
2022  double precision x, y
2023 
2024  double precision hypotx, r
2025  r = hypotx(x, y)
2026  x = x/r
2027  y = y/r
2028 
2029  return
2030  end
2031 
2032  double precision function log1px(x)
2033 * input
2034  double precision x
2035 
2036  double precision y, z
2037  y = 1 + x
2038  z = y - 1
2039  if (z .eq. 0) then
2040  log1px = x
2041  else
2042  log1px = x * log(y) / z
2043  end if
2044 
2045  return
2046  end
2047 
2048  double precision function atanhx(x)
2049 * input
2050  double precision x
2051 
2052 * With Fortran 2008, this becomes: atanhx = atanh(x)
2053  double precision log1px, y
2054  y = abs(x)
2055  y = log1px(2 * y/(1 - y))/2
2056  atanhx = sign(y, x)
2057 
2058  return
2059  end
2060 
2061  double precision function cbrt(x)
2062 * input
2063  double precision x
2064 
2065  cbrt = sign(abs(x)**(1/3d0), x)
2066 
2067  return
2068  end
2069 
2070  double precision function trgsum(sinp, sinx, cosx, c, n)
2071 * Evaluate
2072 * y = sinp ? sum(c[i] * sin( 2*i * x), i, 1, n) :
2073 * sum(c[i] * cos((2*i-1) * x), i, 1, n)
2074 * using Clenshaw summation.
2075 * Approx operation count = (n + 5) mult and (2 * n + 2) add
2076 * input
2077  logical sinp
2078  integer n
2079  double precision sinx, cosx, c(n)
2080 
2081  double precision ar, y0, y1
2082  integer n2, k
2083 
2084 * 2 * cos(2 * x)
2085  ar = 2 * (cosx - sinx) * (cosx + sinx)
2086 * accumulators for sum
2087  if (mod(n, 2) .eq. 1) then
2088  y0 = c(n)
2089  n2 = n - 1
2090  else
2091  y0 = 0
2092  n2 = n
2093  end if
2094  y1 = 0
2095 * Now n2 is even
2096  do 10 k = n2, 2, -2
2097 * Unroll loop x 2, so accumulators return to their original role
2098  y1 = ar * y0 - y1 + c(k)
2099  y0 = ar * y1 - y0 + c(k-1)
2100  10 continue
2101  if (sinp) then
2102 * sin(2 * x) * y0
2103  trgsum = 2 * sinx * cosx * y0
2104  else
2105 * cos(x) * (y0 - y1)
2106  trgsum = cosx * (y0 - y1)
2107  end if
2108 
2109  return
2110  end
2111 
2112  integer function trnsit(lon1, lon2)
2113 * input
2114  double precision lon1, lon2
2115 
2116  double precision lon1x, lon2x, lon12, AngNm, AngDif, e
2117  lon1x = angnm(lon1)
2118  lon2x = angnm(lon2)
2119  lon12 = angdif(lon1x, lon2x, e)
2120  trnsit = 0
2121  if (lon1x .le. 0 .and. lon2x .gt. 0 .and. lon12 .gt. 0) then
2122  trnsit = 1
2123  else if (lon2x .le. 0 .and. lon1x .gt. 0 .and. lon12 .lt. 0) then
2124  trnsit = -1
2125  end if
2126 
2127  return
2128  end
2129 
2130  subroutine accini(s)
2131 * Initialize an accumulator; this is an array with two elements.
2132 * input/output
2133  double precision s(2)
2134 
2135  s(1) = 0
2136  s(2) = 0
2137 
2138  return
2139  end
2140 
2141  subroutine accadd(s, y)
2142 * Add y to an accumulator.
2143 * input
2144  double precision y
2145 * input/output
2146  double precision s(2)
2147 
2148  double precision z, u, sumx
2149  z = sumx(y, s(2), u)
2150  s(1) = sumx(z, s(1), s(2))
2151  if (s(1) .eq. 0) then
2152  s(1) = u
2153  else
2154  s(2) = s(2) + u
2155  end if
2156 
2157  return
2158  end
2159 
2160  subroutine accrem(s, y)
2161 * Reduce s to [-y/2, y/2].
2162 * input
2163  double precision y
2164 * input/output
2165  double precision s(2)
2166 
2167  double precision remx
2168  s(1) = remx(s(1), y)
2169  call accadd(s, 0d0)
2170 
2171  return
2172  end
2173 
2174  subroutine sncsdx(x, sinx, cosx)
2175 * Compute sin(x) and cos(x) with x in degrees
2176 * input
2177  double precision x
2178 * input/output
2179  double precision sinx, cosx
2180 
2181  double precision dblmin, dbleps, pi, degree, tiny,
2182  + tol0, tol1, tol2, tolb, xthrsh
2183  integer digits, maxit1, maxit2
2184  logical init
2185  common /geocom/ dblmin, dbleps, pi, degree, tiny,
2186  + tol0, tol1, tol2, tolb, xthrsh, digits, maxit1, maxit2, init
2187 
2188  double precision r, s, c
2189  integer q
2190  r = mod(x, 360d0)
2191  q = nint(r / 90)
2192  r = (r - 90 * q) * degree
2193  s = sin(r)
2194 * sin(-0) isn't reliably -0, so...
2195  if (x .eq. 0) s = x
2196  c = cos(r)
2197  q = mod(q + 4, 4)
2198  if (q .eq. 0) then
2199  sinx = s
2200  cosx = c
2201  else if (q .eq. 1) then
2202  sinx = c
2203  cosx = -s
2204  else if (q .eq. 2) then
2205  sinx = -s
2206  cosx = -c
2207  else
2208 * q.eq.3
2209  sinx = -c
2210  cosx = s
2211  end if
2212 
2213  if (x .ne. 0) then
2214  sinx = 0d0 + sinx
2215  cosx = 0d0 + cosx
2216  end if
2217 
2218  return
2219  end
2220 
2221  double precision function atn2dx(y, x)
2222 * input
2223  double precision x, y
2224 
2225  double precision dblmin, dbleps, pi, degree, tiny,
2226  + tol0, tol1, tol2, tolb, xthrsh
2227  integer digits, maxit1, maxit2
2228  logical init
2229  common /geocom/ dblmin, dbleps, pi, degree, tiny,
2230  + tol0, tol1, tol2, tolb, xthrsh, digits, maxit1, maxit2, init
2231 
2232  double precision xx, yy
2233  integer q
2234  if (abs(y) .gt. abs(x)) then
2235  xx = y
2236  yy = x
2237  q = 2
2238  else
2239  xx = x
2240  yy = y
2241  q = 0
2242  end if
2243  if (xx .lt. 0) then
2244  xx = -xx
2245  q = q + 1
2246  end if
2247  atn2dx = atan2(yy, xx) / degree
2248  if (q .eq. 1) then
2249  if (yy .ge. 0) then
2250  atn2dx = 180 - atn2dx
2251  else
2252  atn2dx = -180 - atn2dx
2253  end if
2254  else if (q .eq. 2) then
2255  atn2dx = 90 - atn2dx
2256  else if (q .eq. 3) then
2257  atn2dx = -90 + atn2dx
2258  end if
2259 
2260  return
2261  end
2262 
2263  double precision function polval(N, p, x)
2264 * input
2265  integer N
2266  double precision p(0:N), x
2267 
2268  integer i
2269  if (N .lt. 0) then
2270  polval = 0
2271  else
2272  polval = p(0)
2273  end if
2274  do 10 i = 1, n
2275  polval = polval * x + p(i)
2276  10 continue
2277 
2278  return
2279  end
2280 
2281 * Table of name abbreviations to conform to the 6-char limit and
2282 * potential name conflicts.
2283 * A3coeff A3cof
2284 * C3coeff C3cof
2285 * C4coeff C4cof
2286 * AngNormalize AngNm
2287 * AngDiff AngDif
2288 * AngRound AngRnd
2289 * arcmode arcmod
2290 * Astroid Astrd
2291 * betscale betscl
2292 * lamscale lamscl
2293 * cbet12a cbt12a
2294 * sbet12a sbt12a
2295 * epsilon dbleps
2296 * realmin dblmin
2297 * geodesic geod
2298 * inverse invers
2299 * InverseStart InvSta
2300 * Lambda12 Lam12f
2301 * latsign latsgn
2302 * lonsign lonsgn
2303 * Lengths Lengs
2304 * meridian merid
2305 * outmask omask
2306 * shortline shortp
2307 * norm norm2x
2308 * SinCosSeries TrgSum
2309 * xthresh xthrsh
2310 * transit trnsit
2311 * polyval polval
2312 * LONG_UNROLL unroll
2313 * sincosd sncsdx
2314 * atan2d atn2dx
2315 *> @endcond SKIP