russell_lab/math/
bessel_1.rs

1use super::{PI, SQRT_PI, TWO_129, TWO_M27};
2
3// This implementation is based on j1.go file from Go (1.22.1),
4// which, in turn, is based on the FreeBSD code as explained below.
5//
6// Copyright 2010 The Go Authors. All rights reserved.
7// Use of this source code is governed by a BSD-style
8// license that can be found in the LICENSE file.
9//
10// Bessel function of the first and second kinds of order one.
11//
12// The original C code and the long comment below are
13// from FreeBSD's /usr/src/lib/msun/src/e_j1.c and
14// came with this notice. The go code is a simplified
15// version of the original C.
16//
17// ====================================================
18// Copyright (C) 1993 by Sun Microsystems, Inc. All rights reserved.
19//
20// Developed at SunPro, a Sun Microsystems, Inc. business.
21// Permission to use, copy, modify, and distribute this
22// software is freely granted, provided that this notice
23// is preserved.
24// ====================================================
25
26// R0/S0 on [0, 2]
27const R00: f64 = -6.25000000000000000000e-02; // 0xBFB0000000000000
28const R01: f64 = 1.40705666955189706048e-03; // 0x3F570D9F98472C61
29const R02: f64 = -1.59955631084035597520e-05; // 0xBEF0C5C6BA169668
30const R03: f64 = 4.96727999609584448412e-08; // 0x3E6AAAFA46CA0BD9
31const S01: f64 = 1.91537599538363460805e-02; // 0x3F939D0B12637E53
32const S02: f64 = 1.85946785588630915560e-04; // 0x3F285F56B9CDF664
33const S03: f64 = 1.17718464042623683263e-06; // 0x3EB3BFF8333F8498
34const S04: f64 = 5.04636257076217042715e-09; // 0x3E35AC88C97DFF2C
35const S05: f64 = 1.23542274426137913908e-11; // 0x3DAB2ACFCFB97ED8
36
37/// Evaluates the Bessel function J1(x) for any real x
38///
39/// See: <https://mathworld.wolfram.com/BesselFunctionoftheFirstKind.html>
40///
41/// See also: <https://en.wikipedia.org/wiki/Bessel_function>
42///
43/// # Special cases
44///
45/// * `J1(NaN)  = NaN`
46/// * `J1(±Inf) = 0.0`
47/// * `J1(0.0)  = 0.0`
48///
49/// # Examples
50///
51/// ![Bessel J1](https://raw.githubusercontent.com/cpmech/russell/main/russell_lab/data/figures/math_bessel_functions_j1.svg)
52///
53/// ```
54/// use russell_lab::{approx_eq, math};
55///
56/// approx_eq(math::bessel_j1(2.0), 0.57672480775687339, 1e-15);
57/// ```
58pub fn bessel_j1(x: f64) -> f64 {
59    //
60    // For tiny x, use j1(x) = x/2 - x**3/16 + x**5/384 - ...
61    //
62    // Reduce x to |x| since j1(x)=-j1(-x), and:
63    //
64    // for x in (0,2)
65    //
66    // j1(x) = x/2 + x*z*R0/S0,  where z = x*x;
67    //
68    // (precision:  |j1/x - 1/2 - R0/S0 |<2**-61.51 )
69    //
70    // for x in (2,inf)
71    //
72    // j1(x) = sqrt(2/(pi*x))*(p1(x)*cos(x1)-q1(x)*sin(x1))
73    // y1(x) = sqrt(2/(pi*x))*(p1(x)*sin(x1)+q1(x)*cos(x1))
74    //
75    // where x1 = x-3*pi/4.
76    //
77    // Compute sin(x1),cos(x1) as follows:
78    //
79    // cos(x1) =  cos(x)cos(3pi/4)+sin(x)sin(3pi/4)
80    //         =  1/sqrt(2) * (sin(x) - cos(x))
81    // sin(x1) =  sin(x)cos(3pi/4)-cos(x)sin(3pi/4)
82    //         = -1/sqrt(2) * (sin(x) + cos(x))
83    //
84    // To avoid cancellation, use:
85    //
86    // sin(x) +- cos(x) = -cos(2x)/(sin(x) -+ cos(x))
87
88    if f64::is_nan(x) {
89        return f64::NAN;
90    } else if f64::is_infinite(x) {
91        return 0.0;
92    } else if x == 0.0 {
93        return 0.0;
94    }
95
96    let (xx, negative) = if x < 0.0 { (-x, true) } else { (x, false) };
97
98    if xx >= 2.0 {
99        let (s, c) = f64::sin_cos(xx);
100        let mut ss = -s - c;
101        let mut cc = s - c;
102
103        // make sure x+x does not overflow
104        if xx < f64::MAX / 2.0 {
105            let z = f64::cos(xx + xx);
106            if s * c > 0.0 {
107                cc = z / ss;
108            } else {
109                ss = z / cc;
110            }
111        }
112
113        let z = if xx > TWO_129 {
114            (1.0 / SQRT_PI) * cc / f64::sqrt(xx)
115        } else {
116            let u = pone(xx);
117            let v = qone(xx);
118            (1.0 / SQRT_PI) * (u * cc - v * ss) / f64::sqrt(xx)
119        };
120
121        if negative {
122            return -z;
123        } else {
124            return z;
125        }
126    }
127
128    if xx < TWO_M27 {
129        return 0.5 * xx;
130    }
131
132    let mut z = xx * xx;
133    let mut r = z * (R00 + z * (R01 + z * (R02 + z * R03)));
134    let s = 1.0 + z * (S01 + z * (S02 + z * (S03 + z * (S04 + z * S05))));
135    r *= xx;
136    z = 0.5 * xx + r / s;
137
138    if negative {
139        return -z;
140    } else {
141        return z;
142    }
143}
144
145// constant computed with Mathematica
146const TWO_M54: f64 = 5.5511151231257827021181583404541015625000000000000e-17; // 2**-54 0x3c90000000000000 Mathematica: N[2^-54, 50]
147
148const U00: f64 = -1.96057090646238940668e-01; // 0xBFC91866143CBC8A
149const U01: f64 = 5.04438716639811282616e-02; // 0x3FA9D3C776292CD1
150const U02: f64 = -1.91256895875763547298e-03; // 0xBF5F55E54844F50F
151const U03: f64 = 2.35252600561610495928e-05; // 0x3EF8AB038FA6B88E
152const U04: f64 = -9.19099158039878874504e-08; // 0xBE78AC00569105B8
153const V00: f64 = 1.99167318236649903973e-02; // 0x3F94650D3F4DA9F0
154const V01: f64 = 2.02552581025135171496e-04; // 0x3F2A8C896C257764
155const V02: f64 = 1.35608801097516229404e-06; // 0x3EB6C05A894E8CA6
156const V03: f64 = 6.22741452364621501295e-09; // 0x3E3ABF1D5BA69A86
157const V04: f64 = 1.66559246207992079114e-11; // 0x3DB25039DACA772A
158
159/// Evaluates the Bessel function Y1(x) for positive real x
160///
161/// See: <https://mathworld.wolfram.com/BesselFunctionoftheSecondKind.html>
162///
163/// See also: <https://en.wikipedia.org/wiki/Bessel_function>
164///
165/// # Special cases
166///
167/// * `Y1(x < 0.0) = NaN`
168/// * `Y1(NaN)     = NaN`
169/// * `Y1(+Inf)    = 0.0`
170/// * `Y1(0.0)     = -Inf`
171///
172/// # Examples
173///
174/// ![Bessel Y1](https://raw.githubusercontent.com/cpmech/russell/main/russell_lab/data/figures/math_bessel_functions_y1.svg)
175///
176/// ```
177/// use russell_lab::{approx_eq, math};
178///
179/// approx_eq(math::bessel_y1(2.0), -0.10703243154093755, 1e-15);
180/// ```
181pub fn bessel_y1(x: f64) -> f64 {
182    //
183    // Screen out x<=0 cases: y1(0)=-inf, y1(x<0)=NaN
184    //
185    // For x<2. Since
186    //
187    // y1(x) = 2/pi*(j1(x)*(ln(x/2)+Euler)-1/x-x/2+5/64*x**3-...)
188    //
189    // therefore y1(x)-2/pi*j1(x)*ln(x)-1/x is an odd function.
190    //
191    // Use the following function to approximate y1,
192    //
193    // y1(x) = x*U(z)/V(z) + (2/pi)*(j1(x)*ln(x)-1/x), z= x**2
194    //
195    // where for x in [0,2] (abs err less than 2**-65.89)
196    //
197    // U(z) = U0[0] + U0[1]*z + ... + U0[4]*z**4
198    // V(z) = 1  + v0[0]*z + ... + v0[4]*z**5
199    //
200    // Note: For tiny x, 1/x dominate y1 and hence
201    // y1(tiny) = -2/pi/tiny, (choose tiny<2**-54)
202    //
203    // For x>=2.
204    //
205    // y1(x) = sqrt(2/(pi*x))*(p1(x)*sin(x1)+q1(x)*cos(x1))
206    //
207    // where x1 = x-3*pi/4.
208    //
209    // Compute sin(x1),cos(x1) by method mentioned in bessel_j1
210
211    if x < 0.0 || f64::is_nan(x) {
212        return f64::NAN;
213    } else if f64::is_infinite(x) {
214        // this handles +Inf since x < 0.0 was handled already
215        return 0.0;
216    } else if x == 0.0 {
217        return f64::NEG_INFINITY;
218    }
219
220    if x >= 2.0 {
221        let (s, c) = f64::sin_cos(x);
222        let mut ss = -s - c;
223        let mut cc = s - c;
224
225        // make sure x+x does not overflow
226        if x < f64::MAX / 2.0 {
227            let z = f64::cos(x + x);
228            if s * c > 0.0 {
229                cc = z / ss;
230            } else {
231                ss = z / cc;
232            }
233        }
234
235        let z = if x > TWO_129 {
236            (1.0 / SQRT_PI) * ss / f64::sqrt(x)
237        } else {
238            let u = pone(x);
239            let v = qone(x);
240            (1.0 / SQRT_PI) * (u * ss + v * cc) / f64::sqrt(x)
241        };
242
243        return z;
244    }
245
246    if x <= TWO_M54 {
247        return -(2.0 / PI) / x;
248    }
249
250    let z = x * x;
251    let u = U00 + z * (U01 + z * (U02 + z * (U03 + z * U04)));
252    let v = 1.0 + z * (V00 + z * (V01 + z * (V02 + z * (V03 + z * V04))));
253    x * (u / v) + (2.0 / PI) * (bessel_j1(x) * f64::ln(x) - 1.0 / x)
254}
255
256// For x >= 8, the asymptotic expansions of pone is
257//
258// 1 + 15/128 s**2 - 4725/2**15 s**4 - ..., where s = 1/x.
259//
260// We approximate pone by
261//
262// pone(x) = 1 + (R/S)
263//
264// where R = pr0 + pr1*s**2 + pr2*s**4 + ... + pr5*s**10
265//
266// S = 1 + ps0*s**2 + ... + ps4*s**10
267//
268// and
269//
270// | pone(x)-1-R/S | <= 2**(-60.06)
271
272// for x in [inf, 8]=1/[0,0.125]
273const P1R8: [f64; 6] = [
274    0.00000000000000000000e+00, // 0x0000000000000000
275    1.17187499999988647970e-01, // 0x3FBDFFFFFFFFFCCE
276    1.32394806593073575129e+01, // 0x402A7A9D357F7FCE
277    4.12051854307378562225e+02, // 0x4079C0D4652EA590
278    3.87474538913960532227e+03, // 0x40AE457DA3A532CC
279    7.91447954031891731574e+03, // 0x40BEEA7AC32782DD
280];
281const P1S8: [f64; 5] = [
282    1.14207370375678408436e+02, // 0x405C8D458E656CAC
283    3.65093083420853463394e+03, // 0x40AC85DC964D274F
284    3.69562060269033463555e+04, // 0x40E20B8697C5BB7F
285    9.76027935934950801311e+04, // 0x40F7D42CB28F17BB
286    3.08042720627888811578e+04, // 0x40DE1511697A0B2D
287];
288
289// for x in [8,4.5454] = 1/[0.125,0.22001]
290const P1R5: [f64; 6] = [
291    1.31990519556243522749e-11, // 0x3DAD0667DAE1CA7D
292    1.17187493190614097638e-01, // 0x3FBDFFFFE2C10043
293    6.80275127868432871736e+00, // 0x401B36046E6315E3
294    1.08308182990189109773e+02, // 0x405B13B9452602ED
295    5.17636139533199752805e+02, // 0x40802D16D052D649
296    5.28715201363337541807e+02, // 0x408085B8BB7E0CB7
297];
298const P1S5: [f64; 5] = [
299    5.92805987221131331921e+01, // 0x404DA3EAA8AF633D
300    9.91401418733614377743e+02, // 0x408EFB361B066701
301    5.35326695291487976647e+03, // 0x40B4E9445706B6FB
302    7.84469031749551231769e+03, // 0x40BEA4B0B8A5BB15
303    1.50404688810361062679e+03, // 0x40978030036F5E51
304];
305
306// for x in[4.5453,2.8571] = 1/[0.2199,0.35001]
307const P1R3: [f64; 6] = [
308    3.02503916137373618024e-09, // 0x3E29FC21A7AD9EDD
309    1.17186865567253592491e-01, // 0x3FBDFFF55B21D17B
310    3.93297750033315640650e+00, // 0x400F76BCE85EAD8A
311    3.51194035591636932736e+01, // 0x40418F489DA6D129
312    9.10550110750781271918e+01, // 0x4056C3854D2C1837
313    4.85590685197364919645e+01, // 0x4048478F8EA83EE5
314];
315const P1S3: [f64; 5] = [
316    3.47913095001251519989e+01, // 0x40416549A134069C
317    3.36762458747825746741e+02, // 0x40750C3307F1A75F
318    1.04687139975775130551e+03, // 0x40905B7C5037D523
319    8.90811346398256432622e+02, // 0x408BD67DA32E31E9
320    1.03787932439639277504e+02, // 0x4059F26D7C2EED53
321];
322
323// for x in [2.8570,2] = 1/[0.3499,0.5]
324const P1R2: [f64; 6] = [
325    1.07710830106873743082e-07, // 0x3E7CE9D4F65544F4
326    1.17176219462683348094e-01, // 0x3FBDFF42BE760D83
327    2.36851496667608785174e+00, // 0x4002F2B7F98FAEC0
328    1.22426109148261232917e+01, // 0x40287C377F71A964
329    1.76939711271687727390e+01, // 0x4031B1A8177F8EE2
330    5.07352312588818499250e+00, // 0x40144B49A574C1FE
331];
332const P1S2: [f64; 5] = [
333    2.14364859363821409488e+01, // 0x40356FBD8AD5ECDC
334    1.25290227168402751090e+02, // 0x405F529314F92CD5
335    2.32276469057162813669e+02, // 0x406D08D8D5A2DBD9
336    1.17679373287147100768e+02, // 0x405D6B7ADA1884A9
337    8.36463893371618283368e+00, // 0x4020BAB1F44E5192
338];
339
340fn pone(x: f64) -> f64 {
341    let (p, q) = if x >= 8.0 {
342        (&P1R8, &P1S8)
343    } else if x >= 4.5454 {
344        (&P1R5, &P1S5)
345    } else if x >= 2.8571 {
346        (&P1R3, &P1S3)
347    } else if x >= 2.0 {
348        (&P1R2, &P1S2)
349    } else {
350        panic!("INTERNAL ERROR: x must be ≥ 2.0 for pone");
351    };
352    let z = 1.0 / (x * x);
353    let r = p[0] + z * (p[1] + z * (p[2] + z * (p[3] + z * (p[4] + z * p[5]))));
354    let s = 1.0 + z * (q[0] + z * (q[1] + z * (q[2] + z * (q[3] + z * q[4]))));
355    1.0 + r / s
356}
357
358// For x >= 8, the asymptotic expansions of qone is
359//
360// 3/8 s - 105/1024 s**3 - ..., where s = 1/x.
361//
362// Approximate qone by
363//
364// qone(x) = s*(0.375 + (R/S))
365//
366// where
367//
368// R = qr1*s**2 + qr2*s**4 + ... + qr5*s**10
369// S = 1 + qs1*s**2 + ... + qs6*s**12
370//
371// and
372//
373// | qone(x)/s -0.375-R/S | <= 2**(-61.13)
374
375// for x in [inf, 8] = 1/[0,0.125]
376const Q1R8: [f64; 6] = [
377    0.00000000000000000000e+00,  // 0x0000000000000000
378    -1.02539062499992714161e-01, // 0xBFBA3FFFFFFFFDF3
379    -1.62717534544589987888e+01, // 0xC0304591A26779F7
380    -7.59601722513950107896e+02, // 0xC087BCD053E4B576
381    -1.18498066702429587167e+04, // 0xC0C724E740F87415
382    -4.84385124285750353010e+04, // 0xC0E7A6D065D09C6A
383];
384const Q1S8: [f64; 6] = [
385    1.61395369700722909556e+02,  // 0x40642CA6DE5BCDE5
386    7.82538599923348465381e+03,  // 0x40BE9162D0D88419
387    1.33875336287249578163e+05,  // 0x4100579AB0B75E98
388    7.19657723683240939863e+05,  // 0x4125F65372869C19
389    6.66601232617776375264e+05,  // 0x412457D27719AD5C
390    -2.94490264303834643215e+05, // 0xC111F9690EA5AA18
391];
392
393// for x in [8,4.5454] = 1/[0.125,0.22001]
394const Q1R5: [f64; 6] = [
395    -2.08979931141764104297e-11, // 0xBDB6FA431AA1A098
396    -1.02539050241375426231e-01, // 0xBFBA3FFFCB597FEF
397    -8.05644828123936029840e+00, // 0xC0201CE6CA03AD4B
398    -1.83669607474888380239e+02, // 0xC066F56D6CA7B9B0
399    -1.37319376065508163265e+03, // 0xC09574C66931734F
400    -2.61244440453215656817e+03, // 0xC0A468E388FDA79D
401];
402const Q1S5: [f64; 6] = [
403    8.12765501384335777857e+01,  // 0x405451B2FF5A11B2
404    1.99179873460485964642e+03,  // 0x409F1F31E77BF839
405    1.74684851924908907677e+04,  // 0x40D10F1F0D64CE29
406    4.98514270910352279316e+04,  // 0x40E8576DAABAD197
407    2.79480751638918118260e+04,  // 0x40DB4B04CF7C364B
408    -4.71918354795128470869e+03, // 0xC0B26F2EFCFFA004
409];
410
411// for x in [4.5454,2.8571] = 1/[0.2199,0.35001] ???
412const Q1R3: [f64; 6] = [
413    -5.07831226461766561369e-09, // 0xBE35CFA9D38FC84F
414    -1.02537829820837089745e-01, // 0xBFBA3FEB51AEED54
415    -4.61011581139473403113e+00, // 0xC01270C23302D9FF
416    -5.78472216562783643212e+01, // 0xC04CEC71C25D16DA
417    -2.28244540737631695038e+02, // 0xC06C87D34718D55F
418    -2.19210128478909325622e+02, // 0xC06B66B95F5C1BF6
419];
420const Q1S3: [f64; 6] = [
421    4.76651550323729509273e+01,  // 0x4047D523CCD367E4
422    6.73865112676699709482e+02,  // 0x40850EEBC031EE3E
423    3.38015286679526343505e+03,  // 0x40AA684E448E7C9A
424    5.54772909720722782367e+03,  // 0x40B5ABBAA61D54A6
425    1.90311919338810798763e+03,  // 0x409DBC7A0DD4DF4B
426    -1.35201191444307340817e+02, // 0xC060E670290A311F
427];
428
429// for x in [2.8570,2] = 1/[0.3499,0.5]
430const Q1R2: [f64; 6] = [
431    -1.78381727510958865572e-07, // 0xBE87F12644C626D2
432    -1.02517042607985553460e-01, // 0xBFBA3E8E9148B010
433    -2.75220568278187460720e+00, // 0xC006048469BB4EDA
434    -1.96636162643703720221e+01, // 0xC033A9E2C168907F
435    -4.23253133372830490089e+01, // 0xC04529A3DE104AAA
436    -2.13719211703704061733e+01, // 0xC0355F3639CF6E52
437];
438const Q1S2: [f64; 6] = [
439    2.95333629060523854548e+01,  // 0x403D888A78AE64FF
440    2.52981549982190529136e+02,  // 0x406F9F68DB821CBA
441    7.57502834868645436472e+02,  // 0x4087AC05CE49A0F7
442    7.39393205320467245656e+02,  // 0x40871B2548D4C029
443    1.55949003336666123687e+02,  // 0x40637E5E3C3ED8D4
444    -4.95949898822628210127e+00, // 0xC013D686E71BE86B
445];
446
447fn qone(x: f64) -> f64 {
448    let (p, q) = if x >= 8.0 {
449        (&Q1R8, &Q1S8)
450    } else if x >= 4.5454 {
451        (&Q1R5, &Q1S5)
452    } else if x >= 2.8571 {
453        (&Q1R3, &Q1S3)
454    } else if x >= 2.0 {
455        (&Q1R2, &Q1S2)
456    } else {
457        panic!("INTERNAL ERROR: x must be ≥ 2.0 for qone");
458    };
459    let z = 1.0 / (x * x);
460    let r = p[0] + z * (p[1] + z * (p[2] + z * (p[3] + z * (p[4] + z * p[5]))));
461    let s = 1.0 + z * (q[0] + z * (q[1] + z * (q[2] + z * (q[3] + z * (q[4] + z * q[5])))));
462    (0.375 + r / s) / x
463}
464
465////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
466
467#[cfg(test)]
468mod tests {
469    use super::{bessel_j1, bessel_y1, pone, qone, TWO_129, TWO_M54};
470    use crate::{approx_eq, assert_alike};
471
472    #[test]
473    #[should_panic(expected = "INTERNAL ERROR: x must be ≥ 2.0 for pone")]
474    fn pone_panics_on_wrong_input() {
475        pone(1.99);
476    }
477
478    #[test]
479    #[should_panic(expected = "INTERNAL ERROR: x must be ≥ 2.0 for qone")]
480    fn qone_panics_on_wrong_input() {
481        qone(1.99);
482    }
483
484    #[test]
485    fn bessel_j1_handles_special_cases() {
486        assert!(bessel_j1(f64::NAN).is_nan());
487        assert_eq!(bessel_j1(f64::NEG_INFINITY), 0.0);
488        assert_eq!(bessel_j1(f64::INFINITY), 0.0);
489        assert_eq!(bessel_j1(0.0), 0.0);
490    }
491
492    #[test]
493    fn bessel_y1_handles_special_cases() {
494        assert!(bessel_y1(f64::NEG_INFINITY).is_nan());
495        assert!(bessel_y1(-0.01).is_nan());
496        assert!(bessel_y1(f64::NAN).is_nan());
497        assert_eq!(bessel_y1(f64::INFINITY), 0.0);
498        assert_eq!(bessel_y1(0.0), f64::NEG_INFINITY);
499    }
500
501    #[test]
502    fn bessel_j1_works() {
503        // Mathematica: N[BesselJ[1, -123], 100]
504        approx_eq(
505            bessel_j1(-123.0),
506            -0.02156735149890660940086328069036361009670280158394365333670828669644907165480535175323397489229874718,
507            1e-17,
508        );
509
510        // Mathematica: N[BesselJ[1, -5], 100]
511        assert_eq!(
512            bessel_j1(-5.0),
513            0.3275791375914652220377343219101691327608499046240540186864806450648753089914574086157808718640701942
514        );
515
516        // Mathematica: N[BesselJ[1, -2], 100]
517        approx_eq(
518            bessel_j1(-2.0),
519            -0.5767248077568733872024482422691370869203026897196754401211390207640871162896121849483995433063402923,
520            1e-15,
521        );
522
523        // Mathematica: N[BesselJ[1, -1], 100]
524        assert_eq!(
525            bessel_j1(-1.0),
526            -0.4400505857449335159596822037189149131273723019927652511367581717801382224780155479307965923811982542
527        );
528
529        // Mathematica: N[BesselJ[1, 10^-9], 100]
530        assert_eq!(
531            bessel_j1(1e-9),
532            4.999999999999999999375000000000000000026041666666666666666124131944444444444451226128472222222222166e-10
533        );
534
535        // Mathematica: N[BesselJ[1, 1], 100]
536        assert_eq!(
537            bessel_j1(1.0),
538            0.4400505857449335159596822037189149131273723019927652511367581717801382224780155479307965923811982542
539        );
540
541        // Mathematica: N[BesselJ[1, 2], 100]
542        approx_eq(
543            bessel_j1(2.0),
544            0.5767248077568733872024482422691370869203026897196754401211390207640871162896121849483995433063402923,
545            1e-15,
546        );
547
548        // Mathematica: N[BesselJ[1, 5], 100]
549        assert_eq!(
550            bessel_j1(5.0),
551            -0.3275791375914652220377343219101691327608499046240540186864806450648753089914574086157808718640701942
552        );
553
554        // Mathematica: N[BesselJ[1, 123], 100]
555        approx_eq(
556            bessel_j1(123.0),
557            0.02156735149890660940086328069036361009670280158394365333670828669644907165480535175323397489229874718,
558            1e-17,
559        );
560    }
561
562    #[test]
563    fn bessel_y1_works() {
564        // Mathematica: N[BesselY[0, 10^-9], 100]
565        approx_eq(
566            bessel_y1(1e-9),
567            -6.366197723675813498680125340511460635200799230620073016920009919302213806667149128128951250692787669e8,
568            1e-6,
569        );
570
571        // Mathematica: N[BesselY[0, 1], 100]
572        approx_eq(
573            bessel_y1(1.0),
574            -0.7812128213002887165471500000479648205499063907164446078438332461277843915385602167276292380048056346,
575            1e-17,
576        );
577
578        // Mathematica: N[BesselY[0, 2], 100]
579        approx_eq(
580            bessel_y1(2.0),
581            -0.1070324315409375468883707722774766366874808982350538605257945572313199774994906724960750929874561422,
582            1e-16,
583        );
584
585        // Mathematica: N[BesselY[0, 5], 100]
586        approx_eq(
587            bessel_y1(5.0),
588            0.1478631433912268448010506754880372053734652184104756133006976595086882659262171067506925531074677045,
589            1e-16,
590        );
591
592        // Mathematica: N[BesselY[0, 123], 100]
593        approx_eq(
594            bessel_y1(123.0),
595            0.06863489010254926652096999352892516176822669293542798648532852667873692449570610866356639175413101604,
596            1e-17,
597        );
598    }
599
600    #[test]
601    fn bessel_j1_edge_cases_work() {
602        //
603        // x > f64::MAX / 2.0
604        //
605        // println!("x = {:?}", (f64::MAX / 2.0)); // 8.988465674311579e307
606        // println!("j1(x) = {:?}", bessel_j1(8.988465674311579e307));
607        // Reference value from Go 1.22.1
608        approx_eq(bessel_j1(f64::MAX / 2.0), 5.936112522662019e-155, 1e-310);
609
610        //
611        // x > TWO_129
612        //
613        // Mathematica: N[BesselJ[1, 2  2^129], 100]
614        approx_eq(
615            bessel_j1(2.0 * TWO_129),
616            -2.148812412212001397545607680138144226557017053953123779067461809570040430645106137686812136124982505e-20,
617            1e-100,
618        );
619    }
620
621    #[test]
622    fn bessel_y1_edge_cases_work() {
623        //
624        // x > f64::MAX / 2.0
625        //
626        // println!("x = {:?}", (f64::MAX / 2.0)); // 8.988465674311579e307
627        // println!("y1(x) = {:?}", bessel_y1(8.988465674311579e307));
628        // Reference value from Go 1.22.1
629        approx_eq(bessel_y1(f64::MAX / 2.0), -5.965640685080747e-155, 1e-310);
630
631        //
632        // x <= TWO_M54
633        //
634        approx_eq(
635            bessel_y1(TWO_M54),
636            -1.146832227844531729100095246803011708838780681972379423086390888013625612852380226098068118499203402e16,
637            1e-100,
638        );
639    }
640
641    // The code below is based on all_test.go file from Go (1.22.1)
642
643    const VALUES: [f64; 10] = [
644        4.9790119248836735e+00,
645        7.7388724745781045e+00,
646        -2.7688005719200159e-01,
647        -5.0106036182710749e+00,
648        9.6362937071984173e+00,
649        2.9263772392439646e+00,
650        5.2290834314593066e+00,
651        2.7279399104360102e+00,
652        1.8253080916808550e+00,
653        -8.6859247685756013e+00,
654    ];
655
656    const SOLUTION_J1: [f64; 10] = [
657        -3.251526395295203422162967e-01,
658        1.893581711430515718062564e-01,
659        -1.3711761352467242914491514e-01,
660        3.287486536269617297529617e-01,
661        1.3133899188830978473849215e-01,
662        3.660243417832986825301766e-01,
663        -3.4436769271848174665420672e-01,
664        4.329481396640773768835036e-01,
665        5.8181350531954794639333955e-01,
666        -2.7030574577733036112996607e-01,
667    ];
668
669    const SOLUTION_Y1: [f64; 10] = [
670        0.15494213737457922210218611,
671        -0.2165955142081145245075746,
672        -2.4644949631241895201032829,
673        0.1442740489541836405154505,
674        0.2215379960518984777080163,
675        0.3038800915160754150565448,
676        0.0691107642452362383808547,
677        0.2380116417809914424860165,
678        -0.20849492979459761009678934,
679        0.0242503179793232308250804,
680    ];
681
682    const SC_VALUES: [f64; 4] = [f64::NEG_INFINITY, 0.0, f64::INFINITY, f64::NAN];
683
684    const SC_SOLUTION_J1: [f64; 4] = [0.0, 0.0, 0.0, f64::NAN];
685
686    const SC_SOLUTION_Y1: [f64; 5] = [f64::NAN, f64::NEG_INFINITY, 0.0, f64::NAN, f64::NAN];
687
688    #[test]
689    fn test_bessel_j1() {
690        for (i, v) in VALUES.iter().enumerate() {
691            let f = bessel_j1(*v);
692            approx_eq(SOLUTION_J1[i], f, 1e-15);
693        }
694        for (i, v) in SC_VALUES.iter().enumerate() {
695            let f = bessel_j1(*v);
696            assert_alike(SC_SOLUTION_J1[i], f);
697        }
698    }
699
700    #[test]
701    fn test_bessel_y1() {
702        for (i, v) in VALUES.iter().enumerate() {
703            let a = f64::abs(*v);
704            let f = bessel_y1(a);
705            approx_eq(SOLUTION_Y1[i], f, 1e-15);
706        }
707        for (i, v) in SC_VALUES.iter().enumerate() {
708            let f = bessel_y1(*v);
709            assert_alike(SC_SOLUTION_Y1[i], f);
710        }
711    }
712}