Theoretica
A C++ numerical and automatic mathematical library
integral.h
Go to the documentation of this file.
1 
5 
6 #ifndef THEORETICA_INTEGRATION_H
7 #define THEORETICA_INTEGRATION_H
8 
9 #include "../core/constants.h"
10 #include "../core/function.h"
11 #include "../polynomial/polynomial.h"
12 #include "../polynomial/orthogonal.h"
13 #include "./gauss.h"
14 
15 
16 namespace theoretica {
17 
18 
23  template<typename T = real>
25 
26  polynomial<T> P;
27  P.coeff.resize(p.size() + 1);
28  P[0] = 0;
29 
30  for (unsigned int i = 0; i < p.size(); ++i)
31  P[i + 1] = p[i] / T(i + 1);
32 
33  return P;
34  }
35 
36 
45  template<typename RealFunction>
46  inline real integral_midpoint(RealFunction f, real a, real b,
47  unsigned int steps = CALCULUS_INTEGRAL_STEPS) {
48 
49  if(steps == 0) {
50  TH_MATH_ERROR("integral_midpoint", steps, DIV_BY_ZERO);
51  return nan();
52  }
53 
54  const real dx = (b - a) / steps;
55  real res = 0;
56 
57  for (unsigned int i = 0; i < steps; ++i)
58  res += f(a + (i + 0.5) * dx);
59 
60  return res * dx;
61  }
62 
63 
72  template<typename RealFunction>
73  inline real integral_trapezoid(RealFunction f, real a, real b,
74  unsigned int steps = CALCULUS_INTEGRAL_STEPS) {
75 
76  if(steps == 0) {
77  TH_MATH_ERROR("integral_trapezoid", steps, DIV_BY_ZERO);
78  return nan();
79  }
80 
81  const real dx = (b - a) / steps;
82  real res = 0;
83 
84  res += 0.5 * (f(a) + f(b));
85 
86  for (unsigned int i = 1; i < steps; ++i)
87  res += f(a + i * dx);
88 
89  return res * dx;
90  }
91 
92 
101  template<typename RealFunction>
102  inline real integral_simpson(RealFunction f, real a, real b,
103  unsigned int steps = CALCULUS_INTEGRAL_STEPS) {
104 
105  if(steps == 0) {
106  TH_MATH_ERROR("integral_simpson", steps, DIV_BY_ZERO);
107  return nan();
108  }
109 
110  const real dx = (b - a) / (real) steps;
111  real res = 0;
112 
113  // Sum terms by order of magnitude supposing that
114  // f stays at the same order inside the interval,
115  // to alleviate truncation errors
116 
117  res += f(a) + f(b);
118 
119  for (unsigned int i = 2; i < steps; i += 2)
120  res += 2 * f(a + i * dx);
121 
122  for (unsigned int i = 1; i < steps; i += 2)
123  res += 4 * f(a + i * dx);
124 
125  return res * dx / 3.0;
126  }
127 
128 
137  template<typename RealFunction>
139  RealFunction f,
140  real a, real b,
141  unsigned int iter = 8) {
142 
143  real T[iter][iter];
144 
145  T[0][0] = (f(a) + f(b)) * (b - a) / 2.0;
146 
147  for (unsigned int j = 1; j < iter; ++j) {
148 
149  // Composite Trapezoidal Rule
150  T[j][0] = integral_trapezoid(f, a, b, 1 << j);
151 
152  // Richardson extrapolation
153  for (unsigned int k = 1; k <= j; ++k) {
154  const uint64_t coeff = uint64_t(1) << (2 * k);
155  T[j][k] = (coeff * T[j][k - 1] - T[j - 1][k - 1]) / (coeff - 1);
156  }
157  }
158 
159  // Return the best approximation
160  return T[iter - 1][iter - 1];
161  }
162 
163 
172  template<typename RealFunction>
174  RealFunction f,
175  real a, real b,
176  real tolerance = CALCULUS_INTEGRAL_TOL) {
177 
178  const unsigned int MAX_ROMBERG_ITER = 16;
179  real T[MAX_ROMBERG_ITER][MAX_ROMBERG_ITER];
180 
181  T[0][0] = (f(a) + f(b)) * (b - a) / 2.0;
182 
183  for (unsigned int j = 1; j < 16; ++j) {
184 
185  // Composite Trapezoidal Rule
186  T[j][0] = integral_trapezoid(f, a, b, 1 << j);
187 
188  // Richardson extrapolation
189  for (unsigned int k = 1; k <= j; ++k) {
190  const uint64_t coeff = uint64_t(1) << (2 * k);
191  T[j][k] = (coeff * T[j][k - 1] - T[j - 1][k - 1]) / (coeff - 1);
192  }
193 
194  // Stop the algorithm when the desired precision has been reached
195  if(abs(T[j][j] - T[j - 1][j - 1]) < tolerance)
196  return T[j][j];
197  }
198 
199  // Return the best approximation
200  return T[MAX_ROMBERG_ITER - 1][MAX_ROMBERG_ITER - 1];
201  }
202 
203 
209  template<typename RealFunction>
211  RealFunction f, const std::vector<real>& x, const std::vector<real>& w) {
212 
213  if(x.size() != w.size()) {
214  TH_MATH_ERROR("integral_gauss", x.size(), INVALID_ARGUMENT);
215  return nan();
216  }
217 
218  real res = 0;
219 
220  for (int i = 0; i < x.size(); ++i)
221  res += w[i] * f(x[i]);
222 
223  return res;
224  }
225 
226 
233  template<typename RealFunction>
235  RealFunction f, real* x, real* w, unsigned int n) {
236 
237  real res = 0;
238 
239  for (unsigned int i = 0; i < n; ++i)
240  res += w[i] * f(x[i]);
241 
242  return res;
243  }
244 
245 
253  template<typename RealFunction>
255  RealFunction f, real* x, real* w, unsigned int n,
256  real_function Winv) {
257 
258  real res = 0;
259 
260  for (unsigned int i = 0; i < n; ++i)
261  res += w[i] * f(x[i]) * Winv(x[i]);
262 
263  return res;
264  }
265 
266 
276  template<typename RealFunction>
278  RealFunction f, real a, real b, real* x, real* w, unsigned int n) {
279 
280  const real mean = (b + a) / 2.0;
281  const real halfdiff = (b - a) / 2.0;
282 
283  real res = 0;
284 
285  for (int i = n - 1; i >= 0; --i)
286  res += w[i] * f(halfdiff * x[i] + mean);
287 
288  return res * halfdiff;
289  }
290 
291 
301  template<typename RealFunction>
303  RealFunction f, real a, real b,
304  const std::vector<real>& x, const std::vector<real>& w) {
305 
306  const real mean = (b + a) / 2.0;
307  const real halfdiff = (b - a) / 2.0;
308 
309  real res = 0;
310 
311  for (int i = x.size() - 1; i >= 0; --i)
312  res += w[i] * f(halfdiff * x[i] + mean);
313 
314  return res * halfdiff;
315  }
316 
317 
326  template<typename RealFunction>
328  RealFunction f, real a, real b, const std::vector<real>& x) {
329 
330  return integral_legendre(f, a, b, x, legendre_weights(x));
331  }
332 
333 
344  template<typename RealFunction>
345  inline real integral_legendre(RealFunction f, real a, real b, unsigned int n = 16) {
346 
347  switch(n) {
348  case 2: return integral_legendre(f, a, b,
349  tables::legendre_roots_2, tables::legendre_weights_2, 2); break;
350  case 4: return integral_legendre(f, a, b,
351  tables::legendre_roots_4, tables::legendre_weights_4, 4); break;
352  case 8: return integral_legendre(f, a, b,
353  tables::legendre_roots_8, tables::legendre_weights_8, 8); break;
354  case 16: return integral_legendre(f, a, b,
355  tables::legendre_roots_16, tables::legendre_weights_16, 16); break;
356  // case 32: return integral_legendre(f, a, b,
357  // tables::legendre_roots_32, tables::legendre_weights_32, 32); break;
358  default: return integral_legendre(f, a, b, legendre_roots(n)); break;
359  }
360  }
361 
362 
369  template<typename RealFunction>
370  inline real integral_laguerre(RealFunction f, const std::vector<real>& x) {
371 
372  return integral_gauss(f, x, laguerre_weights(x));
373  }
374 
375 
384  template<typename RealFunction>
386  RealFunction f, real a, real b, const std::vector<real>& x) {
387 
388  const std::vector<real> weights = laguerre_weights(x);
389 
390  const real exp_a = exp(-a);
391  const real exp_b = exp(-b);
392 
393  real res = 0;
394 
395  for (int i = x.size() - 1; i >= 0; --i)
396  res += weights[i] * (exp_a * f(x[i] + a) - exp_b * f(x[i] + b));
397 
398  return res;
399  }
400 
401 
410  template<typename RealFunction>
411  inline real integral_laguerre(RealFunction f, unsigned int n = 16) {
412 
413  switch(n) {
414  case 2: return integral_gauss(f,
415  tables::laguerre_roots_2, tables::laguerre_weights_2, 2); break;
416  case 4: return integral_gauss(f,
417  tables::laguerre_roots_4, tables::laguerre_weights_4, 4); break;
418  case 8: return integral_gauss(f,
419  tables::laguerre_roots_8, tables::laguerre_weights_8, 8); break;
420  case 16: return integral_gauss(f,
421  tables::laguerre_roots_16, tables::laguerre_weights_16, 16); break;
422  // case 32: return integral_gauss(f,
423  // tables::laguerre_roots_32, tables::laguerre_weights_32, 32); break;
424  default: {
425  TH_MATH_ERROR("integral_laguerre", n, INVALID_ARGUMENT);
426  return nan(); break;
427  }
428  }
429  }
430 
431 
438  template<typename RealFunction>
439  inline real integral_hermite(RealFunction f, const std::vector<real>& x) {
440 
441  return integral_gauss(f, x, hermite_weights(x));
442  }
443 
444 
453  template<typename RealFunction>
454  inline real integral_hermite(RealFunction f, unsigned int n = 16) {
455 
456  switch(n) {
457  case 2: return integral_gauss(f,
458  tables::hermite_roots_2, tables::hermite_weights_2, 2); break;
459  case 4: return integral_gauss(f,
460  tables::hermite_roots_4, tables::hermite_weights_4, 4); break;
461  case 8: return integral_gauss(f,
462  tables::hermite_roots_8, tables::hermite_weights_8, 8); break;
463  case 16: return integral_gauss(f,
464  tables::hermite_roots_16, tables::hermite_weights_16, 16); break;
465  default: {
466  TH_MATH_ERROR("integral_hermite", n, INVALID_ARGUMENT);
467  return nan(); break;
468  }
469  }
470  }
471 
472 
478  real_function f, real a, real step_sz = 1,
479  real tol = CALCULUS_INTEGRAL_TOL, unsigned int max_iter = 100) {
480 
481  // Current lower extreme of the interval
482  real x_n = a + step_sz;
483 
484  // Total integral sum
485  real sum = integral_romberg_tol(f, a, x_n, tol);
486 
487  // Variation between steps
488  real delta = inf();
489 
490  // Number of steps performed
491  unsigned int i = 0;
492 
493  while(abs(delta) > tol && i <= max_iter) {
494 
495  delta = integral_romberg_tol(f, x_n, x_n + step_sz, tol);
496  sum += delta;
497  x_n += step_sz;
498  i++;
499  }
500 
501  if(i >= max_iter) {
502  TH_MATH_ERROR("integral_inf_riemann", i, NO_ALGO_CONVERGENCE);
503  return nan();
504  }
505 
506  return sum;
507  }
508 
509 
517  template<typename RealFunction>
518  inline real integral(RealFunction f, real a, real b) {
519  return integral_romberg_tol(f, a, b);
520  }
521 
522 
532  template<typename RealFunction>
533  inline real integral(RealFunction f, real a, real b, real tol) {
534  return integral_romberg_tol(f, a, b, tol);
535  }
536 
537 }
538 
539 
540 #endif
A polynomial of arbitrary order.
Definition: polynomial.h:25
#define TH_MATH_ERROR(F_NAME, VALUE, EXCEPTION)
TH_MATH_ERROR is a macro which throws exceptions or modifies errno (depending on which compiling opti...
Definition: error.h:219
Roots and weights for Gaussian quadrature.
real mean(const histogram &h)
Compute the mean of the values of a histogram.
Definition: histogram.h:296
Main namespace of the library which contains all functions and objects.
Definition: algebra.h:27
double real
A real number, defined as a floating point type.
Definition: constants.h:188
real inf()
Return positive infinity in floating point representation.
Definition: error.h:76
dual2 abs(dual2 x)
Compute the absolute value of a second order dual number.
Definition: dual2_functions.h:183
std::vector< real > hermite_weights(const std::vector< real > &roots)
Hermite weights for Gauss-Hermite quadrature of n-th order.
Definition: orthogonal.h:316
dual2 exp(dual2 x)
Compute the exponential of a second order dual number.
Definition: dual2_functions.h:130
auto sum(const Vector &X)
Compute the sum of a vector of real values using pairwise summation to reduce round-off error.
Definition: dataset.h:219
real integral_inf_riemann(real_function f, real a, real step_sz=1, real tol=CALCULUS_INTEGRAL_TOL, unsigned int max_iter=100)
Integrate a function from a point up to infinity by integrating it by steps, stopping execution when ...
Definition: integral.h:477
std::vector< real > legendre_weights(const std::vector< real > &roots)
Legendre weights for Gauss-Legendre quadrature of n-th order.
Definition: orthogonal.h:274
real integral_romberg(RealFunction f, real a, real b, unsigned int iter=8)
Approximate the definite integral of an arbitrary function using Romberg's method accurate to the giv...
Definition: integral.h:138
constexpr int CALCULUS_INTEGRAL_STEPS
Default number of steps for integral approximation.
Definition: constants.h:272
real integral_gauss(RealFunction f, const std::vector< real > &x, const std::vector< real > &w)
Use Gaussian quadrature using the given points and weights.
Definition: integral.h:210
real integral_laguerre(RealFunction f, const std::vector< real > &x)
Use Gauss-Laguerre quadrature of arbitrary degree to approximate an integral over [0,...
Definition: integral.h:370
real integral_trapezoid(RealFunction f, real a, real b, unsigned int steps=CALCULUS_INTEGRAL_STEPS)
Approximate the definite integral of an arbitrary function using the trapezoid method.
Definition: integral.h:73
real integral_romberg_tol(RealFunction f, real a, real b, real tolerance=CALCULUS_INTEGRAL_TOL)
Approximate the definite integral of an arbitrary function using Romberg's method to the given tolera...
Definition: integral.h:173
real integral_hermite(RealFunction f, const std::vector< real > &x)
Use Gauss-Hermite quadrature of arbitrary degree to approximate an integral over (-inf,...
Definition: integral.h:439
std::vector< real > legendre_roots(unsigned int n)
Roots of the n-th Legendre polynomial.
Definition: orthogonal.h:249
polynomial< T > integral(const polynomial< T > &p)
Compute the indefinite integral of a polynomial.
Definition: integral.h:24
real integral_simpson(RealFunction f, real a, real b, unsigned int steps=CALCULUS_INTEGRAL_STEPS)
Approximate the definite integral of an arbitrary function using Simpson's method.
Definition: integral.h:102
real integral_legendre(RealFunction f, real a, real b, real *x, real *w, unsigned int n)
Use Gauss-Legendre quadrature of arbitrary degree to approximate a definite integral providing the ro...
Definition: integral.h:277
std::function< real(real)> real_function
Function pointer to a real function of real variable.
Definition: function.h:20
real integral_midpoint(RealFunction f, real a, real b, unsigned int steps=CALCULUS_INTEGRAL_STEPS)
Approximate the definite integral of an arbitrary function using the midpoint method.
Definition: integral.h:46
real nan()
Return a quiet NaN number in floating point representation.
Definition: error.h:54
std::vector< real > laguerre_weights(const std::vector< real > &roots)
Laguerre weights for Gauss-Laguerre quadrature of n-th order.
Definition: orthogonal.h:295