GeographicLib 2.6
Loading...
Searching...
No Matches
Trigfun.cpp
Go to the documentation of this file.
1/**
2 * \file Trigfun.cpp
3 * \brief Implementation for GeographicLib::Trigfun class
4 *
5 * Copyright (c) Charles Karney (2024-2025) <karney@alum.mit.edu> and licensed
6 * under the MIT/X11 License. For more information, see
7 * https://geographiclib.sourceforge.io/
8 **********************************************************************/
9
10/// \cond SKIP
11// For to_string
12#include <string>
13#include <iostream>
14#include <iomanip>
16#include "kissfft.hh"
17
18#define USE_ANGLE 0
19
20namespace GeographicLib {
21
22 using namespace std;
23
24 Trigfun::Trigfun(int n, const std::function<real(real)>& f,
25 bool odd, bool sym, real halfp, bool centerp) {
26 if (n > 0) {
27 int M = n + (!(odd || sym || centerp) ? 1 : 0);
28 real p = halfp / (sym ? 2 : 1), d = p / n,
29 o = centerp ? d/2 : ( odd ? d : 0 );
30 vector<real> F(M);
31 for (int i = 0; i < M; ++i)
32 F[i] = f(o + d * i);
33 *this = initbysamples(F, odd, sym, halfp, centerp);
34 } else
35 *this = Trigfun{vector<real>{1,0}, odd,sym, halfp};
36 }
37
38 Trigfun::Trigfun(const function<real(real)>& f, bool odd, bool sym,
39 real halfp, int nmax, real tol, real scale) {
40 int n = 16;
41 Trigfun t(n, f, odd, sym, halfp, false);
42 while (n <= nmax) {
43 int K = chop(t._coeff, tol, scale);
44 if (K < n) {
45 t._m = K;
46 t._n = t._sym ? K : K - 1;
47 t._coeff.resize(K);
48 *this = t;
49 return;
50 }
51 Trigfun tx(n, f, odd, sym, halfp, true);
52 t.refine(tx);
53 n *= 2;
54 }
55 *this = t;
56 }
57
58 Trigfun::Trigfun(const function<real(real, real)>& f, bool odd, bool sym,
59 real halfp, int nmax, real tol, real scale) {
60 // Initialize with 2 samples
61 Trigfun t(2,
62 [&f] (real x) -> real
63 { return f(x, Math::NaN()); },
64 odd, sym, halfp, false);
65 while (t._n <= nmax) {
66 int K = chop(t._coeff, tol, scale);
67 if (K < t._n) {
68 t._m = K;
69 t._n = t._sym ? K : K - 1;
70 t._coeff.resize(K);
71 *this = t;
72 return;
73 }
74 // bool centerp = true;
75 int M = t._n;
76 real p = halfp / (sym ? 2 : 1), d = p / t._n, o = d/2;
77 vector<real> F(M);
78 for (int i = 0; i < M; ++i)
79 F[i] = f(o + d * i, t(o + d * i));
80 t.refine(initbysamples(F, odd, sym, halfp, true));
81 // n *= 2;
82 }
83 *this = t;
84 return;
85 }
86
87 Trigfun Trigfun::initbysamples(const vector<real>& F,
88 bool odd, bool sym, real halfp, bool centerp) {
89 if (!(isfinite(halfp) && halfp > 0))
90 throw GeographicErr("Trigfun::initbysamples halfp not positive");
91 using fft_t = kissfft<real>;
92 int n = int(F.size()) - (!(odd || sym || centerp) ? 1 : 0),
93 M = n * (sym ? 4 : 2); // The size of the sample array over a period
94 vector<real> H(M, Math::NaN());
95 if (!centerp) {
96 if (odd) H[0] = 0;
97 // real slope = (odd & !sym) ? F[n-1] / n : 0;
98 for (int i = 0; i < n; ++i)
99 // H[i + (odd ? 1 : 0)] = F[i] - slope * i;
100 H[i + (odd ? 1 : 0)] = F[i];
101 if (!odd) {
102 H[n] = sym ? 0 : F[n];
103 }
104 // Now H[0:n] is populated
105 if (sym) {
106 for (int i = 0; i < n; ++i)
107 H[2*n - i] = (odd ? 1 : -1) * H[i];
108 }
109 // Now H[0:M/2] is populated
110 for (int i = 1; i < M/2; ++i)
111 H[M - i] = (odd ? -1 : 1) * H[i];
112 // Now H[0:M-1] is populated
113 } else {
114 for (int i = 0; i < n; ++i)
115 H[i] = F[i];
116 // Now H[0:n-1] is populated
117 if (sym) {
118 for (int i = 0; i < n; ++i)
119 H[2*n - i - 1] = (odd ? 1 : -1) * H[i];
120 }
121 // Now H[0:M/2-1] is populated
122 for (int i = 0; i < M/2; ++i)
123 H[M - i - 1] = (odd ? -1 : 1) * H[i];
124 // Now H[0:M-1] is populated
125 }
126 // cout << "FFT size " << M/2 << "\n";
127 fft_t fft(M/2, false);
128 // Leave an extra slot
129 vector<complex<real>> cF(M/2 + 1);
130 fft.transform_real(H.data(), cF.data());
131 cF[M/2] = cF[0].imag(); cF[0] = cF[0].real();
132 if (centerp) {
133 for (int i = 1; i <= M/2; ++i)
134 cF[i] *= exp(complex<real>(0, i * (-Math::pi() / M)));
135 }
136 if (!sym) {
137 H.resize(n+1);
138 if (!odd) {
139 for (int i = 0; i <= n; ++i)
140 H[i] = cF[i].real() / n;
141 H[0] /= 2;
142 H[n] = centerp ? 0 : H[n]/2;
143 } else {
144 for (int i = 0; i <= n; ++i)
145 H[i] = -cF[i].imag() / n;
146 H[0] = 0;
147 H[n] = !centerp ? 0 : H[n]/2;
148 }
149 } else { // sym
150 H.resize(n);
151 if (!odd) {
152 for (int i = 0; i < n; ++i)
153 H[i] = cF[2*i+1].real() / (2*n);
154 } else {
155 for (int i = 0; i < n; ++i)
156 H[i] = -cF[2*i+1].imag() / (2*n);
157 }
158 }
159 // if (centerp) cout << "SIZE " << F.size() << " " << H.size() << "\n";
160 Trigfun t(H, odd, sym, halfp);
161 if constexpr (debug_) {
162 real err = t.check(F, centerp, 0);
163 if (err > 100) {
164 cout << t._n << " " << err << endl;
165 throw GeographicErr("initbysamples error");
166 }
167 }
168 return t;
169 }
170
171 Math::real Trigfun::check(const vector<real>& F, bool centerp, real tol)
172 const {
173 real err = 0, maxval = 0;
174 real d = (_sym ? _h/2 : _h) / _n;
175 for (int i = 0; i < (centerp ? _n : _n + 1); ++i) {
176 real a = centerp ? F[i] :
177 (_odd ? (i == 0 ? 0 : F[i-1]) :
178 (_sym && i == _n ? 0 : F[i])),
179 x = d * i + (centerp ? d/2 : 0),
180 b = (*this)(x);
181 maxval = fmax(maxval, fabs(a));
182 err = err + fabs(a - b);
183 }
184 // cout << "Maxval " << maxval << "\n";
185 return err / (tolerance(tol) *
186 maxval * (centerp ? _n : _n + 1));
187 }
188
189 void Trigfun::refine(const Trigfun& tb) {
190 int m = 2 * _n + (_sym ? 0 : 1);
191 _coeff.resize(m);
192 for (int i = 0; i < _n; ++i)
193 _coeff[2*_n + (_sym ? 0 : 1) - 1 - i] =
194 (_odd ? -1 : 1) * (_coeff[i] - tb._coeff[i])/2;
195 if (_odd && !_sym) _coeff[_n] = tb._coeff[_n];
196 for (int i = 0; i < _n; ++i)
197 _coeff[i] = (_coeff[i] + tb._coeff[i])/2;
198 _max = -1;
199 _n *= 2;
200 _m = m;
201 }
202
203 void Trigfun::setsecular(real f0) {
204 if (!(_odd && !_sym))
205 throw GeographicErr
206 ("Trigfun: cannot set secular term unless odd && !sym");
207 _coeff[0] = f0 / Math::pi();
208 }
209
210 Math::real Trigfun::Max() const {
211 if (_max < 0) {
212 _max = 0;
213 for (int k = _m; k > (_sym ? 0 : 1);)
214 _max += fabs(_coeff[--k]);
215 }
216 return _max;
217 }
218
219#if !USE_ANGLE
221 // Evaluate
222 // y = sum(c[k] * sin((k+1/2) * pi/q * z), k, 0, n - 1) if odd && sym
223 // y = sum(c[k] * cos((k+1/2) * pi/q * z), k, 0, n - 1) if !odd && sym
224 // y = c[0] * pi/h * z +
225 // sum(c[k] * sin(k * pi/h * z), k, 1, n) if odd && !sym
226 // y = c[0] +
227 // sum(c[k] * cos(k * pi/h * z), k, 1, n) if !odd && !sym
228 if (_coeff.empty()) return 0;
229 real y = Math::pi()/(_sym ? _q : _h) * z;
230 int k = _m,
231 k0 = _odd && !_sym ? 1 : 0; // add secular term at the end
232 real u0 = 0, u1 = 0, // accumulators for sum
233 x = 2 * cos(y);
234 for (; k > k0;) {
235 real t = x * u0 - u1 + _coeff[--k];
236 u1 = u0; u0 = t;
237 }
238 // sym
239 // v = u0*f0(zeta) - u1*fm1(zeta)
240 // f0 = odd ? sin(y/2) : cos(y/2)
241 // fm1 = odd ? -sin(y/2) : cos(y/2)
242 // v = odd ? sin(y/2) * (u0 + u1) : cos(y/2) * (u0 - u1)
243 // !sym && !odd
244 // v = u0*f0(zeta) - u1*fm1(zeta)
245 // f0 = 1
246 // fm1 = cos(y)
247 // v = u0 - cos(y) * u1 = u0 - (x/2) * u1
248 // !sym && odd
249 // v = u0*f1(zeta) - u1*f0(zeta)
250 // f1 = sin(y)
251 // f0 = 0
252 // v = sin(y) * u0 + secular term
253 return _sym ? (_odd ? sin(y/2) * (u0 + u1) : cos(y/2) * (u0 - u1)) :
254 !_odd ? u0 - (x/2) * u1 : _coeff[0] * y + sin(y) * u0;
255 }
256
257#else
258 // Implementation in terms of Angle
260 // Evaluate
261 // y = sum(c[k] * sin((k+1/2) * pi/q * z), k, 0, n - 1) if odd && sym
262 // y = sum(c[k] * cos((k+1/2) * pi/q * z), k, 0, n - 1) if !odd && sym
263 // y = c[0] * pi/h * z +
264 // sum(c[k] * sin(k * pi/h * z), k, 1, n) if odd && !sym
265 // y = c[0] +
266 // sum(c[k] * cos(k * pi/h * z), k, 1, n) if !odd && !sym
267 if (_coeff.empty()) return 0;
268 real y = Math::pi()/(_sym ? _q : _h) * z;
269 return (_odd && !_sym ? _coeff[0] * y : 0) + eval(Angle::radians(y/2));
270 }
271
272 Math::real Trigfun::eval(Angle phi) const {
273 // Evaluate
274 // y = sum(c[k] * sin((2*k+1) * phi), k, 0, n - 1) if odd && sym
275 // y = sum(c[k] * cos((2*k+1) * phi), k, 0, n - 1) if !odd && sym
276 // y = c[0] * 0 +
277 // sum(c[k] * sin(2*k * phi), k, 1, n) if odd && !sym
278 // y = c[0] +
279 // sum(c[k] * cos(2*k * phi), k, 1, n) if !odd && !sym
280 // Note that secular term c[0] is ignored for odd && !sym.
281 if (_coeff.empty()) return 0;
282 int k = _m,
283 k0 = _odd && !_sym ? 1 : 0; // secular term excluded
284 real u0 = 0, u1 = 0, // accumulators for sum
285 x = 2 * (phi.c() - phi.s()) * (phi.c() + phi.s());
286 for (; k > k0;) {
287 real t = x * u0 - u1 + _coeff[--k];
288 u1 = u0; u0 = t;
289 }
290 // sym
291 // v = u0*f0(zeta) - u1*fm1(zeta)
292 // f0 = odd ? sin(phi) : cos(phi)
293 // fm1 = odd ? -sin(phi) : cos(phi)
294 // v = odd ? sin(phi) * (u0 + u1) : cos(phi) * (u0 - u1)
295 // !sym && !odd
296 // v = u0*f0(zeta) - u1*fm1(zeta)
297 // f0 = 1
298 // fm1 = cos(2*phi)
299 // v = u0 - cos(2*phi) * u1 = u2 - (x/2) * u1
300 // !sym && odd
301 // v = u0*f1(zeta) - u1*f0(zeta)
302 // f1 = sin(2*phi)
303 // f0 = 0
304 // v = sin(2*phi) * u0 + no secular term
305 return _sym ? (_odd ? phi.s() * (u0 + u1) : phi.c() * (u0 - u1)) :
306 !_odd ? u0 - (x/2) * u1 : 2 * phi.s() * phi.c() * u0;
307 }
308#endif
309
310 Trigfun Trigfun::integral() const {
311 if (_odd && !_sym && _coeff[0] != 0)
312 throw GeographicErr("Trigfun: cannot integrate a secular term");
313 vector <real> c(_coeff);
314 real mult = (_odd ? -1 : 1) * (_sym ? _q : _h) / Math::pi();
315 for (int i = 0; i < _m; ++i)
316 c[i] *= mult / (i + (_sym ? real(0.5) : 0));
317 if (!_sym)
318 c[0] = _odd ? 0 : _coeff[0] * mult;
319 return Trigfun(c, !_odd, _sym, _h);
320 }
321
322 // root sig 2
323 Math::real Trigfun::root(ind indicator,
324 real z, const function<real(real)>& fp,
325 real x0,
326 int* countn, int* countb, real tol) const {
327 if (!(_odd && !_sym))
328 throw GeographicErr
329 ("Trigfun: cannot take root unless odd && !sym");
330 // y = pi/h * x
331 // f(x) = c[0] * y + sum(c[k] * sin(k * y), k, 1, n)
332 real hr = Math::pi() * _coeff[0], s = _h / hr,
333 x00 = s * z, dx = fabs(s) * Max();
334 x0 = isfinite(x0) ? fmin(x00 + dx, fmax(x00 - dx, x0)) : x00;
335 // cout << "QQG " << dx << "\n";
336 if (dx == 0) return x0;
337 auto ffp = [this, &fp]
338 (real x) -> pair<real, real>
339 { return pair<real, real>(this->operator()(x), fp(x)); };
340 return root(indicator,ffp, z, x0, x00 - dx, x00 + dx, _h, fabs(hr),
341 s > 0 ? 1 : -1, countn, countb, tol);
342 }
343
344 // root sig 4
345 Math::real Trigfun::root(ind indicator,
346 const function<pair<real, real>(real)>& ffp,
347 real z,
348 real x0, real xa, real xb,
349 real xscale, real zscale, int s,
350 int* countn, int* countb,
351 real tol) {
352 // Solve v = f(x) - z = 0
353 if (!(xa <= x0 && x0 <= xb)) return Math::NaN();
354 if (x0 == xa && x0 == xb)
355 return x0;
356 tol = tolerance(tol);
357 real vtol = tol * zscale/100,
358 xtol = pow(tol, real(0.75)) * xscale,
359 x = x0, oldx = Math::infinity(), oldv = oldx, olddx = oldx;
360 int k = 0, b = 0;
361 real p = Math::pi()/2 * 0;
362 if constexpr (debug_) {
363 cout << "SCALE " << xscale << " " << zscale << "\n";
364 pair<real, real> vala = ffp(xa);
365 pair<real, real> val0 = ffp(x0);
366 pair<real, real> valb = ffp(xb);
367 cout << "DAT " << s << " " << x0-xa << " " << xb-x0 << " " << z << "\n";
368 cout << "DAT "
369 << xa << " " << vala.first - z << " " << vala.second << "\n";
370 cout << "DAT "
371 << x0 << " " << val0.first - z << " " << val0.second << "\n";
372 cout << "DAT "
373 << xb << " " << valb.first - z << " " << valb.second << "\n";
374 if ((vala.first - z) * (valb.first - z) > 0)
375 cout << "DATBAD\n";
376 }
377 for (; k < maxit_ ||
378 (throw_ && (throw GeographicLib::GeographicErr
379 ("Convergence failure Trigfun::root case=" +
380 to_string(indicator)), false));) {
381 // TODO: This inverse problem uses lots of iterations
382 // 20 60 -90 180 127.4974 24.6254 2.4377
383 // Need to figure out why. (Probably fixed by now.)
384 ++k;
385 pair<real, real> val = ffp(x);
386 real v = val.first - z,
387 vp = val.second,
388 dx = - v/vp;
389 if constexpr (debug_)
390 cout << "XX " << k << " " << xa-p << " " << x-p << " " << xb-p << " "
391 << dx << " " << x + dx-p << " " << v << " " << vp << endl;
392 if (!(fabs(v) > (k < 2 ? 0 : vtol))) {
393 if constexpr (debug_) cout << "break1 " << k << " " << fabs(v) << endl;
394 break;
395 } else if (s*v > 0)
396 xb = fmin(xb, x);
397 else
398 xa = fmax(xa, x);
399 x += dx;
400 if (!(xa <= x && x <= xb) || fabs(v) > oldv ||
401 (k > 2 && 2 * fabs(dx) > olddx)) {
402 if constexpr (debug_)
403 cout << "bis " << k << " " << xa-x << " " << x-xb << " ";
404 x = (xa + xb)/2;
405 ++b;
406 if (x == oldx) {
407 if constexpr (debug_)
408 cout << "break3 " << k << " " << x << " " << dx << "\n";
409 break;
410 }
411 } else if (!(fabs(dx) > xtol)) {
412 if constexpr (debug_)
413 cout << "break2 " << k << " " << dx << " " << xtol << endl;
414 break;
415 }
416 if constexpr (debug_)
417 cout << "GAPS " << k << " " << dx << " " << x-xa << " " << xb-x << " "
418 << oldx << " " << x << " " << (oldx - x) << "\n";
419 oldx = x;
420 oldv = fabs(v);
421 olddx = fabs(dx);
422 }
423 if (countn) *countn += k;
424 if (countb) *countb += b;
425 if constexpr (debug_)
426 cout << "return " << x << "\n";
427 return x;
428 }
429
430 Math::real Trigfun::inversep(real z,
431 const function<real(real)>& fp,
432 real dx0,
433 int* countn, int* countb, real tol) const {
434 real hr = Math::pi() * _coeff[0], nslope = _h / hr;
435 return root(INVERSEP, z, fp, z * nslope + dx0, countn, countb, tol) -
436 nslope * z;
437 }
438
439 Trigfun Trigfun::inverse(const function<real(real)>& fp,
440 int* countn, int* countb,
441 int nmax, real tol, real scale) const {
442 if (!(_odd && !_sym && isfinite(_coeff[0]) && _coeff[0] != 0))
443 throw GeographicErr("Can only invert Trigfun with a secular term");
444 int s = _coeff[0] > 0 ? 1 : -1;
445 real hp = _h, hr = Math::pi() * _coeff[0],
446 nhp = hr * s, nhr = hp * s, c0p = nhr / Math::pi();
447 Trigfun t(
448 [this, &fp, countn, countb, tol]
449 (real z, real dx0) -> real
450 { return inversep(z, fp, dx0, countn, countb, tol); },
451 _odd, _sym, nhp, nmax, tol, scale);
452 t._coeff[0] = c0p;
453 return t;
454 }
455
456 int Trigfun::chop(const vector<real>& c, real tol, real scale) {
457 // This is a clone of Chebfun's standardChop function. For C++, the return
458 // value is number of terms to retain. Index of last term is one less than
459 // this.
460 //
461 // See J. L. Aurentz and L. N. Trefethen, "Chopping a Chebyshev series",
462 // https://doi.org/10.1145/2998442 (2017) and
463 // https://arxiv.org/abs/1512.01803 (2015).
464 //
465 // Input:
466 //
467 // COEFFS A nonempty row or column vector of real or complex numbers
468 // which typically will be Chebyshev or Fourier coefficients.
469 //
470 // TOL A number in (0,1) representing a target relative accuracy.
471 // TOL will typically will be set to the Chebfun EPS parameter,
472 // sometimes multiplied by a factor such as vglobal/vlocal in
473 // construction of local pieces of global chebfuns.
474 // Default value: machine epsilon (MATLAB EPS).
475 //
476 // Output:
477 //
478 // CUTOFF A positive integer.
479 // If CUTOFF == length(COEFFS), then we are "not happy":
480 // a satisfactory chopping point has not been found.
481 // If CUTOFF < length(COEFFS), we are "happy" and CUTOFF
482 // represents the last index of COEFFS that should be retained.
483 //
484 // Examples:
485 //
486 // coeffs = 10.^-(1:50); random = cos((1:50).^2);
487 // standardChop(coeffs) // = 18
488 // standardChop(coeffs + 1e-16*random) // = 15
489 // standardChop(coeffs + 1e-13*random) // = 13
490 // standardChop(coeffs + 1e-10*random) // = 50
491 // standardChop(coeffs + 1e-10*random, 1e-10) // = 10
492
493 // Jared Aurentz and Nick Trefethen, July 2015.
494 //
495 // Copyright 2017 by The University of Oxford and The Chebfun Developers.
496 // See http://www.chebfun.org/ for Chebfun information.
497
498 // STANDARDCHOP normally chops COEFFS at a point beyond which it is smaller
499 // than TOL^(2/3). COEFFS will never be chopped unless it is of length at
500 // least 17 and falls at least below TOL^(1/3). It will always be chopped
501 // if it has a long enough final segment below TOL, and the final entry
502 // COEFFS(CUTOFF) will never be smaller than TOL^(7/6). All these
503 // statements are relative to MAX(ABS(COEFFS)) and assume CUTOFF > 1.
504 // These parameters result from extensive experimentation involving
505 // functions such as those presented in the paper cited above. They are
506 // not derived from first principles and there is no claim that they are
507 // optimal.
508
509 // Check magnitude of TOL:
510 tol = tolerance(tol);
511 if (tol >= 1) return 1;
512
513 // Make sure c has length at least 17:
514 int n = int(c.size());
515 // Change 17 in original code to 16 to accommodate trig expansions which
516 // may only have 2^n terms.
517 if (n < 16) return n;
518
519 // Step 1: Convert c to a new monotonically nonincreasing
520 // vector ENVELOPE normalized to begin with the value 1.
521
522 vector<real> m(n);
523 int j = n;
524 m[--j] = fabs(c[n - 1]);
525 for (; j;) {
526 --j;
527 m[j] = fmax(fabs(c[j]), m[j + 1]);
528 }
529 if (m[0] == 0) return 1;
530 if (scale >= 0) m[0] = fmax(scale, m[0]);
531 for (j = n; j;)
532 m[--j] /= m[0];
533
534 // Step 2: Scan ENVELOPE for a value PLATEAUPOINT, the first point J-1, if
535 // any, that is followed by a plateau. A plateau is a stretch of
536 // coefficients ENVELOPE(J),...,ENVELOPE(J2), J2 = round(1.25*J+5) <= n,
537 // with the property that ENVELOPE(J2)/ENVELOPE(J) > R. The number R
538 // ranges from R = 0 if ENVELOPE(J) = TOL up to R = 1 if ENVELOPE(J) =
539 // TOL^(2/3). Thus a potential plateau whose starting value is ENVELOPE(J)
540 // ~ TOL^(2/3) has to be perfectly flat to count, whereas with ENVELOPE(J)
541 // ~ TOL it doesn't have to be flat at all. If a plateau point is found,
542 // then we know we are going to chop the vector, but the precise chopping
543 // point CUTOFF still remains to be determined in Step 3.
544
545 int j2 = 0, plateauPoint = n;
546 real logtol = log(tol);
547 for (j = 2; j <= n; ++j) { // j is a MATLAB index (starts at 1)
548 j2 = int(round(1.25*j + 5));
549 if (j2 > n) return n;
550 real e1 = m[j-1],
551 e2 = m[j2-1],
552 r = 3 * (1 - log(e1)/logtol);
553 if ( e1 == 0 || e2/e1 > r ) {
554 // a plateau has been found: go to Step 3
555 plateauPoint = j - 1;
556 break;
557 }
558 }
559
560 // Step 3: fix CUTOFF at a point where ENVELOPE, plus a linear function
561 // included to bias the result towards the left end, is minimal.
562 //
563 // Some explanation is needed here. One might imagine that if a plateau is
564 // found, then one should simply set CUTOFF = PLATEAUPOINT and be done,
565 // without the need for a Step 3. However, sometimes CUTOFF should be
566 // smaller or larger than PLATEAUPOINT, and that is what Step 3 achieves.
567 //
568 // CUTOFF should be smaller than PLATEAUPOINT if the last few coefficients
569 // made negligible improvement but just managed to bring the vector
570 // ENVELOPE below the level TOL^(2/3), above which no plateau will ever be
571 // detected. This part of the code is important for avoiding situations
572 // where a coefficient vector is chopped at a point that looks "obviously
573 // wrong" with PLOTCOEFFS.
574 //
575 // CUTOFF should be larger than PLATEAUPOINT if, although a plateau has
576 // been found, one can nevertheless reduce the amplitude of the
577 // coefficients a good deal further by taking more of them. This will
578 // happen most often when a plateau is detected at an amplitude close to
579 // TOL, because in this case, the "plateau" need not be very flat. This
580 // part of the code is important to getting an extra digit or two beyond
581 // the minimal prescribed accuracy when it is easy to do so.
582
583 if ( m[plateauPoint - 1] == 0 ) return plateauPoint;
584 real tol76 = tol * sqrt(cbrt(tol)); // tol^(7/6)
585 int j3 = 0;
586 for (j = 0; j < n; ++j) {
587 if (m[j] >= tol76) ++j3;
588 }
589 if ( j3 < j2 ) {
590 j2 = j3 + 1;
591 m[j2 - 1] = tol76;
592 }
593 vector<real> cc(j2);
594 // Replace log10 by log. This involved no change in the logic.
595 real tol3 = logtol/(3 * (j2 - 1));
596 int d = 0;
597 for (j = 0; j < j2; ++j) {
598 cc[j] = log(m[j]) - tol3 * j;
599 if (j > 0 && cc[j] < cc[d]) d = j;
600 }
601 return max(d, 1);
602 }
603
604 TrigfunExt::TrigfunExt(const function<real(real)>& fp, real halfp,
605 bool sym, real scale)
606 : _fp(fp)
607 , _sym(sym)
608 // N.B. tol defaults to epsilon() here. We need to compute the
609 // integral accurately.
610 , _f(Trigfun(_fp, false, _sym, halfp, 1 << 16, 0, scale).integral())
611 , _tol(sqrt(numeric_limits<real>::epsilon()))
612 , _nmax(int(ceil(real(1.5) * _f.NCoeffs())))
613 , _invp(false)
614 {}
615
616} // namespace GeographicLib
617/// \endcond
GeographicLib::Math::real real
Header for GeographicLib::Trigfun class.
static AngleT radians(Math::real rad)
static T infinity()
Definition Math.cpp:308
static T pi()
Definition Math.hpp:187
static T NaN()
Definition Math.cpp:301
Representing a function by a Fourier series.
Definition Trigfun.hpp:65
real operator()(real x) const
Trigfun integral() const
void setsecular(real f0)
Namespace for GeographicLib.
AngleT< Math::real > Angle
Definition Angle.hpp:760