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