ndelement/
polynomials.rs

1//! Orthonormal polynomials.
2//!
3//! The orthonormal polynomials in ndelement span the degree $k$ natural polynomial (or rationomial
4//! for a pyramid) space on a cell. As given at <https://defelement.org/ciarlet.html#The+degree+of+a+finite+element>,
5//! these natural spaces are defined for an interval, triangle, quadrilateral, tetrahedron,
6//! hexahedron, triangular prism and square-based pyramid by
7//!
8//! - $\mathbb{P}^{\text{interval}}_k=\operatorname{span}\left\{x^{p_0}\,\middle|\,p_0\in\mathbb{N},\,p_0\leqslant k\right\},$
9//! - $\mathbb{P}^{\text{triangle}}_k=\operatorname{span}\left\{x^{p_0}y^{p_1}\,\middle|\,p_0,p_1\in\mathbb{N}_0,\,p_0+p_1\leqslant k\right\},$
10//! - $\mathbb{P}^{\text{quadrilateral}}_k=\operatorname{span}\left\{x^{p_0}y^{p_1}\,\middle|\,p_0,p_1\in\mathbb{N}_0,\,p_0\leqslant k,\,p_1\leqslant k\right\},$
11//! - $\mathbb{P}^{\text{tetrahedron}}_k=\operatorname{span}\left\{x^{p_0}y^{p_1}z^{p_2}\,\middle|\,p_0,p_1,p_2\in\mathbb{N}_0,\,p_0+p_1+p_2\leqslant k\right\},$
12//! - $\mathbb{P}^{\text{hexahedron}}_k=\operatorname{span}\left\{x^{p_0}y^{p_1}z^{p_2}\,\middle|\,p_0,p_1,p_2\in\mathbb{N}_0,\,p_0\leqslant k,\,p_1\leqslant k,\,p_2\leqslant k\right\},$
13//! - $\mathbb{P}^{\text{prism}}_k=\operatorname{span}\left\{x^{p_0}y^{p_1}z^{p_2}\,\middle|\,p_0,p_1,p_2\in\mathbb{N}_0,\,p_0+p_1\leqslant k,\,p_2\leqslant k\right\},$
14//! - $\mathbb{P}^{\text{pyramid}}_k=\operatorname{span}\left\{\frac{x^{p_0}y^{p_1}z^{p_2}}{(1-z)^{p_0+p_1}}\,\middle|\,p_0,p_1,p_2\in\mathbb{N}_0,\,p_0\leqslant k,\,p_1\leqslant k,\,p_2\leqslant k\right\}.$
15//!
16//! Note that for non-pyramid cells, these coincide with the polynomial space for the degree $k$
17//! Lagrange element on the cell.
18
19mod legendre;
20pub use legendre::{shape as legendre_shape, tabulate as tabulate_legendre_polynomials};
21
22use crate::reference_cell;
23use crate::types::ReferenceCellType;
24
25/// Return the number of polynomials in the natural polynomial set for a given cell type and degree.
26pub fn polynomial_count(cell_type: ReferenceCellType, degree: usize) -> usize {
27    match cell_type {
28        ReferenceCellType::Interval => degree + 1,
29        ReferenceCellType::Triangle => (degree + 1) * (degree + 2) / 2,
30        ReferenceCellType::Quadrilateral => (degree + 1).pow(2),
31        ReferenceCellType::Tetrahedron => (degree + 1) * (degree + 2) * (degree + 3) / 6,
32        ReferenceCellType::Hexahedron => (degree + 1).pow(3),
33        _ => {
34            panic!("Unsupported cell type: {cell_type:?}");
35        }
36    }
37}
38
39/// Return the total number of partial derivatives up to a given degree.
40pub fn derivative_count(cell_type: ReferenceCellType, degree: usize) -> usize {
41    let mut num = 1;
42    let mut denom = 1;
43    for i in 0..reference_cell::dim(cell_type) {
44        num *= degree + i + 1;
45        denom *= i + 1;
46    }
47    num / denom
48}
49
50#[cfg(test)]
51mod test {
52    use super::*;
53    use crate::quadrature::gauss_jacobi_rule;
54    use crate::types::ReferenceCellType;
55    use approx::*;
56    use paste::paste;
57    use rlst::{DynArray, SliceArray, rlst_dynamic_array};
58
59    macro_rules! test_orthogonal {
60        ($cell:ident, $degree:expr) => {
61            paste! {
62                #[test]
63                fn [<test_orthogonal_ $cell:lower _ $degree>]() {
64                    let rule = gauss_jacobi_rule(
65                        ReferenceCellType::[<$cell>],
66                        2 * [<$degree>],
67                    ).unwrap();
68
69
70                    let points = SliceArray::<f64, 2>::from_shape(&rule.points, [rule.dim, rule.npoints]);
71
72                    let mut data = DynArray::<f64, 3>::from_shape(
73                        legendre_shape(ReferenceCellType::[<$cell>], &points, [<$degree>], 0,)
74                    );
75                    tabulate_legendre_polynomials(ReferenceCellType::[<$cell>], &points, [<$degree>], 0, &mut data);
76
77                    for i in 0..[<$degree>] + 1 {
78                        for j in 0..[<$degree>] + 1 {
79                            let mut product = 0.0;
80                            for k in 0..rule.npoints {
81                                product += data.get([0, i, k]).unwrap()
82                                    * data.get([0, j, k]).unwrap()
83                                    * rule.weights[k];
84                            }
85                            if i == j {
86                                assert_relative_eq!(product, 1.0, epsilon = 1e-12);
87                            } else {
88                                assert_relative_eq!(product, 0.0, epsilon = 1e-12);
89                            }
90                        }
91                    }
92                }
93            }
94        };
95    }
96
97    test_orthogonal!(Interval, 2);
98    test_orthogonal!(Interval, 3);
99    test_orthogonal!(Interval, 4);
100    test_orthogonal!(Interval, 5);
101    test_orthogonal!(Interval, 6);
102    test_orthogonal!(Triangle, 2);
103    test_orthogonal!(Triangle, 3);
104    test_orthogonal!(Triangle, 4);
105    test_orthogonal!(Triangle, 5);
106    test_orthogonal!(Triangle, 6);
107    test_orthogonal!(Quadrilateral, 2);
108    test_orthogonal!(Quadrilateral, 3);
109    test_orthogonal!(Quadrilateral, 4);
110    test_orthogonal!(Quadrilateral, 5);
111    test_orthogonal!(Quadrilateral, 6);
112    test_orthogonal!(Tetrahedron, 2);
113    test_orthogonal!(Tetrahedron, 3);
114    test_orthogonal!(Tetrahedron, 4);
115    test_orthogonal!(Tetrahedron, 5);
116    test_orthogonal!(Tetrahedron, 6);
117    test_orthogonal!(Hexahedron, 2);
118    test_orthogonal!(Hexahedron, 3);
119    test_orthogonal!(Hexahedron, 4);
120    test_orthogonal!(Hexahedron, 5);
121    test_orthogonal!(Hexahedron, 6);
122
123    fn generate_points(cell: ReferenceCellType, epsilon: f64) -> DynArray<f64, 2> {
124        let mut points = match cell {
125            ReferenceCellType::Interval => {
126                let mut points = rlst_dynamic_array!(f64, [1, 20]);
127                for i in 0..10 {
128                    *points.get_mut([0, 2 * i]).unwrap() = i as f64 / 10.0;
129                }
130                points
131            }
132            ReferenceCellType::Triangle => {
133                let mut points = rlst_dynamic_array!(f64, [2, 165]);
134                let mut index = 0;
135                for i in 0..10 {
136                    for j in 0..10 - i {
137                        *points.get_mut([0, 3 * index]).unwrap() = i as f64 / 10.0;
138                        *points.get_mut([1, 3 * index]).unwrap() = j as f64 / 10.0;
139                        index += 1;
140                    }
141                }
142                points
143            }
144            ReferenceCellType::Quadrilateral => {
145                let mut points = rlst_dynamic_array!(f64, [2, 300]);
146                for i in 0..10 {
147                    for j in 0..10 {
148                        let index = 10 * i + j;
149                        *points.get_mut([0, 3 * index]).unwrap() = i as f64 / 10.0;
150                        *points.get_mut([1, 3 * index]).unwrap() = j as f64 / 10.0;
151                    }
152                }
153                points
154            }
155            ReferenceCellType::Tetrahedron => {
156                let mut points = rlst_dynamic_array!(f64, [3, 140]);
157                let mut index = 0;
158                for i in 0..5 {
159                    for j in 0..5 - i {
160                        for k in 0..5 - i - j {
161                            *points.get_mut([0, 4 * index]).unwrap() = i as f64 / 5.0;
162                            *points.get_mut([1, 4 * index]).unwrap() = j as f64 / 5.0;
163                            *points.get_mut([2, 4 * index]).unwrap() = k as f64 / 5.0;
164                            index += 1;
165                        }
166                    }
167                }
168                points
169            }
170            ReferenceCellType::Hexahedron => {
171                let mut points = rlst_dynamic_array!(f64, [3, 500]);
172                for i in 0..5 {
173                    for j in 0..5 {
174                        for k in 0..5 {
175                            let index = 25 * i + 5 * j + k;
176                            *points.get_mut([0, 4 * index]).unwrap() = i as f64 / 5.0;
177                            *points.get_mut([1, 4 * index]).unwrap() = j as f64 / 5.0;
178                            *points.get_mut([2, 4 * index]).unwrap() = k as f64 / 5.0;
179                        }
180                    }
181                }
182                points
183            }
184            _ => {
185                panic!("Unsupported cell type: {cell:?}");
186            }
187        };
188        let dim = reference_cell::dim(cell);
189        for index in 0..points.shape()[1] / (dim + 1) {
190            for d in 0..dim {
191                for i in 0..dim {
192                    *points.get_mut([i, (dim + 1) * index + d + 1]).unwrap() =
193                        *points.get([i, (dim + 1) * index]).unwrap()
194                            + if i == d { epsilon } else { 0.0 };
195                }
196            }
197        }
198        points
199    }
200    macro_rules! test_derivatives {
201        ($cell:ident, $degree:expr) => {
202            paste! {
203                #[test]
204                fn [<test_legendre_derivatives_ $cell:lower _ $degree>]() {
205                    let dim = reference_cell::dim(ReferenceCellType::[<$cell>]);
206                    let epsilon = 1e-10;
207                    let points = generate_points(ReferenceCellType::[<$cell>], epsilon);
208
209                    let mut data = DynArray::<f64, 3>::from_shape(
210                        legendre_shape(ReferenceCellType::[<$cell>], &points, [<$degree>], 1,)
211                    );
212                    tabulate_legendre_polynomials(ReferenceCellType::[<$cell>], &points, [<$degree>], 1, &mut data);
213
214                    for i in 0..data.shape()[1] {
215                        for k in 0..points.shape()[1] / (dim + 1) {
216                            for d in 1..dim + 1 {
217                                assert_relative_eq!(
218                                    *data.get([d, i, (dim + 1) * k]).unwrap(),
219                                    (data.get([0, i, (dim + 1) * k + d]).unwrap() - data.get([0, i, (dim + 1) * k]).unwrap())
220                                        / epsilon,
221                                    epsilon = 1e-3
222                                );
223                            }
224                        }
225                    }
226                }
227            }
228        };
229    }
230
231    test_derivatives!(Interval, 1);
232    test_derivatives!(Interval, 2);
233    test_derivatives!(Interval, 3);
234    test_derivatives!(Interval, 4);
235    test_derivatives!(Interval, 5);
236    test_derivatives!(Interval, 6);
237    test_derivatives!(Triangle, 1);
238    test_derivatives!(Triangle, 2);
239    test_derivatives!(Triangle, 3);
240    test_derivatives!(Triangle, 4);
241    test_derivatives!(Triangle, 5);
242    test_derivatives!(Triangle, 6);
243    test_derivatives!(Quadrilateral, 1);
244    test_derivatives!(Quadrilateral, 2);
245    test_derivatives!(Quadrilateral, 3);
246    test_derivatives!(Quadrilateral, 4);
247    test_derivatives!(Quadrilateral, 5);
248    test_derivatives!(Quadrilateral, 6);
249    test_derivatives!(Tetrahedron, 1);
250    test_derivatives!(Tetrahedron, 2);
251    test_derivatives!(Tetrahedron, 3);
252    test_derivatives!(Tetrahedron, 4);
253    test_derivatives!(Tetrahedron, 5);
254    test_derivatives!(Tetrahedron, 6);
255    test_derivatives!(Hexahedron, 1);
256    test_derivatives!(Hexahedron, 2);
257    test_derivatives!(Hexahedron, 3);
258    test_derivatives!(Hexahedron, 4);
259    test_derivatives!(Hexahedron, 5);
260    test_derivatives!(Hexahedron, 6);
261
262    #[test]
263    fn test_legendre_interval_against_known_polynomials() {
264        let degree = 3;
265
266        let mut points = rlst_dynamic_array!(f64, [1, 11]);
267        for i in 0..11 {
268            *points.get_mut([0, i]).unwrap() = i as f64 / 10.0;
269        }
270
271        let mut data = DynArray::<f64, 3>::from_shape(legendre_shape(
272            ReferenceCellType::Interval,
273            &points,
274            degree,
275            3,
276        ));
277        tabulate_legendre_polynomials(ReferenceCellType::Interval, &points, degree, 3, &mut data);
278
279        for k in 0..points.shape()[0] {
280            let x = *points.get([0, k]).unwrap();
281
282            // 0 => 1
283            assert_relative_eq!(*data.get([0, 0, k]).unwrap(), 1.0, epsilon = 1e-12);
284            assert_relative_eq!(*data.get([1, 0, k]).unwrap(), 0.0, epsilon = 1e-12);
285            assert_relative_eq!(*data.get([2, 0, k]).unwrap(), 0.0, epsilon = 1e-12);
286            assert_relative_eq!(*data.get([3, 0, k]).unwrap(), 0.0, epsilon = 1e-12);
287
288            // 1 => sqrt(3)*(2x - 1)
289            assert_relative_eq!(
290                *data.get([0, 1, k]).unwrap(),
291                f64::sqrt(3.0) * (2.0 * x - 1.0),
292                epsilon = 1e-12
293            );
294            assert_relative_eq!(
295                *data.get([1, 1, k]).unwrap(),
296                2.0 * f64::sqrt(3.0),
297                epsilon = 1e-12
298            );
299            assert_relative_eq!(*data.get([2, 1, k]).unwrap(), 0.0, epsilon = 1e-12);
300            assert_relative_eq!(*data.get([3, 1, k]).unwrap(), 0.0, epsilon = 1e-12);
301
302            // 2 => sqrt(5)*(6x^2 - 6x + 1)
303            assert_relative_eq!(
304                *data.get([0, 2, k]).unwrap(),
305                f64::sqrt(5.0) * (6.0 * x * x - 6.0 * x + 1.0),
306                epsilon = 1e-12
307            );
308            assert_relative_eq!(
309                *data.get([1, 2, k]).unwrap(),
310                f64::sqrt(5.0) * (12.0 * x - 6.0),
311                epsilon = 1e-12
312            );
313            assert_relative_eq!(
314                *data.get([2, 2, k]).unwrap(),
315                f64::sqrt(5.0) * 12.0,
316                epsilon = 1e-12
317            );
318            assert_relative_eq!(*data.get([3, 2, k]).unwrap(), 0.0, epsilon = 1e-12);
319
320            // 3 => sqrt(7)*(20x^3 - 30x^2 + 12x - 1)
321            assert_relative_eq!(
322                *data.get([0, 3, k]).unwrap(),
323                f64::sqrt(7.0) * (20.0 * x * x * x - 30.0 * x * x + 12.0 * x - 1.0),
324                epsilon = 1e-12
325            );
326            assert_relative_eq!(
327                *data.get([1, 3, k]).unwrap(),
328                f64::sqrt(7.0) * (60.0 * x * x - 60.0 * x + 12.0),
329                epsilon = 1e-12
330            );
331            assert_relative_eq!(
332                *data.get([2, 3, k]).unwrap(),
333                f64::sqrt(7.0) * (120.0 * x - 60.0),
334                epsilon = 1e-12
335            );
336            assert_relative_eq!(
337                *data.get([3, 3, k]).unwrap(),
338                f64::sqrt(7.0) * 120.0,
339                epsilon = 1e-12
340            );
341        }
342    }
343
344    #[test]
345    fn test_legendre_quadrilateral_against_known_polynomials() {
346        let degree = 2;
347
348        let mut points = rlst_dynamic_array!(f64, [2, 121]);
349        for i in 0..11 {
350            for j in 0..11 {
351                *points.get_mut([0, 11 * i + j]).unwrap() = i as f64 / 10.0;
352                *points.get_mut([1, 11 * i + j]).unwrap() = j as f64 / 10.0;
353            }
354        }
355
356        let mut data = DynArray::<f64, 3>::from_shape(legendre_shape(
357            ReferenceCellType::Quadrilateral,
358            &points,
359            degree,
360            1,
361        ));
362        tabulate_legendre_polynomials(
363            ReferenceCellType::Quadrilateral,
364            &points,
365            degree,
366            1,
367            &mut data,
368        );
369
370        for k in 0..points.shape()[1] {
371            let x = *points.get([0, k]).unwrap();
372            let y = *points.get([1, k]).unwrap();
373
374            // 0 => 1
375            assert_relative_eq!(*data.get([0, 0, k]).unwrap(), 1.0, epsilon = 1e-12);
376            assert_relative_eq!(*data.get([1, 0, k]).unwrap(), 0.0, epsilon = 1e-12);
377            assert_relative_eq!(*data.get([2, 0, k]).unwrap(), 0.0, epsilon = 1e-12);
378
379            // 3 => sqrt(3)*(2x - 1)
380            assert_relative_eq!(
381                *data.get([0, 3, k]).unwrap(),
382                f64::sqrt(3.0) * (2.0 * x - 1.0),
383                epsilon = 1e-12
384            );
385            assert_relative_eq!(
386                *data.get([1, 3, k]).unwrap(),
387                2.0 * f64::sqrt(3.0),
388                epsilon = 1e-12
389            );
390            assert_relative_eq!(*data.get([2, 3, k]).unwrap(), 0.0, epsilon = 1e-12);
391
392            // 6 => sqrt(5)*(6x^2 - 6x + 1)
393            assert_relative_eq!(
394                *data.get([0, 6, k]).unwrap(),
395                f64::sqrt(5.0) * (6.0 * x * x - 6.0 * x + 1.0),
396                epsilon = 1e-12
397            );
398            assert_relative_eq!(
399                *data.get([1, 6, k]).unwrap(),
400                f64::sqrt(5.0) * (12.0 * x - 6.0),
401                epsilon = 1e-12
402            );
403            assert_relative_eq!(*data.get([2, 6, k]).unwrap(), 0.0, epsilon = 1e-12);
404
405            // 1 => sqrt(3)*(2y - 1)
406            assert_relative_eq!(
407                *data.get([0, 1, k]).unwrap(),
408                f64::sqrt(3.0) * (2.0 * y - 1.0),
409                epsilon = 1e-12
410            );
411
412            assert_relative_eq!(*data.get([1, 1, k]).unwrap(), 0.0, epsilon = 1e-12);
413            assert_relative_eq!(
414                *data.get([2, 1, k]).unwrap(),
415                2.0 * f64::sqrt(3.0),
416                epsilon = 1e-12
417            );
418
419            // 4 => 3*(2x - 1)*(2y - 1)
420            assert_relative_eq!(
421                *data.get([0, 4, k]).unwrap(),
422                3.0 * (2.0 * x - 1.0) * (2.0 * y - 1.0),
423                epsilon = 1e-12
424            );
425            assert_relative_eq!(
426                *data.get([1, 4, k]).unwrap(),
427                6.0 * (2.0 * y - 1.0),
428                epsilon = 1e-12
429            );
430            assert_relative_eq!(
431                *data.get([2, 4, k]).unwrap(),
432                6.0 * (2.0 * x - 1.0),
433                epsilon = 1e-12
434            );
435
436            // 7 => sqrt(15)*(6x^2 - 6x + 1)*(2y - 1)
437            assert_relative_eq!(
438                *data.get([0, 7, k]).unwrap(),
439                f64::sqrt(15.0) * (6.0 * x * x - 6.0 * x + 1.0) * (2.0 * y - 1.0),
440                epsilon = 1e-12
441            );
442            assert_relative_eq!(
443                *data.get([1, 7, k]).unwrap(),
444                f64::sqrt(15.0) * (12.0 * x - 6.0) * (2.0 * y - 1.0),
445                epsilon = 1e-12
446            );
447            assert_relative_eq!(
448                *data.get([2, 7, k]).unwrap(),
449                2.0 * f64::sqrt(15.0) * (6.0 * x * x - 6.0 * x + 1.0),
450                epsilon = 1e-12
451            );
452
453            // 2 => sqrt(5)*(6y^2 - 6y + 1)
454            assert_relative_eq!(
455                *data.get([0, 2, k]).unwrap(),
456                f64::sqrt(5.0) * (6.0 * y * y - 6.0 * y + 1.0),
457                epsilon = 1e-12
458            );
459            assert_relative_eq!(*data.get([1, 2, k]).unwrap(), 0.0, epsilon = 1e-12);
460            assert_relative_eq!(
461                *data.get([2, 2, k]).unwrap(),
462                f64::sqrt(5.0) * (12.0 * y - 6.0),
463                epsilon = 1e-12
464            );
465
466            // 5 => sqrt(15)*(2x - 1)*(6y^2 - 6y + 1)
467            assert_relative_eq!(
468                *data.get([0, 5, k]).unwrap(),
469                f64::sqrt(15.0) * (2.0 * x - 1.0) * (6.0 * y * y - 6.0 * y + 1.0),
470                epsilon = 1e-12
471            );
472            assert_relative_eq!(
473                *data.get([1, 5, k]).unwrap(),
474                2.0 * f64::sqrt(15.0) * (6.0 * y * y - 6.0 * y + 1.0),
475                epsilon = 1e-12
476            );
477            assert_relative_eq!(
478                *data.get([2, 5, k]).unwrap(),
479                f64::sqrt(15.0) * (2.0 * x - 1.0) * (12.0 * y - 6.0),
480                epsilon = 1e-12
481            );
482
483            // 8 => 5*(6x^2 - 6x + 1)*(6y^2 - 6y + 1)
484            assert_relative_eq!(
485                *data.get([0, 8, k]).unwrap(),
486                5.0 * (6.0 * x * x - 6.0 * x + 1.0) * (6.0 * y * y - 6.0 * y + 1.0),
487                epsilon = 1e-12
488            );
489            assert_relative_eq!(
490                *data.get([1, 8, k]).unwrap(),
491                5.0 * (12.0 * x - 6.0) * (6.0 * y * y - 6.0 * y + 1.0),
492                epsilon = 1e-12
493            );
494            assert_relative_eq!(
495                *data.get([2, 8, k]).unwrap(),
496                5.0 * (12.0 * y - 6.0) * (6.0 * x * x - 6.0 * x + 1.0),
497                epsilon = 1e-12
498            );
499        }
500    }
501}