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