13 #ifndef HEX_CLENSHAWCURTIS
14 #define HEX_CLENSHAWCURTIS
22 #include <gsl/gsl_sf.h>
48 Limit(true), Recurrence(true), NNest(5), NStack(5), L(1.0), Verbose(false),
49 vName(
"[ClenshawCurtis_ff]"), Throw(true) {}
52 inline double eps()
const {
return EpsRel; }
55 inline void setEps(
double epsrel) { EpsRel = epsrel; }
58 inline double tol()
const {
return EpsAbs; }
61 inline void setTol(
double epsabs) { EpsAbs = epsabs; }
70 inline int subdiv()
const {
return NNest; }
73 inline void setSubdiv(
int nlevel) { NNest = nlevel; }
76 inline int stack()
const {
return NStack; }
79 inline void setStack(
int nlevel) { NStack = nlevel; }
82 inline bool lim()
const {
return Limit; }
85 inline void setLim(
bool limit) { Limit = limit; }
88 inline bool Rec()
const {
return Recurrence; }
91 inline void setRec(
bool recurrence) { Recurrence = recurrence; }
94 inline bool verbose()
const {
return Verbose; }
116 FType
integrate (
double x1,
double x2,
int * n =
nullptr)
const
122 if (not std::isfinite(x1) and not std::isfinite(x2))
126 if (not std::isfinite(x1))
146 if (not std::isfinite(x2))
193 auto f = [&](
double x) -> FType
195 return F(x1 + 0.5 * (1.0 + x) * (x2 - x1));
199 std::vector<FType> coefs;
202 std::vector<Complex> fvals_prev, fvals, ftraf;
208 int maxN = gsl_sf_pow_int(2,NNest);
211 for (
int N = 4; N <= maxN ; N *= 2)
217 fvals.resize(2*N + 1);
221 fftw_plan plan = fftw_plan_dft_1d
224 reinterpret_cast<fftw_complex*>(&fvals[0]),
225 reinterpret_cast<fftw_complex*>(&ftraf[0]),
234 double pi_over_N = M_PI / N;
235 for (
int k = 0; k < N; k++)
237 fvals[k] = fvals[2*N-k] = f(
cos(k * pi_over_N));
239 if (not std::isfinite(
std::abs(fvals[k])))
240 throw exception(
"%s \"%g\" when evaluating function at %g", vName.c_str(), fvals[k],
cos(k * pi_over_N));
244 if (not std::isfinite(
std::abs(fvals[N])))
245 throw exception(
"%s \"%g\" when evaluating function at -1.", vName.c_str(), fvals[N]);
250 double pi_over_N = M_PI / N;
251 for (
int k = 0; k < N; k++)
253 fvals[k] = fvals[2*N-k] = (k % 2 == 0) ? fvals_prev[k/2] : f(
cos(k * pi_over_N));
255 if (not std::isfinite(
std::abs(fvals[k])))
256 throw exception(
"%s \"%g\" when evaluating function.", vName.c_str(), fvals[k]);
258 fvals[N] = fvals_prev[N/2];
265 fftw_destroy_plan(plan);
268 FType
const * ftraf_ptr =
reinterpret_cast<FType*
>(&ftraf[0]);
271 if (
typeid(FType) ==
typeid(
Complex))
274 for (
int i = 0; i <= N; i++)
275 coefs[i] = 0.5 * (*(ftraf_ptr + i));
277 else if (
typeid(FType) ==
typeid(double))
280 for (
int i = 0; i <= N; i++)
281 coefs[i] = 0.5 * (*(ftraf_ptr + 2*i));
285 throw exception(
"%s Can't handle datatype \"%s\".", vName.c_str(),
typeid(FType).name());
289 sum = 0.5 * (coefs[0] - coefs[N] / (N*N - 1.));
290 for (
int twok = 2; twok < N; twok += 2)
291 sum -= coefs[twok] / (twok*twok - 1.);
295 std::cout << vName <<
" N = " << N <<
", Sum = " << FType(2.*(x2-x1)/N)*sum <<
"\n";
301 std::cout << vName <<
" Convergence for N = " << N <<
", sum = " << FType(2. * (x2 - x1) / N) * sum <<
"\n";
303 return FType(2. * (x2 - x1) / N) *
sum;
305 else if ( std::isfinite(
std::abs(sum))
306 and std::isfinite(
std::abs(sum_prev))
310 std::cout << vName <<
" EpsAbs limit matched, " << EpsAbs <<
" on (" << x1 <<
"," << x2 <<
").\n";
320 fvals_prev = std::move(fvals);
331 throw exception(
"%s Insufficient evaluation limit %d", vName.c_str(), maxN);
335 std::cout << vName <<
" WARNING: Insufficient evaluation limit " << maxN <<
".\n";
336 return FType(2. * (x2 - x1) / maxN) *
sum;
348 std::cout << vName <<
" Interval smaller than " << EpsAbs <<
"\n";
356 std::cout << vName <<
" Bisection inhibited due to internal stack limit.\n";
357 return FType(2. * (x2 - x1) / maxN) *
sum;
362 std::cout << vName <<
" Bisecting to ("
363 << x1 <<
"," << (x2+x1)/2 <<
") and ("
364 << (x2+x1)/2 <<
"," << x2 <<
")\n";
370 double this_EpsAbs = EpsAbs;
373 EpsAbs = 0.5 * (2. * (x2 - x1) / maxN) *
std::abs(sum) * EpsRel;
379 EpsAbs = this_EpsAbs;
381 if (n !=
nullptr) *n = n1 + n2;
395 std::vector<FType> fvals, fvals_prev;
398 std::vector<double> weights, weights_prev;
404 for (
int N = 2; ; N *= 2)
410 if (fvals_prev.empty())
413 for (
int i = 1; i <= N - 1; i++)
415 double x = i * M_PI / N;
416 fvals[i] = F(L / tan(x));
417 if (not std::isfinite(
std::abs(fvals[i])))
419 weights[i] = 1. / gsl_sf_pow_int(
sin(x), 2);
425 for (
int i = 1; i < N - 1; i++)
429 fvals[i] = fvals_prev[i/2];
430 weights[i] = weights_prev[i/2];
434 double x = i * M_PI / N;
435 fvals[i] = F(L / tan(x));
436 if (not std::isfinite(
std::abs(fvals[i])))
438 weights[i] = 1. / gsl_sf_pow_int(
sin(x), 2);
445 for (
int i = 1; i <= N - 1; i++)
446 sum += weights[i] * fvals[i];
449 if (n !=
nullptr) *n = N;
450 return FType(L * M_PI / N) *
sum;
456 weights_prev = weights;
469 mutable double EpsAbs;
ClenshawCurtis(Functor const &f)
Definition: clenshawcurtis.h:47
bool lim() const
Get limit flag.
Definition: clenshawcurtis.h:82
void setEps(double epsrel)
Set relative tolerance.
Definition: clenshawcurtis.h:55
void setStack(int nlevel)
Set subdivision level limit.
Definition: clenshawcurtis.h:79
NumberArray< T > cos(NumberArray< T > const &A)
Return per-element cosine.
Definition: arrays.h:1747
FType integrate(double x1, double x2, int *n=nullptr) const
General integration.
Definition: clenshawcurtis.h:116
double getRange() const
Get compacification parameter.
Definition: clenshawcurtis.h:64
Clenshaw-Curtis quadrature.
Definition: clenshawcurtis.h:39
bool verbose() const
Get verbose flag.
Definition: clenshawcurtis.h:94
const double e
Definition: special.h:44
void setLim(bool limit)
Set limit flag.
Definition: clenshawcurtis.h:85
void setThrowAll(bool t)
Set warn flag.
Definition: clenshawcurtis.h:100
FType integrate_ii(int *n=nullptr) const
Improper integration.
Definition: clenshawcurtis.h:392
void setVerbose(bool verbose, std::string name="ClenshawCurtis_ff")
Set verbose flag.
Definition: clenshawcurtis.h:97
class NodeIntegrator exception
T sum(const ArrayView< T > v)
Sum elements in array.
Definition: arrays.h:1770
double tol() const
Get absolute tolerance.
Definition: clenshawcurtis.h:58
NumberArray< T > sin(NumberArray< T > const &A)
Return per-element sine.
Definition: arrays.h:1723
void setRange(double l)
Set compacification parameter.
Definition: clenshawcurtis.h:67
void setSubdiv(int nlevel)
Set 2-log of maximal subdivision interval count.
Definition: clenshawcurtis.h:73
int stack() const
Get subdivision level limit.
Definition: clenshawcurtis.h:76
void setRec(bool recurrence)
Set recurrence flag.
Definition: clenshawcurtis.h:91
void setTol(double epsabs)
Set absolute tolerance.
Definition: clenshawcurtis.h:61
rArray abs(const cArrayView u)
Definition: arrays.cpp:19
std::complex< double > Complex
Definition: complex.h:20
int subdiv() const
Get 2-log of maximal subdivision interval count.
Definition: clenshawcurtis.h:70
bool throwall() const
Get warn flag.
Definition: clenshawcurtis.h:103
T max(const ArrayView< T > a)
Maximal element of array.
Definition: arrays.h:1673
const double Nan
Definition: special.h:68
double eps() const
Get relative tolerance.
Definition: clenshawcurtis.h:52
Compactification multiplied by its jacobian.
Definition: compact.h:311
bool Rec() const
Get recurrence flag.
Definition: clenshawcurtis.h:88
FType integrate_ff(double x1, double x2, int *n=nullptr) const
Finite integration.
Definition: clenshawcurtis.h:178