mathhook_core/core/polynomial/
special_families.rs1pub use crate::functions::polynomials::{
33 chebyshev, evaluation, hermite, laguerre, legendre, symbolic, PolynomialIntelligence,
34};
35
36pub use crate::functions::polynomials::evaluation::{
38 evaluate_chebyshev_first_numerical, evaluate_chebyshev_second_numerical,
39 evaluate_hermite_numerical, evaluate_laguerre_numerical, evaluate_legendre_numerical,
40};
41pub use crate::functions::polynomials::symbolic::{
42 expand_chebyshev_first_symbolic, expand_chebyshev_second_symbolic, expand_hermite_symbolic,
43 expand_laguerre_symbolic, expand_legendre_symbolic,
44};
45
46use crate::core::{Expression, Symbol};
47use crate::pattern::Substitutable;
48
49pub trait OrthogonalPolynomial {
54 fn polynomial(n: usize, var: &Symbol) -> Expression;
65
66 fn evaluate(n: usize, x: f64) -> f64;
77
78 fn recurrence_coefficients(n: usize) -> (f64, f64, f64);
86}
87
88pub struct Legendre;
90
91impl OrthogonalPolynomial for Legendre {
92 fn polynomial(n: usize, var: &Symbol) -> Expression {
93 let expr = expand_legendre_symbolic(n);
95 if var.name() == "x" {
96 expr
97 } else {
98 expr.subs(&Expression::symbol("x"), &Expression::symbol(var.clone()))
99 }
100 }
101
102 fn evaluate(n: usize, x: f64) -> f64 {
103 let result = evaluate_legendre_numerical(&[n as f64, x]);
104 result.first().copied().unwrap_or(0.0)
105 }
106
107 fn recurrence_coefficients(n: usize) -> (f64, f64, f64) {
108 let n_f64 = n as f64;
109 let a_n = (2.0 * n_f64 + 1.0) / (n_f64 + 1.0);
110 let b_n = 0.0;
111 let c_n = n_f64 / (n_f64 + 1.0);
112 (a_n, b_n, c_n)
113 }
114}
115
116pub struct ChebyshevT;
118
119impl OrthogonalPolynomial for ChebyshevT {
120 fn polynomial(n: usize, var: &Symbol) -> Expression {
121 let expr = expand_chebyshev_first_symbolic(n);
122 if var.name() == "x" {
123 expr
124 } else {
125 expr.subs(&Expression::symbol("x"), &Expression::symbol(var.clone()))
126 }
127 }
128
129 fn evaluate(n: usize, x: f64) -> f64 {
130 let result = evaluate_chebyshev_first_numerical(&[n as f64, x]);
131 result.first().copied().unwrap_or(0.0)
132 }
133
134 fn recurrence_coefficients(n: usize) -> (f64, f64, f64) {
135 if n == 0 {
136 (1.0, 0.0, 0.0)
137 } else {
138 (2.0, 0.0, 1.0)
139 }
140 }
141}
142
143pub struct ChebyshevU;
145
146impl OrthogonalPolynomial for ChebyshevU {
147 fn polynomial(n: usize, var: &Symbol) -> Expression {
148 let expr = expand_chebyshev_second_symbolic(n);
149 if var.name() == "x" {
150 expr
151 } else {
152 expr.subs(&Expression::symbol("x"), &Expression::symbol(var.clone()))
153 }
154 }
155
156 fn evaluate(n: usize, x: f64) -> f64 {
157 let result = evaluate_chebyshev_second_numerical(&[n as f64, x]);
158 result.first().copied().unwrap_or(0.0)
159 }
160
161 fn recurrence_coefficients(n: usize) -> (f64, f64, f64) {
162 if n == 0 {
163 (2.0, 0.0, 0.0)
164 } else {
165 (2.0, 0.0, 1.0)
166 }
167 }
168}
169
170pub struct Hermite;
172
173impl OrthogonalPolynomial for Hermite {
174 fn polynomial(n: usize, var: &Symbol) -> Expression {
175 let expr = expand_hermite_symbolic(n);
176 if var.name() == "x" {
177 expr
178 } else {
179 expr.subs(&Expression::symbol("x"), &Expression::symbol(var.clone()))
180 }
181 }
182
183 fn evaluate(n: usize, x: f64) -> f64 {
184 let result = evaluate_hermite_numerical(&[n as f64, x]);
185 result.first().copied().unwrap_or(0.0)
186 }
187
188 fn recurrence_coefficients(n: usize) -> (f64, f64, f64) {
189 let n_f64 = n as f64;
190 (2.0, 0.0, 2.0 * n_f64)
191 }
192}
193
194pub struct Laguerre;
196
197impl OrthogonalPolynomial for Laguerre {
198 fn polynomial(n: usize, var: &Symbol) -> Expression {
199 let expr = expand_laguerre_symbolic(n);
200 if var.name() == "x" {
201 expr
202 } else {
203 expr.subs(&Expression::symbol("x"), &Expression::symbol(var.clone()))
204 }
205 }
206
207 fn evaluate(n: usize, x: f64) -> f64 {
208 let result = evaluate_laguerre_numerical(&[n as f64, x]);
209 result.first().copied().unwrap_or(0.0)
210 }
211
212 fn recurrence_coefficients(n: usize) -> (f64, f64, f64) {
213 let n_f64 = n as f64;
214 let a_n = -1.0 / (n_f64 + 1.0);
215 let b_n = (2.0 * n_f64 + 1.0) / (n_f64 + 1.0);
216 let c_n = n_f64 / (n_f64 + 1.0);
217 (a_n, b_n, c_n)
218 }
219}
220
221#[cfg(test)]
222mod tests {
223 use super::*;
224 use crate::symbol;
225
226 #[test]
227 fn test_legendre_base_cases() {
228 let x = symbol!(x);
229 let p0 = Legendre::polynomial(0, &x);
230 let p1 = Legendre::polynomial(1, &x);
231
232 assert_eq!(p0, Expression::integer(1));
233 assert_eq!(p1, Expression::symbol(x));
234 }
235
236 #[test]
237 fn test_legendre_eval() {
238 assert!((Legendre::evaluate(0, 0.5) - 1.0).abs() < 1e-10);
240 assert!((Legendre::evaluate(1, 0.5) - 0.5).abs() < 1e-10);
242 assert!((Legendre::evaluate(2, 0.5) - (-0.125)).abs() < 1e-10);
244 }
245
246 #[test]
247 fn test_chebyshev_t_base_cases() {
248 let x = symbol!(x);
249 let t0 = ChebyshevT::polynomial(0, &x);
250 let t1 = ChebyshevT::polynomial(1, &x);
251
252 assert_eq!(t0, Expression::integer(1));
253 assert_eq!(t1, Expression::symbol(x));
254 }
255
256 #[test]
257 fn test_chebyshev_t_eval() {
258 assert!((ChebyshevT::evaluate(0, 0.5) - 1.0).abs() < 1e-10);
260 assert!((ChebyshevT::evaluate(1, 0.5) - 0.5).abs() < 1e-10);
262 assert!((ChebyshevT::evaluate(2, 0.5) - (-0.5)).abs() < 1e-10);
264 }
265
266 #[test]
267 fn test_chebyshev_u_eval() {
268 assert!((ChebyshevU::evaluate(0, 0.5) - 1.0).abs() < 1e-10);
270 assert!((ChebyshevU::evaluate(1, 0.5) - 1.0).abs() < 1e-10);
272 }
273
274 #[test]
275 fn test_hermite_base_cases() {
276 let x = symbol!(x);
277 let h0 = Hermite::polynomial(0, &x);
278 let _h1 = Hermite::polynomial(1, &x);
279
280 assert_eq!(h0, Expression::integer(1));
281 }
283
284 #[test]
285 fn test_hermite_eval() {
286 assert!((Hermite::evaluate(0, 0.5) - 1.0).abs() < 1e-10);
288 assert!((Hermite::evaluate(1, 0.5) - 1.0).abs() < 1e-10);
290 }
291
292 #[test]
293 fn test_laguerre_base_cases() {
294 let x = symbol!(x);
295 let l0 = Laguerre::polynomial(0, &x);
296
297 assert_eq!(l0, Expression::integer(1));
298 }
299
300 #[test]
301 fn test_laguerre_eval() {
302 assert!((Laguerre::evaluate(0, 0.5) - 1.0).abs() < 1e-10);
304 assert!((Laguerre::evaluate(1, 0.5) - 0.5).abs() < 1e-10);
306 }
307
308 #[test]
309 fn test_recurrence_coefficients() {
310 let (a, b, c) = Legendre::recurrence_coefficients(2);
312 assert!(a > 0.0);
313 assert_eq!(b, 0.0);
314 assert!(c > 0.0);
315
316 let (a, _, c) = ChebyshevT::recurrence_coefficients(2);
317 assert_eq!(a, 2.0);
318 assert_eq!(c, 1.0);
319 }
320
321 #[test]
322 fn test_variable_substitution() {
323 let t = symbol!(t);
324 let p1 = Legendre::polynomial(1, &t);
325 assert_eq!(p1, Expression::symbol(t));
327 }
328}