1#![allow(dead_code)]
21#![allow(clippy::too_many_arguments)]
22
23use std::f64::consts::PI;
24
25#[derive(Debug, Clone, Copy, PartialEq)]
31pub struct Complex {
32 pub re: f64,
34 pub im: f64,
36}
37
38impl Complex {
39 pub fn new(re: f64, im: f64) -> Self {
41 Self { re, im }
42 }
43
44 pub fn from_real(re: f64) -> Self {
46 Self { re, im: 0.0 }
47 }
48
49 pub fn from_imag(im: f64) -> Self {
51 Self { re: 0.0, im }
52 }
53
54 pub fn zero() -> Self {
56 Self::new(0.0, 0.0)
57 }
58
59 pub fn one() -> Self {
61 Self::new(1.0, 0.0)
62 }
63
64 pub fn i() -> Self {
66 Self::new(0.0, 1.0)
67 }
68
69 pub fn abs(&self) -> f64 {
71 self.re.hypot(self.im)
72 }
73
74 pub fn arg(&self) -> f64 {
76 self.im.atan2(self.re)
77 }
78
79 pub fn conj(&self) -> Self {
81 Self::new(self.re, -self.im)
82 }
83
84 pub fn norm_sq(&self) -> f64 {
86 self.re * self.re + self.im * self.im
87 }
88
89 pub fn inv(&self) -> Self {
93 let n = self.norm_sq();
94 Self::new(self.re / n, -self.im / n)
95 }
96
97 pub fn exp(&self) -> Self {
99 let r = self.re.exp();
100 Self::new(r * self.im.cos(), r * self.im.sin())
101 }
102
103 pub fn ln(&self) -> Self {
105 Self::new(self.abs().ln(), self.arg())
106 }
107
108 pub fn sqrt(&self) -> Self {
110 let r = self.abs().sqrt();
111 let theta = self.arg() / 2.0;
112 Self::new(r * theta.cos(), r * theta.sin())
113 }
114
115 pub fn pow(&self, n: f64) -> Self {
117 if self.re == 0.0 && self.im == 0.0 {
118 return Self::zero();
119 }
120 let r = self.abs().powf(n);
121 let theta = self.arg() * n;
122 Self::new(r * theta.cos(), r * theta.sin())
123 }
124
125 pub fn cpow(&self, w: Self) -> Self {
127 (w * self.ln()).exp()
128 }
129
130 pub fn sin(&self) -> Self {
132 Self::new(
133 self.re.sin() * self.im.cosh(),
134 self.re.cos() * self.im.sinh(),
135 )
136 }
137
138 pub fn cos(&self) -> Self {
140 Self::new(
141 self.re.cos() * self.im.cosh(),
142 -self.re.sin() * self.im.sinh(),
143 )
144 }
145
146 pub fn sinh(&self) -> Self {
148 Self::new(
149 self.re.sinh() * self.im.cos(),
150 self.re.cosh() * self.im.sin(),
151 )
152 }
153
154 pub fn cosh(&self) -> Self {
156 Self::new(
157 self.re.cosh() * self.im.cos(),
158 self.re.sinh() * self.im.sin(),
159 )
160 }
161
162 pub fn tan(&self) -> Self {
164 self.sin() / self.cos()
165 }
166
167 pub fn approx_eq(&self, other: Self, tol: f64) -> bool {
169 (self.re - other.re).abs() < tol && (self.im - other.im).abs() < tol
170 }
171}
172
173impl std::ops::Add for Complex {
176 type Output = Self;
177 fn add(self, rhs: Self) -> Self {
178 Self::new(self.re + rhs.re, self.im + rhs.im)
179 }
180}
181
182impl std::ops::Sub for Complex {
183 type Output = Self;
184 fn sub(self, rhs: Self) -> Self {
185 Self::new(self.re - rhs.re, self.im - rhs.im)
186 }
187}
188
189impl std::ops::Mul for Complex {
190 type Output = Self;
191 fn mul(self, rhs: Self) -> Self {
192 Self::new(
193 self.re * rhs.re - self.im * rhs.im,
194 self.re * rhs.im + self.im * rhs.re,
195 )
196 }
197}
198
199impl std::ops::Div for Complex {
200 type Output = Self;
201 #[allow(clippy::suspicious_arithmetic_impl)]
202 fn div(self, rhs: Self) -> Self {
203 self * rhs.inv()
204 }
205}
206
207impl std::ops::Neg for Complex {
208 type Output = Self;
209 fn neg(self) -> Self {
210 Self::new(-self.re, -self.im)
211 }
212}
213
214impl std::ops::Add<f64> for Complex {
215 type Output = Self;
216 fn add(self, rhs: f64) -> Self {
217 Self::new(self.re + rhs, self.im)
218 }
219}
220
221impl std::ops::Mul<f64> for Complex {
222 type Output = Self;
223 fn mul(self, rhs: f64) -> Self {
224 Self::new(self.re * rhs, self.im * rhs)
225 }
226}
227
228impl std::ops::Div<f64> for Complex {
229 type Output = Self;
230 fn div(self, rhs: f64) -> Self {
231 Self::new(self.re / rhs, self.im / rhs)
232 }
233}
234
235impl std::fmt::Display for Complex {
236 fn fmt(&self, f: &mut std::fmt::Formatter<'_>) -> std::fmt::Result {
237 if self.im >= 0.0 {
238 write!(f, "{} + {}i", self.re, self.im)
239 } else {
240 write!(f, "{} - {}i", self.re, -self.im)
241 }
242 }
243}
244
245pub fn cauchy_integral(f: impl Fn(Complex) -> Complex, contour: &[Complex]) -> Complex {
262 if contour.len() < 2 {
263 return Complex::zero();
264 }
265 let n = contour.len();
266 let mut result = Complex::zero();
267 for i in 0..n {
268 let z0 = contour[i];
269 let z1 = contour[(i + 1) % n];
270 let dz = z1 - z0;
271 let mid = (f(z0) + f(z1)) / 2.0;
273 result = result + mid * dz;
274 }
275 result
276}
277
278pub fn residue(f: impl Fn(Complex) -> Complex, pole: Complex, order: usize) -> Complex {
293 let eps = 1e-5_f64;
294 let n_points = 512_usize;
295 let contour: Vec<Complex> = (0..n_points)
301 .map(|k| {
302 let theta = 2.0 * PI * (k as f64) / (n_points as f64);
303 pole + Complex::new(eps * theta.cos(), eps * theta.sin())
304 })
305 .collect();
306 let integral = cauchy_integral(
307 |z| {
308 let zp = z - pole;
309 let factor = zp.pow((order as f64) - 1.0);
310 f(z) * factor
311 },
312 &contour,
313 );
314 integral / Complex::new(0.0, 2.0 * PI)
316}
317
318pub fn conformal_map_joukowski(z: Complex, c: f64) -> Complex {
333 let c2 = Complex::from_real(c * c);
334 z + c2 / z
335}
336
337pub fn conformal_map_mobius(z: Complex, a: Complex, b: Complex, c: Complex, d: Complex) -> Complex {
348 (a * z + b) / (c * z + d)
349}
350
351#[derive(Debug, Clone)]
361pub struct ContourIntegral {
362 pub path: Vec<Complex>,
364}
365
366impl ContourIntegral {
367 pub fn new(path: Vec<Complex>) -> Self {
369 Self { path }
370 }
371
372 pub fn integrate(&self, f: impl Fn(Complex) -> Complex) -> Complex {
376 cauchy_integral(f, &self.path)
377 }
378
379 pub fn winding_number(&self, p: Complex) -> i32 {
390 if self.path.len() < 2 {
391 return 0;
392 }
393 let n = self.path.len();
394 let mut winding = 0.0_f64;
395 for i in 0..n {
396 let z0 = self.path[i] - p;
397 let z1 = self.path[(i + 1) % n] - p;
398 let dtheta = (z1 / z0).arg();
400 winding += dtheta;
401 }
402 (winding / (2.0 * PI)).round() as i32
403 }
404}
405
406const LANCZOS_G: f64 = 7.0;
412const LANCZOS_C: [f64; 9] = [
413 0.999_999_999_999_809_3,
414 676.520_368_121_885_1,
415 -1_259.139_216_722_402_8,
416 771.323_428_777_653_1,
417 -176.615_029_162_140_6,
418 12.507_343_278_686_905,
419 -0.138_571_095_265_720_12,
420 9.984_369_578_019_572e-6,
421 1.505_632_735_149_311_6e-7,
422];
423
424pub fn gamma_lanczos(z: Complex) -> Complex {
435 if z.re < 0.5 {
436 let pi_z = Complex::new(PI * z.re, PI * z.im);
438 let sin_piz = pi_z.sin();
439 let g1mz = gamma_lanczos(Complex::new(1.0 - z.re, -z.im));
440 Complex::from_real(PI) / (sin_piz * g1mz)
441 } else {
442 let z1 = z - Complex::from_real(1.0);
443 let mut x = Complex::from_real(LANCZOS_C[0]);
444 for (i, &c) in LANCZOS_C[1..].iter().enumerate() {
445 x = x + Complex::from_real(c) / (z1 + Complex::from_real((i + 1) as f64));
446 }
447 let t = z1 + Complex::from_real(LANCZOS_G + 0.5);
448 let sqrt2pi = (2.0 * PI).sqrt();
449 Complex::from_real(sqrt2pi) * t.cpow(z1 + Complex::from_real(0.5)) * (-t).exp() * x
450 }
451}
452
453pub fn beta_function(a: f64, b: f64) -> f64 {
468 let ga = gamma_lanczos(Complex::from_real(a));
469 let gb = gamma_lanczos(Complex::from_real(b));
470 let gab = gamma_lanczos(Complex::from_real(a + b));
471 (ga * gb / gab).re
472}
473
474pub fn zeta_euler_maclaurin(s: Complex, n_terms: usize) -> Complex {
490 let n = n_terms.max(2);
491 let mut sum = Complex::zero();
493 for k in 1..n {
494 sum = sum + Complex::from_real(k as f64).cpow(-s);
495 }
496 let nf = Complex::from_real(n as f64);
498 let half_term = nf.cpow(-s) / 2.0;
499 let integral_term = nf.cpow(Complex::from_real(1.0) - s) / (s - Complex::from_real(1.0));
500 let b2_term = Complex::from_real(1.0 / 12.0) * s * nf.cpow(-s - Complex::from_real(1.0));
502 sum + half_term + integral_term + b2_term
503}
504
505pub fn elliptic_k(k: f64) -> f64 {
521 debug_assert!(k.abs() < 1.0, "modulus |k| must be < 1");
522 let mut a = 1.0_f64;
523 let mut b = (1.0 - k * k).sqrt();
524 for _ in 0..50 {
525 let a_new = (a + b) / 2.0;
526 let b_new = (a * b).sqrt();
527 a = a_new;
528 b = b_new;
529 if (a - b).abs() < 1e-15 {
530 break;
531 }
532 }
533 PI / (2.0 * a)
534}
535
536pub fn elliptic_e(k: f64) -> f64 {
548 debug_assert!(k.abs() < 1.0, "modulus |k| must be < 1");
549 let mut a = 1.0_f64;
550 let mut b = (1.0 - k * k).sqrt();
551 let mut c = k;
552 let mut two_n = 1.0_f64;
553 let mut sum = 1.0 - 0.5 * c * c;
554 for _ in 0..50 {
555 let a_new = (a + b) / 2.0;
556 let b_new = (a * b).sqrt();
557 c = (a - b) / 2.0;
558 a = a_new;
559 b = b_new;
560 two_n *= 2.0;
561 sum -= two_n * c * c;
562 if c.abs() < 1e-15 {
563 break;
564 }
565 }
566 sum * PI / (2.0 * a)
567}
568
569pub fn bessel_j(n: i32, x: f64) -> f64 {
585 if x == 0.0 {
586 return if n == 0 { 1.0 } else { 0.0 };
587 }
588 if n < 0 {
589 let sign = if (-n) % 2 == 0 { 1.0 } else { -1.0 };
591 return sign * bessel_j(-n, x);
592 }
593 let j0 = bessel_j0(x);
595 let j1 = bessel_j1(x);
596 if n == 0 {
597 return j0;
598 }
599 if n == 1 {
600 return j1;
601 }
602 let mut jm1 = j0;
604 let mut j_curr = j1;
605 for k in 1..n {
606 let jp1 = (2.0 * (k as f64) / x) * j_curr - jm1;
607 jm1 = j_curr;
608 j_curr = jp1;
609 }
610 j_curr
611}
612
613fn bessel_j0(x: f64) -> f64 {
615 let x2 = x * x;
617 let mut term = 1.0_f64;
618 let mut sum = 1.0_f64;
619 for m in 1..=40_usize {
620 term *= -x2 / (4.0 * (m as f64) * (m as f64));
621 sum += term;
622 if term.abs() < 1e-16 * sum.abs() {
623 break;
624 }
625 }
626 sum
627}
628
629fn bessel_j1(x: f64) -> f64 {
631 let x2 = x * x;
633 let mut term = x / 2.0;
634 let mut sum = term;
635 for m in 1..=40_usize {
636 term *= -x2 / (4.0 * (m as f64) * ((m + 1) as f64));
637 sum += term;
638 if term.abs() < 1e-16 * sum.abs() {
639 break;
640 }
641 }
642 sum
643}
644
645pub fn bessel_y(n: i32, x: f64) -> f64 {
656 debug_assert!(x > 0.0, "Y_n(x) requires x > 0");
657 if n < 0 {
658 let sign = if (-n) % 2 == 0 { 1.0 } else { -1.0 };
659 return sign * bessel_y(-n, x);
660 }
661 let y0 = bessel_y0(x);
662 let y1 = bessel_y1(x);
663 if n == 0 {
664 return y0;
665 }
666 if n == 1 {
667 return y1;
668 }
669 let mut ym1 = y0;
670 let mut y_curr = y1;
671 for k in 1..n {
672 let yp1 = (2.0 * (k as f64) / x) * y_curr - ym1;
673 ym1 = y_curr;
674 y_curr = yp1;
675 }
676 y_curr
677}
678
679fn bessel_y0(x: f64) -> f64 {
681 let gamma_euler = 0.577_215_664_901_532_9;
684 let j0 = bessel_j0(x);
685 let x2 = x * x;
686 let mut term = 1.0_f64;
688 let mut h = 0.0_f64;
689 let mut correction = 0.0_f64;
690 for m in 1..=40_usize {
691 term *= -x2 / (4.0 * (m as f64) * (m as f64));
692 h += 1.0 / (m as f64);
693 correction += term * h;
694 if (term * h).abs() < 1e-16 {
695 break;
696 }
697 }
698 (2.0 / PI) * ((x / 2.0).ln() + gamma_euler) * j0 + (2.0 / PI) * correction
699}
700
701fn bessel_y1(x: f64) -> f64 {
703 let gamma_euler = 0.577_215_664_901_532_9;
704 let j1 = bessel_j1(x);
705 let x2 = x * x;
706 let mut t2 = x / 2.0;
708 let mut h = 1.0_f64;
709 let mut correction = 0.0_f64;
710 for m in 1..=40_usize {
711 t2 *= -x2 / (4.0 * (m as f64) * ((m + 1) as f64));
712 h += 1.0 / ((m + 1) as f64);
713 correction += t2 * (h + h - 1.0 / (m as f64));
714 if (t2 * h).abs() < 1e-16 * correction.abs().max(1e-30) {
715 break;
716 }
717 }
718 (2.0 / PI) * ((x / 2.0).ln() + gamma_euler - 0.5) * j1 - (1.0 / (PI * x))
719 + (2.0 / PI) * correction
720}
721
722pub fn legendre_p(n: usize, x: f64) -> f64 {
738 if n == 0 {
739 return 1.0;
740 }
741 if n == 1 {
742 return x;
743 }
744 let mut p_prev = 1.0_f64;
745 let mut p_curr = x;
746 for k in 1..n {
747 let p_next = ((2 * k + 1) as f64 * x * p_curr - k as f64 * p_prev) / ((k + 1) as f64);
748 p_prev = p_curr;
749 p_curr = p_next;
750 }
751 p_curr
752}
753
754pub fn hermite_h(n: usize, x: f64) -> f64 {
765 if n == 0 {
766 return 1.0;
767 }
768 if n == 1 {
769 return 2.0 * x;
770 }
771 let mut h_prev = 1.0_f64;
772 let mut h_curr = 2.0 * x;
773 for k in 1..n {
774 let h_next = 2.0 * x * h_curr - 2.0 * k as f64 * h_prev;
775 h_prev = h_curr;
776 h_curr = h_next;
777 }
778 h_curr
779}
780
781#[cfg(test)]
786mod tests {
787 use super::*;
788
789 const TOL: f64 = 1e-9;
790 const MED_TOL: f64 = 1e-6;
791
792 #[test]
795 fn test_complex_add() {
796 let a = Complex::new(1.0, 2.0);
797 let b = Complex::new(3.0, 4.0);
798 let c = a + b;
799 assert!((c.re - 4.0).abs() < TOL);
800 assert!((c.im - 6.0).abs() < TOL);
801 }
802
803 #[test]
804 fn test_complex_mul() {
805 let a = Complex::new(1.0, 2.0);
807 let b = Complex::new(3.0, 4.0);
808 let c = a * b;
809 assert!((c.re - (-5.0)).abs() < TOL);
810 assert!((c.im - 10.0).abs() < TOL);
811 }
812
813 #[test]
814 fn test_complex_div() {
815 let a = Complex::new(1.0, 2.0);
817 let b = Complex::new(1.0, 1.0);
818 let c = a / b;
819 assert!((c.re - 1.5).abs() < TOL);
820 assert!((c.im - 0.5).abs() < TOL);
821 }
822
823 #[test]
824 fn test_complex_abs_arg() {
825 let z = Complex::new(0.0, 1.0);
826 assert!((z.abs() - 1.0).abs() < TOL);
827 assert!((z.arg() - PI / 2.0).abs() < TOL);
828 }
829
830 #[test]
831 fn test_complex_exp() {
832 let z = Complex::new(0.0, PI);
834 let ez = z.exp();
835 assert!((ez.re - (-1.0)).abs() < 1e-14);
836 assert!(ez.im.abs() < 1e-14);
837 }
838
839 #[test]
840 fn test_complex_ln() {
841 let z = Complex::from_real(std::f64::consts::E);
843 let lnz = z.ln();
844 assert!((lnz.re - 1.0).abs() < TOL);
845 assert!(lnz.im.abs() < TOL);
846 }
847
848 #[test]
849 fn test_complex_sqrt() {
850 let z = Complex::from_real(-1.0);
852 let s = z.sqrt();
853 assert!(s.re.abs() < 1e-14);
854 assert!((s.im - 1.0).abs() < 1e-14);
855 }
856
857 #[test]
858 fn test_complex_pow() {
859 let z = Complex::new(0.0, 1.0);
861 let z2 = z.pow(2.0);
862 assert!((z2.re - (-1.0)).abs() < 1e-14);
863 assert!(z2.im.abs() < 1e-14);
864 }
865
866 #[test]
867 fn test_complex_sin_cos() {
868 let z = Complex::zero();
870 assert!(z.sin().abs() < TOL);
871 assert!((z.cos().re - 1.0).abs() < TOL);
872 }
873
874 #[test]
875 fn test_complex_sinh_cosh() {
876 let z = Complex::new(1.0, 0.5);
878 let ch = z.cosh();
879 let sh = z.sinh();
880 let diff = ch * ch - sh * sh;
881 assert!((diff.re - 1.0).abs() < 1e-12, "re diff = {}", diff.re);
882 assert!(diff.im.abs() < 1e-12, "im = {}", diff.im);
883 }
884
885 #[test]
886 fn test_complex_conj() {
887 let z = Complex::new(3.0, 4.0);
888 let zc = z.conj();
889 assert!((zc.re - 3.0).abs() < TOL);
890 assert!((zc.im - (-4.0)).abs() < TOL);
891 }
892
893 #[test]
894 fn test_complex_inv() {
895 let z = Complex::new(1.0, 1.0);
897 let inv = z.inv();
898 assert!((inv.re - 0.5).abs() < TOL);
899 assert!((inv.im - (-0.5)).abs() < TOL);
900 }
901
902 #[test]
905 fn test_cauchy_integral_constant() {
906 let n = 64_usize;
908 let contour: Vec<Complex> = (0..n)
909 .map(|k| {
910 let t = 2.0 * PI * k as f64 / n as f64;
911 Complex::new(t.cos(), t.sin())
912 })
913 .collect();
914 let result = cauchy_integral(|_z| Complex::one(), &contour);
915 assert!(result.abs() < 1e-10);
916 }
917
918 #[test]
919 fn test_cauchy_integral_z() {
920 let n = 256_usize;
922 let contour: Vec<Complex> = (0..n)
923 .map(|k| {
924 let t = 2.0 * PI * k as f64 / n as f64;
925 Complex::new(t.cos(), t.sin())
926 })
927 .collect();
928 let result = cauchy_integral(|z| z, &contour);
929 assert!(result.abs() < 1e-10);
930 }
931
932 #[test]
933 fn test_cauchy_integral_one_over_z() {
934 let n = 512_usize;
936 let contour: Vec<Complex> = (0..n)
937 .map(|k| {
938 let t = 2.0 * PI * k as f64 / n as f64;
939 Complex::new(t.cos(), t.sin())
940 })
941 .collect();
942 let result = cauchy_integral(|z| z.inv(), &contour);
943 assert!(result.re.abs() < 1e-3);
944 assert!((result.im - 2.0 * PI).abs() < 1e-2);
945 }
946
947 #[test]
950 fn test_residue_simple_pole() {
951 let r = residue(|z| z.inv(), Complex::zero(), 1);
953 assert!((r.re - 1.0).abs() < 1e-4);
954 assert!(r.im.abs() < 1e-4);
955 }
956
957 #[test]
960 fn test_joukowski_identity_at_two() {
961 let w = conformal_map_joukowski(Complex::from_real(2.0), 1.0);
963 assert!((w.re - 2.5).abs() < TOL);
964 assert!(w.im.abs() < TOL);
965 }
966
967 #[test]
968 fn test_mobius_identity() {
969 let z = Complex::new(2.0, 3.0);
971 let w = conformal_map_mobius(
972 z,
973 Complex::one(),
974 Complex::zero(),
975 Complex::zero(),
976 Complex::one(),
977 );
978 assert!(z.approx_eq(w, TOL));
979 }
980
981 #[test]
984 fn test_winding_number_inside() {
985 let n = 64_usize;
986 let circle: Vec<Complex> = (0..n)
987 .map(|k| {
988 let t = 2.0 * PI * k as f64 / n as f64;
989 Complex::new(t.cos(), t.sin())
990 })
991 .collect();
992 let ci = ContourIntegral::new(circle);
993 assert_eq!(ci.winding_number(Complex::zero()), 1);
994 }
995
996 #[test]
997 fn test_winding_number_outside() {
998 let n = 64_usize;
999 let circle: Vec<Complex> = (0..n)
1000 .map(|k| {
1001 let t = 2.0 * PI * k as f64 / n as f64;
1002 Complex::new(t.cos(), t.sin())
1003 })
1004 .collect();
1005 let ci = ContourIntegral::new(circle);
1006 assert_eq!(ci.winding_number(Complex::new(10.0, 0.0)), 0);
1008 }
1009
1010 #[test]
1013 fn test_gamma_factorial() {
1014 for (n, expected) in [
1016 (1.0_f64, 1.0),
1017 (2.0, 1.0),
1018 (3.0, 2.0),
1019 (4.0, 6.0),
1020 (5.0, 24.0),
1021 ] {
1022 let g = gamma_lanczos(Complex::from_real(n)).re;
1023 assert!(
1024 (g - expected).abs() < MED_TOL,
1025 "Γ({n}) = {g}, expected {expected}"
1026 );
1027 }
1028 }
1029
1030 #[test]
1031 fn test_gamma_half() {
1032 let g = gamma_lanczos(Complex::from_real(0.5)).re;
1034 assert!((g - PI.sqrt()).abs() < MED_TOL);
1035 }
1036
1037 #[test]
1040 fn test_beta_symmetric() {
1041 let a = 2.5;
1043 let b = 3.5;
1044 assert!((beta_function(a, b) - beta_function(b, a)).abs() < MED_TOL);
1045 }
1046
1047 #[test]
1048 fn test_beta_known_value() {
1049 assert!((beta_function(1.0, 1.0) - 1.0).abs() < MED_TOL);
1051 assert!((beta_function(2.0, 2.0) - 1.0 / 6.0).abs() < MED_TOL);
1053 }
1054
1055 #[test]
1058 fn test_zeta_known_value() {
1059 let z = zeta_euler_maclaurin(Complex::from_real(2.0), 50);
1061 assert!((z.re - PI * PI / 6.0).abs() < 1e-4);
1062 }
1063
1064 #[test]
1067 fn test_elliptic_k_zero() {
1068 assert!((elliptic_k(0.0) - PI / 2.0).abs() < TOL);
1070 }
1071
1072 #[test]
1073 fn test_elliptic_e_zero() {
1074 assert!((elliptic_e(0.0) - PI / 2.0).abs() < TOL);
1076 }
1077
1078 #[test]
1079 fn test_elliptic_ke_relation() {
1080 let k = 0.7;
1082 assert!(elliptic_e(k) <= elliptic_k(k) + 1e-10);
1083 }
1084
1085 #[test]
1088 fn test_bessel_j0_zero() {
1089 assert!((bessel_j(0, 0.0) - 1.0).abs() < TOL);
1091 }
1092
1093 #[test]
1094 fn test_bessel_j1_zero() {
1095 assert!(bessel_j(1, 0.0).abs() < TOL);
1097 }
1098
1099 #[test]
1100 fn test_bessel_j0_known() {
1101 let v = bessel_j(0, 2.40483).abs();
1103 assert!(v < 1e-4);
1104 }
1105
1106 #[test]
1107 fn test_bessel_recurrence() {
1108 let x = 3.0;
1110 let n = 2;
1111 let lhs = bessel_j(n - 1, x) + bessel_j(n + 1, x);
1112 let rhs = (2.0 * n as f64 / x) * bessel_j(n, x);
1113 assert!((lhs - rhs).abs() < 1e-10);
1114 }
1115
1116 #[test]
1119 fn test_legendre_p0_p1() {
1120 assert!((legendre_p(0, 0.5) - 1.0).abs() < TOL);
1121 assert!((legendre_p(1, 0.5) - 0.5).abs() < TOL);
1122 }
1123
1124 #[test]
1125 fn test_legendre_p2() {
1126 let x = 0.7;
1128 let expected = (3.0 * x * x - 1.0) / 2.0;
1129 assert!((legendre_p(2, x) - expected).abs() < TOL);
1130 }
1131
1132 #[test]
1133 fn test_legendre_orthogonality_approx() {
1134 let n = 200_usize;
1136 let mut sum = 0.0_f64;
1137 for i in 0..n {
1138 let x = -1.0 + 2.0 * (i as f64 + 0.5) / n as f64;
1139 sum += legendre_p(2, x) * legendre_p(0, x) * (2.0 / n as f64);
1140 }
1141 assert!(sum.abs() < 1e-4);
1142 }
1143
1144 #[test]
1147 fn test_hermite_h0_h1() {
1148 assert!((hermite_h(0, 1.0) - 1.0).abs() < TOL);
1149 assert!((hermite_h(1, 1.0) - 2.0).abs() < TOL);
1150 }
1151
1152 #[test]
1153 fn test_hermite_h2() {
1154 let x = 1.5;
1156 let expected = 4.0 * x * x - 2.0;
1157 assert!((hermite_h(2, x) - expected).abs() < TOL);
1158 }
1159
1160 #[test]
1161 fn test_hermite_h3() {
1162 let x = 2.0;
1164 let expected = 8.0 * x * x * x - 12.0 * x;
1165 assert!((hermite_h(3, x) - expected).abs() < TOL);
1166 }
1167}