Theoretica
A C++ numerical and automatic mathematical library
orthogonal.h
Go to the documentation of this file.
1 
5 
6 #ifndef THEORETICA_POLYN_ORTHOGONAL_H
7 #define THEORETICA_POLYN_ORTHOGONAL_H
8 
9 #include "./polynomial.h"
10 #include "../optimization/roots.h"
11 #include <functional>
12 
13 
14 namespace theoretica {
15 
16 
20  = std::function<polynomial<real>(polynomial<real>, polynomial<real>, unsigned int)>;
21 
22 
33  unsigned int n) {
34 
35  if(n == 0)
36  return P0;
37 
38  if(n == 1)
39  return P1;
40 
42  for (unsigned int l = 2; l <= n; ++l) {
43 
44  Pl = f(P0, P1, l);
45 
46  P0 = P1;
47  P1 = Pl;
48  }
49 
50  return Pl;
51  }
52 
53 
54  // Legendre polynomials
55 
56 
59  polynomial<real> P0, polynomial<real> P1, unsigned int l) {
60 
61  return ((2 * l - 1) * P1 * polynomial<real>({0, 1}) - (l - 1) * P0) / l;
62  }
63 
64 
67  inline polynomial<real> legendre_polynomial(unsigned int n) {
68 
69  // P0 = 1
70  // P1 = x
71 
72  return gen_polyn_recurr({1}, {0, 1}, legendre_polyn_recurr, n);
73  }
74 
75 
77  inline real legendre_polyn_normalization(unsigned int n) {
78 
79  return sqrt((2 * n + 1) / 2.0);
80  }
81 
82 
83  // Associated Legendre polynomial
84  inline std::function<real(real)> assoc_legendre_polynomial(unsigned int l, int m) {
85 
86  real K = pow(-1, m);
87  polynomial<real> L = legendre_polynomial(l);
88 
89 
90  if(m < 0) {
91  m = -m;
92  K *= pow(-1, m) * fact(l + m) / fact(l - m);
93  }
94 
95  for (int i = 0; i < m; ++i)
96  L = deriv(L);
97 
98 
99  if(m % 2 == 0) {
100 
101  return [=](real x) {
102  return pow(1 - x * x, m / 2) * L(x) / K;
103  };
104 
105  } else {
106 
107  return [=](real x) {
108  return sqrt(pow(1 - x * x, m)) * L(x) / K;
109  };
110  }
111  }
112 
113 
114  // Associated Legendre polynomial
115  inline polynomial<real> assoc_legendre_polynomial_even(unsigned int l, int m) {
116 
117  if(m % 2 != 0)
118  TH_MATH_ERROR("assoc_legendre_polynomial_even", m, IMPOSSIBLE_OPERATION);
119 
120  real K = 1;
121 
122  if(m < 0) {
123  m = -m;
124  K *= fact(l + m) / fact(l - m);
125  }
126 
127  polynomial<real> L = legendre_polynomial(l);
128  polynomial<real> P = ipow(polynomial<real>({1, 0, -1}), m / 2);
129 
130  for (int i = 0; i < m; ++i)
131  L = deriv(L);
132 
133  return L * P / K;
134  }
135 
136 
137  // Laguerre polynomials
138 
139 
142  polynomial<real> L0, polynomial<real> L1, unsigned int i) {
143 
144  return (polynomial<real>({2 * (real) i - 1, -1}) * L1 - (i - 1) * L0) / i;
145  }
146 
147 
149  inline polynomial<real> laguerre_polynomial(unsigned int n) {
150 
151  // L0 = 1
152  // L1 = 1 - x
153 
154  return gen_polyn_recurr({1}, {1, -1}, laguerre_polyn_recurr, n);
155  }
156 
157 
158  // Generalized Laguerre polynomials
159 
160 
163  polynomial<real> L0, polynomial<real> L1, real alpha, unsigned int i) {
164 
165  return (polynomial<real>({2 * (real) i + alpha - 1, -1}) * L1 - (i + alpha - 1) * L0) / i;
166  }
167 
168 
170  inline polynomial<real> general_laguerre_polynomial(real alpha, unsigned int n) {
171 
172  // L0 = 1
173  // L1 = 1 + alpha - x
174 
175  return gen_polyn_recurr(
176  {1}, {1 + alpha, -1},
177  [alpha](polynomial<real> L0, polynomial<real> L1, unsigned int i) {
178  return general_laguerre_polyn_recurr(L0, L1, alpha, i);
179  }, n);
180  }
181 
182 
183  // Hermite polynomials
184 
185 
188  polynomial<real> H0, polynomial<real> H1, unsigned int i) {
189 
190  return polynomial<real>({0, 2}) * H1 - 2 * (i - 1) * H0;
191  }
192 
193 
196  inline polynomial<real> hermite_polynomial(unsigned int n) {
197 
198  // H0 = 1
199  // H1 = 2x
200 
201  return gen_polyn_recurr({1}, {0, 2}, hermite_polyn_recurr, n);
202  }
203 
204 
206  inline real hermite_polyn_normalization(unsigned int n) {
207 
208  return 1.0 / sqrt((1 << n) * fact(n) * SQRTPI);
209  }
210 
211 
212  // Chebyshev polynomials
213 
214 
218  polynomial<real> T0, polynomial<real> T1, unsigned int i) {
219 
220  return polynomial<real>({0, 2}) * T1 - T0;
221  }
222 
223 
226  inline polynomial<real> chebyshev1_polynomial(unsigned int n) {
227 
228  // T0 = 1
229  // T1 = x
230 
231  return gen_polyn_recurr({1}, {0, 1}, chebyshev_polyn_recurr, n);
232  }
233 
234 
237  inline polynomial<real> chebyshev2_polynomial(unsigned int n) {
238 
239  // U0 = 1
240  // U1 = 2x
241 
242  return gen_polyn_recurr({1}, {0, 2}, chebyshev_polyn_recurr, n);
243  }
244 
245 
249  inline std::vector<real> legendre_roots(unsigned int n) {
250 
251  if(n == 0)
252  return {};
253 
254  if(n == 1)
255  return {0};
256 
258  std::vector<real> roots;
259  roots.reserve(n);
260 
261  for (unsigned int i = 1; i <= n; ++i) {
262  roots.push_back(
263  root_newton(P, (2.0 / (n + 1)) * i - 1.0)
264  );
265  }
266 
267  return roots;
268  }
269 
270 
274  inline std::vector<real> legendre_weights(const std::vector<real>& roots) {
275 
276  const unsigned int n = roots.size();
278 
279  std::vector<real> weights;
280  weights.reserve(n);
281 
282  for (unsigned int i = 0; i < n; ++i) {
283  weights.push_back(
284  2.0 / ((1 - square(roots[i])) * square(dP(roots[i])))
285  );
286  }
287 
288  return weights;
289  }
290 
291 
295  inline std::vector<real> laguerre_weights(const std::vector<real>& roots) {
296 
297  const unsigned int n = roots.size();
298  const polynomial<real> L = laguerre_polynomial(n + 1);
299 
300  std::vector<real> weights;
301  weights.reserve(n);
302 
303  for (unsigned int i = 0; i < n; ++i) {
304  weights.push_back(
305  roots[i] / square((n + 1) * L(roots[i]))
306  );
307  }
308 
309  return weights;
310  }
311 
312 
316  inline std::vector<real> hermite_weights(const std::vector<real>& roots) {
317 
318  const unsigned int n = roots.size();
319  const polynomial<real> H = hermite_polynomial(n - 1);
320 
321  std::vector<real> weights;
322  weights.reserve(n);
323 
324  for (unsigned int i = 0; i < n; ++i) {
325  weights.push_back(
326  ((1 << (n - 1)) * fact(n) * SQRTPI) / square(n * H(roots[i]))
327  );
328  }
329 
330  return weights;
331  }
332 
333 }
334 
335 
336 #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
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
polynomial< real > legendre_polyn_recurr(polynomial< real > P0, polynomial< real > P1, unsigned int l)
Recursion formula for Legendre polynomials.
Definition: orthogonal.h:58
dual2 sqrt(dual2 x)
Compute the square root of a second order dual number.
Definition: dual2_functions.h:54
polynomial< real > chebyshev_polyn_recurr(polynomial< real > T0, polynomial< real > T1, unsigned int i)
Recursion formula for Chebyshev polynomials The formula is the same for first and second kind polynom...
Definition: orthogonal.h:217
std::vector< real > hermite_weights(const std::vector< real > &roots)
Hermite weights for Gauss-Hermite quadrature of n-th order.
Definition: orthogonal.h:316
polynomial< real > legendre_polynomial(unsigned int n)
Compute the nth Legendre polynomial.
Definition: orthogonal.h:67
polynomial< real > laguerre_polyn_recurr(polynomial< real > L0, polynomial< real > L1, unsigned int i)
Recursion formula for Laguerre polynomials.
Definition: orthogonal.h:141
constexpr real SQRTPI
The square root of Pi.
Definition: constants.h:224
polynomial< real > hermite_polynomial(unsigned int n)
Compute the nth Hermite polynomial.
Definition: orthogonal.h:196
TH_CONSTEXPR T ipow(T x, unsigned int n, T neutral_element=T(1))
Compute the n-th positive power of x (where n is natural)
Definition: real_analysis.h:644
polynomial< real > general_laguerre_polynomial(real alpha, unsigned int n)
Compute the nth Laguerre polynomial.
Definition: orthogonal.h:170
std::vector< real > legendre_weights(const std::vector< real > &roots)
Legendre weights for Gauss-Legendre quadrature of n-th order.
Definition: orthogonal.h:274
polynomial< real > gen_polyn_recurr(polynomial< real > P0, polynomial< real > P1, polyn_recurr_formula f, unsigned int n)
Generate a polynomial basis using a recursion formula.
Definition: orthogonal.h:29
real hermite_polyn_normalization(unsigned int n)
Normalization constant for the nth Hermite polynomial.
Definition: orthogonal.h:206
TH_CONSTEXPR IntType fact(unsigned int n)
Compute the factorial of n.
Definition: real_analysis.h:671
polynomial< real > general_laguerre_polyn_recurr(polynomial< real > L0, polynomial< real > L1, real alpha, unsigned int i)
Recursion formula for Generalized Laguerre polynomials.
Definition: orthogonal.h:162
std::vector< real > roots(RealFunction f, real a, real b, real tolerance=OPTIMIZATION_TOL, real steps=10)
Find the roots of a function inside a given interval.
Definition: roots.h:462
real root_newton(RealFunction f, RealFunction Df, real guess=0)
Approximate a root of an arbitrary function using Newton's method.
Definition: roots.h:100
polynomial< real > chebyshev2_polynomial(unsigned int n)
Compute the nth Chebyshev polynomial of the second kind.
Definition: orthogonal.h:237
std::vector< real > legendre_roots(unsigned int n)
Roots of the n-th Legendre polynomial.
Definition: orthogonal.h:249
polynomial< real > hermite_polyn_recurr(polynomial< real > H0, polynomial< real > H1, unsigned int i)
Recursion formula for Hermite polynomials.
Definition: orthogonal.h:187
real legendre_polyn_normalization(unsigned int n)
Normalization constant for the nth Legendre polynomial.
Definition: orthogonal.h:77
polynomial< Field > deriv(const polynomial< Field > &p)
Compute the exact derivative of a polynomial function.
Definition: deriv.h:21
polynomial< real > laguerre_polynomial(unsigned int n)
Compute the nth Laguerre polynomial.
Definition: orthogonal.h:149
std::function< polynomial< real >(polynomial< real >, polynomial< real >, unsigned int)> polyn_recurr_formula
Polynomial sequence recurrence formula type Used for computing orthogonal polynomial basis elements.
Definition: orthogonal.h:20
dual2 square(dual2 x)
Return the square of a second order dual number.
Definition: dual2_functions.h:23
std::vector< real > laguerre_weights(const std::vector< real > &roots)
Laguerre weights for Gauss-Laguerre quadrature of n-th order.
Definition: orthogonal.h:295
polynomial< real > chebyshev1_polynomial(unsigned int n)
Compute the nth Chebyshev polynomial of the first kind.
Definition: orthogonal.h:226
dual2 pow(dual2 x, int n)
Compute the n-th power of a second order dual number.
Definition: dual2_functions.h:41
Polynomial storage and manipulation.