GeographicLib 2.3
AuxLatitude.cpp
Go to the documentation of this file.
1/**
2 * \file AuxLatitude.cpp
3 * \brief Implementation for the GeographicLib::AuxLatitude class.
4 *
5 * \note This is just sample code. It is not part of GeographicLib itself.
6 *
7 * This file is an implementation of the methods described in
8 * - C. F. F. Karney,
9 * <a href="https://doi.org/10.1080/00396265.2023.2217604">
10 * On auxiliary latitudes,</a>
11 * Survey Review (2023);
12 * preprint
13 * <a href="https://arxiv.org/abs/2212.05818">arXiv:2212.05818</a>.
14 * .
15 * Copyright (c) Charles Karney (2022-2023) <karney@alum.mit.edu> and licensed
16 * under the MIT/X11 License. For more information, see
17 * https://geographiclib.sourceforge.io/
18 **********************************************************************/
19
22
23#if defined(_MSC_VER)
24// Squelch warnings about constant conditional and enum-float expressions
25# pragma warning (disable: 4127 5055)
26#endif
27
28namespace GeographicLib {
29
30 using namespace std;
31
32 AuxLatitude::AuxLatitude(real a, real f)
33 : tol_( sqrt(numeric_limits<real>::epsilon()) )
34 , bmin_( log2(numeric_limits<real>::min()) )
35 , bmax_( log2(numeric_limits<real>::max()) )
36 , _a(a)
37 , _b(_a * (1 - f))
38 , _f( f )
39 , _fm1( 1 - _f )
40 , _e2( _f * (2 - _f) )
41 , _e2m1( _fm1 * _fm1 )
42 , _e12( _e2/(1 - _e2) )
43 , _e12p1( 1 / _e2m1 )
44 , _n( _f/(2 - _f) )
45 , _e( sqrt(fabs(_e2)) )
46 , _e1( sqrt(fabs(_e12)) )
47 , _n2( _n * _n )
48 , _q( _e12p1 + (_f == 0 ? 1 : (_f > 0 ? asinh(_e1) : atan(_e)) / _e) )
49 {
50 if (!(isfinite(_a) && _a > 0))
51 throw GeographicErr("Equatorial radius is not positive");
52 if (!(isfinite(_b) && _b > 0))
53 throw GeographicErr("Polar semi-axis is not positive");
54 fill(_c, _c + Lmax * AUXNUMBER * AUXNUMBER,
55 numeric_limits<real>::quiet_NaN());
56 }
57
58 /// \cond SKIP
59 AuxLatitude::AuxLatitude(const pair<real, real>& axes)
60 : tol_( sqrt(numeric_limits<real>::epsilon()) )
61 , bmin_( log2(numeric_limits<real>::min()) )
62 , bmax_( log2(numeric_limits<real>::max()) )
63 , _a(axes.first)
64 , _b(axes.second)
65 , _f( (_a - _b) / _a )
66 , _fm1( _b / _a )
67 , _e2( ((_a - _b) * (_a + _b)) / (_a * _a) )
68 , _e2m1( (_b * _b) / (_a * _a) )
69 , _e12( ((_a - _b) * (_a + _b)) / (_b * _b) )
70 , _e12p1( (_a * _a) / (_b * _b) )
71 , _n( (_a - _b) / (_a + _b) )
72 , _e( sqrt(fabs(_a - _b) * (_a + _b)) / _a )
73 , _e1( sqrt(fabs(_a - _b) * (_a + _b)) / _b )
74 , _n2( _n * _n )
75 , _q( _e12p1 + (_f == 0 ? 1 : (_f > 0 ? asinh(_e1) : atan(_e)) / _e) )
76 {
77 if (!(isfinite(_a) && _a > 0))
78 throw GeographicErr("Equatorial radius is not positive");
79 if (!(isfinite(_b) && _b > 0))
80 throw GeographicErr("Polar semi-axis is not positive");
81 fill(_c, _c + Lmax * AUXNUMBER * AUXNUMBER,
82 numeric_limits<real>::quiet_NaN());
83 }
84 /// \endcond
85
87 static const AuxLatitude wgs84(Constants::WGS84_a(), Constants::WGS84_f());
88 return wgs84;
89 }
90
91 AuxAngle AuxLatitude::Parametric(const AuxAngle& phi, real* diff) const {
92 if (diff) *diff = _fm1;
93 return AuxAngle(phi.y() * _fm1, phi.x());
94 }
95
96 AuxAngle AuxLatitude::Geocentric(const AuxAngle& phi, real* diff) const {
97 if (diff) *diff = _e2m1;
98 return AuxAngle(phi.y() * _e2m1, phi.x());
99 }
100
101 AuxAngle AuxLatitude::Rectifying(const AuxAngle& phi, real* diff) const {
102 using std::isinf; // Needed for Centos 7, ubuntu 14
103 AuxAngle beta(Parametric(phi).normalized());
104 real sbeta = fabs(beta.y()), cbeta = fabs(beta.x());
105 real a = 1, b = _fm1, ka = _e2, kb = -_e12, ka1 = _e2m1, kb1 = _e12p1,
106 smu, cmu, mr;
107 if (_f < 0) {
108 swap(a, b); swap(ka, kb); swap(ka1, kb1); swap(sbeta, cbeta);
109 }
110 // now a,b = larger/smaller semiaxis
111 // beta now measured from larger semiaxis
112 // kb,ka = modulus-squared for distance from beta = 0,pi/2
113 // NB kb <= 0; 0 <= ka <= 1
114 // sa = b*E(beta,sqrt(kb)), sb = a*E(beta',sqrt(ka))
115 // 1 - ka * (1 - sb2) = 1 -ka + ka*sb2
116 real
117 sb2 = sbeta * sbeta,
118 cb2 = cbeta * cbeta,
119 db2 = 1 - kb * sb2,
120 da2 = ka1 + ka * sb2,
121 // DLMF Eq. 19.25.9
122 sa = b * sbeta * ( EllipticFunction::RF(cb2, db2, 1) -
123 kb * sb2 * EllipticFunction::RD(cb2, db2, 1)/3 ),
124 // DLMF Eq. 19.25.10 with complementary angles
125 sb = a * cbeta * ( ka1 * EllipticFunction::RF(sb2, da2, 1)
126 + ka * ka1 * cb2 * EllipticFunction::RD(sb2, 1, da2)/3
127 + ka * sbeta / sqrt(da2) );
128 // sa + sb = 2*EllipticFunction::RG(a*a, b*b) = a*E(e) = b*E(i*e')
129 // mr = a*E(e)*(2/pi) = b*E(i*e')*(2/pi)
130 mr = (2 * (sa + sb)) / Math::pi();
131 smu = sin(sa / mr);
132 cmu = sin(sb / mr);
133 if (_f < 0) { swap(smu, cmu); swap(a, b); }
134 // mu is normalized
135 AuxAngle mu(AuxAngle(smu, cmu).copyquadrant(phi));
136 if (diff) {
137 real cphi = phi.normalized().x(), tphi = phi.tan();
138 if (!isinf(tphi)) {
139 cmu = mu.x(); cbeta = beta.x();
140 *diff = _fm1 * b/mr * Math::sq(cbeta / cmu) * (cbeta / cphi);
141 } else
142 *diff = _fm1 * mr/a;
143 }
144 return mu;
145 }
146
147 AuxAngle AuxLatitude::Conformal(const AuxAngle& phi, real* diff) const {
148 using std::isinf; // Needed for Centos 7, ubuntu 14
149 real tphi = fabs(phi.tan()), tchi = tphi;
150 if ( !( !isfinite(tphi) || tphi == 0 || _f == 0 ) ) {
151 real scphi = sc(tphi),
152 sig = sinh(_e2 * atanhee(tphi) ),
153 scsig = sc(sig);
154 if (_f <= 0) {
155 tchi = tphi * scsig - sig * scphi;
156 } else {
157 // The general expression for tchi is
158 // tphi * scsig - sig * scphi
159 // This involves cancellation if f > 0, so change to
160 // (tphi - sig) * (tphi + sig) / (tphi * scsig + sig * scphi)
161 // To control overflow, write as (sigtphi = sig / tphi)
162 // (tphi - sig) * (1 + sigtphi) / (scsig + sigtphi * scphi)
163 real sigtphi = sig / tphi, tphimsig;
164 if (sig < tphi / 2)
165 tphimsig = tphi - sig;
166 else {
167 // Still have possibly dangerous cancellation in tphi - sig.
168 //
169 // Write tphi - sig = (1 - e) * Dg(1, e)
170 // Dg(x, y) = (g(x) - g(y)) / (x - y)
171 // g(x) = sinh(x * atanh(sphi * x))
172 // Note sinh(atanh(sphi)) = tphi
173 // Turn the crank on divided differences, substitute
174 // sphi = tphi/sc(tphi)
175 // atanh(x) = asinh(x/sqrt(1-x^2))
176 real em1 = _e2m1 / (1 + _e), // 1 - e
177 atanhs = asinh(tphi), // atanh(sphi)
178 scbeta = sc(_fm1 * tphi), // sec(beta)
179 scphibeta = sc(tphi) / scbeta, // sec(phi)/sec(beta)
180 atanhes = asinh(_e * tphi / scbeta), // atanh(e * sphi)
181 t1 = (atanhs - _e * atanhes)/2,
182 t2 = asinh(em1 * (tphi * scphibeta)) / em1,
183 Dg = cosh((atanhs + _e * atanhes)/2) * (sinh(t1) / t1)
184 * ((atanhs + atanhes)/2 + (1 + _e)/2 * t2);
185 tphimsig = em1 * Dg; // tphi - sig
186 }
187 tchi = tphimsig * (1 + sigtphi) / (scsig + sigtphi * scphi);
188 }
189 }
190 AuxAngle chi(AuxAngle(tchi).copyquadrant(phi));
191 if (diff) {
192 if (!isinf(tphi)) {
193 real cchi = chi.normalized().x(),
194 cphi = phi.normalized().x(),
195 cbeta = Parametric(phi).normalized().x();
196 *diff = _e2m1 * (cbeta / cchi) * (cbeta / cphi);
197 } else {
198 real ss = _f > 0 ? sinh(_e * asinh(_e1)) : sinh(-_e * atan(_e));
199 *diff = _f > 0 ? 1/( sc(ss) + ss ) : sc(ss) - ss;
200 }
201 }
202 return chi;
203 }
204
205 AuxAngle AuxLatitude::Authalic(const AuxAngle& phi, real* diff) const {
206 using std::isnan; // Needed for Centos 7, ubuntu 14
207 real tphi = fabs(phi.tan());
208 AuxAngle xi(phi), phin(phi.normalized());
209 if ( !( !isfinite(tphi) || tphi == 0 || _f == 0 ) ) {
210 real qv = q(tphi),
211 Dqp = Dq(tphi),
212 Dqm = (_q + qv) / (1 + fabs(phin.y())); // Dq(-tphi)
213 xi = AuxAngle( copysign(qv, phi.y()), phin.x() * sqrt(Dqp * Dqm) );
214 }
215 if (diff) {
216 if (!isnan(tphi)) {
217 real cbeta = Parametric(phi).normalized().x(),
218 cxi = xi.normalized().x();
219 *diff =
220 (2/_q) * Math::sq(cbeta / cxi) * (cbeta / cxi) * (cbeta / phin.x());
221 } else
222 *diff = _e2m1 * sqrt(_q/2);
223 }
224 return xi;
225 }
226
228 real* diff) const {
229 switch (auxout) {
230 case GEOGRAPHIC: if (diff) *diff = 1; return phi; break;
231 case PARAMETRIC: return Parametric(phi, diff); break;
232 case GEOCENTRIC: return Geocentric(phi, diff); break;
233 case RECTIFYING: return Rectifying(phi, diff); break;
234 case CONFORMAL : return Conformal (phi, diff); break;
235 case AUTHALIC : return Authalic (phi, diff); break;
236 default:
237 if (diff) *diff = numeric_limits<real>::quiet_NaN();
238 return AuxAngle::NaN();
239 break;
240 }
241 }
242
244 int* niter) const {
245 int n = 0; if (niter) *niter = n;
246 real tphi = _fm1;
247 switch (auxin) {
248 case GEOGRAPHIC: return zeta; break;
249 // case PARAMETRIC: break;
250 case PARAMETRIC: return AuxAngle(zeta.y() / _fm1, zeta.x()); break;
251 // case GEOCENTRIC: tphi *= _fm1 ; break;
252 case GEOCENTRIC: return AuxAngle(zeta.y() / _e2m1, zeta.x()); break;
253 case RECTIFYING: tphi *= sqrt(_fm1); break;
254 case CONFORMAL : tphi *= _fm1 ; break;
255 case AUTHALIC : tphi *= cbrt(_fm1); break;
256 default: return AuxAngle::NaN(); break;
257 }
258
259 // Drop through to solution by Newton's method
260 real tzeta = fabs(zeta.tan()), ltzeta = log2(tzeta);
261 if (!isfinite(ltzeta)) return zeta;
262 tphi = tzeta / tphi;
263 real ltphi = log2(tphi),
264 bmin = fmin(ltphi, bmin_), bmax = fmax(ltphi, bmax_);
265 for (int sign = 0, osign = 0, ntrip = 0; n < numit_;) {
266 ++n;
267 real diff;
268 AuxAngle zeta1(ToAuxiliary(auxin, AuxAngle(tphi), &diff));
269 real tzeta1 = zeta1.tan(), ltzeta1 = log2(tzeta1);
270 // Convert derivative from dtan(zeta)/dtan(phi) to
271 // dlog(tan(zeta))/dlog(tan(phi))
272 diff *= tphi/tzeta1;
273 osign = sign;
274 if (tzeta1 == tzeta)
275 break;
276 else if (tzeta1 > tzeta) {
277 sign = 1;
278 bmax = ltphi;
279 } else {
280 sign = -1;
281 bmin = ltphi;
282 }
283 real dltphi = -(ltzeta1 - ltzeta) / diff;
284 ltphi += dltphi;
285 tphi = exp2(ltphi);
286 if (!(fabs(dltphi) >= tol_)) {
287 ++n;
288 // Final Newton iteration without the logs
289 zeta1 = ToAuxiliary(auxin, AuxAngle(tphi), &diff);
290 tphi -= (zeta1.tan() - tzeta) / diff;
291 break;
292 }
293 if ((sign * osign < 0 && n - ntrip > 2) ||
294 ltphi >= bmax || ltphi <= bmin) {
295 sign = 0; ntrip = n;
296 ltphi = (bmin + bmax) / 2;
297 tphi = exp2(ltphi);
298 }
299 }
300 if (niter) *niter = n;
301 return AuxAngle(tphi).copyquadrant(zeta);
302 }
303
304 AuxAngle AuxLatitude::Convert(int auxin, int auxout, const AuxAngle& zeta,
305 bool exact) const {
306 using std::isnan; // Needed for Centos 7, ubuntu 14
307 int k = ind(auxout, auxin);
308 if (k < 0) return AuxAngle::NaN();
309 if (auxin == auxout) return zeta;
310 if (exact) {
311 if (auxin < 3 && auxout < 3)
312 // Need extra real because, since C++11, pow(float, int) returns double
313 return AuxAngle(zeta.y() * real(pow(_fm1, auxout - auxin)), zeta.x());
314 else
315 return ToAuxiliary(auxout, FromAuxiliary(auxin, zeta));
316 } else {
317 if ( isnan(_c[Lmax * (k + 1) - 1]) ) fillcoeff(auxin, auxout, k);
318 AuxAngle zetan(zeta.normalized());
319 real d = Clenshaw(true, zetan.y(), zetan.x(), _c + Lmax * k, Lmax);
320 zetan += AuxAngle::radians(d);
321 return zetan;
322 }
323 }
324
325 Math::real AuxLatitude::Convert(int auxin, int auxout, real zeta,
326 bool exact) const {
327 AuxAngle zetaa(AuxAngle::degrees(zeta));
328 real m = round((zeta - zetaa.degrees()) / Math::td);
329 return Math::td * m + Convert(auxin, auxout, zetaa, exact).degrees();
330 }
331
333 if (exact) {
334 return EllipticFunction::RG(Math::sq(_a), Math::sq(_b)) * 4 / Math::pi();
335 } else {
336 // Maxima code for these coefficients:
337 // df[i]:=if i<0 then df[i+2]/(i+2) else i!!$
338 // R(Lmax):=sum((df[2*j-3]/df[2*j])^2*n^(2*j),j,0,floor(Lmax/2))$
339 // cf(Lmax):=block([t:R(Lmax)],
340 // t:makelist(coeff(t,n,2*(floor(Lmax/2)-j)),j,0,floor(Lmax/2)),
341 // map(lambda([x],num(x)/
342 // (if denom(x) = 1 then 1 else real(denom(x)))),t))$
343#if GEOGRAPHICLIB_AUXLATITUDE_ORDER == 4
344 static const real coeff[] = {1/real(64), 1/real(4), 1};
345#elif GEOGRAPHICLIB_AUXLATITUDE_ORDER == 6
346 static const real coeff[] = {1/real(256), 1/real(64), 1/real(4), 1};
347#elif GEOGRAPHICLIB_AUXLATITUDE_ORDER == 8
348 static const real coeff[] = {
349 25/real(16384), 1/real(256), 1/real(64), 1/real(4), 1
350 };
351#else
352#error "Unsupported value for GEOGRAPHICLIB_AUXLATITUDE_ORDER"
353#endif
354 int m = Lmax/2;
355 return (_a + _b) / 2 * Math::polyval(m, coeff, _n2);
356 }
357 }
358
360 if (exact) {
361 return Math::sq(_b) * _q / 2;
362 } else {
363 // Using a * (a + b) / 2 as the multiplying factor leads to a rapidly
364 // converging series in n. Of course, using this series isn't really
365 // necessary, since the exact expression is simple to evaluate. However,
366 // we do it for consistency with RectifyingRadius; and, presumably, the
367 // roundoff error is smaller compared to that for the exact expression.
368 //
369 // Maxima code for these coefficients:
370 // c2:subst(e=2*sqrt(n)/(1+n),
371 // (atanh(e)/e * (1-n)^2 + (1+n)^2)/(2*(1+n)))$
372 // cf(Lmax):=block([t:expand(ratdisrep(taylor(c2,n,0,Lmax)))],
373 // t:makelist(coeff(t,n,Lmax-j),j,0,Lmax),
374 // map(lambda([x],num(x)/
375 // (if denom(x) = 1 then 1 else real(denom(x)))),t))$
376 // N.B. Coeff of n^j = 1 for j = 0
377 // -1/3 for j = 1
378 // 4*(2*j-5)!!/(2*j+1)!! for j > 1
379#if GEOGRAPHICLIB_AUXLATITUDE_ORDER == 4
380 static const real coeff[] = {
381 4/real(315), 4/real(105), 4/real(15), -1/real(3), 1
382 };
383#elif GEOGRAPHICLIB_AUXLATITUDE_ORDER == 6
384 static const real coeff[] = {
385 4/real(1287), 4/real(693), 4/real(315), 4/real(105), 4/real(15),
386 -1/real(3), 1
387 };
388#elif GEOGRAPHICLIB_AUXLATITUDE_ORDER == 8
389 static const real coeff[] = {
390 4/real(3315), 4/real(2145), 4/real(1287), 4/real(693), 4/real(315),
391 4/real(105), 4/real(15), -1/real(3), 1
392 };
393#else
394#error "Unsupported value for GEOGRAPHICLIB_AUXLATITUDE_ORDER"
395#endif
396 int m = Lmax;
397 return _a * (_a + _b) / 2 * Math::polyval(m, coeff, _n);
398 }
399 }
400
401 /// \cond SKIP
402 Math::real AuxLatitude::atanhee(real tphi) const {
403 real s = _f <= 0 ? sn(tphi) : sn(_fm1 * tphi);
404 return _f == 0 ? s :
405 // atanh(e * sphi) = asinh(e' * sbeta)
406 (_f < 0 ? atan( _e * s ) : asinh( _e1 * s )) / _e;
407 }
408 /// \endcond
409
410 Math::real AuxLatitude::q(real tphi) const {
411 real scbeta = sc(_fm1 * tphi);
412 return atanhee(tphi) + (tphi / scbeta) * (sc(tphi) / scbeta);
413 }
414
415 Math::real AuxLatitude::Dq(real tphi) const {
416 real scphi = sc(tphi), sphi = sn(tphi),
417 // d = (1 - sphi) can underflow to zero for large tphi
418 d = tphi > 0 ? 1 / (scphi * scphi * (1 + sphi)) : 1 - sphi;
419 if (tphi <= 0)
420 // This branch is not reached; this case is open-coded in Authalic.
421 return (_q - q(tphi)) / d;
422 else if (d == 0)
423 return 2 / Math::sq(_e2m1);
424 else {
425 // General expression for Dq(1, sphi) is
426 // atanh(e * d / (1 - e2 * sphi)) / (e * d) +
427 // (1 + e2 * sphi) / ((1 - e2 * sphi * sphi) * e2m1);
428 // atanh( e * d / (1 - e2 * sphi))
429 // = atanh( e * d * scphi/(scphi - e2 * tphi))
430 // =
431 real scbeta = sc(_fm1 * tphi);
432 return (_f == 0 ? 1 :
433 (_f > 0 ? asinh(_e1 * d * scphi / scbeta) :
434 atan(_e * d / (1 - _e2 * sphi))) / (_e * d) ) +
435 (_f > 0 ?
436 ((scphi + _e2 * tphi) / (_e2m1 * scbeta)) * (scphi / scbeta) :
437 (1 + _e2 * sphi) / ((1 - _e2 * sphi*sphi) * _e2m1) );
438 }
439 }
440
441 /// \cond SKIP
442 void AuxLatitude::fillcoeff(int auxin, int auxout, int k) const {
443#if GEOGRAPHICLIB_AUXLATITUDE_ORDER == 4
444 static const real coeffs[] = {
445 // C[phi,phi] skipped
446 // C[phi,beta]; even coeffs only
447 0, 1,
448 0, 1/real(2),
449 1/real(3),
450 1/real(4),
451 // C[phi,theta]; even coeffs only
452 -2, 2,
453 -4, 2,
454 8/real(3),
455 4,
456 // C[phi,mu]; even coeffs only
457 -27/real(32), 3/real(2),
458 -55/real(32), 21/real(16),
459 151/real(96),
460 1097/real(512),
461 // C[phi,chi]
462 116/real(45), -2, -2/real(3), 2,
463 -227/real(45), -8/real(5), 7/real(3),
464 -136/real(35), 56/real(15),
465 4279/real(630),
466 // C[phi,xi]
467 -2582/real(14175), -16/real(35), 4/real(45), 4/real(3),
468 -11966/real(14175), 152/real(945), 46/real(45),
469 3802/real(14175), 3044/real(2835),
470 6059/real(4725),
471 // C[beta,phi]; even coeffs only
472 0, -1,
473 0, 1/real(2),
474 -1/real(3),
475 1/real(4),
476 // C[beta,beta] skipped
477 // C[beta,theta]; even coeffs only
478 0, 1,
479 0, 1/real(2),
480 1/real(3),
481 1/real(4),
482 // C[beta,mu]; even coeffs only
483 -9/real(32), 1/real(2),
484 -37/real(96), 5/real(16),
485 29/real(96),
486 539/real(1536),
487 // C[beta,chi]
488 38/real(45), -1/real(3), -2/real(3), 1,
489 -7/real(9), -14/real(15), 5/real(6),
490 -34/real(21), 16/real(15),
491 2069/real(1260),
492 // C[beta,xi]
493 -1082/real(14175), -46/real(315), 4/real(45), 1/real(3),
494 -338/real(2025), 68/real(945), 17/real(90),
495 1102/real(14175), 461/real(2835),
496 3161/real(18900),
497 // C[theta,phi]; even coeffs only
498 2, -2,
499 -4, 2,
500 -8/real(3),
501 4,
502 // C[theta,beta]; even coeffs only
503 0, -1,
504 0, 1/real(2),
505 -1/real(3),
506 1/real(4),
507 // C[theta,theta] skipped
508 // C[theta,mu]; even coeffs only
509 -23/real(32), -1/real(2),
510 -5/real(96), 5/real(16),
511 1/real(32),
512 283/real(1536),
513 // C[theta,chi]
514 4/real(9), -2/real(3), -2/real(3), 0,
515 -23/real(45), -4/real(15), 1/real(3),
516 -24/real(35), 2/real(5),
517 83/real(126),
518 // C[theta,xi]
519 -2102/real(14175), -158/real(315), 4/real(45), -2/real(3),
520 934/real(14175), -16/real(945), 16/real(45),
521 922/real(14175), -232/real(2835),
522 719/real(4725),
523 // C[mu,phi]; even coeffs only
524 9/real(16), -3/real(2),
525 -15/real(32), 15/real(16),
526 -35/real(48),
527 315/real(512),
528 // C[mu,beta]; even coeffs only
529 3/real(16), -1/real(2),
530 1/real(32), -1/real(16),
531 -1/real(48),
532 -5/real(512),
533 // C[mu,theta]; even coeffs only
534 13/real(16), 1/real(2),
535 33/real(32), -1/real(16),
536 -5/real(16),
537 -261/real(512),
538 // C[mu,mu] skipped
539 // C[mu,chi]
540 41/real(180), 5/real(16), -2/real(3), 1/real(2),
541 557/real(1440), -3/real(5), 13/real(48),
542 -103/real(140), 61/real(240),
543 49561/real(161280),
544 // C[mu,xi]
545 -1609/real(28350), 121/real(1680), 4/real(45), -1/real(6),
546 16463/real(453600), 26/real(945), -29/real(720),
547 449/real(28350), -1003/real(45360),
548 -40457/real(2419200),
549 // C[chi,phi]
550 -82/real(45), 4/real(3), 2/real(3), -2,
551 -13/real(9), -16/real(15), 5/real(3),
552 34/real(21), -26/real(15),
553 1237/real(630),
554 // C[chi,beta]
555 -16/real(45), 0, 2/real(3), -1,
556 19/real(45), -2/real(5), 1/real(6),
557 16/real(105), -1/real(15),
558 17/real(1260),
559 // C[chi,theta]
560 -2/real(9), 2/real(3), 2/real(3), 0,
561 43/real(45), 4/real(15), -1/real(3),
562 2/real(105), -2/real(5),
563 -55/real(126),
564 // C[chi,mu]
565 1/real(360), -37/real(96), 2/real(3), -1/real(2),
566 437/real(1440), -1/real(15), -1/real(48),
567 37/real(840), -17/real(480),
568 -4397/real(161280),
569 // C[chi,chi] skipped
570 // C[chi,xi]
571 -2312/real(14175), -88/real(315), 34/real(45), -2/real(3),
572 6079/real(14175), -184/real(945), 1/real(45),
573 772/real(14175), -106/real(2835),
574 -167/real(9450),
575 // C[xi,phi]
576 538/real(4725), 88/real(315), -4/real(45), -4/real(3),
577 -2482/real(14175), 8/real(105), 34/real(45),
578 -898/real(14175), -1532/real(2835),
579 6007/real(14175),
580 // C[xi,beta]
581 34/real(675), 32/real(315), -4/real(45), -1/real(3),
582 74/real(2025), -4/real(315), -7/real(90),
583 2/real(14175), -83/real(2835),
584 -797/real(56700),
585 // C[xi,theta]
586 778/real(4725), 62/real(105), -4/real(45), 2/real(3),
587 12338/real(14175), -32/real(315), 4/real(45),
588 -1618/real(14175), -524/real(2835),
589 -5933/real(14175),
590 // C[xi,mu]
591 1297/real(18900), -817/real(10080), -4/real(45), 1/real(6),
592 -29609/real(453600), -2/real(35), 49/real(720),
593 -2917/real(56700), 4463/real(90720),
594 331799/real(7257600),
595 // C[xi,chi]
596 2458/real(4725), 46/real(315), -34/real(45), 2/real(3),
597 3413/real(14175), -256/real(315), 19/real(45),
598 -15958/real(14175), 248/real(567),
599 16049/real(28350),
600 // C[xi,xi] skipped
601 };
602 static const int ptrs[] = {
603 0, 0, 6, 12, 18, 28, 38, 44, 44, 50, 56, 66, 76, 82, 88, 88, 94, 104,
604 114, 120, 126, 132, 132, 142, 152, 162, 172, 182, 192, 192, 202, 212,
605 222, 232, 242, 252, 252,
606 };
607#elif GEOGRAPHICLIB_AUXLATITUDE_ORDER == 6
608 static const real coeffs[] = {
609 // C[phi,phi] skipped
610 // C[phi,beta]; even coeffs only
611 0, 0, 1,
612 0, 0, 1/real(2),
613 0, 1/real(3),
614 0, 1/real(4),
615 1/real(5),
616 1/real(6),
617 // C[phi,theta]; even coeffs only
618 2, -2, 2,
619 6, -4, 2,
620 -8, 8/real(3),
621 -16, 4,
622 32/real(5),
623 32/real(3),
624 // C[phi,mu]; even coeffs only
625 269/real(512), -27/real(32), 3/real(2),
626 6759/real(4096), -55/real(32), 21/real(16),
627 -417/real(128), 151/real(96),
628 -15543/real(2560), 1097/real(512),
629 8011/real(2560),
630 293393/real(61440),
631 // C[phi,chi]
632 -2854/real(675), 26/real(45), 116/real(45), -2, -2/real(3), 2,
633 2323/real(945), 2704/real(315), -227/real(45), -8/real(5), 7/real(3),
634 73814/real(2835), -1262/real(105), -136/real(35), 56/real(15),
635 -399572/real(14175), -332/real(35), 4279/real(630),
636 -144838/real(6237), 4174/real(315),
637 601676/real(22275),
638 // C[phi,xi]
639 28112932/real(212837625), 60136/real(467775), -2582/real(14175),
640 -16/real(35), 4/real(45), 4/real(3),
641 251310128/real(638512875), -21016/real(51975), -11966/real(14175),
642 152/real(945), 46/real(45),
643 -8797648/real(10945935), -94388/real(66825), 3802/real(14175),
644 3044/real(2835),
645 -1472637812/real(638512875), 41072/real(93555), 6059/real(4725),
646 455935736/real(638512875), 768272/real(467775),
647 4210684958LL/real(1915538625),
648 // C[beta,phi]; even coeffs only
649 0, 0, -1,
650 0, 0, 1/real(2),
651 0, -1/real(3),
652 0, 1/real(4),
653 -1/real(5),
654 1/real(6),
655 // C[beta,beta] skipped
656 // C[beta,theta]; even coeffs only
657 0, 0, 1,
658 0, 0, 1/real(2),
659 0, 1/real(3),
660 0, 1/real(4),
661 1/real(5),
662 1/real(6),
663 // C[beta,mu]; even coeffs only
664 205/real(1536), -9/real(32), 1/real(2),
665 1335/real(4096), -37/real(96), 5/real(16),
666 -75/real(128), 29/real(96),
667 -2391/real(2560), 539/real(1536),
668 3467/real(7680),
669 38081/real(61440),
670 // C[beta,chi]
671 -3118/real(4725), -1/real(3), 38/real(45), -1/real(3), -2/real(3), 1,
672 -247/real(270), 50/real(21), -7/real(9), -14/real(15), 5/real(6),
673 17564/real(2835), -5/real(3), -34/real(21), 16/real(15),
674 -49877/real(14175), -28/real(9), 2069/real(1260),
675 -28244/real(4455), 883/real(315),
676 797222/real(155925),
677 // C[beta,xi]
678 7947332/real(212837625), 11824/real(467775), -1082/real(14175),
679 -46/real(315), 4/real(45), 1/real(3),
680 39946703/real(638512875), -16672/real(155925), -338/real(2025),
681 68/real(945), 17/real(90),
682 -255454/real(1563705), -101069/real(467775), 1102/real(14175),
683 461/real(2835),
684 -189032762/real(638512875), 1786/real(18711), 3161/real(18900),
685 80274086/real(638512875), 88868/real(467775),
686 880980241/real(3831077250LL),
687 // C[theta,phi]; even coeffs only
688 -2, 2, -2,
689 6, -4, 2,
690 8, -8/real(3),
691 -16, 4,
692 -32/real(5),
693 32/real(3),
694 // C[theta,beta]; even coeffs only
695 0, 0, -1,
696 0, 0, 1/real(2),
697 0, -1/real(3),
698 0, 1/real(4),
699 -1/real(5),
700 1/real(6),
701 // C[theta,theta] skipped
702 // C[theta,mu]; even coeffs only
703 499/real(1536), -23/real(32), -1/real(2),
704 6565/real(12288), -5/real(96), 5/real(16),
705 -77/real(128), 1/real(32),
706 -4037/real(7680), 283/real(1536),
707 1301/real(7680),
708 17089/real(61440),
709 // C[theta,chi]
710 -3658/real(4725), 2/real(9), 4/real(9), -2/real(3), -2/real(3), 0,
711 61/real(135), 68/real(45), -23/real(45), -4/real(15), 1/real(3),
712 9446/real(2835), -46/real(35), -24/real(35), 2/real(5),
713 -34712/real(14175), -80/real(63), 83/real(126),
714 -2362/real(891), 52/real(45),
715 335882/real(155925),
716 // C[theta,xi]
717 216932/real(2627625), 109042/real(467775), -2102/real(14175),
718 -158/real(315), 4/real(45), -2/real(3),
719 117952358/real(638512875), -7256/real(155925), 934/real(14175),
720 -16/real(945), 16/real(45),
721 -7391576/real(54729675), -25286/real(66825), 922/real(14175),
722 -232/real(2835),
723 -67048172/real(638512875), 268/real(18711), 719/real(4725),
724 46774256/real(638512875), 14354/real(467775),
725 253129538/real(1915538625),
726 // C[mu,phi]; even coeffs only
727 -3/real(32), 9/real(16), -3/real(2),
728 135/real(2048), -15/real(32), 15/real(16),
729 105/real(256), -35/real(48),
730 -189/real(512), 315/real(512),
731 -693/real(1280),
732 1001/real(2048),
733 // C[mu,beta]; even coeffs only
734 -1/real(32), 3/real(16), -1/real(2),
735 -9/real(2048), 1/real(32), -1/real(16),
736 3/real(256), -1/real(48),
737 3/real(512), -5/real(512),
738 -7/real(1280),
739 -7/real(2048),
740 // C[mu,theta]; even coeffs only
741 -15/real(32), 13/real(16), 1/real(2),
742 -1673/real(2048), 33/real(32), -1/real(16),
743 349/real(256), -5/real(16),
744 963/real(512), -261/real(512),
745 -921/real(1280),
746 -6037/real(6144),
747 // C[mu,mu] skipped
748 // C[mu,chi]
749 7891/real(37800), -127/real(288), 41/real(180), 5/real(16), -2/real(3),
750 1/real(2),
751 -1983433/real(1935360), 281/real(630), 557/real(1440), -3/real(5),
752 13/real(48),
753 167603/real(181440), 15061/real(26880), -103/real(140), 61/real(240),
754 6601661/real(7257600), -179/real(168), 49561/real(161280),
755 -3418889/real(1995840), 34729/real(80640),
756 212378941/real(319334400),
757 // C[mu,xi]
758 12674323/real(851350500), -384229/real(14968800), -1609/real(28350),
759 121/real(1680), 4/real(45), -1/real(6),
760 -31621753811LL/real(1307674368000LL), -431/real(17325),
761 16463/real(453600), 26/real(945), -29/real(720),
762 -32844781/real(1751349600), 3746047/real(119750400), 449/real(28350),
763 -1003/real(45360),
764 10650637121LL/real(326918592000LL), 629/real(53460),
765 -40457/real(2419200),
766 205072597/real(20432412000LL), -1800439/real(119750400),
767 -59109051671LL/real(3923023104000LL),
768 // C[chi,phi]
769 4642/real(4725), 32/real(45), -82/real(45), 4/real(3), 2/real(3), -2,
770 -1522/real(945), 904/real(315), -13/real(9), -16/real(15), 5/real(3),
771 -12686/real(2835), 8/real(5), 34/real(21), -26/real(15),
772 -24832/real(14175), -12/real(5), 1237/real(630),
773 109598/real(31185), -734/real(315),
774 444337/real(155925),
775 // C[chi,beta]
776 -998/real(4725), 2/real(5), -16/real(45), 0, 2/real(3), -1,
777 -2/real(27), -22/real(105), 19/real(45), -2/real(5), 1/real(6),
778 116/real(567), -22/real(105), 16/real(105), -1/real(15),
779 2123/real(14175), -8/real(105), 17/real(1260),
780 128/real(4455), -1/real(105),
781 149/real(311850),
782 // C[chi,theta]
783 1042/real(4725), -14/real(45), -2/real(9), 2/real(3), 2/real(3), 0,
784 -712/real(945), -4/real(45), 43/real(45), 4/real(15), -1/real(3),
785 274/real(2835), 124/real(105), 2/real(105), -2/real(5),
786 21068/real(14175), -16/real(105), -55/real(126),
787 -9202/real(31185), -22/real(45),
788 -90263/real(155925),
789 // C[chi,mu]
790 -96199/real(604800), 81/real(512), 1/real(360), -37/real(96), 2/real(3),
791 -1/real(2),
792 1118711/real(3870720), -46/real(105), 437/real(1440), -1/real(15),
793 -1/real(48),
794 -5569/real(90720), 209/real(4480), 37/real(840), -17/real(480),
795 830251/real(7257600), 11/real(504), -4397/real(161280),
796 108847/real(3991680), -4583/real(161280),
797 -20648693/real(638668800),
798 // C[chi,chi] skipped
799 // C[chi,xi]
800 -55271278/real(212837625), 27128/real(93555), -2312/real(14175),
801 -88/real(315), 34/real(45), -2/real(3),
802 106691108/real(638512875), -65864/real(155925), 6079/real(14175),
803 -184/real(945), 1/real(45),
804 5921152/real(54729675), -14246/real(467775), 772/real(14175),
805 -106/real(2835),
806 75594328/real(638512875), -5312/real(467775), -167/real(9450),
807 2837636/real(638512875), -248/real(13365),
808 -34761247/real(1915538625),
809 // C[xi,phi]
810 -44732/real(2837835), 20824/real(467775), 538/real(4725), 88/real(315),
811 -4/real(45), -4/real(3),
812 -12467764/real(212837625), -37192/real(467775), -2482/real(14175),
813 8/real(105), 34/real(45),
814 100320856/real(1915538625), 54968/real(467775), -898/real(14175),
815 -1532/real(2835),
816 -5884124/real(70945875), 24496/real(467775), 6007/real(14175),
817 -839792/real(19348875), -23356/real(66825),
818 570284222/real(1915538625),
819 // C[xi,beta]
820 -70496/real(8513505), 2476/real(467775), 34/real(675), 32/real(315),
821 -4/real(45), -1/real(3),
822 53836/real(212837625), 3992/real(467775), 74/real(2025), -4/real(315),
823 -7/real(90),
824 -661844/real(1915538625), 7052/real(467775), 2/real(14175),
825 -83/real(2835),
826 1425778/real(212837625), 934/real(467775), -797/real(56700),
827 390088/real(212837625), -3673/real(467775),
828 -18623681/real(3831077250LL),
829 // C[xi,theta]
830 -4286228/real(42567525), -193082/real(467775), 778/real(4725),
831 62/real(105), -4/real(45), 2/real(3),
832 -61623938/real(70945875), 92696/real(467775), 12338/real(14175),
833 -32/real(315), 4/real(45),
834 427003576/real(1915538625), 612536/real(467775), -1618/real(14175),
835 -524/real(2835),
836 427770788/real(212837625), -8324/real(66825), -5933/real(14175),
837 -9153184/real(70945875), -320044/real(467775),
838 -1978771378/real(1915538625),
839 // C[xi,mu]
840 -9292991/real(302702400), 7764059/real(239500800), 1297/real(18900),
841 -817/real(10080), -4/real(45), 1/real(6),
842 36019108271LL/real(871782912000LL), 35474/real(467775),
843 -29609/real(453600), -2/real(35), 49/real(720),
844 3026004511LL/real(30648618000LL), -4306823/real(59875200),
845 -2917/real(56700), 4463/real(90720),
846 -368661577/real(4036032000LL), -102293/real(1871100),
847 331799/real(7257600),
848 -875457073/real(13621608000LL), 11744233/real(239500800),
849 453002260127LL/real(7846046208000LL),
850 // C[xi,chi]
851 2706758/real(42567525), -55222/real(93555), 2458/real(4725),
852 46/real(315), -34/real(45), 2/real(3),
853 -340492279/real(212837625), 516944/real(467775), 3413/real(14175),
854 -256/real(315), 19/real(45),
855 4430783356LL/real(1915538625), 206834/real(467775), -15958/real(14175),
856 248/real(567),
857 62016436/real(70945875), -832976/real(467775), 16049/real(28350),
858 -651151712/real(212837625), 15602/real(18711),
859 2561772812LL/real(1915538625),
860 // C[xi,xi] skipped
861 };
862 static const int ptrs[] = {
863 0, 0, 12, 24, 36, 57, 78, 90, 90, 102, 114, 135, 156, 168, 180, 180, 192,
864 213, 234, 246, 258, 270, 270, 291, 312, 333, 354, 375, 396, 396, 417,
865 438, 459, 480, 501, 522, 522,
866 };
867#elif GEOGRAPHICLIB_AUXLATITUDE_ORDER == 8
868 static const real coeffs[] = {
869 // C[phi,phi] skipped
870 // C[phi,beta]; even coeffs only
871 0, 0, 0, 1,
872 0, 0, 0, 1/real(2),
873 0, 0, 1/real(3),
874 0, 0, 1/real(4),
875 0, 1/real(5),
876 0, 1/real(6),
877 1/real(7),
878 1/real(8),
879 // C[phi,theta]; even coeffs only
880 -2, 2, -2, 2,
881 -8, 6, -4, 2,
882 16, -8, 8/real(3),
883 40, -16, 4,
884 -32, 32/real(5),
885 -64, 32/real(3),
886 128/real(7),
887 32,
888 // C[phi,mu]; even coeffs only
889 -6607/real(24576), 269/real(512), -27/real(32), 3/real(2),
890 -155113/real(122880), 6759/real(4096), -55/real(32), 21/real(16),
891 87963/real(20480), -417/real(128), 151/real(96),
892 2514467/real(245760), -15543/real(2560), 1097/real(512),
893 -69119/real(6144), 8011/real(2560),
894 -5962461/real(286720), 293393/real(61440),
895 6459601/real(860160),
896 332287993/real(27525120),
897 // C[phi,chi]
898 189416/real(99225), 16822/real(4725), -2854/real(675), 26/real(45),
899 116/real(45), -2, -2/real(3), 2,
900 141514/real(8505), -31256/real(1575), 2323/real(945), 2704/real(315),
901 -227/real(45), -8/real(5), 7/real(3),
902 -2363828/real(31185), 98738/real(14175), 73814/real(2835),
903 -1262/real(105), -136/real(35), 56/real(15),
904 14416399/real(935550), 11763988/real(155925), -399572/real(14175),
905 -332/real(35), 4279/real(630),
906 258316372/real(1216215), -2046082/real(31185), -144838/real(6237),
907 4174/real(315),
908 -2155215124LL/real(14189175), -115444544/real(2027025),
909 601676/real(22275),
910 -170079376/real(1216215), 38341552/real(675675),
911 1383243703/real(11351340),
912 // C[phi,xi]
913 -1683291094/real(37574026875LL), 22947844/real(1915538625),
914 28112932/real(212837625), 60136/real(467775), -2582/real(14175),
915 -16/real(35), 4/real(45), 4/real(3),
916 -14351220203LL/real(488462349375LL), 1228352/real(3007125),
917 251310128/real(638512875), -21016/real(51975), -11966/real(14175),
918 152/real(945), 46/real(45),
919 505559334506LL/real(488462349375LL), 138128272/real(147349125),
920 -8797648/real(10945935), -94388/real(66825), 3802/real(14175),
921 3044/real(2835),
922 973080708361LL/real(488462349375LL), -45079184/real(29469825),
923 -1472637812/real(638512875), 41072/real(93555), 6059/real(4725),
924 -1385645336626LL/real(488462349375LL), -550000184/real(147349125),
925 455935736/real(638512875), 768272/real(467775),
926 -2939205114427LL/real(488462349375LL), 443810768/real(383107725),
927 4210684958LL/real(1915538625),
928 101885255158LL/real(54273594375LL), 387227992/real(127702575),
929 1392441148867LL/real(325641566250LL),
930 // C[beta,phi]; even coeffs only
931 0, 0, 0, -1,
932 0, 0, 0, 1/real(2),
933 0, 0, -1/real(3),
934 0, 0, 1/real(4),
935 0, -1/real(5),
936 0, 1/real(6),
937 -1/real(7),
938 1/real(8),
939 // C[beta,beta] skipped
940 // C[beta,theta]; even coeffs only
941 0, 0, 0, 1,
942 0, 0, 0, 1/real(2),
943 0, 0, 1/real(3),
944 0, 0, 1/real(4),
945 0, 1/real(5),
946 0, 1/real(6),
947 1/real(7),
948 1/real(8),
949 // C[beta,mu]; even coeffs only
950 -4879/real(73728), 205/real(1536), -9/real(32), 1/real(2),
951 -86171/real(368640), 1335/real(4096), -37/real(96), 5/real(16),
952 2901/real(4096), -75/real(128), 29/real(96),
953 1082857/real(737280), -2391/real(2560), 539/real(1536),
954 -28223/real(18432), 3467/real(7680),
955 -733437/real(286720), 38081/real(61440),
956 459485/real(516096),
957 109167851/real(82575360),
958 // C[beta,chi]
959 -25666/real(99225), 4769/real(4725), -3118/real(4725), -1/real(3),
960 38/real(45), -1/real(3), -2/real(3), 1,
961 193931/real(42525), -14404/real(4725), -247/real(270), 50/real(21),
962 -7/real(9), -14/real(15), 5/real(6),
963 -1709614/real(155925), -36521/real(14175), 17564/real(2835), -5/real(3),
964 -34/real(21), 16/real(15),
965 -637699/real(85050), 2454416/real(155925), -49877/real(14175),
966 -28/real(9), 2069/real(1260),
967 48124558/real(1216215), -20989/real(2835), -28244/real(4455),
968 883/real(315),
969 -16969807/real(1091475), -2471888/real(184275), 797222/real(155925),
970 -1238578/real(42525), 2199332/real(225225),
971 87600385/real(4540536),
972 // C[beta,xi]
973 -5946082372LL/real(488462349375LL), 9708931/real(1915538625),
974 7947332/real(212837625), 11824/real(467775), -1082/real(14175),
975 -46/real(315), 4/real(45), 1/real(3),
976 190673521/real(69780335625LL), 164328266/real(1915538625),
977 39946703/real(638512875), -16672/real(155925), -338/real(2025),
978 68/real(945), 17/real(90),
979 86402898356LL/real(488462349375LL), 236067184/real(1915538625),
980 -255454/real(1563705), -101069/real(467775), 1102/real(14175),
981 461/real(2835),
982 110123070361LL/real(488462349375LL), -98401826/real(383107725),
983 -189032762/real(638512875), 1786/real(18711), 3161/real(18900),
984 -200020620676LL/real(488462349375LL), -802887278/real(1915538625),
985 80274086/real(638512875), 88868/real(467775),
986 -296107325077LL/real(488462349375LL), 66263486/real(383107725),
987 880980241/real(3831077250LL),
988 4433064236LL/real(18091198125LL), 37151038/real(127702575),
989 495248998393LL/real(1302566265000LL),
990 // C[theta,phi]; even coeffs only
991 2, -2, 2, -2,
992 -8, 6, -4, 2,
993 -16, 8, -8/real(3),
994 40, -16, 4,
995 32, -32/real(5),
996 -64, 32/real(3),
997 -128/real(7),
998 32,
999 // C[theta,beta]; even coeffs only
1000 0, 0, 0, -1,
1001 0, 0, 0, 1/real(2),
1002 0, 0, -1/real(3),
1003 0, 0, 1/real(4),
1004 0, -1/real(5),
1005 0, 1/real(6),
1006 -1/real(7),
1007 1/real(8),
1008 // C[theta,theta] skipped
1009 // C[theta,mu]; even coeffs only
1010 -14321/real(73728), 499/real(1536), -23/real(32), -1/real(2),
1011 -201467/real(368640), 6565/real(12288), -5/real(96), 5/real(16),
1012 2939/real(4096), -77/real(128), 1/real(32),
1013 1155049/real(737280), -4037/real(7680), 283/real(1536),
1014 -19465/real(18432), 1301/real(7680),
1015 -442269/real(286720), 17089/real(61440),
1016 198115/real(516096),
1017 48689387/real(82575360),
1018 // C[theta,chi]
1019 64424/real(99225), 76/real(225), -3658/real(4725), 2/real(9), 4/real(9),
1020 -2/real(3), -2/real(3), 0,
1021 2146/real(1215), -2728/real(945), 61/real(135), 68/real(45),
1022 -23/real(45), -4/real(15), 1/real(3),
1023 -95948/real(10395), 428/real(945), 9446/real(2835), -46/real(35),
1024 -24/real(35), 2/real(5),
1025 29741/real(85050), 4472/real(525), -34712/real(14175), -80/real(63),
1026 83/real(126),
1027 280108/real(13365), -17432/real(3465), -2362/real(891), 52/real(45),
1028 -48965632/real(4729725), -548752/real(96525), 335882/real(155925),
1029 -197456/real(15795), 51368/real(12285),
1030 1461335/real(174636),
1031 // C[theta,xi]
1032 -230886326/real(6343666875LL), -189115382/real(1915538625),
1033 216932/real(2627625), 109042/real(467775), -2102/real(14175),
1034 -158/real(315), 4/real(45), -2/real(3),
1035 -11696145869LL/real(69780335625LL), 288456008/real(1915538625),
1036 117952358/real(638512875), -7256/real(155925), 934/real(14175),
1037 -16/real(945), 16/real(45),
1038 91546732346LL/real(488462349375LL), 478700902/real(1915538625),
1039 -7391576/real(54729675), -25286/real(66825), 922/real(14175),
1040 -232/real(2835),
1041 218929662961LL/real(488462349375LL), -67330724/real(383107725),
1042 -67048172/real(638512875), 268/real(18711), 719/real(4725),
1043 -129039188386LL/real(488462349375LL), -117954842/real(273648375),
1044 46774256/real(638512875), 14354/real(467775),
1045 -178084928947LL/real(488462349375LL), 2114368/real(34827975),
1046 253129538/real(1915538625),
1047 6489189398LL/real(54273594375LL), 13805944/real(127702575),
1048 59983985827LL/real(325641566250LL),
1049 // C[mu,phi]; even coeffs only
1050 57/real(2048), -3/real(32), 9/real(16), -3/real(2),
1051 -105/real(4096), 135/real(2048), -15/real(32), 15/real(16),
1052 -105/real(2048), 105/real(256), -35/real(48),
1053 693/real(16384), -189/real(512), 315/real(512),
1054 693/real(2048), -693/real(1280),
1055 -1287/real(4096), 1001/real(2048),
1056 -6435/real(14336),
1057 109395/real(262144),
1058 // C[mu,beta]; even coeffs only
1059 19/real(2048), -1/real(32), 3/real(16), -1/real(2),
1060 7/real(4096), -9/real(2048), 1/real(32), -1/real(16),
1061 -3/real(2048), 3/real(256), -1/real(48),
1062 -11/real(16384), 3/real(512), -5/real(512),
1063 7/real(2048), -7/real(1280),
1064 9/real(4096), -7/real(2048),
1065 -33/real(14336),
1066 -429/real(262144),
1067 // C[mu,theta]; even coeffs only
1068 509/real(2048), -15/real(32), 13/real(16), 1/real(2),
1069 2599/real(4096), -1673/real(2048), 33/real(32), -1/real(16),
1070 -2989/real(2048), 349/real(256), -5/real(16),
1071 -43531/real(16384), 963/real(512), -261/real(512),
1072 5545/real(2048), -921/real(1280),
1073 16617/real(4096), -6037/real(6144),
1074 -19279/real(14336),
1075 -490925/real(262144),
1076 // C[mu,mu] skipped
1077 // C[mu,chi]
1078 -18975107/real(50803200), 72161/real(387072), 7891/real(37800),
1079 -127/real(288), 41/real(180), 5/real(16), -2/real(3), 1/real(2),
1080 148003883/real(174182400), 13769/real(28800), -1983433/real(1935360),
1081 281/real(630), 557/real(1440), -3/real(5), 13/real(48),
1082 79682431/real(79833600), -67102379/real(29030400), 167603/real(181440),
1083 15061/real(26880), -103/real(140), 61/real(240),
1084 -40176129013LL/real(7664025600LL), 97445/real(49896),
1085 6601661/real(7257600), -179/real(168), 49561/real(161280),
1086 2605413599LL/real(622702080), 14644087/real(9123840),
1087 -3418889/real(1995840), 34729/real(80640),
1088 175214326799LL/real(58118860800LL), -30705481/real(10378368),
1089 212378941/real(319334400),
1090 -16759934899LL/real(3113510400LL), 1522256789/real(1383782400),
1091 1424729850961LL/real(743921418240LL),
1092 // C[mu,xi]
1093 -375027460897LL/real(125046361440000LL),
1094 7183403063LL/real(560431872000LL), 12674323/real(851350500),
1095 -384229/real(14968800), -1609/real(28350), 121/real(1680), 4/real(45),
1096 -1/real(6),
1097 30410873385097LL/real(2000741783040000LL),
1098 1117820213/real(122594472000LL), -31621753811LL/real(1307674368000LL),
1099 -431/real(17325), 16463/real(453600), 26/real(945), -29/real(720),
1100 151567502183LL/real(17863765920000LL),
1101 -116359346641LL/real(3923023104000LL), -32844781/real(1751349600),
1102 3746047/real(119750400), 449/real(28350), -1003/real(45360),
1103 -317251099510901LL/real(8002967132160000LL), -13060303/real(766215450),
1104 10650637121LL/real(326918592000LL), 629/real(53460),
1105 -40457/real(2419200),
1106 -2105440822861LL/real(125046361440000LL),
1107 146875240637LL/real(3923023104000LL), 205072597/real(20432412000LL),
1108 -1800439/real(119750400),
1109 91496147778023LL/real(2000741783040000LL), 228253559/real(24518894400LL),
1110 -59109051671LL/real(3923023104000LL),
1111 126430355893LL/real(13894040160000LL),
1112 -4255034947LL/real(261534873600LL),
1113 -791820407649841LL/real(42682491371520000LL),
1114 // C[chi,phi]
1115 1514/real(1323), -8384/real(4725), 4642/real(4725), 32/real(45),
1116 -82/real(45), 4/real(3), 2/real(3), -2,
1117 142607/real(42525), -2288/real(1575), -1522/real(945), 904/real(315),
1118 -13/real(9), -16/real(15), 5/real(3),
1119 120202/real(51975), 44644/real(14175), -12686/real(2835), 8/real(5),
1120 34/real(21), -26/real(15),
1121 -1097407/real(187110), 1077964/real(155925), -24832/real(14175),
1122 -12/real(5), 1237/real(630),
1123 -12870194/real(1216215), 1040/real(567), 109598/real(31185),
1124 -734/real(315),
1125 -126463/real(72765), -941912/real(184275), 444337/real(155925),
1126 3463678/real(467775), -2405834/real(675675),
1127 256663081/real(56756700),
1128 // C[chi,beta]
1129 1384/real(11025), -34/real(4725), -998/real(4725), 2/real(5),
1130 -16/real(45), 0, 2/real(3), -1,
1131 -12616/real(42525), 1268/real(4725), -2/real(27), -22/real(105),
1132 19/real(45), -2/real(5), 1/real(6),
1133 1724/real(51975), -1858/real(14175), 116/real(567), -22/real(105),
1134 16/real(105), -1/real(15),
1135 115249/real(935550), -26836/real(155925), 2123/real(14175), -8/real(105),
1136 17/real(1260),
1137 140836/real(1216215), -424/real(6237), 128/real(4455), -1/real(105),
1138 210152/real(4729725), -31232/real(2027025), 149/real(311850),
1139 30208/real(6081075), -499/real(225225),
1140 -68251/real(113513400),
1141 // C[chi,theta]
1142 -1738/real(11025), 18/real(175), 1042/real(4725), -14/real(45),
1143 -2/real(9), 2/real(3), 2/real(3), 0,
1144 23159/real(42525), 332/real(945), -712/real(945), -4/real(45),
1145 43/real(45), 4/real(15), -1/real(3),
1146 13102/real(31185), -1352/real(945), 274/real(2835), 124/real(105),
1147 2/real(105), -2/real(5),
1148 -2414843/real(935550), 1528/real(4725), 21068/real(14175), -16/real(105),
1149 -55/real(126),
1150 60334/real(93555), 20704/real(10395), -9202/real(31185), -22/real(45),
1151 40458083/real(14189175), -299444/real(675675), -90263/real(155925),
1152 -3818498/real(6081075), -8962/real(12285),
1153 -4259027/real(4365900),
1154 // C[chi,mu]
1155 -7944359/real(67737600), 5406467/real(38707200), -96199/real(604800),
1156 81/real(512), 1/real(360), -37/real(96), 2/real(3), -1/real(2),
1157 -24749483/real(348364800), -51841/real(1209600), 1118711/real(3870720),
1158 -46/real(105), 437/real(1440), -1/real(15), -1/real(48),
1159 6457463/real(17740800), -9261899/real(58060800), -5569/real(90720),
1160 209/real(4480), 37/real(840), -17/real(480),
1161 -324154477/real(7664025600LL), -466511/real(2494800),
1162 830251/real(7257600), 11/real(504), -4397/real(161280),
1163 -22894433/real(124540416), 8005831/real(63866880), 108847/real(3991680),
1164 -4583/real(161280),
1165 2204645983LL/real(12915302400LL), 16363163/real(518918400),
1166 -20648693/real(638668800),
1167 497323811/real(12454041600LL), -219941297/real(5535129600LL),
1168 -191773887257LL/real(3719607091200LL),
1169 // C[chi,chi] skipped
1170 // C[chi,xi]
1171 -17451293242LL/real(488462349375LL), 308365186/real(1915538625),
1172 -55271278/real(212837625), 27128/real(93555), -2312/real(14175),
1173 -88/real(315), 34/real(45), -2/real(3),
1174 -101520127208LL/real(488462349375LL), 149984636/real(1915538625),
1175 106691108/real(638512875), -65864/real(155925), 6079/real(14175),
1176 -184/real(945), 1/real(45),
1177 10010741462LL/real(37574026875LL), -99534832/real(383107725),
1178 5921152/real(54729675), -14246/real(467775), 772/real(14175),
1179 -106/real(2835),
1180 1615002539/real(75148053750LL), -35573728/real(273648375),
1181 75594328/real(638512875), -5312/real(467775), -167/real(9450),
1182 -3358119706LL/real(488462349375LL), 130601488/real(1915538625),
1183 2837636/real(638512875), -248/real(13365),
1184 46771947158LL/real(488462349375LL), -3196/real(3553875),
1185 -34761247/real(1915538625),
1186 -18696014/real(18091198125LL), -2530364/real(127702575),
1187 -14744861191LL/real(651283132500LL),
1188 // C[xi,phi]
1189 -88002076/real(13956067125LL), -86728/real(16372125),
1190 -44732/real(2837835), 20824/real(467775), 538/real(4725), 88/real(315),
1191 -4/real(45), -4/real(3),
1192 -2641983469LL/real(488462349375LL), -895712/real(147349125),
1193 -12467764/real(212837625), -37192/real(467775), -2482/real(14175),
1194 8/real(105), 34/real(45),
1195 8457703444LL/real(488462349375LL), 240616/real(4209975),
1196 100320856/real(1915538625), 54968/real(467775), -898/real(14175),
1197 -1532/real(2835),
1198 -4910552477LL/real(97692469875LL), -4832848/real(147349125),
1199 -5884124/real(70945875), 24496/real(467775), 6007/real(14175),
1200 9393713176LL/real(488462349375LL), 816824/real(13395375),
1201 -839792/real(19348875), -23356/real(66825),
1202 -4532926649LL/real(97692469875LL), 1980656/real(54729675),
1203 570284222/real(1915538625),
1204 -14848113968LL/real(488462349375LL), -496894276/real(1915538625),
1205 224557742191LL/real(976924698750LL),
1206 // C[xi,beta]
1207 29232878/real(97692469875LL), -18484/real(4343625), -70496/real(8513505),
1208 2476/real(467775), 34/real(675), 32/real(315), -4/real(45), -1/real(3),
1209 -324943819/real(488462349375LL), -4160804/real(1915538625),
1210 53836/real(212837625), 3992/real(467775), 74/real(2025), -4/real(315),
1211 -7/real(90),
1212 -168643106/real(488462349375LL), 237052/real(383107725),
1213 -661844/real(1915538625), 7052/real(467775), 2/real(14175),
1214 -83/real(2835),
1215 113042383/real(97692469875LL), -2915326/real(1915538625),
1216 1425778/real(212837625), 934/real(467775), -797/real(56700),
1217 -558526274/real(488462349375LL), 6064888/real(1915538625),
1218 390088/real(212837625), -3673/real(467775),
1219 155665021/real(97692469875LL), 41288/real(29469825),
1220 -18623681/real(3831077250LL),
1221 504234982/real(488462349375LL), -6205669/real(1915538625),
1222 -8913001661LL/real(3907698795000LL),
1223 // C[xi,theta]
1224 182466964/real(8881133625LL), 53702182/real(212837625),
1225 -4286228/real(42567525), -193082/real(467775), 778/real(4725),
1226 62/real(105), -4/real(45), 2/real(3),
1227 367082779691LL/real(488462349375LL), -32500616/real(273648375),
1228 -61623938/real(70945875), 92696/real(467775), 12338/real(14175),
1229 -32/real(315), 4/real(45),
1230 -42668482796LL/real(488462349375LL), -663111728/real(383107725),
1231 427003576/real(1915538625), 612536/real(467775), -1618/real(14175),
1232 -524/real(2835),
1233 -327791986997LL/real(97692469875LL), 421877252/real(1915538625),
1234 427770788/real(212837625), -8324/real(66825), -5933/real(14175),
1235 74612072536LL/real(488462349375LL), 6024982024LL/real(1915538625),
1236 -9153184/real(70945875), -320044/real(467775),
1237 489898512247LL/real(97692469875LL), -46140784/real(383107725),
1238 -1978771378/real(1915538625),
1239 -42056042768LL/real(488462349375LL), -2926201612LL/real(1915538625),
1240 -2209250801969LL/real(976924698750LL),
1241 // C[xi,mu]
1242 39534358147LL/real(2858202547200LL),
1243 -25359310709LL/real(1743565824000LL), -9292991/real(302702400),
1244 7764059/real(239500800), 1297/real(18900), -817/real(10080), -4/real(45),
1245 1/real(6),
1246 -13216941177599LL/real(571640509440000LL),
1247 -14814966289LL/real(245188944000LL), 36019108271LL/real(871782912000LL),
1248 35474/real(467775), -29609/real(453600), -2/real(35), 49/real(720),
1249 -27782109847927LL/real(250092722880000LL),
1250 99871724539LL/real(1569209241600LL), 3026004511LL/real(30648618000LL),
1251 -4306823/real(59875200), -2917/real(56700), 4463/real(90720),
1252 168979300892599LL/real(1600593426432000LL),
1253 2123926699/real(15324309000LL), -368661577/real(4036032000LL),
1254 -102293/real(1871100), 331799/real(7257600),
1255 1959350112697LL/real(9618950880000LL),
1256 -493031379277LL/real(3923023104000LL), -875457073/real(13621608000LL),
1257 11744233/real(239500800),
1258 -145659994071373LL/real(800296713216000LL),
1259 -793693009/real(9807557760LL), 453002260127LL/real(7846046208000LL),
1260 -53583096419057LL/real(500185445760000LL),
1261 103558761539LL/real(1426553856000LL),
1262 real(12272105438887727LL)/real(128047474114560000LL),
1263 // C[xi,chi]
1264 -64724382148LL/real(97692469875LL), 16676974/real(30405375),
1265 2706758/real(42567525), -55222/real(93555), 2458/real(4725),
1266 46/real(315), -34/real(45), 2/real(3),
1267 85904355287LL/real(37574026875LL), 158999572/real(1915538625),
1268 -340492279/real(212837625), 516944/real(467775), 3413/real(14175),
1269 -256/real(315), 19/real(45),
1270 2986003168LL/real(37574026875LL), -7597644214LL/real(1915538625),
1271 4430783356LL/real(1915538625), 206834/real(467775), -15958/real(14175),
1272 248/real(567),
1273 -375566203/real(39037950), 851209552/real(174139875),
1274 62016436/real(70945875), -832976/real(467775), 16049/real(28350),
1275 5106181018156LL/real(488462349375LL), 3475643362LL/real(1915538625),
1276 -651151712/real(212837625), 15602/real(18711),
1277 34581190223LL/real(8881133625LL), -10656173804LL/real(1915538625),
1278 2561772812LL/real(1915538625),
1279 -5150169424688LL/real(488462349375LL), 873037408/real(383107725),
1280 7939103697617LL/real(1953849397500LL),
1281 // C[xi,xi] skipped
1282 };
1283 static const int ptrs[] = {
1284 0, 0, 20, 40, 60, 96, 132, 152, 152, 172, 192, 228, 264, 284, 304, 304,
1285 324, 360, 396, 416, 436, 456, 456, 492, 528, 564, 600, 636, 672, 672,
1286 708, 744, 780, 816, 852, 888, 888,
1287 };
1288#else
1289#error "Unsupported value for GEOGRAPHICLIB_AUXLATITUDE_ORDER"
1290#endif
1291
1292 static_assert(sizeof(ptrs) / sizeof(int) == AUXNUMBER*AUXNUMBER+1,
1293 "Mismatch in size of ptrs array");
1294 static_assert(sizeof(coeffs) / sizeof(real) ==
1296 (Lmax * (Lmax + 3) - 2*(Lmax/2))/4 // Even only arrays
1298 (Lmax * (Lmax + 1))/2,
1299 "Mismatch in size of coeffs array");
1300
1301 if (k < 0) return; // auxout or auxin out of range
1302 if (auxout == auxin)
1303 fill(_c + Lmax * k, _c + Lmax * (k + 1), 0);
1304 else {
1305 int o = ptrs[k];
1306 real d = _n;
1307 if (auxin <= RECTIFYING && auxout <= RECTIFYING) {
1308 for (int l = 0; l < Lmax; ++l) {
1309 int m = (Lmax - l - 1) / 2; // order of polynomial in n^2
1310 _c[Lmax * k + l] = d * Math::polyval(m, coeffs + o, _n2);
1311 o += m + 1;
1312 d *= _n;
1313 }
1314 } else {
1315 for (int l = 0; l < Lmax; ++l) {
1316 int m = (Lmax - l - 1); // order of polynomial in n
1317 _c[Lmax * k + l] = d * Math::polyval(m, coeffs + o, _n);
1318 o += m + 1;
1319 d *= _n;
1320 }
1321 }
1322 // assert (o == ptrs[AUXNUMBER * auxout + auxin + 1])
1323 }
1324 }
1325
1326 Math::real AuxLatitude::Clenshaw(bool sinp, real szeta, real czeta,
1327 const real c[], int K) {
1328 // Evaluate
1329 // y = sum(c[k] * sin( (2*k+2) * zeta), i, 0, K-1) if sinp
1330 // y = sum(c[k] * cos( (2*k+2) * zeta), i, 0, K-1) if !sinp
1331 // Approx operation count = (K + 5) mult and (2 * K + 2) add
1332 int k = K;
1333 real u0 = 0, u1 = 0, // accumulators for sum
1334 x = 2 * (czeta - szeta) * (czeta + szeta); // 2 * cos(2*zeta)
1335 for (; k > 0;) {
1336 real t = x * u0 - u1 + c[--k];
1337 u1 = u0; u0 = t;
1338 }
1339 // u0*f0(zeta) - u1*fm1(zeta)
1340 // f0 = sinp ? sin(2*zeta) : cos(2*zeta)
1341 // fm1 = sinp ? 0 : 1
1342 real f0 = sinp ? 2 * szeta * czeta : x / 2, fm1 = sinp ? 0 : 1;
1343 return f0 * u0 - fm1 * u1;
1344 }
1345 /// \endcond
1346
1347} // namespace GeographicLib
Header for the GeographicLib::AuxLatitude class.
Header for GeographicLib::EllipticFunction class.
GeographicLib::Math::real real
Definition: GeodSolve.cpp:29
An accurate representation of angles.
Definition: AuxAngle.hpp:47
Math::real y() const
Definition: AuxAngle.hpp:70
Math::real x() const
Definition: AuxAngle.hpp:75
Math::real radians() const
Definition: AuxAngle.hpp:228
AuxAngle normalized() const
Definition: AuxAngle.cpp:31
Math::real degrees() const
Definition: AuxAngle.hpp:224
static AuxAngle NaN()
Definition: AuxAngle.cpp:26
Math::real tan() const
Definition: AuxAngle.hpp:113
AuxAngle copyquadrant(const AuxAngle &p) const
Definition: AuxAngle.cpp:47
Conversions between auxiliary latitudes.
Definition: AuxLatitude.hpp:75
AuxAngle Conformal(const AuxAngle &phi, real *diff=nullptr) const
AuxAngle Convert(int auxin, int auxout, const AuxAngle &zeta, bool exact=false) const
Math::real AuthalicRadiusSquared(bool exact=false) const
AuxAngle FromAuxiliary(int auxin, const AuxAngle &zeta, int *niter=nullptr) const
AuxAngle ToAuxiliary(int auxout, const AuxAngle &phi, real *diff=nullptr) const
Math::real RectifyingRadius(bool exact=false) const
AuxAngle Authalic(const AuxAngle &phi, real *diff=nullptr) const
AuxAngle Geocentric(const AuxAngle &phi, real *diff=nullptr) const
Definition: AuxLatitude.cpp:96
AuxAngle Parametric(const AuxAngle &phi, real *diff=nullptr) const
Definition: AuxLatitude.cpp:91
static Math::real Clenshaw(bool sinp, real szeta, real czeta, const real c[], int K)
static const AuxLatitude & WGS84()
Definition: AuxLatitude.cpp:86
AuxAngle Rectifying(const AuxAngle &phi, real *diff=nullptr) const
static real RG(real x, real y, real z)
static real RD(real x, real y, real z)
static real RF(real x, real y, real z)
Exception handling for GeographicLib.
Definition: Constants.hpp:316
static T sq(T x)
Definition: Math.hpp:205
static T pi()
Definition: Math.hpp:183
static T polyval(int N, const T p[], T x)
Definition: Math.hpp:264
@ td
degrees per turn
Definition: Math.hpp:142
Namespace for GeographicLib.
Definition: Accumulator.cpp:12
void swap(GeographicLib::NearestNeighbor< dist_t, pos_t, distfun_t > &a, GeographicLib::NearestNeighbor< dist_t, pos_t, distfun_t > &b)