GeographicLib 2.5
Geodesic.cpp
Go to the documentation of this file.
1/**
2 * \file Geodesic.cpp
3 * \brief Implementation for GeographicLib::Geodesic class
4 *
5 * Copyright (c) Charles Karney (2009-2023) <karney@alum.mit.edu> and licensed
6 * under the MIT/X11 License. For more information, see
7 * https://geographiclib.sourceforge.io/
8 *
9 * This is a reformulation of the geodesic problem. The notation is as
10 * follows:
11 * - at a general point (no suffix or 1 or 2 as suffix)
12 * - phi = latitude
13 * - beta = latitude on auxiliary sphere
14 * - omega = longitude on auxiliary sphere
15 * - lambda = longitude
16 * - alpha = azimuth of great circle
17 * - sigma = arc length along great circle
18 * - s = distance
19 * - tau = scaled distance (= sigma at multiples of pi/2)
20 * - at northwards equator crossing
21 * - beta = phi = 0
22 * - omega = lambda = 0
23 * - alpha = alpha0
24 * - sigma = s = 0
25 * - a 12 suffix means a difference, e.g., s12 = s2 - s1.
26 * - s and c prefixes mean sin and cos
27 **********************************************************************/
28
31
32#if defined(_MSC_VER)
33// Squelch warnings about potentially uninitialized local variables
34# pragma warning (disable: 4701)
35#endif
36
37namespace GeographicLib {
38
39 using namespace std;
40
41 Geodesic::Geodesic(real a, real f, bool exact)
42 : maxit2_(maxit1_ + Math::digits() + 10)
43 // Underflow guard. We require
44 // tiny_ * epsilon() > 0
45 // tiny_ + epsilon() == epsilon()
46 , tiny_(sqrt(numeric_limits<real>::min()))
47 , tol0_(numeric_limits<real>::epsilon())
48 // Increase multiplier in defn of tol1_ from 100 to 200 to fix inverse
49 // case 52.784459512564 0 -52.784459512563990912 179.634407464943777557
50 // which otherwise failed for Visual Studio 10 (Release and Debug)
51 , tol1_(200 * tol0_)
52 , tol2_(sqrt(tol0_))
53 , tolb_(tol0_) // Check on bisection interval
54 , xthresh_(1000 * tol2_)
55 , _a(a)
56 , _f(f)
57 , _exact(exact)
58 , _f1(1 - _f)
59 , _e2(_f * (2 - _f))
60 , _ep2(_e2 / Math::sq(_f1)) // e2 / (1 - e2)
61 , _n(_f / ( 2 - _f))
62 , _b(_a * _f1)
63 , _c2((Math::sq(_a) + Math::sq(_b) *
64 (_e2 == 0 ? 1 :
65 Math::eatanhe(real(1), (_f < 0 ? -1 : 1) * sqrt(fabs(_e2))) / _e2))
66 / 2) // authalic radius squared
67 // The sig12 threshold for "really short". Using the auxiliary sphere
68 // solution with dnm computed at (bet1 + bet2) / 2, the relative error in
69 // the azimuth consistency check is sig12^2 * abs(f) * min(1, 1-f/2) / 2.
70 // (Error measured for 1/100 < b/a < 100 and abs(f) >= 1/1000. For a
71 // given f and sig12, the max error occurs for lines near the pole. If
72 // the old rule for computing dnm = (dn1 + dn2)/2 is used, then the error
73 // increases by a factor of 2.) Setting this equal to epsilon gives
74 // sig12 = etol2. Here 0.1 is a safety factor (error decreased by 100)
75 // and max(0.001, abs(f)) stops etol2 getting too large in the nearly
76 // spherical case.
77 , _etol2(real(0.1) * tol2_ /
78 sqrt( fmax(real(0.001), fabs(_f)) * fmin(real(1), 1 - _f/2) / 2 ))
79 , _geodexact(_exact ? GeodesicExact(a, f) : GeodesicExact())
80 {
81 if (_exact)
82 _c2 = _geodexact._c2;
83 else {
84 if (!(isfinite(_a) && _a > 0))
85 throw GeographicErr("Equatorial radius is not positive");
86 if (!(isfinite(_b) && _b > 0))
87 throw GeographicErr("Polar semi-axis is not positive");
88 A3coeff();
89 C3coeff();
90 C4coeff();
91 }
92 }
93
95 static const Geodesic wgs84(Constants::WGS84_a(), Constants::WGS84_f());
96 return wgs84;
97 }
98
99 Math::real Geodesic::SinCosSeries(bool sinp,
100 real sinx, real cosx,
101 const real c[], int n) {
102 // Evaluate
103 // y = sinp ? sum(c[i] * sin( 2*i * x), i, 1, n) :
104 // sum(c[i] * cos((2*i+1) * x), i, 0, n-1)
105 // using Clenshaw summation. N.B. c[0] is unused for sin series
106 // Approx operation count = (n + 5) mult and (2 * n + 2) add
107 c += (n + sinp); // Point to one beyond last element
108 real
109 ar = 2 * (cosx - sinx) * (cosx + sinx), // 2 * cos(2 * x)
110 y0 = n & 1 ? *--c : 0, y1 = 0; // accumulators for sum
111 // Now n is even
112 n /= 2;
113 while (n--) {
114 // Unroll loop x 2, so accumulators return to their original role
115 y1 = ar * y0 - y1 + *--c;
116 y0 = ar * y1 - y0 + *--c;
117 }
118 return sinp
119 ? 2 * sinx * cosx * y0 // sin(2 * x) * y0
120 : cosx * (y0 - y1); // cos(x) * (y0 - y1)
121 }
122
123 GeodesicLine Geodesic::Line(real lat1, real lon1, real azi1,
124 unsigned caps) const {
125 return GeodesicLine(*this, lat1, lon1, azi1, caps);
126 }
127
128 Math::real Geodesic::GenDirect(real lat1, real lon1, real azi1,
129 bool arcmode, real s12_a12, unsigned outmask,
130 real& lat2, real& lon2, real& azi2,
131 real& s12, real& m12, real& M12, real& M21,
132 real& S12) const {
133 if (_exact)
134 return _geodexact.GenDirect(lat1, lon1, azi1, arcmode, s12_a12, outmask,
135 lat2, lon2, azi2,
136 s12, m12, M12, M21, S12);
137 // Automatically supply DISTANCE_IN if necessary
138 if (!arcmode) outmask |= DISTANCE_IN;
139 return GeodesicLine(*this, lat1, lon1, azi1, outmask)
140 . // Note the dot!
141 GenPosition(arcmode, s12_a12, outmask,
142 lat2, lon2, azi2, s12, m12, M12, M21, S12);
143 }
144
145 GeodesicLine Geodesic::GenDirectLine(real lat1, real lon1, real azi1,
146 bool arcmode, real s12_a12,
147 unsigned caps) const {
148 azi1 = Math::AngNormalize(azi1);
149 real salp1, calp1;
150 // Guard against underflow in salp0. Also -0 is converted to +0.
151 Math::sincosd(Math::AngRound(azi1), salp1, calp1);
152 // Automatically supply DISTANCE_IN if necessary
153 if (!arcmode) caps |= DISTANCE_IN;
154 return GeodesicLine(*this, lat1, lon1, azi1, salp1, calp1,
155 caps, arcmode, s12_a12);
156 }
157
158 GeodesicLine Geodesic::DirectLine(real lat1, real lon1, real azi1, real s12,
159 unsigned caps) const {
160 return GenDirectLine(lat1, lon1, azi1, false, s12, caps);
161 }
162
163 GeodesicLine Geodesic::ArcDirectLine(real lat1, real lon1, real azi1,
164 real a12, unsigned caps) const {
165 return GenDirectLine(lat1, lon1, azi1, true, a12, caps);
166 }
167
168 Math::real Geodesic::GenInverse(real lat1, real lon1, real lat2, real lon2,
169 unsigned outmask, real& s12,
170 real& salp1, real& calp1,
171 real& salp2, real& calp2,
172 real& m12, real& M12, real& M21,
173 real& S12) const {
174 if (_exact)
175 return _geodexact.GenInverse(lat1, lon1, lat2, lon2,
176 outmask, s12,
177 salp1, calp1, salp2, calp2,
178 m12, M12, M21, S12);
179 // Compute longitude difference (AngDiff does this carefully).
180 using std::isnan; // Needed for Centos 7, ubuntu 14
181 real lon12s, lon12 = Math::AngDiff(lon1, lon2, lon12s);
182 // Make longitude difference positive.
183 int lonsign = signbit(lon12) ? -1 : 1;
184 lon12 *= lonsign; lon12s *= lonsign;
185 real
186 lam12 = lon12 * Math::degree(),
187 slam12, clam12;
188 // Calculate sincos of lon12 + error (this applies AngRound internally).
189 Math::sincosde(lon12, lon12s, slam12, clam12);
190 // the supplementary longitude difference
191 lon12s = (Math::hd - lon12) - lon12s;
192
193 // If really close to the equator, treat as on equator.
194 lat1 = Math::AngRound(Math::LatFix(lat1));
195 lat2 = Math::AngRound(Math::LatFix(lat2));
196 // Swap points so that point with higher (abs) latitude is point 1.
197 // If one latitude is a nan, then it becomes lat1.
198 int swapp = fabs(lat1) < fabs(lat2) || isnan(lat2) ? -1 : 1;
199 if (swapp < 0) {
200 lonsign *= -1;
201 swap(lat1, lat2);
202 }
203 // Make lat1 <= -0
204 int latsign = signbit(lat1) ? 1 : -1;
205 lat1 *= latsign;
206 lat2 *= latsign;
207 // Now we have
208 //
209 // 0 <= lon12 <= 180
210 // -90 <= lat1 <= -0
211 // lat1 <= lat2 <= -lat1
212 //
213 // longsign, swapp, latsign register the transformation to bring the
214 // coordinates to this canonical form. In all cases, 1 means no change was
215 // made. We make these transformations so that there are few cases to
216 // check, e.g., on verifying quadrants in atan2. In addition, this
217 // enforces some symmetries in the results returned.
218
219 real sbet1, cbet1, sbet2, cbet2, s12x, m12x = Math::NaN();
220
221 Math::sincosd(lat1, sbet1, cbet1); sbet1 *= _f1;
222 // Ensure cbet1 = +epsilon at poles; doing the fix on beta means that sig12
223 // will be <= 2*tiny for two points at the same pole.
224 Math::norm(sbet1, cbet1); cbet1 = fmax(tiny_, cbet1);
225
226 Math::sincosd(lat2, sbet2, cbet2); sbet2 *= _f1;
227 // Ensure cbet2 = +epsilon at poles
228 Math::norm(sbet2, cbet2); cbet2 = fmax(tiny_, cbet2);
229
230 // If cbet1 < -sbet1, then cbet2 - cbet1 is a sensitive measure of the
231 // |bet1| - |bet2|. Alternatively (cbet1 >= -sbet1), abs(sbet2) + sbet1 is
232 // a better measure. This logic is used in assigning calp2 in Lambda12.
233 // Sometimes these quantities vanish and in that case we force bet2 = +/-
234 // bet1 exactly. An example where is is necessary is the inverse problem
235 // 48.522876735459 0 -48.52287673545898293 179.599720456223079643
236 // which failed with Visual Studio 10 (Release and Debug)
237
238 if (cbet1 < -sbet1) {
239 if (cbet2 == cbet1)
240 sbet2 = copysign(sbet1, sbet2);
241 } else {
242 if (fabs(sbet2) == -sbet1)
243 cbet2 = cbet1;
244 }
245
246 real
247 dn1 = sqrt(1 + _ep2 * Math::sq(sbet1)),
248 dn2 = sqrt(1 + _ep2 * Math::sq(sbet2));
249
250 real a12, sig12;
251 // index zero element of this array is unused
252 real Ca[nC_];
253
254 bool meridian = lat1 == -Math::qd || slam12 == 0;
255
256 if (meridian) {
257
258 // Endpoints are on a single full meridian, so the geodesic might lie on
259 // a meridian.
260
261 calp1 = clam12; salp1 = slam12; // Head to the target longitude
262 calp2 = 1; salp2 = 0; // At the target we're heading north
263
264 real
265 // tan(bet) = tan(sig) * cos(alp)
266 ssig1 = sbet1, csig1 = calp1 * cbet1,
267 ssig2 = sbet2, csig2 = calp2 * cbet2;
268
269 // sig12 = sig2 - sig1
270 sig12 = atan2(fmax(real(0), csig1 * ssig2 - ssig1 * csig2) + real(0),
271 csig1 * csig2 + ssig1 * ssig2);
272 {
273 real dummy;
274 Lengths(_n, sig12, ssig1, csig1, dn1, ssig2, csig2, dn2, cbet1, cbet2,
275 outmask | DISTANCE | REDUCEDLENGTH,
276 s12x, m12x, dummy, M12, M21, Ca);
277 }
278 // Add the check for sig12 since zero length geodesics might yield m12 <
279 // 0. Test case was
280 //
281 // echo 20.001 0 20.001 0 | GeodSolve -i
282 //
283 // In fact, we will have sig12 > pi/2 for meridional geodesic which is
284 // not a shortest path.
285 // TODO: investigate m12 < 0 result for aarch/ppc (with -f -p 20)
286 // 20.001000000000001 0.000000000000000 180.000000000000000
287 // 20.001000000000001 0.000000000000000 180.000000000000000
288 // 0.0000000002 0.000000000000001 -0.0000000001
289 // 0.99999999999999989 0.99999999999999989 0.000
290 if (sig12 < 1 || m12x >= 0) {
291 // Need at least 2, to handle 90 0 90 180
292 if (sig12 < 3 * tiny_ ||
293 // Prevent negative s12 or m12 for short lines
294 (sig12 < tol0_ && (s12x < 0 || m12x < 0)))
295 sig12 = m12x = s12x = 0;
296 m12x *= _b;
297 s12x *= _b;
298 a12 = sig12 / Math::degree();
299 } else
300 // m12 < 0, i.e., prolate and too close to anti-podal
301 meridian = false;
302 }
303
304 // somg12 == 2 marks that it needs to be calculated
305 real omg12 = 0, somg12 = 2, comg12 = 0;
306 if (!meridian &&
307 sbet1 == 0 && // and sbet2 == 0
308 (_f <= 0 || lon12s >= _f * Math::hd)) {
309
310 // Geodesic runs along equator
311 calp1 = calp2 = 0; salp1 = salp2 = 1;
312 s12x = _a * lam12;
313 sig12 = omg12 = lam12 / _f1;
314 m12x = _b * sin(sig12);
315 if (outmask & GEODESICSCALE)
316 M12 = M21 = cos(sig12);
317 a12 = lon12 / _f1;
318
319 } else if (!meridian) {
320
321 // Now point1 and point2 belong within a hemisphere bounded by a
322 // meridian and geodesic is neither meridional or equatorial.
323
324 // Figure a starting point for Newton's method
325 real dnm;
326 sig12 = InverseStart(sbet1, cbet1, dn1, sbet2, cbet2, dn2,
327 lam12, slam12, clam12,
328 salp1, calp1, salp2, calp2, dnm,
329 Ca);
330
331 if (sig12 >= 0) {
332 // Short lines (InverseStart sets salp2, calp2, dnm)
333 s12x = sig12 * _b * dnm;
334 m12x = Math::sq(dnm) * _b * sin(sig12 / dnm);
335 if (outmask & GEODESICSCALE)
336 M12 = M21 = cos(sig12 / dnm);
337 a12 = sig12 / Math::degree();
338 omg12 = lam12 / (_f1 * dnm);
339 } else {
340
341 // Newton's method. This is a straightforward solution of f(alp1) =
342 // lambda12(alp1) - lam12 = 0 with one wrinkle. f(alp) has exactly one
343 // root in the interval (0, pi) and its derivative is positive at the
344 // root. Thus f(alp) is positive for alp > alp1 and negative for alp <
345 // alp1. During the course of the iteration, a range (alp1a, alp1b) is
346 // maintained which brackets the root and with each evaluation of
347 // f(alp) the range is shrunk, if possible. Newton's method is
348 // restarted whenever the derivative of f is negative (because the new
349 // value of alp1 is then further from the solution) or if the new
350 // estimate of alp1 lies outside (0,pi); in this case, the new starting
351 // guess is taken to be (alp1a + alp1b) / 2.
352 //
353 // initial values to suppress warnings (if loop is executed 0 times)
354 real ssig1 = 0, csig1 = 0, ssig2 = 0, csig2 = 0, eps = 0, domg12 = 0;
355 unsigned numit = 0;
356 // Bracketing range
357 real salp1a = tiny_, calp1a = 1, salp1b = tiny_, calp1b = -1;
358 for (bool tripn = false, tripb = false;; ++numit) {
359 // the WGS84 test set: mean = 1.47, sd = 1.25, max = 16
360 // WGS84 and random input: mean = 2.85, sd = 0.60
361 real dv;
362 real v = Lambda12(sbet1, cbet1, dn1, sbet2, cbet2, dn2, salp1, calp1,
363 slam12, clam12,
364 salp2, calp2, sig12, ssig1, csig1, ssig2, csig2,
365 eps, domg12, numit < maxit1_, dv, Ca);
366 if (tripb ||
367 // Reversed test to allow escape with NaNs
368 !(fabs(v) >= (tripn ? 8 : 1) * tol0_) ||
369 // Enough bisections to get accurate result
370 numit == maxit2_)
371 break;
372 // Update bracketing values
373 if (v > 0 && (numit > maxit1_ || calp1/salp1 > calp1b/salp1b))
374 { salp1b = salp1; calp1b = calp1; }
375 else if (v < 0 && (numit > maxit1_ || calp1/salp1 < calp1a/salp1a))
376 { salp1a = salp1; calp1a = calp1; }
377 if (numit < maxit1_ && dv > 0) {
378 real
379 dalp1 = -v/dv;
380 // |dalp1| < pi test moved earlier because GEOGRAPHICLIB_PRECISION
381 // = 5 can result in dalp1 = 10^(10^8). Then sin(dalp1) takes ages
382 // (because of the need to do accurate range reduction).
383 if (fabs(dalp1) < Math::pi()) {
384 real
385 sdalp1 = sin(dalp1), cdalp1 = cos(dalp1),
386 nsalp1 = salp1 * cdalp1 + calp1 * sdalp1;
387 if (nsalp1 > 0) {
388 calp1 = calp1 * cdalp1 - salp1 * sdalp1;
389 salp1 = nsalp1;
390 Math::norm(salp1, calp1);
391 // In some regimes we don't get quadratic convergence because
392 // slope -> 0. So use convergence conditions based on epsilon
393 // instead of sqrt(epsilon).
394 tripn = fabs(v) <= 16 * tol0_;
395 continue;
396 }
397 }
398 }
399 // Either dv was not positive or updated value was outside legal
400 // range. Use the midpoint of the bracket as the next estimate.
401 // This mechanism is not needed for the WGS84 ellipsoid, but it does
402 // catch problems with more eccentric ellipsoids. Its efficacy is
403 // such for the WGS84 test set with the starting guess set to alp1 =
404 // 90deg:
405 // the WGS84 test set: mean = 5.21, sd = 3.93, max = 24
406 // WGS84 and random input: mean = 4.74, sd = 0.99
407 salp1 = (salp1a + salp1b)/2;
408 calp1 = (calp1a + calp1b)/2;
409 Math::norm(salp1, calp1);
410 tripn = false;
411 tripb = (fabs(salp1a - salp1) + (calp1a - calp1) < tolb_ ||
412 fabs(salp1 - salp1b) + (calp1 - calp1b) < tolb_);
413 }
414 {
415 real dummy;
416 // Ensure that the reduced length and geodesic scale are computed in
417 // a "canonical" way, with the I2 integral.
418 unsigned lengthmask = outmask |
419 (outmask & (REDUCEDLENGTH | GEODESICSCALE) ? DISTANCE : NONE);
420 Lengths(eps, sig12, ssig1, csig1, dn1, ssig2, csig2, dn2,
421 cbet1, cbet2, lengthmask, s12x, m12x, dummy, M12, M21, Ca);
422 }
423 m12x *= _b;
424 s12x *= _b;
425 a12 = sig12 / Math::degree();
426 if (outmask & AREA) {
427 // omg12 = lam12 - domg12
428 real sdomg12 = sin(domg12), cdomg12 = cos(domg12);
429 somg12 = slam12 * cdomg12 - clam12 * sdomg12;
430 comg12 = clam12 * cdomg12 + slam12 * sdomg12;
431 }
432 }
433 }
434
435 if (outmask & DISTANCE)
436 s12 = real(0) + s12x; // Convert -0 to 0
437
438 if (outmask & REDUCEDLENGTH)
439 m12 = real(0) + m12x; // Convert -0 to 0
440
441 if (outmask & AREA) {
442 real
443 // From Lambda12: sin(alp1) * cos(bet1) = sin(alp0)
444 salp0 = salp1 * cbet1,
445 calp0 = hypot(calp1, salp1 * sbet1); // calp0 > 0
446 real alp12;
447 if (calp0 != 0 && salp0 != 0) {
448 real
449 // From Lambda12: tan(bet) = tan(sig) * cos(alp)
450 ssig1 = sbet1, csig1 = calp1 * cbet1,
451 ssig2 = sbet2, csig2 = calp2 * cbet2,
452 k2 = Math::sq(calp0) * _ep2,
453 eps = k2 / (2 * (1 + sqrt(1 + k2)) + k2),
454 // Multiplier = a^2 * e^2 * cos(alpha0) * sin(alpha0).
455 A4 = Math::sq(_a) * calp0 * salp0 * _e2;
456 Math::norm(ssig1, csig1);
457 Math::norm(ssig2, csig2);
458 C4f(eps, Ca);
459 real
460 B41 = SinCosSeries(false, ssig1, csig1, Ca, nC4_),
461 B42 = SinCosSeries(false, ssig2, csig2, Ca, nC4_);
462 S12 = A4 * (B42 - B41);
463 } else
464 // Avoid problems with indeterminate sig1, sig2 on equator
465 S12 = 0;
466 if (!meridian && somg12 == 2) {
467 somg12 = sin(omg12); comg12 = cos(omg12);
468 }
469
470 if (!meridian &&
471 // omg12 < 3/4 * pi
472 comg12 > -real(0.7071) && // Long difference not too big
473 sbet2 - sbet1 < real(1.75)) { // Lat difference not too big
474 // Use tan(Gamma/2) = tan(omg12/2)
475 // * (tan(bet1/2)+tan(bet2/2))/(1+tan(bet1/2)*tan(bet2/2))
476 // with tan(x/2) = sin(x)/(1+cos(x))
477 real domg12 = 1 + comg12, dbet1 = 1 + cbet1, dbet2 = 1 + cbet2;
478 alp12 = 2 * atan2( somg12 * ( sbet1 * dbet2 + sbet2 * dbet1 ),
479 domg12 * ( sbet1 * sbet2 + dbet1 * dbet2 ) );
480 } else {
481 // alp12 = alp2 - alp1, used in atan2 so no need to normalize
482 real
483 salp12 = salp2 * calp1 - calp2 * salp1,
484 calp12 = calp2 * calp1 + salp2 * salp1;
485 // The right thing appears to happen if alp1 = +/-180 and alp2 = 0, viz
486 // salp12 = -0 and alp12 = -180. However this depends on the sign
487 // being attached to 0 correctly. The following ensures the correct
488 // behavior.
489 if (salp12 == 0 && calp12 < 0) {
490 salp12 = tiny_ * calp1;
491 calp12 = -1;
492 }
493 alp12 = atan2(salp12, calp12);
494 }
495 S12 += _c2 * alp12;
496 S12 *= swapp * lonsign * latsign;
497 // Convert -0 to 0
498 S12 += 0;
499 }
500
501 // Convert calp, salp to azimuth accounting for lonsign, swapp, latsign.
502 if (swapp < 0) {
503 swap(salp1, salp2);
504 swap(calp1, calp2);
505 if (outmask & GEODESICSCALE)
506 swap(M12, M21);
507 }
508
509 salp1 *= swapp * lonsign; calp1 *= swapp * latsign;
510 salp2 *= swapp * lonsign; calp2 *= swapp * latsign;
511 // Returned value in [0, 180]
512 return a12;
513 }
514
515 Math::real Geodesic::GenInverse(real lat1, real lon1, real lat2, real lon2,
516 unsigned outmask,
517 real& s12, real& azi1, real& azi2,
518 real& m12, real& M12, real& M21,
519 real& S12) const {
520 outmask &= OUT_MASK;
521 real salp1, calp1, salp2, calp2,
522 a12 = GenInverse(lat1, lon1, lat2, lon2,
523 outmask, s12, salp1, calp1, salp2, calp2,
524 m12, M12, M21, S12);
525 if (outmask & AZIMUTH) {
526 azi1 = Math::atan2d(salp1, calp1);
527 azi2 = Math::atan2d(salp2, calp2);
528 }
529 return a12;
530 }
531
533 real lat2, real lon2,
534 unsigned caps) const {
535 real t, salp1, calp1, salp2, calp2,
536 a12 = GenInverse(lat1, lon1, lat2, lon2,
537 // No need to specify AZIMUTH here
538 0u, t, salp1, calp1, salp2, calp2,
539 t, t, t, t),
540 azi1 = Math::atan2d(salp1, calp1);
541 // Ensure that a12 can be converted to a distance
542 if (caps & (OUT_MASK & DISTANCE_IN)) caps |= DISTANCE;
543 return
544 GeodesicLine(*this, lat1, lon1, azi1, salp1, calp1, caps, true, a12);
545 }
546
547 void Geodesic::Lengths(real eps, real sig12,
548 real ssig1, real csig1, real dn1,
549 real ssig2, real csig2, real dn2,
550 real cbet1, real cbet2, unsigned outmask,
551 real& s12b, real& m12b, real& m0,
552 real& M12, real& M21,
553 // Scratch area of the right size
554 real Ca[]) const {
555 // Return m12b = (reduced length)/_b; also calculate s12b = distance/_b,
556 // and m0 = coefficient of secular term in expression for reduced length.
557
558 outmask &= OUT_MASK;
559 // outmask & DISTANCE: set s12b
560 // outmask & REDUCEDLENGTH: set m12b & m0
561 // outmask & GEODESICSCALE: set M12 & M21
562
563 real m0x = 0, J12 = 0, A1 = 0, A2 = 0;
564 real Cb[nC2_ + 1];
565 if (outmask & (DISTANCE | REDUCEDLENGTH | GEODESICSCALE)) {
566 A1 = A1m1f(eps);
567 C1f(eps, Ca);
568 if (outmask & (REDUCEDLENGTH | GEODESICSCALE)) {
569 A2 = A2m1f(eps);
570 C2f(eps, Cb);
571 m0x = A1 - A2;
572 A2 = 1 + A2;
573 }
574 A1 = 1 + A1;
575 }
576 if (outmask & DISTANCE) {
577 real B1 = SinCosSeries(true, ssig2, csig2, Ca, nC1_) -
578 SinCosSeries(true, ssig1, csig1, Ca, nC1_);
579 // Missing a factor of _b
580 s12b = A1 * (sig12 + B1);
581 if (outmask & (REDUCEDLENGTH | GEODESICSCALE)) {
582 real B2 = SinCosSeries(true, ssig2, csig2, Cb, nC2_) -
583 SinCosSeries(true, ssig1, csig1, Cb, nC2_);
584 J12 = m0x * sig12 + (A1 * B1 - A2 * B2);
585 }
586 } else if (outmask & (REDUCEDLENGTH | GEODESICSCALE)) {
587 // Assume here that nC1_ >= nC2_
588 for (int l = 1; l <= nC2_; ++l)
589 Cb[l] = A1 * Ca[l] - A2 * Cb[l];
590 J12 = m0x * sig12 + (SinCosSeries(true, ssig2, csig2, Cb, nC2_) -
591 SinCosSeries(true, ssig1, csig1, Cb, nC2_));
592 }
593 if (outmask & REDUCEDLENGTH) {
594 m0 = m0x;
595 // Missing a factor of _b.
596 // Add parens around (csig1 * ssig2) and (ssig1 * csig2) to ensure
597 // accurate cancellation in the case of coincident points.
598 m12b = dn2 * (csig1 * ssig2) - dn1 * (ssig1 * csig2) -
599 csig1 * csig2 * J12;
600 }
601 if (outmask & GEODESICSCALE) {
602 real csig12 = csig1 * csig2 + ssig1 * ssig2;
603 real t = _ep2 * (cbet1 - cbet2) * (cbet1 + cbet2) / (dn1 + dn2);
604 M12 = csig12 + (t * ssig2 - csig2 * J12) * ssig1 / dn1;
605 M21 = csig12 - (t * ssig1 - csig1 * J12) * ssig2 / dn2;
606 }
607 }
608
609 Math::real Geodesic::Astroid(real x, real y) {
610 // Solve k^4+2*k^3-(x^2+y^2-1)*k^2-2*y^2*k-y^2 = 0 for positive root k.
611 // This solution is adapted from Geocentric::Reverse.
612 real k;
613 real
614 p = Math::sq(x),
615 q = Math::sq(y),
616 r = (p + q - 1) / 6;
617 if ( !(q == 0 && r <= 0) ) {
618 real
619 // Avoid possible division by zero when r = 0 by multiplying equations
620 // for s and t by r^3 and r, resp.
621 S = p * q / 4, // S = r^3 * s
622 r2 = Math::sq(r),
623 r3 = r * r2,
624 // The discriminant of the quadratic equation for T3. This is zero on
625 // the evolute curve p^(1/3)+q^(1/3) = 1
626 disc = S * (S + 2 * r3);
627 real u = r;
628 if (disc >= 0) {
629 real T3 = S + r3;
630 // Pick the sign on the sqrt to maximize abs(T3). This minimizes loss
631 // of precision due to cancellation. The result is unchanged because
632 // of the way the T is used in definition of u.
633 T3 += T3 < 0 ? -sqrt(disc) : sqrt(disc); // T3 = (r * t)^3
634 // N.B. cbrt always returns the real root. cbrt(-8) = -2.
635 real T = cbrt(T3); // T = r * t
636 // T can be zero; but then r2 / T -> 0.
637 u += T + (T != 0 ? r2 / T : 0);
638 } else {
639 // T is complex, but the way u is defined the result is real.
640 real ang = atan2(sqrt(-disc), -(S + r3));
641 // There are three possible cube roots. We choose the root which
642 // avoids cancellation. Note that disc < 0 implies that r < 0.
643 u += 2 * r * cos(ang / 3);
644 }
645 real
646 v = sqrt(Math::sq(u) + q), // guaranteed positive
647 // Avoid loss of accuracy when u < 0.
648 uv = u < 0 ? q / (v - u) : u + v, // u+v, guaranteed positive
649 w = (uv - q) / (2 * v); // positive?
650 // Rearrange expression for k to avoid loss of accuracy due to
651 // subtraction. Division by 0 not possible because uv > 0, w >= 0.
652 k = uv / (sqrt(uv + Math::sq(w)) + w); // guaranteed positive
653 } else { // q == 0 && r <= 0
654 // y = 0 with |x| <= 1. Handle this case directly.
655 // for y small, positive root is k = abs(y)/sqrt(1-x^2)
656 k = 0;
657 }
658 return k;
659 }
660
661 Math::real Geodesic::InverseStart(real sbet1, real cbet1, real dn1,
662 real sbet2, real cbet2, real dn2,
663 real lam12, real slam12, real clam12,
664 real& salp1, real& calp1,
665 // Only updated if return val >= 0
666 real& salp2, real& calp2,
667 // Only updated for short lines
668 real& dnm,
669 // Scratch area of the right size
670 real Ca[]) const {
671 // Return a starting point for Newton's method in salp1 and calp1 (function
672 // value is -1). If Newton's method doesn't need to be used, return also
673 // salp2 and calp2 and function value is sig12.
674 real
675 sig12 = -1, // Return value
676 // bet12 = bet2 - bet1 in [0, pi); bet12a = bet2 + bet1 in (-pi, 0]
677 sbet12 = sbet2 * cbet1 - cbet2 * sbet1,
678 cbet12 = cbet2 * cbet1 + sbet2 * sbet1;
679 real sbet12a = sbet2 * cbet1 + cbet2 * sbet1;
680 bool shortline = cbet12 >= 0 && sbet12 < real(0.5) &&
681 cbet2 * lam12 < real(0.5);
682 real somg12, comg12;
683 if (shortline) {
684 real sbetm2 = Math::sq(sbet1 + sbet2);
685 // sin((bet1+bet2)/2)^2
686 // = (sbet1 + sbet2)^2 / ((sbet1 + sbet2)^2 + (cbet1 + cbet2)^2)
687 sbetm2 /= sbetm2 + Math::sq(cbet1 + cbet2);
688 dnm = sqrt(1 + _ep2 * sbetm2);
689 real omg12 = lam12 / (_f1 * dnm);
690 somg12 = sin(omg12); comg12 = cos(omg12);
691 } else {
692 somg12 = slam12; comg12 = clam12;
693 }
694
695 salp1 = cbet2 * somg12;
696 calp1 = comg12 >= 0 ?
697 sbet12 + cbet2 * sbet1 * Math::sq(somg12) / (1 + comg12) :
698 sbet12a - cbet2 * sbet1 * Math::sq(somg12) / (1 - comg12);
699
700 real
701 ssig12 = hypot(salp1, calp1),
702 csig12 = sbet1 * sbet2 + cbet1 * cbet2 * comg12;
703
704 if (shortline && ssig12 < _etol2) {
705 // really short lines
706 salp2 = cbet1 * somg12;
707 calp2 = sbet12 - cbet1 * sbet2 *
708 (comg12 >= 0 ? Math::sq(somg12) / (1 + comg12) : 1 - comg12);
709 Math::norm(salp2, calp2);
710 // Set return value
711 sig12 = atan2(ssig12, csig12);
712 } else if (fabs(_n) > real(0.1) || // Skip astroid calc if too eccentric
713 csig12 >= 0 ||
714 ssig12 >= 6 * fabs(_n) * Math::pi() * Math::sq(cbet1)) {
715 // Nothing to do, zeroth order spherical approximation is OK
716 } else {
717 // Scale lam12 and bet2 to x, y coordinate system where antipodal point
718 // is at origin and singular point is at y = 0, x = -1.
719 real x, y, lamscale, betscale;
720 real lam12x = atan2(-slam12, -clam12); // lam12 - pi
721 if (_f >= 0) { // In fact f == 0 does not get here
722 // x = dlong, y = dlat
723 {
724 real
725 k2 = Math::sq(sbet1) * _ep2,
726 eps = k2 / (2 * (1 + sqrt(1 + k2)) + k2);
727 lamscale = _f * cbet1 * A3f(eps) * Math::pi();
728 }
729 betscale = lamscale * cbet1;
730
731 x = lam12x / lamscale;
732 y = sbet12a / betscale;
733 } else { // _f < 0
734 // x = dlat, y = dlong
735 real
736 cbet12a = cbet2 * cbet1 - sbet2 * sbet1,
737 bet12a = atan2(sbet12a, cbet12a);
738 real m12b, m0, dummy;
739 // In the case of lon12 = 180, this repeats a calculation made in
740 // Inverse.
741 Lengths(_n, Math::pi() + bet12a,
742 sbet1, -cbet1, dn1, sbet2, cbet2, dn2,
743 cbet1, cbet2,
744 REDUCEDLENGTH, dummy, m12b, m0, dummy, dummy, Ca);
745 x = -1 + m12b / (cbet1 * cbet2 * m0 * Math::pi());
746 betscale = x < -real(0.01) ? sbet12a / x :
747 -_f * Math::sq(cbet1) * Math::pi();
748 lamscale = betscale / cbet1;
749 y = lam12x / lamscale;
750 }
751
752 if (y > -tol1_ && x > -1 - xthresh_) {
753 // strip near cut
754 // Need real(x) here to cast away the volatility of x for min/max
755 if (_f >= 0) {
756 salp1 = fmin(real(1), -x); calp1 = - sqrt(1 - Math::sq(salp1));
757 } else {
758 calp1 = fmax(real(x > -tol1_ ? 0 : -1), x);
759 salp1 = sqrt(1 - Math::sq(calp1));
760 }
761 } else {
762 // Estimate alp1, by solving the astroid problem.
763 //
764 // Could estimate alpha1 = theta + pi/2, directly, i.e.,
765 // calp1 = y/k; salp1 = -x/(1+k); for _f >= 0
766 // calp1 = x/(1+k); salp1 = -y/k; for _f < 0 (need to check)
767 //
768 // However, it's better to estimate omg12 from astroid and use
769 // spherical formula to compute alp1. This reduces the mean number of
770 // Newton iterations for astroid cases from 2.24 (min 0, max 6) to 2.12
771 // (min 0 max 5). The changes in the number of iterations are as
772 // follows:
773 //
774 // change percent
775 // 1 5
776 // 0 78
777 // -1 16
778 // -2 0.6
779 // -3 0.04
780 // -4 0.002
781 //
782 // The histogram of iterations is (m = number of iterations estimating
783 // alp1 directly, n = number of iterations estimating via omg12, total
784 // number of trials = 148605):
785 //
786 // iter m n
787 // 0 148 186
788 // 1 13046 13845
789 // 2 93315 102225
790 // 3 36189 32341
791 // 4 5396 7
792 // 5 455 1
793 // 6 56 0
794 //
795 // Because omg12 is near pi, estimate work with omg12a = pi - omg12
796 real k = Astroid(x, y);
797 real
798 omg12a = lamscale * ( _f >= 0 ? -x * k/(1 + k) : -y * (1 + k)/k );
799 somg12 = sin(omg12a); comg12 = -cos(omg12a);
800 // Update spherical estimate of alp1 using omg12 instead of lam12
801 salp1 = cbet2 * somg12;
802 calp1 = sbet12a - cbet2 * sbet1 * Math::sq(somg12) / (1 - comg12);
803 }
804 }
805 // Sanity check on starting guess. Backwards check allows NaN through.
806 if (!(salp1 <= 0))
807 Math::norm(salp1, calp1);
808 else {
809 salp1 = 1; calp1 = 0;
810 }
811 return sig12;
812 }
813
814 Math::real Geodesic::Lambda12(real sbet1, real cbet1, real dn1,
815 real sbet2, real cbet2, real dn2,
816 real salp1, real calp1,
817 real slam120, real clam120,
818 real& salp2, real& calp2,
819 real& sig12,
820 real& ssig1, real& csig1,
821 real& ssig2, real& csig2,
822 real& eps, real& domg12,
823 bool diffp, real& dlam12,
824 // Scratch area of the right size
825 real Ca[]) const {
826
827 if (sbet1 == 0 && calp1 == 0)
828 // Break degeneracy of equatorial line. This case has already been
829 // handled.
830 calp1 = -tiny_;
831
832 real
833 // sin(alp1) * cos(bet1) = sin(alp0)
834 salp0 = salp1 * cbet1,
835 calp0 = hypot(calp1, salp1 * sbet1); // calp0 > 0
836
837 real somg1, comg1, somg2, comg2, somg12, comg12, lam12;
838 // tan(bet1) = tan(sig1) * cos(alp1)
839 // tan(omg1) = sin(alp0) * tan(sig1) = tan(omg1)=tan(alp1)*sin(bet1)
840 ssig1 = sbet1; somg1 = salp0 * sbet1;
841 csig1 = comg1 = calp1 * cbet1;
842 Math::norm(ssig1, csig1);
843 // Math::norm(somg1, comg1); -- don't need to normalize!
844
845 // Enforce symmetries in the case abs(bet2) = -bet1. Need to be careful
846 // about this case, since this can yield singularities in the Newton
847 // iteration.
848 // sin(alp2) * cos(bet2) = sin(alp0)
849 salp2 = cbet2 != cbet1 ? salp0 / cbet2 : salp1;
850 // calp2 = sqrt(1 - sq(salp2))
851 // = sqrt(sq(calp0) - sq(sbet2)) / cbet2
852 // and subst for calp0 and rearrange to give (choose positive sqrt
853 // to give alp2 in [0, pi/2]).
854 calp2 = cbet2 != cbet1 || fabs(sbet2) != -sbet1 ?
855 sqrt(Math::sq(calp1 * cbet1) +
856 (cbet1 < -sbet1 ?
857 (cbet2 - cbet1) * (cbet1 + cbet2) :
858 (sbet1 - sbet2) * (sbet1 + sbet2))) / cbet2 :
859 fabs(calp1);
860 // tan(bet2) = tan(sig2) * cos(alp2)
861 // tan(omg2) = sin(alp0) * tan(sig2).
862 ssig2 = sbet2; somg2 = salp0 * sbet2;
863 csig2 = comg2 = calp2 * cbet2;
864 Math::norm(ssig2, csig2);
865 // Math::norm(somg2, comg2); -- don't need to normalize!
866
867 // sig12 = sig2 - sig1, limit to [0, pi]
868 sig12 = atan2(fmax(real(0), csig1 * ssig2 - ssig1 * csig2) + real(0),
869 csig1 * csig2 + ssig1 * ssig2);
870
871 // omg12 = omg2 - omg1, limit to [0, pi]
872 somg12 = fmax(real(0), comg1 * somg2 - somg1 * comg2) + real(0);
873 comg12 = comg1 * comg2 + somg1 * somg2;
874 // eta = omg12 - lam120
875 real eta = atan2(somg12 * clam120 - comg12 * slam120,
876 comg12 * clam120 + somg12 * slam120);
877 real B312;
878 real k2 = Math::sq(calp0) * _ep2;
879 eps = k2 / (2 * (1 + sqrt(1 + k2)) + k2);
880 C3f(eps, Ca);
881 B312 = (SinCosSeries(true, ssig2, csig2, Ca, nC3_-1) -
882 SinCosSeries(true, ssig1, csig1, Ca, nC3_-1));
883 domg12 = -_f * A3f(eps) * salp0 * (sig12 + B312);
884 lam12 = eta + domg12;
885
886 if (diffp) {
887 if (calp2 == 0)
888 dlam12 = - 2 * _f1 * dn1 / sbet1;
889 else {
890 real dummy;
891 Lengths(eps, sig12, ssig1, csig1, dn1, ssig2, csig2, dn2,
892 cbet1, cbet2, REDUCEDLENGTH,
893 dummy, dlam12, dummy, dummy, dummy, Ca);
894 dlam12 *= _f1 / (calp2 * cbet2);
895 }
896 }
897
898 return lam12;
899 }
900
901 Math::real Geodesic::A3f(real eps) const {
902 // Evaluate A3
903 return Math::polyval(nA3_ - 1, _aA3x, eps);
904 }
905
906 void Geodesic::C3f(real eps, real c[]) const {
907 // Evaluate C3 coeffs
908 // Elements c[1] thru c[nC3_ - 1] are set
909 real mult = 1;
910 int o = 0;
911 for (int l = 1; l < nC3_; ++l) { // l is index of C3[l]
912 int m = nC3_ - l - 1; // order of polynomial in eps
913 mult *= eps;
914 c[l] = mult * Math::polyval(m, _cC3x + o, eps);
915 o += m + 1;
916 }
917 // Post condition: o == nC3x_
918 }
919
920 void Geodesic::C4f(real eps, real c[]) const {
921 // Evaluate C4 coeffs
922 // Elements c[0] thru c[nC4_ - 1] are set
923 real mult = 1;
924 int o = 0;
925 for (int l = 0; l < nC4_; ++l) { // l is index of C4[l]
926 int m = nC4_ - l - 1; // order of polynomial in eps
927 c[l] = mult * Math::polyval(m, _cC4x + o, eps);
928 o += m + 1;
929 mult *= eps;
930 }
931 // Post condition: o == nC4x_
932 }
933
934 // The static const coefficient arrays in the following functions are
935 // generated by Maxima and give the coefficients of the Taylor expansions for
936 // the geodesics. The convention on the order of these coefficients is as
937 // follows:
938 //
939 // ascending order in the trigonometric expansion,
940 // then powers of eps in descending order,
941 // finally powers of n in descending order.
942 //
943 // (For some expansions, only a subset of levels occur.) For each polynomial
944 // of order n at the lowest level, the (n+1) coefficients of the polynomial
945 // are followed by a divisor which is applied to the whole polynomial. In
946 // this way, the coefficients are expressible with no round off error. The
947 // sizes of the coefficient arrays are:
948 //
949 // A1m1f, A2m1f = floor(N/2) + 2
950 // C1f, C1pf, C2f, A3coeff = (N^2 + 7*N - 2*floor(N/2)) / 4
951 // C3coeff = (N - 1) * (N^2 + 7*N - 2*floor(N/2)) / 8
952 // C4coeff = N * (N + 1) * (N + 5) / 6
953 //
954 // where N = GEOGRAPHICLIB_GEODESIC_ORDER
955 // = nA1 = nA2 = nC1 = nC1p = nA3 = nC4
956
957 // The scale factor A1-1 = mean value of (d/dsigma)I1 - 1
958 Math::real Geodesic::A1m1f(real eps) {
959 // Generated by Maxima on 2015-05-05 18:08:12-04:00
960#if GEOGRAPHICLIB_GEODESIC_ORDER/2 == 1
961 static const real coeff[] = {
962 // (1-eps)*A1-1, polynomial in eps2 of order 1
963 1, 0, 4,
964 };
965#elif GEOGRAPHICLIB_GEODESIC_ORDER/2 == 2
966 static const real coeff[] = {
967 // (1-eps)*A1-1, polynomial in eps2 of order 2
968 1, 16, 0, 64,
969 };
970#elif GEOGRAPHICLIB_GEODESIC_ORDER/2 == 3
971 static const real coeff[] = {
972 // (1-eps)*A1-1, polynomial in eps2 of order 3
973 1, 4, 64, 0, 256,
974 };
975#elif GEOGRAPHICLIB_GEODESIC_ORDER/2 == 4
976 static const real coeff[] = {
977 // (1-eps)*A1-1, polynomial in eps2 of order 4
978 25, 64, 256, 4096, 0, 16384,
979 };
980#else
981#error "Bad value for GEOGRAPHICLIB_GEODESIC_ORDER"
982#endif
983 static_assert(sizeof(coeff) / sizeof(real) == nA1_/2 + 2,
984 "Coefficient array size mismatch in A1m1f");
985 int m = nA1_/2;
986 real t = Math::polyval(m, coeff, Math::sq(eps)) / coeff[m + 1];
987 return (t + eps) / (1 - eps);
988 }
989
990 // The coefficients C1[l] in the Fourier expansion of B1
991 void Geodesic::C1f(real eps, real c[]) {
992 // Generated by Maxima on 2015-05-05 18:08:12-04:00
993#if GEOGRAPHICLIB_GEODESIC_ORDER == 3
994 static const real coeff[] = {
995 // C1[1]/eps^1, polynomial in eps2 of order 1
996 3, -8, 16,
997 // C1[2]/eps^2, polynomial in eps2 of order 0
998 -1, 16,
999 // C1[3]/eps^3, polynomial in eps2 of order 0
1000 -1, 48,
1001 };
1002#elif GEOGRAPHICLIB_GEODESIC_ORDER == 4
1003 static const real coeff[] = {
1004 // C1[1]/eps^1, polynomial in eps2 of order 1
1005 3, -8, 16,
1006 // C1[2]/eps^2, polynomial in eps2 of order 1
1007 1, -2, 32,
1008 // C1[3]/eps^3, polynomial in eps2 of order 0
1009 -1, 48,
1010 // C1[4]/eps^4, polynomial in eps2 of order 0
1011 -5, 512,
1012 };
1013#elif GEOGRAPHICLIB_GEODESIC_ORDER == 5
1014 static const real coeff[] = {
1015 // C1[1]/eps^1, polynomial in eps2 of order 2
1016 -1, 6, -16, 32,
1017 // C1[2]/eps^2, polynomial in eps2 of order 1
1018 1, -2, 32,
1019 // C1[3]/eps^3, polynomial in eps2 of order 1
1020 9, -16, 768,
1021 // C1[4]/eps^4, polynomial in eps2 of order 0
1022 -5, 512,
1023 // C1[5]/eps^5, polynomial in eps2 of order 0
1024 -7, 1280,
1025 };
1026#elif GEOGRAPHICLIB_GEODESIC_ORDER == 6
1027 static const real coeff[] = {
1028 // C1[1]/eps^1, polynomial in eps2 of order 2
1029 -1, 6, -16, 32,
1030 // C1[2]/eps^2, polynomial in eps2 of order 2
1031 -9, 64, -128, 2048,
1032 // C1[3]/eps^3, polynomial in eps2 of order 1
1033 9, -16, 768,
1034 // C1[4]/eps^4, polynomial in eps2 of order 1
1035 3, -5, 512,
1036 // C1[5]/eps^5, polynomial in eps2 of order 0
1037 -7, 1280,
1038 // C1[6]/eps^6, polynomial in eps2 of order 0
1039 -7, 2048,
1040 };
1041#elif GEOGRAPHICLIB_GEODESIC_ORDER == 7
1042 static const real coeff[] = {
1043 // C1[1]/eps^1, polynomial in eps2 of order 3
1044 19, -64, 384, -1024, 2048,
1045 // C1[2]/eps^2, polynomial in eps2 of order 2
1046 -9, 64, -128, 2048,
1047 // C1[3]/eps^3, polynomial in eps2 of order 2
1048 -9, 72, -128, 6144,
1049 // C1[4]/eps^4, polynomial in eps2 of order 1
1050 3, -5, 512,
1051 // C1[5]/eps^5, polynomial in eps2 of order 1
1052 35, -56, 10240,
1053 // C1[6]/eps^6, polynomial in eps2 of order 0
1054 -7, 2048,
1055 // C1[7]/eps^7, polynomial in eps2 of order 0
1056 -33, 14336,
1057 };
1058#elif GEOGRAPHICLIB_GEODESIC_ORDER == 8
1059 static const real coeff[] = {
1060 // C1[1]/eps^1, polynomial in eps2 of order 3
1061 19, -64, 384, -1024, 2048,
1062 // C1[2]/eps^2, polynomial in eps2 of order 3
1063 7, -18, 128, -256, 4096,
1064 // C1[3]/eps^3, polynomial in eps2 of order 2
1065 -9, 72, -128, 6144,
1066 // C1[4]/eps^4, polynomial in eps2 of order 2
1067 -11, 96, -160, 16384,
1068 // C1[5]/eps^5, polynomial in eps2 of order 1
1069 35, -56, 10240,
1070 // C1[6]/eps^6, polynomial in eps2 of order 1
1071 9, -14, 4096,
1072 // C1[7]/eps^7, polynomial in eps2 of order 0
1073 -33, 14336,
1074 // C1[8]/eps^8, polynomial in eps2 of order 0
1075 -429, 262144,
1076 };
1077#else
1078#error "Bad value for GEOGRAPHICLIB_GEODESIC_ORDER"
1079#endif
1080 static_assert(sizeof(coeff) / sizeof(real) ==
1081 (nC1_*nC1_ + 7*nC1_ - 2*(nC1_/2)) / 4,
1082 "Coefficient array size mismatch in C1f");
1083 real
1084 eps2 = Math::sq(eps),
1085 d = eps;
1086 int o = 0;
1087 for (int l = 1; l <= nC1_; ++l) { // l is index of C1p[l]
1088 int m = (nC1_ - l) / 2; // order of polynomial in eps^2
1089 c[l] = d * Math::polyval(m, coeff + o, eps2) / coeff[o + m + 1];
1090 o += m + 2;
1091 d *= eps;
1092 }
1093 // Post condition: o == sizeof(coeff) / sizeof(real)
1094 }
1095
1096 // The coefficients C1p[l] in the Fourier expansion of B1p
1097 void Geodesic::C1pf(real eps, real c[]) {
1098 // Generated by Maxima on 2015-05-05 18:08:12-04:00
1099#if GEOGRAPHICLIB_GEODESIC_ORDER == 3
1100 static const real coeff[] = {
1101 // C1p[1]/eps^1, polynomial in eps2 of order 1
1102 -9, 16, 32,
1103 // C1p[2]/eps^2, polynomial in eps2 of order 0
1104 5, 16,
1105 // C1p[3]/eps^3, polynomial in eps2 of order 0
1106 29, 96,
1107 };
1108#elif GEOGRAPHICLIB_GEODESIC_ORDER == 4
1109 static const real coeff[] = {
1110 // C1p[1]/eps^1, polynomial in eps2 of order 1
1111 -9, 16, 32,
1112 // C1p[2]/eps^2, polynomial in eps2 of order 1
1113 -37, 30, 96,
1114 // C1p[3]/eps^3, polynomial in eps2 of order 0
1115 29, 96,
1116 // C1p[4]/eps^4, polynomial in eps2 of order 0
1117 539, 1536,
1118 };
1119#elif GEOGRAPHICLIB_GEODESIC_ORDER == 5
1120 static const real coeff[] = {
1121 // C1p[1]/eps^1, polynomial in eps2 of order 2
1122 205, -432, 768, 1536,
1123 // C1p[2]/eps^2, polynomial in eps2 of order 1
1124 -37, 30, 96,
1125 // C1p[3]/eps^3, polynomial in eps2 of order 1
1126 -225, 116, 384,
1127 // C1p[4]/eps^4, polynomial in eps2 of order 0
1128 539, 1536,
1129 // C1p[5]/eps^5, polynomial in eps2 of order 0
1130 3467, 7680,
1131 };
1132#elif GEOGRAPHICLIB_GEODESIC_ORDER == 6
1133 static const real coeff[] = {
1134 // C1p[1]/eps^1, polynomial in eps2 of order 2
1135 205, -432, 768, 1536,
1136 // C1p[2]/eps^2, polynomial in eps2 of order 2
1137 4005, -4736, 3840, 12288,
1138 // C1p[3]/eps^3, polynomial in eps2 of order 1
1139 -225, 116, 384,
1140 // C1p[4]/eps^4, polynomial in eps2 of order 1
1141 -7173, 2695, 7680,
1142 // C1p[5]/eps^5, polynomial in eps2 of order 0
1143 3467, 7680,
1144 // C1p[6]/eps^6, polynomial in eps2 of order 0
1145 38081, 61440,
1146 };
1147#elif GEOGRAPHICLIB_GEODESIC_ORDER == 7
1148 static const real coeff[] = {
1149 // C1p[1]/eps^1, polynomial in eps2 of order 3
1150 -4879, 9840, -20736, 36864, 73728,
1151 // C1p[2]/eps^2, polynomial in eps2 of order 2
1152 4005, -4736, 3840, 12288,
1153 // C1p[3]/eps^3, polynomial in eps2 of order 2
1154 8703, -7200, 3712, 12288,
1155 // C1p[4]/eps^4, polynomial in eps2 of order 1
1156 -7173, 2695, 7680,
1157 // C1p[5]/eps^5, polynomial in eps2 of order 1
1158 -141115, 41604, 92160,
1159 // C1p[6]/eps^6, polynomial in eps2 of order 0
1160 38081, 61440,
1161 // C1p[7]/eps^7, polynomial in eps2 of order 0
1162 459485, 516096,
1163 };
1164#elif GEOGRAPHICLIB_GEODESIC_ORDER == 8
1165 static const real coeff[] = {
1166 // C1p[1]/eps^1, polynomial in eps2 of order 3
1167 -4879, 9840, -20736, 36864, 73728,
1168 // C1p[2]/eps^2, polynomial in eps2 of order 3
1169 -86171, 120150, -142080, 115200, 368640,
1170 // C1p[3]/eps^3, polynomial in eps2 of order 2
1171 8703, -7200, 3712, 12288,
1172 // C1p[4]/eps^4, polynomial in eps2 of order 2
1173 1082857, -688608, 258720, 737280,
1174 // C1p[5]/eps^5, polynomial in eps2 of order 1
1175 -141115, 41604, 92160,
1176 // C1p[6]/eps^6, polynomial in eps2 of order 1
1177 -2200311, 533134, 860160,
1178 // C1p[7]/eps^7, polynomial in eps2 of order 0
1179 459485, 516096,
1180 // C1p[8]/eps^8, polynomial in eps2 of order 0
1181 109167851, 82575360,
1182 };
1183#else
1184#error "Bad value for GEOGRAPHICLIB_GEODESIC_ORDER"
1185#endif
1186 static_assert(sizeof(coeff) / sizeof(real) ==
1187 (nC1p_*nC1p_ + 7*nC1p_ - 2*(nC1p_/2)) / 4,
1188 "Coefficient array size mismatch in C1pf");
1189 real
1190 eps2 = Math::sq(eps),
1191 d = eps;
1192 int o = 0;
1193 for (int l = 1; l <= nC1p_; ++l) { // l is index of C1p[l]
1194 int m = (nC1p_ - l) / 2; // order of polynomial in eps^2
1195 c[l] = d * Math::polyval(m, coeff + o, eps2) / coeff[o + m + 1];
1196 o += m + 2;
1197 d *= eps;
1198 }
1199 // Post condition: o == sizeof(coeff) / sizeof(real)
1200 }
1201
1202 // The scale factor A2-1 = mean value of (d/dsigma)I2 - 1
1203 Math::real Geodesic::A2m1f(real eps) {
1204 // Generated by Maxima on 2015-05-29 08:09:47-04:00
1205#if GEOGRAPHICLIB_GEODESIC_ORDER/2 == 1
1206 static const real coeff[] = {
1207 // (eps+1)*A2-1, polynomial in eps2 of order 1
1208 -3, 0, 4,
1209 }; // count = 3
1210#elif GEOGRAPHICLIB_GEODESIC_ORDER/2 == 2
1211 static const real coeff[] = {
1212 // (eps+1)*A2-1, polynomial in eps2 of order 2
1213 -7, -48, 0, 64,
1214 }; // count = 4
1215#elif GEOGRAPHICLIB_GEODESIC_ORDER/2 == 3
1216 static const real coeff[] = {
1217 // (eps+1)*A2-1, polynomial in eps2 of order 3
1218 -11, -28, -192, 0, 256,
1219 }; // count = 5
1220#elif GEOGRAPHICLIB_GEODESIC_ORDER/2 == 4
1221 static const real coeff[] = {
1222 // (eps+1)*A2-1, polynomial in eps2 of order 4
1223 -375, -704, -1792, -12288, 0, 16384,
1224 }; // count = 6
1225#else
1226#error "Bad value for GEOGRAPHICLIB_GEODESIC_ORDER"
1227#endif
1228 static_assert(sizeof(coeff) / sizeof(real) == nA2_/2 + 2,
1229 "Coefficient array size mismatch in A2m1f");
1230 int m = nA2_/2;
1231 real t = Math::polyval(m, coeff, Math::sq(eps)) / coeff[m + 1];
1232 return (t - eps) / (1 + eps);
1233 }
1234
1235 // The coefficients C2[l] in the Fourier expansion of B2
1236 void Geodesic::C2f(real eps, real c[]) {
1237 // Generated by Maxima on 2015-05-05 18:08:12-04:00
1238#if GEOGRAPHICLIB_GEODESIC_ORDER == 3
1239 static const real coeff[] = {
1240 // C2[1]/eps^1, polynomial in eps2 of order 1
1241 1, 8, 16,
1242 // C2[2]/eps^2, polynomial in eps2 of order 0
1243 3, 16,
1244 // C2[3]/eps^3, polynomial in eps2 of order 0
1245 5, 48,
1246 };
1247#elif GEOGRAPHICLIB_GEODESIC_ORDER == 4
1248 static const real coeff[] = {
1249 // C2[1]/eps^1, polynomial in eps2 of order 1
1250 1, 8, 16,
1251 // C2[2]/eps^2, polynomial in eps2 of order 1
1252 1, 6, 32,
1253 // C2[3]/eps^3, polynomial in eps2 of order 0
1254 5, 48,
1255 // C2[4]/eps^4, polynomial in eps2 of order 0
1256 35, 512,
1257 };
1258#elif GEOGRAPHICLIB_GEODESIC_ORDER == 5
1259 static const real coeff[] = {
1260 // C2[1]/eps^1, polynomial in eps2 of order 2
1261 1, 2, 16, 32,
1262 // C2[2]/eps^2, polynomial in eps2 of order 1
1263 1, 6, 32,
1264 // C2[3]/eps^3, polynomial in eps2 of order 1
1265 15, 80, 768,
1266 // C2[4]/eps^4, polynomial in eps2 of order 0
1267 35, 512,
1268 // C2[5]/eps^5, polynomial in eps2 of order 0
1269 63, 1280,
1270 };
1271#elif GEOGRAPHICLIB_GEODESIC_ORDER == 6
1272 static const real coeff[] = {
1273 // C2[1]/eps^1, polynomial in eps2 of order 2
1274 1, 2, 16, 32,
1275 // C2[2]/eps^2, polynomial in eps2 of order 2
1276 35, 64, 384, 2048,
1277 // C2[3]/eps^3, polynomial in eps2 of order 1
1278 15, 80, 768,
1279 // C2[4]/eps^4, polynomial in eps2 of order 1
1280 7, 35, 512,
1281 // C2[5]/eps^5, polynomial in eps2 of order 0
1282 63, 1280,
1283 // C2[6]/eps^6, polynomial in eps2 of order 0
1284 77, 2048,
1285 };
1286#elif GEOGRAPHICLIB_GEODESIC_ORDER == 7
1287 static const real coeff[] = {
1288 // C2[1]/eps^1, polynomial in eps2 of order 3
1289 41, 64, 128, 1024, 2048,
1290 // C2[2]/eps^2, polynomial in eps2 of order 2
1291 35, 64, 384, 2048,
1292 // C2[3]/eps^3, polynomial in eps2 of order 2
1293 69, 120, 640, 6144,
1294 // C2[4]/eps^4, polynomial in eps2 of order 1
1295 7, 35, 512,
1296 // C2[5]/eps^5, polynomial in eps2 of order 1
1297 105, 504, 10240,
1298 // C2[6]/eps^6, polynomial in eps2 of order 0
1299 77, 2048,
1300 // C2[7]/eps^7, polynomial in eps2 of order 0
1301 429, 14336,
1302 };
1303#elif GEOGRAPHICLIB_GEODESIC_ORDER == 8
1304 static const real coeff[] = {
1305 // C2[1]/eps^1, polynomial in eps2 of order 3
1306 41, 64, 128, 1024, 2048,
1307 // C2[2]/eps^2, polynomial in eps2 of order 3
1308 47, 70, 128, 768, 4096,
1309 // C2[3]/eps^3, polynomial in eps2 of order 2
1310 69, 120, 640, 6144,
1311 // C2[4]/eps^4, polynomial in eps2 of order 2
1312 133, 224, 1120, 16384,
1313 // C2[5]/eps^5, polynomial in eps2 of order 1
1314 105, 504, 10240,
1315 // C2[6]/eps^6, polynomial in eps2 of order 1
1316 33, 154, 4096,
1317 // C2[7]/eps^7, polynomial in eps2 of order 0
1318 429, 14336,
1319 // C2[8]/eps^8, polynomial in eps2 of order 0
1320 6435, 262144,
1321 };
1322#else
1323#error "Bad value for GEOGRAPHICLIB_GEODESIC_ORDER"
1324#endif
1325 static_assert(sizeof(coeff) / sizeof(real) ==
1326 (nC2_*nC2_ + 7*nC2_ - 2*(nC2_/2)) / 4,
1327 "Coefficient array size mismatch in C2f");
1328 real
1329 eps2 = Math::sq(eps),
1330 d = eps;
1331 int o = 0;
1332 for (int l = 1; l <= nC2_; ++l) { // l is index of C2[l]
1333 int m = (nC2_ - l) / 2; // order of polynomial in eps^2
1334 c[l] = d * Math::polyval(m, coeff + o, eps2) / coeff[o + m + 1];
1335 o += m + 2;
1336 d *= eps;
1337 }
1338 // Post condition: o == sizeof(coeff) / sizeof(real)
1339 }
1340
1341 // The scale factor A3 = mean value of (d/dsigma)I3
1342 void Geodesic::A3coeff() {
1343 // Generated by Maxima on 2015-05-05 18:08:13-04:00
1344#if GEOGRAPHICLIB_GEODESIC_ORDER == 3
1345 static const real coeff[] = {
1346 // A3, coeff of eps^2, polynomial in n of order 0
1347 -1, 4,
1348 // A3, coeff of eps^1, polynomial in n of order 1
1349 1, -1, 2,
1350 // A3, coeff of eps^0, polynomial in n of order 0
1351 1, 1,
1352 };
1353#elif GEOGRAPHICLIB_GEODESIC_ORDER == 4
1354 static const real coeff[] = {
1355 // A3, coeff of eps^3, polynomial in n of order 0
1356 -1, 16,
1357 // A3, coeff of eps^2, polynomial in n of order 1
1358 -1, -2, 8,
1359 // A3, coeff of eps^1, polynomial in n of order 1
1360 1, -1, 2,
1361 // A3, coeff of eps^0, polynomial in n of order 0
1362 1, 1,
1363 };
1364#elif GEOGRAPHICLIB_GEODESIC_ORDER == 5
1365 static const real coeff[] = {
1366 // A3, coeff of eps^4, polynomial in n of order 0
1367 -3, 64,
1368 // A3, coeff of eps^3, polynomial in n of order 1
1369 -3, -1, 16,
1370 // A3, coeff of eps^2, polynomial in n of order 2
1371 3, -1, -2, 8,
1372 // A3, coeff of eps^1, polynomial in n of order 1
1373 1, -1, 2,
1374 // A3, coeff of eps^0, polynomial in n of order 0
1375 1, 1,
1376 };
1377#elif GEOGRAPHICLIB_GEODESIC_ORDER == 6
1378 static const real coeff[] = {
1379 // A3, coeff of eps^5, polynomial in n of order 0
1380 -3, 128,
1381 // A3, coeff of eps^4, polynomial in n of order 1
1382 -2, -3, 64,
1383 // A3, coeff of eps^3, polynomial in n of order 2
1384 -1, -3, -1, 16,
1385 // A3, coeff of eps^2, polynomial in n of order 2
1386 3, -1, -2, 8,
1387 // A3, coeff of eps^1, polynomial in n of order 1
1388 1, -1, 2,
1389 // A3, coeff of eps^0, polynomial in n of order 0
1390 1, 1,
1391 };
1392#elif GEOGRAPHICLIB_GEODESIC_ORDER == 7
1393 static const real coeff[] = {
1394 // A3, coeff of eps^6, polynomial in n of order 0
1395 -5, 256,
1396 // A3, coeff of eps^5, polynomial in n of order 1
1397 -5, -3, 128,
1398 // A3, coeff of eps^4, polynomial in n of order 2
1399 -10, -2, -3, 64,
1400 // A3, coeff of eps^3, polynomial in n of order 3
1401 5, -1, -3, -1, 16,
1402 // A3, coeff of eps^2, polynomial in n of order 2
1403 3, -1, -2, 8,
1404 // A3, coeff of eps^1, polynomial in n of order 1
1405 1, -1, 2,
1406 // A3, coeff of eps^0, polynomial in n of order 0
1407 1, 1,
1408 };
1409#elif GEOGRAPHICLIB_GEODESIC_ORDER == 8
1410 static const real coeff[] = {
1411 // A3, coeff of eps^7, polynomial in n of order 0
1412 -25, 2048,
1413 // A3, coeff of eps^6, polynomial in n of order 1
1414 -15, -20, 1024,
1415 // A3, coeff of eps^5, polynomial in n of order 2
1416 -5, -10, -6, 256,
1417 // A3, coeff of eps^4, polynomial in n of order 3
1418 -5, -20, -4, -6, 128,
1419 // A3, coeff of eps^3, polynomial in n of order 3
1420 5, -1, -3, -1, 16,
1421 // A3, coeff of eps^2, polynomial in n of order 2
1422 3, -1, -2, 8,
1423 // A3, coeff of eps^1, polynomial in n of order 1
1424 1, -1, 2,
1425 // A3, coeff of eps^0, polynomial in n of order 0
1426 1, 1,
1427 };
1428#else
1429#error "Bad value for GEOGRAPHICLIB_GEODESIC_ORDER"
1430#endif
1431 static_assert(sizeof(coeff) / sizeof(real) ==
1432 (nA3_*nA3_ + 7*nA3_ - 2*(nA3_/2)) / 4,
1433 "Coefficient array size mismatch in A3f");
1434 int o = 0, k = 0;
1435 for (int j = nA3_ - 1; j >= 0; --j) { // coeff of eps^j
1436 int m = min(nA3_ - j - 1, j); // order of polynomial in n
1437 _aA3x[k++] = Math::polyval(m, coeff + o, _n) / coeff[o + m + 1];
1438 o += m + 2;
1439 }
1440 // Post condition: o == sizeof(coeff) / sizeof(real) && k == nA3x_
1441 }
1442
1443 // The coefficients C3[l] in the Fourier expansion of B3
1444 void Geodesic::C3coeff() {
1445 // Generated by Maxima on 2015-05-05 18:08:13-04:00
1446#if GEOGRAPHICLIB_GEODESIC_ORDER == 3
1447 static const real coeff[] = {
1448 // C3[1], coeff of eps^2, polynomial in n of order 0
1449 1, 8,
1450 // C3[1], coeff of eps^1, polynomial in n of order 1
1451 -1, 1, 4,
1452 // C3[2], coeff of eps^2, polynomial in n of order 0
1453 1, 16,
1454 };
1455#elif GEOGRAPHICLIB_GEODESIC_ORDER == 4
1456 static const real coeff[] = {
1457 // C3[1], coeff of eps^3, polynomial in n of order 0
1458 3, 64,
1459 // C3[1], coeff of eps^2, polynomial in n of order 1
1460 // This is a case where a leading 0 term has been inserted to maintain the
1461 // pattern in the orders of the polynomials.
1462 0, 1, 8,
1463 // C3[1], coeff of eps^1, polynomial in n of order 1
1464 -1, 1, 4,
1465 // C3[2], coeff of eps^3, polynomial in n of order 0
1466 3, 64,
1467 // C3[2], coeff of eps^2, polynomial in n of order 1
1468 -3, 2, 32,
1469 // C3[3], coeff of eps^3, polynomial in n of order 0
1470 5, 192,
1471 };
1472#elif GEOGRAPHICLIB_GEODESIC_ORDER == 5
1473 static const real coeff[] = {
1474 // C3[1], coeff of eps^4, polynomial in n of order 0
1475 5, 128,
1476 // C3[1], coeff of eps^3, polynomial in n of order 1
1477 3, 3, 64,
1478 // C3[1], coeff of eps^2, polynomial in n of order 2
1479 -1, 0, 1, 8,
1480 // C3[1], coeff of eps^1, polynomial in n of order 1
1481 -1, 1, 4,
1482 // C3[2], coeff of eps^4, polynomial in n of order 0
1483 3, 128,
1484 // C3[2], coeff of eps^3, polynomial in n of order 1
1485 -2, 3, 64,
1486 // C3[2], coeff of eps^2, polynomial in n of order 2
1487 1, -3, 2, 32,
1488 // C3[3], coeff of eps^4, polynomial in n of order 0
1489 3, 128,
1490 // C3[3], coeff of eps^3, polynomial in n of order 1
1491 -9, 5, 192,
1492 // C3[4], coeff of eps^4, polynomial in n of order 0
1493 7, 512,
1494 };
1495#elif GEOGRAPHICLIB_GEODESIC_ORDER == 6
1496 static const real coeff[] = {
1497 // C3[1], coeff of eps^5, polynomial in n of order 0
1498 3, 128,
1499 // C3[1], coeff of eps^4, polynomial in n of order 1
1500 2, 5, 128,
1501 // C3[1], coeff of eps^3, polynomial in n of order 2
1502 -1, 3, 3, 64,
1503 // C3[1], coeff of eps^2, polynomial in n of order 2
1504 -1, 0, 1, 8,
1505 // C3[1], coeff of eps^1, polynomial in n of order 1
1506 -1, 1, 4,
1507 // C3[2], coeff of eps^5, polynomial in n of order 0
1508 5, 256,
1509 // C3[2], coeff of eps^4, polynomial in n of order 1
1510 1, 3, 128,
1511 // C3[2], coeff of eps^3, polynomial in n of order 2
1512 -3, -2, 3, 64,
1513 // C3[2], coeff of eps^2, polynomial in n of order 2
1514 1, -3, 2, 32,
1515 // C3[3], coeff of eps^5, polynomial in n of order 0
1516 7, 512,
1517 // C3[3], coeff of eps^4, polynomial in n of order 1
1518 -10, 9, 384,
1519 // C3[3], coeff of eps^3, polynomial in n of order 2
1520 5, -9, 5, 192,
1521 // C3[4], coeff of eps^5, polynomial in n of order 0
1522 7, 512,
1523 // C3[4], coeff of eps^4, polynomial in n of order 1
1524 -14, 7, 512,
1525 // C3[5], coeff of eps^5, polynomial in n of order 0
1526 21, 2560,
1527 };
1528#elif GEOGRAPHICLIB_GEODESIC_ORDER == 7
1529 static const real coeff[] = {
1530 // C3[1], coeff of eps^6, polynomial in n of order 0
1531 21, 1024,
1532 // C3[1], coeff of eps^5, polynomial in n of order 1
1533 11, 12, 512,
1534 // C3[1], coeff of eps^4, polynomial in n of order 2
1535 2, 2, 5, 128,
1536 // C3[1], coeff of eps^3, polynomial in n of order 3
1537 -5, -1, 3, 3, 64,
1538 // C3[1], coeff of eps^2, polynomial in n of order 2
1539 -1, 0, 1, 8,
1540 // C3[1], coeff of eps^1, polynomial in n of order 1
1541 -1, 1, 4,
1542 // C3[2], coeff of eps^6, polynomial in n of order 0
1543 27, 2048,
1544 // C3[2], coeff of eps^5, polynomial in n of order 1
1545 1, 5, 256,
1546 // C3[2], coeff of eps^4, polynomial in n of order 2
1547 -9, 2, 6, 256,
1548 // C3[2], coeff of eps^3, polynomial in n of order 3
1549 2, -3, -2, 3, 64,
1550 // C3[2], coeff of eps^2, polynomial in n of order 2
1551 1, -3, 2, 32,
1552 // C3[3], coeff of eps^6, polynomial in n of order 0
1553 3, 256,
1554 // C3[3], coeff of eps^5, polynomial in n of order 1
1555 -4, 21, 1536,
1556 // C3[3], coeff of eps^4, polynomial in n of order 2
1557 -6, -10, 9, 384,
1558 // C3[3], coeff of eps^3, polynomial in n of order 3
1559 -1, 5, -9, 5, 192,
1560 // C3[4], coeff of eps^6, polynomial in n of order 0
1561 9, 1024,
1562 // C3[4], coeff of eps^5, polynomial in n of order 1
1563 -10, 7, 512,
1564 // C3[4], coeff of eps^4, polynomial in n of order 2
1565 10, -14, 7, 512,
1566 // C3[5], coeff of eps^6, polynomial in n of order 0
1567 9, 1024,
1568 // C3[5], coeff of eps^5, polynomial in n of order 1
1569 -45, 21, 2560,
1570 // C3[6], coeff of eps^6, polynomial in n of order 0
1571 11, 2048,
1572 };
1573#elif GEOGRAPHICLIB_GEODESIC_ORDER == 8
1574 static const real coeff[] = {
1575 // C3[1], coeff of eps^7, polynomial in n of order 0
1576 243, 16384,
1577 // C3[1], coeff of eps^6, polynomial in n of order 1
1578 10, 21, 1024,
1579 // C3[1], coeff of eps^5, polynomial in n of order 2
1580 3, 11, 12, 512,
1581 // C3[1], coeff of eps^4, polynomial in n of order 3
1582 -2, 2, 2, 5, 128,
1583 // C3[1], coeff of eps^3, polynomial in n of order 3
1584 -5, -1, 3, 3, 64,
1585 // C3[1], coeff of eps^2, polynomial in n of order 2
1586 -1, 0, 1, 8,
1587 // C3[1], coeff of eps^1, polynomial in n of order 1
1588 -1, 1, 4,
1589 // C3[2], coeff of eps^7, polynomial in n of order 0
1590 187, 16384,
1591 // C3[2], coeff of eps^6, polynomial in n of order 1
1592 69, 108, 8192,
1593 // C3[2], coeff of eps^5, polynomial in n of order 2
1594 -2, 1, 5, 256,
1595 // C3[2], coeff of eps^4, polynomial in n of order 3
1596 -6, -9, 2, 6, 256,
1597 // C3[2], coeff of eps^3, polynomial in n of order 3
1598 2, -3, -2, 3, 64,
1599 // C3[2], coeff of eps^2, polynomial in n of order 2
1600 1, -3, 2, 32,
1601 // C3[3], coeff of eps^7, polynomial in n of order 0
1602 139, 16384,
1603 // C3[3], coeff of eps^6, polynomial in n of order 1
1604 -1, 12, 1024,
1605 // C3[3], coeff of eps^5, polynomial in n of order 2
1606 -77, -8, 42, 3072,
1607 // C3[3], coeff of eps^4, polynomial in n of order 3
1608 10, -6, -10, 9, 384,
1609 // C3[3], coeff of eps^3, polynomial in n of order 3
1610 -1, 5, -9, 5, 192,
1611 // C3[4], coeff of eps^7, polynomial in n of order 0
1612 127, 16384,
1613 // C3[4], coeff of eps^6, polynomial in n of order 1
1614 -43, 72, 8192,
1615 // C3[4], coeff of eps^5, polynomial in n of order 2
1616 -7, -40, 28, 2048,
1617 // C3[4], coeff of eps^4, polynomial in n of order 3
1618 -7, 20, -28, 14, 1024,
1619 // C3[5], coeff of eps^7, polynomial in n of order 0
1620 99, 16384,
1621 // C3[5], coeff of eps^6, polynomial in n of order 1
1622 -15, 9, 1024,
1623 // C3[5], coeff of eps^5, polynomial in n of order 2
1624 75, -90, 42, 5120,
1625 // C3[6], coeff of eps^7, polynomial in n of order 0
1626 99, 16384,
1627 // C3[6], coeff of eps^6, polynomial in n of order 1
1628 -99, 44, 8192,
1629 // C3[7], coeff of eps^7, polynomial in n of order 0
1630 429, 114688,
1631 };
1632#else
1633#error "Bad value for GEOGRAPHICLIB_GEODESIC_ORDER"
1634#endif
1635 static_assert(sizeof(coeff) / sizeof(real) ==
1636 ((nC3_-1)*(nC3_*nC3_ + 7*nC3_ - 2*(nC3_/2)))/8,
1637 "Coefficient array size mismatch in C3coeff");
1638 int o = 0, k = 0;
1639 for (int l = 1; l < nC3_; ++l) { // l is index of C3[l]
1640 for (int j = nC3_ - 1; j >= l; --j) { // coeff of eps^j
1641 int m = min(nC3_ - j - 1, j); // order of polynomial in n
1642 _cC3x[k++] = Math::polyval(m, coeff + o, _n) / coeff[o + m + 1];
1643 o += m + 2;
1644 }
1645 }
1646 // Post condition: o == sizeof(coeff) / sizeof(real) && k == nC3x_
1647 }
1648
1649 void Geodesic::C4coeff() {
1650 // Generated by Maxima on 2015-05-05 18:08:13-04:00
1651#if GEOGRAPHICLIB_GEODESIC_ORDER == 3
1652 static const real coeff[] = {
1653 // C4[0], coeff of eps^2, polynomial in n of order 0
1654 -2, 105,
1655 // C4[0], coeff of eps^1, polynomial in n of order 1
1656 16, -7, 35,
1657 // C4[0], coeff of eps^0, polynomial in n of order 2
1658 8, -28, 70, 105,
1659 // C4[1], coeff of eps^2, polynomial in n of order 0
1660 -2, 105,
1661 // C4[1], coeff of eps^1, polynomial in n of order 1
1662 -16, 7, 315,
1663 // C4[2], coeff of eps^2, polynomial in n of order 0
1664 4, 525,
1665 };
1666#elif GEOGRAPHICLIB_GEODESIC_ORDER == 4
1667 static const real coeff[] = {
1668 // C4[0], coeff of eps^3, polynomial in n of order 0
1669 11, 315,
1670 // C4[0], coeff of eps^2, polynomial in n of order 1
1671 -32, -6, 315,
1672 // C4[0], coeff of eps^1, polynomial in n of order 2
1673 -32, 48, -21, 105,
1674 // C4[0], coeff of eps^0, polynomial in n of order 3
1675 4, 24, -84, 210, 315,
1676 // C4[1], coeff of eps^3, polynomial in n of order 0
1677 -1, 105,
1678 // C4[1], coeff of eps^2, polynomial in n of order 1
1679 64, -18, 945,
1680 // C4[1], coeff of eps^1, polynomial in n of order 2
1681 32, -48, 21, 945,
1682 // C4[2], coeff of eps^3, polynomial in n of order 0
1683 -8, 1575,
1684 // C4[2], coeff of eps^2, polynomial in n of order 1
1685 -32, 12, 1575,
1686 // C4[3], coeff of eps^3, polynomial in n of order 0
1687 8, 2205,
1688 };
1689#elif GEOGRAPHICLIB_GEODESIC_ORDER == 5
1690 static const real coeff[] = {
1691 // C4[0], coeff of eps^4, polynomial in n of order 0
1692 4, 1155,
1693 // C4[0], coeff of eps^3, polynomial in n of order 1
1694 -368, 121, 3465,
1695 // C4[0], coeff of eps^2, polynomial in n of order 2
1696 1088, -352, -66, 3465,
1697 // C4[0], coeff of eps^1, polynomial in n of order 3
1698 48, -352, 528, -231, 1155,
1699 // C4[0], coeff of eps^0, polynomial in n of order 4
1700 16, 44, 264, -924, 2310, 3465,
1701 // C4[1], coeff of eps^4, polynomial in n of order 0
1702 4, 1155,
1703 // C4[1], coeff of eps^3, polynomial in n of order 1
1704 80, -99, 10395,
1705 // C4[1], coeff of eps^2, polynomial in n of order 2
1706 -896, 704, -198, 10395,
1707 // C4[1], coeff of eps^1, polynomial in n of order 3
1708 -48, 352, -528, 231, 10395,
1709 // C4[2], coeff of eps^4, polynomial in n of order 0
1710 -8, 1925,
1711 // C4[2], coeff of eps^3, polynomial in n of order 1
1712 384, -88, 17325,
1713 // C4[2], coeff of eps^2, polynomial in n of order 2
1714 320, -352, 132, 17325,
1715 // C4[3], coeff of eps^4, polynomial in n of order 0
1716 -16, 8085,
1717 // C4[3], coeff of eps^3, polynomial in n of order 1
1718 -256, 88, 24255,
1719 // C4[4], coeff of eps^4, polynomial in n of order 0
1720 64, 31185,
1721 };
1722#elif GEOGRAPHICLIB_GEODESIC_ORDER == 6
1723 static const real coeff[] = {
1724 // C4[0], coeff of eps^5, polynomial in n of order 0
1725 97, 15015,
1726 // C4[0], coeff of eps^4, polynomial in n of order 1
1727 1088, 156, 45045,
1728 // C4[0], coeff of eps^3, polynomial in n of order 2
1729 -224, -4784, 1573, 45045,
1730 // C4[0], coeff of eps^2, polynomial in n of order 3
1731 -10656, 14144, -4576, -858, 45045,
1732 // C4[0], coeff of eps^1, polynomial in n of order 4
1733 64, 624, -4576, 6864, -3003, 15015,
1734 // C4[0], coeff of eps^0, polynomial in n of order 5
1735 100, 208, 572, 3432, -12012, 30030, 45045,
1736 // C4[1], coeff of eps^5, polynomial in n of order 0
1737 1, 9009,
1738 // C4[1], coeff of eps^4, polynomial in n of order 1
1739 -2944, 468, 135135,
1740 // C4[1], coeff of eps^3, polynomial in n of order 2
1741 5792, 1040, -1287, 135135,
1742 // C4[1], coeff of eps^2, polynomial in n of order 3
1743 5952, -11648, 9152, -2574, 135135,
1744 // C4[1], coeff of eps^1, polynomial in n of order 4
1745 -64, -624, 4576, -6864, 3003, 135135,
1746 // C4[2], coeff of eps^5, polynomial in n of order 0
1747 8, 10725,
1748 // C4[2], coeff of eps^4, polynomial in n of order 1
1749 1856, -936, 225225,
1750 // C4[2], coeff of eps^3, polynomial in n of order 2
1751 -8448, 4992, -1144, 225225,
1752 // C4[2], coeff of eps^2, polynomial in n of order 3
1753 -1440, 4160, -4576, 1716, 225225,
1754 // C4[3], coeff of eps^5, polynomial in n of order 0
1755 -136, 63063,
1756 // C4[3], coeff of eps^4, polynomial in n of order 1
1757 1024, -208, 105105,
1758 // C4[3], coeff of eps^3, polynomial in n of order 2
1759 3584, -3328, 1144, 315315,
1760 // C4[4], coeff of eps^5, polynomial in n of order 0
1761 -128, 135135,
1762 // C4[4], coeff of eps^4, polynomial in n of order 1
1763 -2560, 832, 405405,
1764 // C4[5], coeff of eps^5, polynomial in n of order 0
1765 128, 99099,
1766 };
1767#elif GEOGRAPHICLIB_GEODESIC_ORDER == 7
1768 static const real coeff[] = {
1769 // C4[0], coeff of eps^6, polynomial in n of order 0
1770 10, 9009,
1771 // C4[0], coeff of eps^5, polynomial in n of order 1
1772 -464, 291, 45045,
1773 // C4[0], coeff of eps^4, polynomial in n of order 2
1774 -4480, 1088, 156, 45045,
1775 // C4[0], coeff of eps^3, polynomial in n of order 3
1776 10736, -224, -4784, 1573, 45045,
1777 // C4[0], coeff of eps^2, polynomial in n of order 4
1778 1664, -10656, 14144, -4576, -858, 45045,
1779 // C4[0], coeff of eps^1, polynomial in n of order 5
1780 16, 64, 624, -4576, 6864, -3003, 15015,
1781 // C4[0], coeff of eps^0, polynomial in n of order 6
1782 56, 100, 208, 572, 3432, -12012, 30030, 45045,
1783 // C4[1], coeff of eps^6, polynomial in n of order 0
1784 10, 9009,
1785 // C4[1], coeff of eps^5, polynomial in n of order 1
1786 112, 15, 135135,
1787 // C4[1], coeff of eps^4, polynomial in n of order 2
1788 3840, -2944, 468, 135135,
1789 // C4[1], coeff of eps^3, polynomial in n of order 3
1790 -10704, 5792, 1040, -1287, 135135,
1791 // C4[1], coeff of eps^2, polynomial in n of order 4
1792 -768, 5952, -11648, 9152, -2574, 135135,
1793 // C4[1], coeff of eps^1, polynomial in n of order 5
1794 -16, -64, -624, 4576, -6864, 3003, 135135,
1795 // C4[2], coeff of eps^6, polynomial in n of order 0
1796 -4, 25025,
1797 // C4[2], coeff of eps^5, polynomial in n of order 1
1798 -1664, 168, 225225,
1799 // C4[2], coeff of eps^4, polynomial in n of order 2
1800 1664, 1856, -936, 225225,
1801 // C4[2], coeff of eps^3, polynomial in n of order 3
1802 6784, -8448, 4992, -1144, 225225,
1803 // C4[2], coeff of eps^2, polynomial in n of order 4
1804 128, -1440, 4160, -4576, 1716, 225225,
1805 // C4[3], coeff of eps^6, polynomial in n of order 0
1806 64, 315315,
1807 // C4[3], coeff of eps^5, polynomial in n of order 1
1808 1792, -680, 315315,
1809 // C4[3], coeff of eps^4, polynomial in n of order 2
1810 -2048, 1024, -208, 105105,
1811 // C4[3], coeff of eps^3, polynomial in n of order 3
1812 -1792, 3584, -3328, 1144, 315315,
1813 // C4[4], coeff of eps^6, polynomial in n of order 0
1814 -512, 405405,
1815 // C4[4], coeff of eps^5, polynomial in n of order 1
1816 2048, -384, 405405,
1817 // C4[4], coeff of eps^4, polynomial in n of order 2
1818 3072, -2560, 832, 405405,
1819 // C4[5], coeff of eps^6, polynomial in n of order 0
1820 -256, 495495,
1821 // C4[5], coeff of eps^5, polynomial in n of order 1
1822 -2048, 640, 495495,
1823 // C4[6], coeff of eps^6, polynomial in n of order 0
1824 512, 585585,
1825 };
1826#elif GEOGRAPHICLIB_GEODESIC_ORDER == 8
1827 static const real coeff[] = {
1828 // C4[0], coeff of eps^7, polynomial in n of order 0
1829 193, 85085,
1830 // C4[0], coeff of eps^6, polynomial in n of order 1
1831 4192, 850, 765765,
1832 // C4[0], coeff of eps^5, polynomial in n of order 2
1833 20960, -7888, 4947, 765765,
1834 // C4[0], coeff of eps^4, polynomial in n of order 3
1835 12480, -76160, 18496, 2652, 765765,
1836 // C4[0], coeff of eps^3, polynomial in n of order 4
1837 -154048, 182512, -3808, -81328, 26741, 765765,
1838 // C4[0], coeff of eps^2, polynomial in n of order 5
1839 3232, 28288, -181152, 240448, -77792, -14586, 765765,
1840 // C4[0], coeff of eps^1, polynomial in n of order 6
1841 96, 272, 1088, 10608, -77792, 116688, -51051, 255255,
1842 // C4[0], coeff of eps^0, polynomial in n of order 7
1843 588, 952, 1700, 3536, 9724, 58344, -204204, 510510, 765765,
1844 // C4[1], coeff of eps^7, polynomial in n of order 0
1845 349, 2297295,
1846 // C4[1], coeff of eps^6, polynomial in n of order 1
1847 -1472, 510, 459459,
1848 // C4[1], coeff of eps^5, polynomial in n of order 2
1849 -39840, 1904, 255, 2297295,
1850 // C4[1], coeff of eps^4, polynomial in n of order 3
1851 52608, 65280, -50048, 7956, 2297295,
1852 // C4[1], coeff of eps^3, polynomial in n of order 4
1853 103744, -181968, 98464, 17680, -21879, 2297295,
1854 // C4[1], coeff of eps^2, polynomial in n of order 5
1855 -1344, -13056, 101184, -198016, 155584, -43758, 2297295,
1856 // C4[1], coeff of eps^1, polynomial in n of order 6
1857 -96, -272, -1088, -10608, 77792, -116688, 51051, 2297295,
1858 // C4[2], coeff of eps^7, polynomial in n of order 0
1859 464, 1276275,
1860 // C4[2], coeff of eps^6, polynomial in n of order 1
1861 -928, -612, 3828825,
1862 // C4[2], coeff of eps^5, polynomial in n of order 2
1863 64256, -28288, 2856, 3828825,
1864 // C4[2], coeff of eps^4, polynomial in n of order 3
1865 -126528, 28288, 31552, -15912, 3828825,
1866 // C4[2], coeff of eps^3, polynomial in n of order 4
1867 -41472, 115328, -143616, 84864, -19448, 3828825,
1868 // C4[2], coeff of eps^2, polynomial in n of order 5
1869 160, 2176, -24480, 70720, -77792, 29172, 3828825,
1870 // C4[3], coeff of eps^7, polynomial in n of order 0
1871 -16, 97461,
1872 // C4[3], coeff of eps^6, polynomial in n of order 1
1873 -16384, 1088, 5360355,
1874 // C4[3], coeff of eps^5, polynomial in n of order 2
1875 -2560, 30464, -11560, 5360355,
1876 // C4[3], coeff of eps^4, polynomial in n of order 3
1877 35840, -34816, 17408, -3536, 1786785,
1878 // C4[3], coeff of eps^3, polynomial in n of order 4
1879 7168, -30464, 60928, -56576, 19448, 5360355,
1880 // C4[4], coeff of eps^7, polynomial in n of order 0
1881 128, 2297295,
1882 // C4[4], coeff of eps^6, polynomial in n of order 1
1883 26624, -8704, 6891885,
1884 // C4[4], coeff of eps^5, polynomial in n of order 2
1885 -77824, 34816, -6528, 6891885,
1886 // C4[4], coeff of eps^4, polynomial in n of order 3
1887 -32256, 52224, -43520, 14144, 6891885,
1888 // C4[5], coeff of eps^7, polynomial in n of order 0
1889 -6784, 8423415,
1890 // C4[5], coeff of eps^6, polynomial in n of order 1
1891 24576, -4352, 8423415,
1892 // C4[5], coeff of eps^5, polynomial in n of order 2
1893 45056, -34816, 10880, 8423415,
1894 // C4[6], coeff of eps^7, polynomial in n of order 0
1895 -1024, 3318315,
1896 // C4[6], coeff of eps^6, polynomial in n of order 1
1897 -28672, 8704, 9954945,
1898 // C4[7], coeff of eps^7, polynomial in n of order 0
1899 1024, 1640925,
1900 };
1901#else
1902#error "Bad value for GEOGRAPHICLIB_GEODESIC_ORDER"
1903#endif
1904 static_assert(sizeof(coeff) / sizeof(real) ==
1905 (nC4_ * (nC4_ + 1) * (nC4_ + 5)) / 6,
1906 "Coefficient array size mismatch in C4coeff");
1907 int o = 0, k = 0;
1908 for (int l = 0; l < nC4_; ++l) { // l is index of C4[l]
1909 for (int j = nC4_ - 1; j >= l; --j) { // coeff of eps^j
1910 int m = nC4_ - j - 1; // order of polynomial in n
1911 _cC4x[k++] = Math::polyval(m, coeff + o, _n) / coeff[o + m + 1];
1912 o += m + 2;
1913 }
1914 }
1915 // Post condition: o == sizeof(coeff) / sizeof(real) && k == nC4x_
1916 }
1917
1918} // namespace GeographicLib
GeographicLib::Math::real real
Definition GeodSolve.cpp:28
Header for GeographicLib::GeodesicLine class.
Header for GeographicLib::Geodesic class.
Exact geodesic calculations.
Math::real GenDirect(real lat1, real lon1, real azi1, bool arcmode, real s12_a12, unsigned outmask, real &lat2, real &lon2, real &azi2, real &s12, real &m12, real &M12, real &M21, real &S12) const
Geodesic calculations
Definition Geodesic.hpp:175
GeodesicLine InverseLine(real lat1, real lon1, real lat2, real lon2, unsigned caps=ALL) const
Definition Geodesic.cpp:532
static const Geodesic & WGS84()
Definition Geodesic.cpp:94
GeodesicLine ArcDirectLine(real lat1, real lon1, real azi1, real a12, unsigned caps=ALL) const
Definition Geodesic.cpp:163
GeodesicLine Line(real lat1, real lon1, real azi1, unsigned caps=ALL) const
Definition Geodesic.cpp:123
GeodesicLine GenDirectLine(real lat1, real lon1, real azi1, bool arcmode, real s12_a12, unsigned caps=ALL) const
Definition Geodesic.cpp:145
friend class GeodesicLine
Definition Geodesic.hpp:178
Math::real GenDirect(real lat1, real lon1, real azi1, bool arcmode, real s12_a12, unsigned outmask, real &lat2, real &lon2, real &azi2, real &s12, real &m12, real &M12, real &M21, real &S12) const
Definition Geodesic.cpp:128
GeodesicLine DirectLine(real lat1, real lon1, real azi1, real s12, unsigned caps=ALL) const
Definition Geodesic.cpp:158
Geodesic(real a, real f, bool exact=false)
Definition Geodesic.cpp:41
Exception handling for GeographicLib.
Mathematical functions needed by GeographicLib.
Definition Math.hpp:77
static T degree()
Definition Math.hpp:209
static T LatFix(T x)
Definition Math.hpp:309
static void sincosd(T x, T &sinx, T &cosx)
Definition Math.cpp:101
static T atan2d(T y, T x)
Definition Math.cpp:199
static void norm(T &x, T &y)
Definition Math.hpp:231
static T AngRound(T x)
Definition Math.cpp:92
static T sq(T x)
Definition Math.hpp:221
static constexpr int qd
degrees per quarter turn
Definition Math.hpp:145
static T AngNormalize(T x)
Definition Math.cpp:66
static void sincosde(T x, T t, T &sinx, T &cosx)
Definition Math.cpp:128
static T pi()
Definition Math.hpp:199
static T NaN()
Definition Math.cpp:277
static T polyval(int N, const T p[], T x)
Definition Math.hpp:280
static T AngDiff(T x, T y, T &e)
Definition Math.cpp:77
static constexpr int hd
degrees per half turn
Definition Math.hpp:148
Namespace for GeographicLib.