Theoretica
Scientific Computing
Loading...
Searching...
No Matches
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
16namespace theoretica {
17
18
23 template<typename Field = real>
25
26 if (p.size() == 0)
27 return polynomial<Field>({});
28
30 P.resize(p.size() + 1);
31
32 if (P.size() != p.size() + 1) {
33 TH_MATH_ERROR("integral(polynomial)", P.size(), MathError::ImpossibleOperation);
34 return polynomial<Field>(nan());
35 }
36
37 P[0] = Field(0.0);
38 for (unsigned int i = 0; i < p.size(); ++i)
39 P[i + 1] = p[i] / Field(i + 1);
40
41 return P;
42 }
43
44
53 inline real integral(const polynomial<real>& p, real a, real b) {
54
55 real P_a = 0.0;
56 real P_b = 0.0;
57
58 for (unsigned int i = 0; i < p.size(); ++i) {
59
60 const unsigned int pos = p.size() - i - 1;
61 P_a = a * (p[pos] / (pos + 1) + P_a);
62 P_b = b * (p[pos] / (pos + 1) + P_b);
63 }
64
65 return P_b - P_a;
66 }
67
68
77 template<typename RealFunction>
79 unsigned int steps = CALCULUS_INTEGRAL_STEPS) {
80
81 if(steps == 0) {
82 TH_MATH_ERROR("integral_midpoint", steps, MathError::DivByZero);
83 return nan();
84 }
85
86 const real dx = (b - a) / steps;
87 real res = 0;
88
89 for (unsigned int i = 0; i < steps; ++i)
90 res += f(a + (i + 0.5) * dx);
91
92 return res * dx;
93 }
94
95
104 template<typename RealFunction>
106 unsigned int steps = CALCULUS_INTEGRAL_STEPS) {
107
108 if(steps == 0) {
109 TH_MATH_ERROR("integral_trapezoid", steps, MathError::DivByZero);
110 return nan();
111 }
112
113 const real dx = (b - a) / steps;
114 real res = 0;
115
116 res += 0.5 * (f(a) + f(b));
117
118 for (unsigned int i = 1; i < steps; ++i)
119 res += f(a + i * dx);
120
121 return res * dx;
122 }
123
124
133 template<typename RealFunction>
135 unsigned int steps = CALCULUS_INTEGRAL_STEPS) {
136
137 if (steps == 0) {
138 TH_MATH_ERROR("integral_simpson", steps, MathError::DivByZero);
139 return nan();
140 }
141
142 if (steps % 2 != 0) {
143 TH_MATH_ERROR("integral_simpson", steps, MathError::InvalidArgument);
144 return nan();
145 }
146
147 // Sum terms by order of magnitude supposing that
148 // f stays at the same order inside the interval,
149 // to alleviate truncation errors
150
151 real res = f(a) + f(b);
152 const real dx = (b - a) / (real) steps;
153
154 for (unsigned int i = 2; i < steps; i += 2)
155 res += 2 * f(a + i * dx);
156
157 for (unsigned int i = 1; i < steps; i += 2)
158 res += 4 * f(a + i * dx);
159
160 return res * dx / 3.0;
161 }
162
163
173 template<typename RealFunction>
175 RealFunction f,
176 real a, real b,
177 real tolerance = CALCULUS_INTEGRAL_TOL,
178 size_t max_iter = 16
179 ) {
180
182
183 real h = (b - a);
184 T[0][0] = (f(a) + f(b)) * (b - a) / 2.0;
185
186 for (unsigned int j = 1; j < max_iter; ++j) {
187
188 h /= 2.0;
189
190 // Reuse previous calculations
191 T[j][0] = T[j - 1][0] / 2.0;
192
193 for (unsigned int i = 0; i < (uint32_t(1) << (j - 1)); i++) {
194 const real x = a + (2 * i + 1) * h;
195 T[j][0] += f(x) * h;
196 }
197
198 // Richardson extrapolation
199 for (unsigned int k = 1; k <= j; ++k) {
200 const uint64_t coeff = uint64_t(1) << (2 * k);
201 T[j][k] = (coeff * T[j][k - 1] - T[j - 1][k - 1]) / (coeff - 1);
202 }
203
204 // Stop the algorithm when the desired precision has been reached
205 if(abs(T[j][j] - T[j - 1][j - 1]) < tolerance)
206 return T[j][j];
207 }
208
209 // Return the best approximation
210 return T[max_iter - 1][max_iter - 1];
211 }
212
213
219 template<typename RealFunction>
221 RealFunction f, const std::vector<real>& x, const std::vector<real>& w) {
222
223 if(x.size() != w.size()) {
224 TH_MATH_ERROR("integral_gauss", x.size(), MathError::InvalidArgument);
225 return nan();
226 }
227
228 real res = 0;
229
230 for (int i = 0; i < x.size(); ++i)
231 res += w[i] * f(x[i]);
232
233 return res;
234 }
235
236
243 template<typename RealFunction>
245 RealFunction f, real* x, real* w, unsigned int n) {
246
247 real res = 0;
248
249 for (unsigned int i = 0; i < n; ++i)
250 res += w[i] * f(x[i]);
251
252 return res;
253 }
254
255
263 template<typename RealFunction>
265 RealFunction f, real* x, real* w, unsigned int n,
267
268 real res = 0;
269
270 for (unsigned int i = 0; i < n; ++i)
271 res += w[i] * f(x[i]) * Winv(x[i]);
272
273 return res;
274 }
275
276
286 template<typename RealFunction>
288 RealFunction f, real a, real b, real* x, real* w, unsigned int n) {
289
290 const real mean = (b + a) / 2.0;
291 const real halfdiff = (b - a) / 2.0;
292
293 real res = 0;
294
295 for (int i = n - 1; i >= 0; --i)
296 res += w[i] * f(halfdiff * x[i] + mean);
297
298 return res * halfdiff;
299 }
300
301
311 template<typename RealFunction>
313 RealFunction f, real a, real b,
314 const std::vector<real>& x, const std::vector<real>& w) {
315
316 const real mean = (b + a) / 2.0;
317 const real halfdiff = (b - a) / 2.0;
318
319 real res = 0;
320
321 for (int i = x.size() - 1; i >= 0; --i)
322 res += w[i] * f(halfdiff * x[i] + mean);
323
324 return res * halfdiff;
325 }
326
327
336 template<typename RealFunction>
338 RealFunction f, real a, real b, const std::vector<real>& x) {
339
340 return integral_legendre(f, a, b, x, legendre_weights(x));
341 }
342
343
356 template<typename RealFunction>
357 inline real integral_legendre(RealFunction f, real a, real b, unsigned int n = 16) {
358
359 switch(n) {
360 case 2: return integral_legendre(f, a, b,
361 tables::legendre_roots_2, tables::legendre_weights_2, 2); break;
362 case 4: return integral_legendre(f, a, b,
363 tables::legendre_roots_4, tables::legendre_weights_4, 4); break;
364 case 8: return integral_legendre(f, a, b,
365 tables::legendre_roots_8, tables::legendre_weights_8, 8); break;
366 case 16: return integral_legendre(f, a, b,
367 tables::legendre_roots_16, tables::legendre_weights_16, 16); break;
368 // case 32: return integral_legendre(f, a, b,
369 // tables::legendre_roots_32, tables::legendre_weights_32, 32); break;
370 default: return integral_legendre(f, a, b, legendre_roots(n)); break;
371 }
372 }
373
374
381 template<typename RealFunction>
382 inline real integral_laguerre(RealFunction f, const std::vector<real>& x) {
383
384 return integral_gauss(f, x, laguerre_weights(x));
385 }
386
387
396 template<typename RealFunction>
398 RealFunction f, real a, real b, const std::vector<real>& x) {
399
400 const std::vector<real> weights = laguerre_weights(x);
401
402 const real exp_a = exp(-a);
403 const real exp_b = exp(-b);
404
405 real res = 0;
406
407 for (int i = x.size() - 1; i >= 0; --i)
408 res += weights[i] * (exp_a * f(x[i] + a) - exp_b * f(x[i] + b));
409
410 return res;
411 }
412
413
422 template<typename RealFunction>
423 inline real integral_laguerre(RealFunction f, unsigned int n = 16) {
424
425 switch(n) {
426 case 2: return integral_gauss(f,
427 tables::laguerre_roots_2, tables::laguerre_weights_2, 2); break;
428 case 4: return integral_gauss(f,
429 tables::laguerre_roots_4, tables::laguerre_weights_4, 4); break;
430 case 8: return integral_gauss(f,
431 tables::laguerre_roots_8, tables::laguerre_weights_8, 8); break;
432 case 16: return integral_gauss(f,
433 tables::laguerre_roots_16, tables::laguerre_weights_16, 16); break;
434 // case 32: return integral_gauss(f,
435 // tables::laguerre_roots_32, tables::laguerre_weights_32, 32); break;
436 default: {
437 TH_MATH_ERROR("integral_laguerre", n, MathError::InvalidArgument);
438 return nan(); break;
439 }
440 }
441 }
442
443
450 template<typename RealFunction>
451 inline real integral_hermite(RealFunction f, const std::vector<real>& x) {
452
453 return integral_gauss(f, x, hermite_weights(x));
454 }
455
456
465 template<typename RealFunction>
466 inline real integral_hermite(RealFunction f, unsigned int n = 16) {
467
468 switch(n) {
469 case 2: return integral_gauss(f,
470 tables::hermite_roots_2, tables::hermite_weights_2, 2); break;
471 case 4: return integral_gauss(f,
472 tables::hermite_roots_4, tables::hermite_weights_4, 4); break;
473 case 8: return integral_gauss(f,
474 tables::hermite_roots_8, tables::hermite_weights_8, 8); break;
475 case 16: return integral_gauss(f,
476 tables::hermite_roots_16, tables::hermite_weights_16, 16); break;
477 default: {
478 TH_MATH_ERROR("integral_hermite", n, MathError::InvalidArgument);
479 return nan(); break;
480 }
481 }
482 }
483
484
490 real_function f, real a, real step_sz = 1,
491 real tol = CALCULUS_INTEGRAL_TOL, unsigned int max_iter = 100) {
492
493 // Current lower extreme of the interval
494 real x_n = a + step_sz;
495
496 // Total integral sum
497 real sum = integral_romberg(f, a, x_n, tol);
498
499 // Variation between steps
500 real delta = inf();
501
502 // Number of steps performed
503 unsigned int i = 0;
504
505 while(abs(delta) > tol && i <= max_iter) {
506
508 sum += delta;
509 x_n += step_sz;
510 i++;
511 }
512
513 if(i >= max_iter) {
514 TH_MATH_ERROR("integral_inf_riemann", i, MathError::NoConvergence);
515 return nan();
516 }
517
518 return sum;
519 }
520
521
529 template<typename RealFunction>
531 return integral_romberg(f, a, b);
532 }
533
534
544 template<typename RealFunction>
546 return integral_romberg(f, a, b, tol);
547 }
548
549}
550
551
552#endif
A polynomial of arbitrary order with coefficients of a specified type.
Definition polynomial.h:28
void resize(size_t sz)
Change the number of coefficients of the polynomial.
Definition polynomial.h:168
#define TH_MATH_ERROR(F_NAME, VALUE, EXCEPTION)
TH_MATH_ERROR is a macro which throws exceptions or modifies errno (depending on which compilation op...
Definition error.h:219
Roots and weights for Gaussian quadrature.
Main namespace of the library which contains all functions and objects.
Definition algebra.h:27
std::vector< real > hermite_weights(const std::vector< real > &roots)
Hermite weights for Gauss-Hermite quadrature of n-th order.
Definition orthogonal.h:316
double real
A real number, defined as a floating point type.
Definition constants.h:207
dual2 abs(dual2 x)
Compute the absolute value of a second order dual number.
Definition dual2_functions.h:242
dual2 exp(dual2 x)
Compute the exponential of a second order dual number.
Definition dual2_functions.h:138
std::vector< real > legendre_weights(const std::vector< real > &roots)
Legendre weights for Gauss-Legendre quadrature of n-th order.
Definition orthogonal.h:274
Vector make_error()
Create a vector representing an error state, with all NaN values.
Definition algebra.h:103
std::vector< real > laguerre_weights(const std::vector< real > &roots)
Laguerre weights for Gauss-Laguerre quadrature of n-th order.
Definition orthogonal.h:295
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:215
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:489
TH_CONSTEXPR real nan()
Return a quiet NaN number in floating point representation.
Definition error.h:74
constexpr int CALCULUS_INTEGRAL_STEPS
Default number of steps for integral approximation.
Definition constants.h:291
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:220
@ InvalidArgument
Invalid argument.
@ ImpossibleOperation
Mathematically impossible operation.
@ NoConvergence
Algorithm did not converge.
@ DivByZero
Division by zero.
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:382
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:105
polynomial< Field > integral(const polynomial< Field > &p)
Compute the indefinite integral of a polynomial.
Definition integral.h:24
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:451
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:134
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:287
std::function< real(real)> real_function
Function pointer to a real function of real variable.
Definition function.h:20
real integral_romberg(RealFunction f, real a, real b, real tolerance=CALCULUS_INTEGRAL_TOL, size_t max_iter=16)
Approximate the definite integral of an arbitrary function using Romberg's method to the given tolera...
Definition integral.h:174
TH_CONSTEXPR real inf()
Get positive infinity in floating point representation.
Definition error.h:96
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:78
std::vector< real > legendre_roots(unsigned int n)
Roots of the n-th Legendre polynomial.
Definition orthogonal.h:249