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


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:42:07