25 bool odd,
bool sym,
real halfp,
bool centerp) {
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 );
31 for (
int i = 0; i < M; ++i)
33 *
this = initbysamples(F, odd, sym, halfp, centerp);
35 *
this = Trigfun{vector<real>{1,0}, odd,sym, halfp};
41 Trigfun t(n, f, odd, sym, halfp,
false);
43 int K = chop(t._coeff, tol, scale);
46 t._n = t._sym ? K : K - 1;
51 Trigfun tx(n, f, odd, sym, halfp,
true);
64 odd, sym, halfp,
false);
65 while (t._n <= nmax) {
66 int K = chop(t._coeff, tol, scale);
69 t._n = t._sym ? K : K - 1;
76 real p = halfp / (sym ? 2 : 1), d = p / t._n, o = d/2;
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));
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);
98 for (
int i = 0; i < n; ++i)
100 H[i + (odd ? 1 : 0)] = F[i];
102 H[n] = sym ? 0 : F[n];
106 for (
int i = 0; i < n; ++i)
107 H[2*n - i] = (odd ? 1 : -1) * H[i];
110 for (
int i = 1; i < M/2; ++i)
111 H[M - i] = (odd ? -1 : 1) * H[i];
114 for (
int i = 0; i < n; ++i)
118 for (
int i = 0; i < n; ++i)
119 H[2*n - i - 1] = (odd ? 1 : -1) * H[i];
122 for (
int i = 0; i < M/2; ++i)
123 H[M - i - 1] = (odd ? -1 : 1) * H[i];
127 fft_t fft(M/2,
false);
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();
133 for (
int i = 1; i <= M/2; ++i)
134 cF[i] *= exp(complex<real>(0, i * (-
Math::pi() / M)));
139 for (
int i = 0; i <= n; ++i)
140 H[i] = cF[i].
real() / n;
142 H[n] = centerp ? 0 : H[n]/2;
144 for (
int i = 0; i <= n; ++i)
145 H[i] = -cF[i].imag() / n;
147 H[n] = !centerp ? 0 : H[n]/2;
152 for (
int i = 0; i < n; ++i)
153 H[i] = cF[2*i+1].
real() / (2*n);
155 for (
int i = 0; i < n; ++i)
156 H[i] = -cF[2*i+1].imag() / (2*n);
160 Trigfun t(H, odd, sym, halfp);
161 if constexpr (debug_) {
162 real err = t.check(F, centerp, 0);
164 cout << t._n <<
" " << err << endl;
165 throw GeographicErr(
"initbysamples error");
171 Math::real Trigfun::check(
const vector<real>& F,
bool centerp,
real tol)
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),
181 maxval = fmax(maxval, fabs(a));
182 err = err + fabs(a - b);
185 return err / (tolerance(tol) *
186 maxval * (centerp ? _n : _n + 1));
189 void Trigfun::refine(
const Trigfun& tb) {
190 int m = 2 * _n + (_sym ? 0 : 1);
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;
204 if (!(_odd && !_sym))
206 (
"Trigfun: cannot set secular term unless odd && !sym");
213 for (
int k = _m; k > (_sym ? 0 : 1);)
214 _max += fabs(_coeff[--k]);
228 if (_coeff.empty())
return 0;
229 real y =
Math::pi()/(_sym ? _q : _h) * z;
231 k0 = _odd && !_sym ? 1 : 0;
235 real t = x * u0 - u1 + _coeff[--k];
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;
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));
281 if (_coeff.empty())
return 0;
283 k0 = _odd && !_sym ? 1 : 0;
285 x = 2 * (phi.c() - phi.s()) * (phi.c() + phi.s());
287 real t = x * u0 - u1 + _coeff[--k];
305 return _sym ? (_odd ? phi.s() * (u0 + u1) : phi.c() * (u0 - u1)) :
306 !_odd ? u0 - (x/2) * u1 : 2 * phi.s() * phi.c() * u0;
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));
318 c[0] = _odd ? 0 : _coeff[0] * mult;
319 return Trigfun(c, !_odd, _sym, _h);
326 int* countn,
int* countb,
real tol)
const {
327 if (!(_odd && !_sym))
329 (
"Trigfun: cannot take root unless odd && !sym");
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;
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);
346 const function<pair<real, real>(
real)>& ffp,
350 int* countn,
int* countb,
353 if (!(xa <= x0 && x0 <= xb))
return Math::NaN();
354 if (x0 == xa && x0 == xb)
356 tol = tolerance(tol);
357 real vtol = tol * zscale/100,
358 xtol = pow(tol,
real(0.75)) * xscale,
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";
369 << xa <<
" " << vala.first - z <<
" " << vala.second <<
"\n";
371 << x0 <<
" " << val0.first - z <<
" " << val0.second <<
"\n";
373 << xb <<
" " << valb.first - z <<
" " << valb.second <<
"\n";
374 if ((vala.first - z) * (valb.first - z) > 0)
378 (throw_ && (
throw GeographicLib::GeographicErr
379 (
"Convergence failure Trigfun::root case=" +
380 to_string(indicator)),
false));) {
385 pair<real, real> val = ffp(x);
386 real v = val.first - z,
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;
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 <<
" ";
407 if constexpr (debug_)
408 cout <<
"break3 " << k <<
" " << x <<
" " << dx <<
"\n";
411 }
else if (!(fabs(dx) > xtol)) {
412 if constexpr (debug_)
413 cout <<
"break2 " << k <<
" " << dx <<
" " << xtol << endl;
416 if constexpr (debug_)
417 cout <<
"GAPS " << k <<
" " << dx <<
" " << x-xa <<
" " << xb-x <<
" "
418 << oldx <<
" " << x <<
" " << (oldx - x) <<
"\n";
423 if (countn) *countn += k;
424 if (countb) *countb += b;
425 if constexpr (debug_)
426 cout <<
"return " << x <<
"\n";
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) -
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();
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);
456 int Trigfun::chop(
const vector<real>& c,
real tol,
real scale) {
510 tol = tolerance(tol);
511 if (tol >= 1)
return 1;
514 int n = int(c.size());
517 if (n < 16)
return n;
524 m[--j] = fabs(c[n - 1]);
527 m[j] = fmax(fabs(c[j]), m[j + 1]);
529 if (m[0] == 0)
return 1;
530 if (scale >= 0) m[0] = fmax(scale, m[0]);
545 int j2 = 0, plateauPoint = n;
546 real logtol = log(tol);
547 for (j = 2; j <= n; ++j) {
548 j2 = int(round(1.25*j + 5));
549 if (j2 > n)
return n;
552 r = 3 * (1 - log(e1)/logtol);
553 if ( e1 == 0 || e2/e1 > r ) {
555 plateauPoint = j - 1;
583 if ( m[plateauPoint - 1] == 0 )
return plateauPoint;
584 real tol76 = tol * sqrt(cbrt(tol));
586 for (j = 0; j < n; ++j) {
587 if (m[j] >= tol76) ++j3;
595 real tol3 = logtol/(3 * (j2 - 1));
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;
605 bool sym,
real scale)
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())))
GeographicLib::Math::real real
Header for GeographicLib::Trigfun class.
static AngleT radians(Math::real rad)
Representing a function by a Fourier series.
real operator()(real x) const
Namespace for GeographicLib.
AngleT< Math::real > Angle