Skip to main content

numra_integrate/
fixed.rs

1//! Fixed-order Gaussian quadrature rules.
2//!
3//! Provides Gauss-Legendre, Gauss-Laguerre, and Gauss-Hermite quadrature
4//! with precomputed nodes and weights.
5//!
6//! Author: Moussa Leblouba
7//! Date: 9 February 2026
8//! Modified: 2 May 2026
9
10// Quadrature tables contain precise mathematical constants that happen to
11// resemble well-known constants like FRAC_1_SQRT_2. These are genuine
12// quadrature node values, not approximations of std constants.
13#![allow(clippy::excessive_precision, clippy::approx_constant)]
14
15use numra_core::Scalar;
16
17use crate::error::IntegrationError;
18
19// ============================================================================
20// Gauss-Legendre nodes and weights for [-1, 1]
21// Stored as (positive_node, weight) pairs; rules are symmetric about 0.
22// For odd n, the zero node and its weight are stored separately.
23// ============================================================================
24
25/// Gauss-Legendre quadrature nodes and weights for n = 1..20.
26/// Each entry is (nodes, weights) where nodes are the positive abscissae
27/// (the rule is symmetric, so negative abscissae have the same weights).
28/// For odd n, the first entry has node = 0.
29///
30/// Source: Abramowitz & Stegun, Table 25.4; verified against DLMF.
31struct GaussLegendreTable {
32    /// For odd n: (weight_at_zero, &[(positive_node, weight)])
33    /// For even n: (0.0, &[(positive_node, weight)])
34    zero_weight: f64,
35    pairs: &'static [(f64, f64)],
36}
37
38// Precomputed tables for n = 1..20
39const GL: [GaussLegendreTable; 20] = [
40    // n=1
41    GaussLegendreTable {
42        zero_weight: 2.0,
43        pairs: &[],
44    },
45    // n=2
46    GaussLegendreTable {
47        zero_weight: 0.0,
48        pairs: &[(0.5773502691896258, 1.0)],
49    },
50    // n=3
51    GaussLegendreTable {
52        zero_weight: 0.8888888888888889,
53        pairs: &[(0.7745966692414834, 0.5555555555555556)],
54    },
55    // n=4
56    GaussLegendreTable {
57        zero_weight: 0.0,
58        pairs: &[
59            (0.3399810435848563, 0.6521451548625461),
60            (0.8611363115940526, 0.3478548451374538),
61        ],
62    },
63    // n=5
64    GaussLegendreTable {
65        zero_weight: 0.5688888888888889,
66        pairs: &[
67            (0.5384693101056831, 0.4786286704993665),
68            (0.9061798459386640, 0.2369268850561891),
69        ],
70    },
71    // n=6
72    GaussLegendreTable {
73        zero_weight: 0.0,
74        pairs: &[
75            (0.2386191860831969, 0.4679139345726910),
76            (0.6612093864662645, 0.3607615730481386),
77            (0.9324695142031521, 0.1713244923791704),
78        ],
79    },
80    // n=7
81    GaussLegendreTable {
82        zero_weight: 0.4179591836734694,
83        pairs: &[
84            (0.4058451513773972, 0.3818300505051189),
85            (0.7415311855993945, 0.2797053914892767),
86            (0.9491079123427585, 0.1294849661688697),
87        ],
88    },
89    // n=8
90    GaussLegendreTable {
91        zero_weight: 0.0,
92        pairs: &[
93            (0.1834346424956498, 0.3626837833783620),
94            (0.5255324099163290, 0.3137066458778873),
95            (0.7966664774136267, 0.2223810344533745),
96            (0.9602898564975363, 0.1012285362903763),
97        ],
98    },
99    // n=9
100    GaussLegendreTable {
101        zero_weight: 0.3302393550012598,
102        pairs: &[
103            (0.3242534234038089, 0.3123470770400029),
104            (0.6133714327005904, 0.2606106964029354),
105            (0.8360311073266358, 0.1806481606948574),
106            (0.9681602395076261, 0.0812743883615744),
107        ],
108    },
109    // n=10
110    GaussLegendreTable {
111        zero_weight: 0.0,
112        pairs: &[
113            (0.1488743389816312, 0.2955242247147529),
114            (0.4333953941292472, 0.2692667193099963),
115            (0.6794095682990244, 0.2190863625159820),
116            (0.8650633666889845, 0.1494513491505806),
117            (0.9739065285171717, 0.0666713443086881),
118        ],
119    },
120    // n=11
121    GaussLegendreTable {
122        zero_weight: 0.2729250867779006,
123        pairs: &[
124            (0.2695431559523450, 0.2628045445102467),
125            (0.5190961292068118, 0.2331937645919905),
126            (0.7301520055740494, 0.1862902109277343),
127            (0.8870625997680953, 0.1255803694649046),
128            (0.9782286581460570, 0.0556685671161737),
129        ],
130    },
131    // n=12
132    GaussLegendreTable {
133        zero_weight: 0.0,
134        pairs: &[
135            (0.1252334085114689, 0.2491470458134028),
136            (0.3678314989981802, 0.2334925365383548),
137            (0.5873179542866175, 0.2031674267230659),
138            (0.7699026741943047, 0.1600783285433462),
139            (0.9041172563704749, 0.1069393259953184),
140            (0.9815606342467192, 0.0471753363865118),
141        ],
142    },
143    // n=13
144    GaussLegendreTable {
145        zero_weight: 0.2325515532308739,
146        pairs: &[
147            (0.2304583159551348, 0.2262831802628972),
148            (0.4484927510364469, 0.2078160475368885),
149            (0.6423493394403402, 0.1781459807619457),
150            (0.8015780907333099, 0.1388735102197872),
151            (0.9175983992229779, 0.0921214998377285),
152            (0.9841830547185881, 0.0404840047653159),
153        ],
154    },
155    // n=14
156    GaussLegendreTable {
157        zero_weight: 0.0,
158        pairs: &[
159            (0.1080549487073437, 0.2152638534631578),
160            (0.3191123689278898, 0.2051984637212956),
161            (0.5152486363581541, 0.1855383974779378),
162            (0.6872929048116855, 0.1572031671581935),
163            (0.8272013150697650, 0.1215185706879032),
164            (0.9284348836635735, 0.0801580871597602),
165            (0.9862838086968123, 0.0351194603317519),
166        ],
167    },
168    // n=15
169    GaussLegendreTable {
170        zero_weight: 0.2025782419255613,
171        pairs: &[
172            (0.2011940939974345, 0.1984314853271116),
173            (0.3941513470775634, 0.1861610000155622),
174            (0.5709721726085388, 0.1662692058169939),
175            (0.7244177313601700, 0.1395706779261543),
176            (0.8482065834104272, 0.1071592204671719),
177            (0.9372733924007060, 0.0703660474881081),
178            (0.9879925180204854, 0.0307532419961173),
179        ],
180    },
181    // n=16
182    GaussLegendreTable {
183        zero_weight: 0.0,
184        pairs: &[
185            (0.0950125098376374, 0.1894506104550685),
186            (0.2816035507792589, 0.1826034150449236),
187            (0.4580167776572274, 0.1691565193950025),
188            (0.6178762444026438, 0.1495959888165767),
189            (0.7554044083550030, 0.1246289712555339),
190            (0.8656312023878318, 0.0951585116824928),
191            (0.9445750230732326, 0.0622535239386479),
192            (0.9894009349916499, 0.0271524594117541),
193        ],
194    },
195    // n=17
196    GaussLegendreTable {
197        zero_weight: 0.1794464703562065,
198        pairs: &[
199            (0.1784841814958479, 0.1765627053669926),
200            (0.3512317634538763, 0.1680041021564500),
201            (0.5126905370864769, 0.1540457610768103),
202            (0.6576711592166907, 0.1351363684685255),
203            (0.7815140038968014, 0.1118838471934040),
204            (0.8802391537269859, 0.0850361483171792),
205            (0.9506755217687678, 0.0554595293739872),
206            (0.9905754753144174, 0.0241483028685479),
207        ],
208    },
209    // n=18
210    GaussLegendreTable {
211        zero_weight: 0.0,
212        pairs: &[
213            (0.0847750130417353, 0.1691423829631436),
214            (0.2518862256915055, 0.1642764837458327),
215            (0.4117511614628426, 0.1546846751262652),
216            (0.5597708310739475, 0.1406429146706507),
217            (0.6916870430603532, 0.1225552067114785),
218            (0.8037049589725231, 0.1009420441062872),
219            (0.8926024664975557, 0.0764043285523263),
220            (0.9558239495713977, 0.0497145488949698),
221            (0.9915651684209309, 0.0216160135264833),
222        ],
223    },
224    // n=19
225    GaussLegendreTable {
226        zero_weight: 0.1610544498487837,
227        pairs: &[
228            (0.1603867852552212, 0.1589688433939543),
229            (0.3165640999636298, 0.1527660420658597),
230            (0.4645707413759609, 0.1426067021736066),
231            (0.6005453046616810, 0.1287539625393362),
232            (0.7209661773352294, 0.1115666455473340),
233            (0.8227146565371428, 0.0914900216224500),
234            (0.9031559036148179, 0.0690445427376412),
235            (0.9602081521348300, 0.0448142267656996),
236            (0.9924068438435844, 0.0194617882297265),
237        ],
238    },
239    // n=20
240    GaussLegendreTable {
241        zero_weight: 0.0,
242        pairs: &[
243            (0.0765265211334973, 0.1527533871307258),
244            (0.2277858511416451, 0.1491729864726037),
245            (0.3737060887154195, 0.1420961093183821),
246            (0.5108670019508271, 0.1316886384491766),
247            (0.6360536807265150, 0.1181945319615184),
248            (0.7463319064601508, 0.1019301198172404),
249            (0.8391169718222188, 0.0832767415767048),
250            (0.9122344282513259, 0.0626720483341091),
251            (0.9639719272779138, 0.0406014298003869),
252            (0.9931285991850949, 0.0176140071391521),
253        ],
254    },
255];
256
257/// Gauss-Legendre quadrature over `[a, b]` with `n` points.
258///
259/// Exact for polynomials of degree `<= 2n - 1`.
260/// Supports `n` from 1 to 20 with precomputed nodes/weights.
261///
262/// # Errors
263///
264/// Returns `IntegrationError::InvalidInput` if `n` is 0 or greater than 20.
265pub fn gauss_legendre<S, F>(f: F, a: S, b: S, n: usize) -> Result<S, IntegrationError>
266where
267    S: Scalar,
268    F: Fn(S) -> S,
269{
270    if n == 0 || n > 20 {
271        return Err(IntegrationError::InvalidInput(format!(
272            "gauss_legendre: n must be 1..20, got {}",
273            n
274        )));
275    }
276
277    let table = &GL[n - 1];
278    let mid = (a + b) * S::HALF;
279    let half_len = (b - a) * S::HALF;
280
281    let mut sum = S::ZERO;
282
283    // Zero node contribution (odd n)
284    if table.zero_weight != 0.0 {
285        sum += S::from_f64(table.zero_weight) * f(mid);
286    }
287
288    // Symmetric pairs
289    for &(node, weight) in table.pairs {
290        let x_pos = mid + half_len * S::from_f64(node);
291        let x_neg = mid - half_len * S::from_f64(node);
292        let w = S::from_f64(weight);
293        sum += w * (f(x_pos) + f(x_neg));
294    }
295
296    Ok(sum * half_len)
297}
298
299// ============================================================================
300// Gauss-Laguerre nodes and weights for [0, inf) with weight e^(-x)
301// ============================================================================
302
303/// Gauss-Laguerre nodes and weights for n = 1..10.
304/// Each entry is a slice of (node, weight) pairs.
305const LAGUERRE: [&[(f64, f64)]; 10] = [
306    // n=1
307    &[(1.0, 1.0)],
308    // n=2
309    &[
310        (0.5857864376269050, 0.8535533905932738),
311        (3.4142135623730950, 0.1464466094067263),
312    ],
313    // n=3
314    &[
315        (0.4157745567834791, 0.7110930099291730),
316        (2.2942803602790417, 0.2785177335692409),
317        (6.2899450829374792, 0.0103892565015861),
318    ],
319    // n=4
320    &[
321        (0.3225476896193923, 0.6031541043416336),
322        (1.7457611011583466, 0.3574186924377997),
323        (4.5366202969681338, 0.0388879085150054),
324        (9.3950709123011331, 0.0005392947055613),
325    ],
326    // n=5
327    &[
328        (0.2635603197181409, 0.5217556105828087),
329        (1.4134030591065168, 0.3986668110831759),
330        (3.5964257710407221, 0.0759424496817076),
331        (7.0858100058588376, 0.0036117586799220),
332        (12.6408008442757826, 0.0000233699723858),
333    ],
334    // n=6
335    &[
336        (0.2228466041792607, 0.4589646739499636),
337        (1.1889321016726231, 0.4170008307721210),
338        (2.9927363260593141, 0.1133733820740450),
339        (5.7751435691045105, 0.0103991974531490),
340        (9.8374674183825899, 0.0002610172028150),
341        (15.9828739806017017, 0.0000008985479064),
342    ],
343    // n=7
344    &[
345        (0.1930436765603623, 0.4093189517012739),
346        (1.0266648953391920, 0.4218312778617198),
347        (2.5678715903161220, 0.1471263486575052),
348        (4.9003530845264845, 0.0206335144687170),
349        (8.1821534445628607, 0.0010740101432807),
350        (12.7341802917978137, 0.0000158654643486),
351        (19.3957278622625403, 0.0000000317031548),
352    ],
353    // n=8
354    &[
355        (0.1702796323051010, 0.3691885893416375),
356        (0.9037017767993799, 0.4187867808143430),
357        (2.2510866298661307, 0.1757949866371718),
358        (4.2667001702876588, 0.0333434922612156),
359        (7.0459054023934657, 0.0027945362352257),
360        (10.7585160101809952, 0.0000907650877336),
361        (15.7406786412780046, 0.0000008485746716),
362        (22.8631317368892641, 0.0000000010480012),
363    ],
364    // n=9
365    &[
366        (0.1523222277318064, 0.3361264114948080),
367        (0.8072049551012728, 0.4112139804309697),
368        (2.0051351556193953, 0.1992875253708959),
369        (3.7834739733309144, 0.0474605627656516),
370        (6.2049567778766197, 0.0055999297143511),
371        (9.3729852516875898, 0.0003054237032950),
372        (13.4662369110920935, 0.0000065921230234),
373        (18.8335977889916966, 0.0000000412268815),
374        (26.3740718909273768, 0.0000000000329188),
375    ],
376    // n=10
377    &[
378        (0.1377934705404924, 0.3084411157650190),
379        (0.7294545495031705, 0.4011199291270945),
380        (1.8083429017403160, 0.2180682876118094),
381        (3.4014336978548995, 0.0620874560986777),
382        (5.5524961400638012, 0.0095015169751811),
383        (8.3301527467644965, 0.0007530083885876),
384        (11.8437858379000656, 0.0000282592334960),
385        (16.2792578313781021, 0.0000000424931398),
386        (21.9965858119807620, 0.0000000000183956),
387        (29.9206970122738766, 0.0000000000000009),
388    ],
389];
390
391/// Gauss-Laguerre quadrature for semi-infinite integrals.
392///
393/// Approximates `integral_0^inf f(x) * exp(-x) dx` using `n` points.
394/// The user provides `f(x)`; the weight function `exp(-x)` is built in.
395///
396/// # Errors
397///
398/// Returns `IntegrationError::InvalidInput` if `n` is 0 or greater than 10.
399pub fn gauss_laguerre<S, F>(f: F, n: usize) -> Result<S, IntegrationError>
400where
401    S: Scalar,
402    F: Fn(S) -> S,
403{
404    if n == 0 || n > 10 {
405        return Err(IntegrationError::InvalidInput(format!(
406            "gauss_laguerre: n must be 1..10, got {}",
407            n
408        )));
409    }
410
411    let table = LAGUERRE[n - 1];
412    let mut sum = S::ZERO;
413    for &(node, weight) in table {
414        sum += S::from_f64(weight) * f(S::from_f64(node));
415    }
416    Ok(sum)
417}
418
419// ============================================================================
420// Gauss-Hermite nodes and weights for (-inf, inf) with weight e^(-x^2)
421// ============================================================================
422
423/// Gauss-Hermite nodes and weights for n = 1..10.
424/// Stored as (positive_node, weight) pairs; symmetric about 0.
425/// For odd n, the zero node weight is stored separately.
426struct GaussHermiteTable {
427    zero_weight: f64,
428    pairs: &'static [(f64, f64)],
429}
430
431const GH: [GaussHermiteTable; 10] = [
432    // n=1
433    GaussHermiteTable {
434        zero_weight: 1.7724538509055159, // sqrt(pi)
435        pairs: &[],
436    },
437    // n=2
438    GaussHermiteTable {
439        zero_weight: 0.0,
440        pairs: &[(0.7071067811865476, 0.8862269254527580)],
441    },
442    // n=3
443    GaussHermiteTable {
444        zero_weight: 1.1816359006036774,
445        pairs: &[(1.2247448713915890, 0.2954089751509193)],
446    },
447    // n=4
448    GaussHermiteTable {
449        zero_weight: 0.0,
450        pairs: &[
451            (0.5246476232752904, 0.8049140900055128),
452            (1.6506801238857846, 0.0813128354472452),
453        ],
454    },
455    // n=5
456    GaussHermiteTable {
457        zero_weight: 0.9453087204829419,
458        pairs: &[
459            (0.9585724646138185, 0.3936193231522412),
460            (2.0201828704560856, 0.0199532420590459),
461        ],
462    },
463    // n=6
464    GaussHermiteTable {
465        zero_weight: 0.0,
466        pairs: &[
467            (0.4360774119276165, 0.7246295952243925),
468            (1.3358490740136970, 0.1570673203228566),
469            (2.3506049736744922, 0.0045300099055088),
470        ],
471    },
472    // n=7
473    GaussHermiteTable {
474        zero_weight: 0.8102646175568073,
475        pairs: &[
476            (0.8162878828589647, 0.4256072526101278),
477            (1.6735516287674714, 0.0545155828191270),
478            (2.6519613568352334, 0.0009717812450995),
479        ],
480    },
481    // n=8
482    GaussHermiteTable {
483        zero_weight: 0.0,
484        pairs: &[
485            (0.3811869902073222, 0.6611470125582413),
486            (1.1571937124467802, 0.2078023258148919),
487            (1.9816567566958429, 0.0170779830074457),
488            (2.9306374202572440, 0.0001996040722114),
489        ],
490    },
491    // n=9
492    GaussHermiteTable {
493        zero_weight: 0.7202352156060510,
494        pairs: &[
495            (0.7235510187528376, 0.4326515590025558),
496            (1.4685532892166679, 0.0886810550706316),
497            (2.2665805845318431, 0.0049436242755369),
498            (3.1909932017815276, 0.0000396069772633),
499        ],
500    },
501    // n=10
502    GaussHermiteTable {
503        zero_weight: 0.0,
504        pairs: &[
505            (0.3429013272237046, 0.6108626337353258),
506            (1.0366108297895137, 0.2401386110823147),
507            (1.7566836492998818, 0.0338743944554811),
508            (2.5327316742327897, 0.0013436457467812),
509            (3.4361591188377376, 0.0000076404328552),
510        ],
511    },
512];
513
514/// Gauss-Hermite quadrature for integrals over `(-inf, inf)`.
515///
516/// Approximates `integral_{-inf}^{inf} f(x) * exp(-x^2) dx` using `n` points.
517/// The user provides `f(x)`; the weight function `exp(-x^2)` is built in.
518///
519/// # Errors
520///
521/// Returns `IntegrationError::InvalidInput` if `n` is 0 or greater than 10.
522pub fn gauss_hermite<S, F>(f: F, n: usize) -> Result<S, IntegrationError>
523where
524    S: Scalar,
525    F: Fn(S) -> S,
526{
527    if n == 0 || n > 10 {
528        return Err(IntegrationError::InvalidInput(format!(
529            "gauss_hermite: n must be 1..10, got {}",
530            n
531        )));
532    }
533
534    let table = &GH[n - 1];
535    let mut sum = S::ZERO;
536
537    // Zero node
538    if table.zero_weight != 0.0 {
539        sum += S::from_f64(table.zero_weight) * f(S::ZERO);
540    }
541
542    // Symmetric pairs
543    for &(node, weight) in table.pairs {
544        let w = S::from_f64(weight);
545        sum += w * (f(S::from_f64(node)) + f(S::from_f64(-node)));
546    }
547
548    Ok(sum)
549}
550
551#[cfg(test)]
552mod tests {
553    use super::*;
554    use approx::assert_relative_eq;
555
556    // --- Gauss-Legendre ---
557
558    #[test]
559    fn test_gl_constant() {
560        // integral of 1 from -1 to 1 = 2
561        let result: f64 = gauss_legendre(|_: f64| 1.0, -1.0, 1.0, 1).unwrap();
562        assert_relative_eq!(result, 2.0, epsilon = 1e-14);
563    }
564
565    #[test]
566    fn test_gl_x_squared() {
567        // n=2 is exact for poly degree <= 3
568        let result: f64 = gauss_legendre(|x: f64| x * x, -1.0, 1.0, 2).unwrap();
569        assert_relative_eq!(result, 2.0 / 3.0, epsilon = 1e-14);
570    }
571
572    #[test]
573    fn test_gl_x_fourth() {
574        // integral of x^4 from -1 to 1 = 2/5, exact for n >= 3
575        let result: f64 = gauss_legendre(|x: f64| x.powi(4), -1.0, 1.0, 3).unwrap();
576        assert_relative_eq!(result, 2.0 / 5.0, epsilon = 1e-14);
577    }
578
579    #[test]
580    fn test_gl_interval_transform() {
581        // integral of x^2 from 0 to 1 = 1/3
582        let result: f64 = gauss_legendre(|x: f64| x * x, 0.0, 1.0, 5).unwrap();
583        assert_relative_eq!(result, 1.0 / 3.0, epsilon = 1e-14);
584    }
585
586    #[test]
587    fn test_gl_high_order() {
588        // integral of x^10 from -1 to 1 = 2/11, needs n >= 6 for exactness
589        let result: f64 = gauss_legendre(|x: f64| x.powi(10), -1.0, 1.0, 6).unwrap();
590        assert_relative_eq!(result, 2.0 / 11.0, epsilon = 1e-12);
591    }
592
593    #[test]
594    fn test_gl_n20_exp() {
595        // integral of exp(x) from -1 to 1 = e - 1/e
596        let result: f64 = gauss_legendre(|x: f64| x.exp(), -1.0, 1.0, 20).unwrap();
597        let expected = core::f64::consts::E - 1.0 / core::f64::consts::E;
598        assert_relative_eq!(result, expected, epsilon = 1e-14);
599    }
600
601    #[test]
602    fn test_gl_invalid_n() {
603        let result: Result<f64, _> = gauss_legendre(|x: f64| x, 0.0, 1.0, 0);
604        assert!(result.is_err());
605        let result: Result<f64, _> = gauss_legendre(|x: f64| x, 0.0, 1.0, 21);
606        assert!(result.is_err());
607    }
608
609    // --- Gauss-Laguerre ---
610
611    #[test]
612    fn test_laguerre_constant() {
613        // integral_0^inf 1 * e^(-x) dx = 1
614        let result: f64 = gauss_laguerre(|_: f64| 1.0, 5).unwrap();
615        assert_relative_eq!(result, 1.0, epsilon = 1e-14);
616    }
617
618    #[test]
619    fn test_laguerre_x() {
620        // integral_0^inf x * e^(-x) dx = 1 (Gamma(2) = 1!)
621        let result: f64 = gauss_laguerre(|x: f64| x, 5).unwrap();
622        assert_relative_eq!(result, 1.0, epsilon = 1e-12);
623    }
624
625    #[test]
626    fn test_laguerre_x_squared() {
627        // integral_0^inf x^2 * e^(-x) dx = 2 (Gamma(3) = 2!)
628        let result: f64 = gauss_laguerre(|x: f64| x * x, 5).unwrap();
629        assert_relative_eq!(result, 2.0, epsilon = 1e-10);
630    }
631
632    #[test]
633    fn test_laguerre_invalid_n() {
634        let result: Result<f64, _> = gauss_laguerre(|x: f64| x, 0);
635        assert!(result.is_err());
636        let result: Result<f64, _> = gauss_laguerre(|x: f64| x, 11);
637        assert!(result.is_err());
638    }
639
640    // --- Gauss-Hermite ---
641
642    #[test]
643    fn test_hermite_constant() {
644        // integral_{-inf}^{inf} 1 * e^(-x^2) dx = sqrt(pi)
645        let result: f64 = gauss_hermite(|_: f64| 1.0, 5).unwrap();
646        assert_relative_eq!(result, core::f64::consts::PI.sqrt(), epsilon = 1e-12);
647    }
648
649    #[test]
650    fn test_hermite_x_squared() {
651        // integral_{-inf}^{inf} x^2 * e^(-x^2) dx = sqrt(pi)/2
652        let result: f64 = gauss_hermite(|x: f64| x * x, 5).unwrap();
653        assert_relative_eq!(result, core::f64::consts::PI.sqrt() / 2.0, epsilon = 1e-12);
654    }
655
656    #[test]
657    fn test_hermite_invalid_n() {
658        let result: Result<f64, _> = gauss_hermite(|x: f64| x, 0);
659        assert!(result.is_err());
660        let result: Result<f64, _> = gauss_hermite(|x: f64| x, 11);
661        assert!(result.is_err());
662    }
663
664    // --- f32 tests ---
665
666    #[test]
667    fn test_gl_f32() {
668        let result: f32 = gauss_legendre(|x: f32| x * x, 0.0f32, 1.0f32, 5).unwrap();
669        assert!((result - 1.0 / 3.0).abs() < 1e-6);
670    }
671
672    #[test]
673    fn test_laguerre_f32() {
674        let result: f32 = gauss_laguerre(|_: f32| 1.0f32, 5).unwrap();
675        assert!((result - 1.0).abs() < 1e-5);
676    }
677}