6#ifndef THEORETICA_INTEGRATION_H
7#define THEORETICA_INTEGRATION_H
9#include "../core/constants.h"
10#include "../core/function.h"
11#include "../polynomial/polynomial.h"
12#include "../polynomial/orthogonal.h"
23 template<
typename Field = real>
32 if (
P.size() !=
p.size() + 1) {
38 for (
unsigned int i = 0; i <
p.size(); ++i)
58 for (
unsigned int i = 0; i <
p.size(); ++i) {
60 const unsigned int pos =
p.size() - i - 1;
77 template<
typename RealFunction>
89 for (
unsigned int i = 0; i <
steps; ++i)
90 res += f(a + (i + 0.5) *
dx);
104 template<
typename RealFunction>
116 res += 0.5 * (f(a) + f(b));
118 for (
unsigned int i = 1; i <
steps; ++i)
119 res += f(a + i *
dx);
133 template<
typename RealFunction>
142 if (
steps % 2 != 0) {
154 for (
unsigned int i = 2; i <
steps; i += 2)
155 res += 2 * f(a + i *
dx);
157 for (
unsigned int i = 1; i <
steps; i += 2)
158 res += 4 * f(a + i *
dx);
160 return res *
dx / 3.0;
173 template<
typename RealFunction>
184 T[0][0] = (f(a) + f(b)) * (b - a) / 2.0;
191 T[
j][0] =
T[
j - 1][0] / 2.0;
193 for (
unsigned int i = 0; i < (
uint32_t(1) << (
j - 1)); i++) {
194 const real x = a + (2 * i + 1) *
h;
199 for (
unsigned int k = 1;
k <=
j; ++
k) {
201 T[
j][
k] = (coeff *
T[
j][
k - 1] -
T[
j - 1][
k - 1]) / (coeff - 1);
219 template<
typename RealFunction>
221 RealFunction f,
const std::vector<real>& x,
const std::vector<real>&
w) {
223 if(x.size() !=
w.size()) {
230 for (
int i = 0; i < x.size(); ++i)
231 res +=
w[i] * f(x[i]);
243 template<
typename RealFunction>
249 for (
unsigned int i = 0; i <
n; ++i)
250 res +=
w[i] * f(x[i]);
263 template<
typename RealFunction>
270 for (
unsigned int i = 0; i <
n; ++i)
286 template<
typename RealFunction>
290 const real mean = (b + a) / 2.0;
295 for (
int i =
n - 1; i >= 0; --i)
311 template<
typename RealFunction>
314 const std::vector<real>& x,
const std::vector<real>&
w) {
316 const real mean = (b + a) / 2.0;
321 for (
int i = x.size() - 1; i >= 0; --i)
336 template<
typename RealFunction>
356 template<
typename RealFunction>
361 tables::legendre_roots_2, tables::legendre_weights_2, 2);
break;
363 tables::legendre_roots_4, tables::legendre_weights_4, 4);
break;
365 tables::legendre_roots_8, tables::legendre_weights_8, 8);
break;
367 tables::legendre_roots_16, tables::legendre_weights_16, 16);
break;
381 template<
typename RealFunction>
396 template<
typename RealFunction>
407 for (
int i = x.size() - 1; i >= 0; --i)
422 template<
typename RealFunction>
427 tables::laguerre_roots_2, tables::laguerre_weights_2, 2);
break;
429 tables::laguerre_roots_4, tables::laguerre_weights_4, 4);
break;
431 tables::laguerre_roots_8, tables::laguerre_weights_8, 8);
break;
433 tables::laguerre_roots_16, tables::laguerre_weights_16, 16);
break;
450 template<
typename RealFunction>
465 template<
typename RealFunction>
470 tables::hermite_roots_2, tables::hermite_weights_2, 2);
break;
472 tables::hermite_roots_4, tables::hermite_weights_4, 4);
break;
474 tables::hermite_roots_8, tables::hermite_weights_8, 8);
break;
476 tables::hermite_roots_16, tables::hermite_weights_16, 16);
break;
529 template<
typename RealFunction>
544 template<
typename RealFunction>
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