]> Git Repo - uclibc-ng.git/blame - libm/double/jv.c
Seems we were lacking an acos() implementation
[uclibc-ng.git] / libm / double / jv.c
CommitLineData
1077fa4d
EA
1/* jv.c
2 *
3 * Bessel function of noninteger order
4 *
5 *
6 *
7 * SYNOPSIS:
8 *
9 * double v, x, y, jv();
10 *
11 * y = jv( v, x );
12 *
13 *
14 *
15 * DESCRIPTION:
16 *
17 * Returns Bessel function of order v of the argument,
18 * where v is real. Negative x is allowed if v is an integer.
19 *
20 * Several expansions are included: the ascending power
21 * series, the Hankel expansion, and two transitional
22 * expansions for large v. If v is not too large, it
23 * is reduced by recurrence to a region of best accuracy.
24 * The transitional expansions give 12D accuracy for v > 500.
25 *
26 *
27 *
28 * ACCURACY:
29 * Results for integer v are indicated by *, where x and v
30 * both vary from -125 to +125. Otherwise,
31 * x ranges from 0 to 125, v ranges as indicated by "domain."
32 * Error criterion is absolute, except relative when |jv()| > 1.
33 *
34 * arithmetic v domain x domain # trials peak rms
35 * IEEE 0,125 0,125 100000 4.6e-15 2.2e-16
36 * IEEE -125,0 0,125 40000 5.4e-11 3.7e-13
37 * IEEE 0,500 0,500 20000 4.4e-15 4.0e-16
38 * Integer v:
39 * IEEE -125,125 -125,125 50000 3.5e-15* 1.9e-16*
40 *
41 */
42\f
43
44/*
45Cephes Math Library Release 2.8: June, 2000
46Copyright 1984, 1987, 1989, 1992, 2000 by Stephen L. Moshier
47*/
48
49
50#include <math.h>
51#define DEBUG 0
52
53#ifdef DEC
54#define MAXGAM 34.84425627277176174
55#else
56#define MAXGAM 171.624376956302725
57#endif
58
59#ifdef ANSIPROT
60extern int airy ( double, double *, double *, double *, double * );
61extern double fabs ( double );
62extern double floor ( double );
63extern double frexp ( double, int * );
64extern double polevl ( double, void *, int );
65extern double j0 ( double );
66extern double j1 ( double );
67extern double sqrt ( double );
68extern double cbrt ( double );
69extern double exp ( double );
70extern double log ( double );
71extern double sin ( double );
72extern double cos ( double );
73extern double acos ( double );
74extern double pow ( double, double );
75extern double gamma ( double );
76extern double lgam ( double );
77static double recur(double *, double, double *, int);
78static double jvs(double, double);
79static double hankel(double, double);
80static double jnx(double, double);
81static double jnt(double, double);
82#else
83int airy();
84double fabs(), floor(), frexp(), polevl(), j0(), j1(), sqrt(), cbrt();
85double exp(), log(), sin(), cos(), acos(), pow(), gamma(), lgam();
86static double recur(), jvs(), hankel(), jnx(), jnt();
87#endif
88
89extern double MAXNUM, MACHEP, MINLOG, MAXLOG;
90#define BIG 1.44115188075855872E+17
91
92double jv( n, x )
93double n, x;
94{
95double k, q, t, y, an;
96int i, sign, nint;
97
98nint = 0; /* Flag for integer n */
99sign = 1; /* Flag for sign inversion */
100an = fabs( n );
101y = floor( an );
102if( y == an )
103 {
104 nint = 1;
105 i = an - 16384.0 * floor( an/16384.0 );
106 if( n < 0.0 )
107 {
108 if( i & 1 )
109 sign = -sign;
110 n = an;
111 }
112 if( x < 0.0 )
113 {
114 if( i & 1 )
115 sign = -sign;
116 x = -x;
117 }
118 if( n == 0.0 )
119 return( j0(x) );
120 if( n == 1.0 )
121 return( sign * j1(x) );
122 }
123
124if( (x < 0.0) && (y != an) )
125 {
126 mtherr( "Jv", DOMAIN );
127 y = 0.0;
128 goto done;
129 }
130
131y = fabs(x);
132
133if( y < MACHEP )
134 goto underf;
135
136k = 3.6 * sqrt(y);
137t = 3.6 * sqrt(an);
138if( (y < t) && (an > 21.0) )
139 return( sign * jvs(n,x) );
140if( (an < k) && (y > 21.0) )
141 return( sign * hankel(n,x) );
142
143if( an < 500.0 )
144 {
145/* Note: if x is too large, the continued
146 * fraction will fail; but then the
147 * Hankel expansion can be used.
148 */
149 if( nint != 0 )
150 {
151 k = 0.0;
152 q = recur( &n, x, &k, 1 );
153 if( k == 0.0 )
154 {
155 y = j0(x)/q;
156 goto done;
157 }
158 if( k == 1.0 )
159 {
160 y = j1(x)/q;
161 goto done;
162 }
163 }
164
165if( an > 2.0 * y )
166 goto rlarger;
167
168 if( (n >= 0.0) && (n < 20.0)
169 && (y > 6.0) && (y < 20.0) )
170 {
171/* Recur backwards from a larger value of n
172 */
173rlarger:
174 k = n;
175
176 y = y + an + 1.0;
177 if( y < 30.0 )
178 y = 30.0;
179 y = n + floor(y-n);
180 q = recur( &y, x, &k, 0 );
181 y = jvs(y,x) * q;
182 goto done;
183 }
184
185 if( k <= 30.0 )
186 {
187 k = 2.0;
188 }
189 else if( k < 90.0 )
190 {
191 k = (3*k)/4;
192 }
193 if( an > (k + 3.0) )
194 {
195 if( n < 0.0 )
196 k = -k;
197 q = n - floor(n);
198 k = floor(k) + q;
199 if( n > 0.0 )
200 q = recur( &n, x, &k, 1 );
201 else
202 {
203 t = k;
204 k = n;
205 q = recur( &t, x, &k, 1 );
206 k = t;
207 }
208 if( q == 0.0 )
209 {
210underf:
211 y = 0.0;
212 goto done;
213 }
214 }
215 else
216 {
217 k = n;
218 q = 1.0;
219 }
220
221/* boundary between convergence of
222 * power series and Hankel expansion
223 */
224 y = fabs(k);
225 if( y < 26.0 )
226 t = (0.0083*y + 0.09)*y + 12.9;
227 else
228 t = 0.9 * y;
229
230 if( x > t )
231 y = hankel(k,x);
232 else
233 y = jvs(k,x);
234#if DEBUG
235printf( "y = %.16e, recur q = %.16e\n", y, q );
236#endif
237 if( n > 0.0 )
238 y /= q;
239 else
240 y *= q;
241 }
242
243else
244 {
245/* For large n, use the uniform expansion
246 * or the transitional expansion.
247 * But if x is of the order of n**2,
248 * these may blow up, whereas the
249 * Hankel expansion will then work.
250 */
251 if( n < 0.0 )
252 {
253 mtherr( "Jv", TLOSS );
254 y = 0.0;
255 goto done;
256 }
257 t = x/n;
258 t /= n;
259 if( t > 0.3 )
260 y = hankel(n,x);
261 else
262 y = jnx(n,x);
263 }
264
265done: return( sign * y);
266}
267\f
268/* Reduce the order by backward recurrence.
269 * AMS55 #9.1.27 and 9.1.73.
270 */
271
272static double recur( n, x, newn, cancel )
273double *n;
274double x;
275double *newn;
276int cancel;
277{
278double pkm2, pkm1, pk, qkm2, qkm1;
279/* double pkp1; */
280double k, ans, qk, xk, yk, r, t, kf;
281static double big = BIG;
282int nflag, ctr;
283
284/* continued fraction for Jn(x)/Jn-1(x) */
285if( *n < 0.0 )
286 nflag = 1;
287else
288 nflag = 0;
289
290fstart:
291
292#if DEBUG
293printf( "recur: n = %.6e, newn = %.6e, cfrac = ", *n, *newn );
294#endif
295
296pkm2 = 0.0;
297qkm2 = 1.0;
298pkm1 = x;
299qkm1 = *n + *n;
300xk = -x * x;
301yk = qkm1;
302ans = 1.0;
303ctr = 0;
304do
305 {
306 yk += 2.0;
307 pk = pkm1 * yk + pkm2 * xk;
308 qk = qkm1 * yk + qkm2 * xk;
309 pkm2 = pkm1;
310 pkm1 = pk;
311 qkm2 = qkm1;
312 qkm1 = qk;
313 if( qk != 0 )
314 r = pk/qk;
315 else
316 r = 0.0;
317 if( r != 0 )
318 {
319 t = fabs( (ans - r)/r );
320 ans = r;
321 }
322 else
323 t = 1.0;
324
325 if( ++ctr > 1000 )
326 {
327 mtherr( "jv", UNDERFLOW );
328 goto done;
329 }
330 if( t < MACHEP )
331 goto done;
332
333 if( fabs(pk) > big )
334 {
335 pkm2 /= big;
336 pkm1 /= big;
337 qkm2 /= big;
338 qkm1 /= big;
339 }
340 }
341while( t > MACHEP );
342
343done:
344
345#if DEBUG
346printf( "%.6e\n", ans );
347#endif
348
349/* Change n to n-1 if n < 0 and the continued fraction is small
350 */
351if( nflag > 0 )
352 {
353 if( fabs(ans) < 0.125 )
354 {
355 nflag = -1;
356 *n = *n - 1.0;
357 goto fstart;
358 }
359 }
360
361
362kf = *newn;
363
364/* backward recurrence
365 * 2k
366 * J (x) = --- J (x) - J (x)
367 * k-1 x k k+1
368 */
369
370pk = 1.0;
371pkm1 = 1.0/ans;
372k = *n - 1.0;
373r = 2 * k;
374do
375 {
376 pkm2 = (pkm1 * r - pk * x) / x;
377 /* pkp1 = pk; */
378 pk = pkm1;
379 pkm1 = pkm2;
380 r -= 2.0;
381/*
382 t = fabs(pkp1) + fabs(pk);
383 if( (k > (kf + 2.5)) && (fabs(pkm1) < 0.25*t) )
384 {
385 k -= 1.0;
386 t = x*x;
387 pkm2 = ( (r*(r+2.0)-t)*pk - r*x*pkp1 )/t;
388 pkp1 = pk;
389 pk = pkm1;
390 pkm1 = pkm2;
391 r -= 2.0;
392 }
393*/
394 k -= 1.0;
395 }
396while( k > (kf + 0.5) );
397
398/* Take the larger of the last two iterates
399 * on the theory that it may have less cancellation error.
400 */
401
402if( cancel )
403 {
404 if( (kf >= 0.0) && (fabs(pk) > fabs(pkm1)) )
405 {
406 k += 1.0;
407 pkm2 = pk;
408 }
409 }
410*newn = k;
411#if DEBUG
412printf( "newn %.6e rans %.6e\n", k, pkm2 );
413#endif
414return( pkm2 );
415}
416
417
418
419/* Ascending power series for Jv(x).
420 * AMS55 #9.1.10.
421 */
422
423extern double PI;
424extern int sgngam;
425
426static double jvs( n, x )
427double n, x;
428{
429double t, u, y, z, k;
430int ex;
431
432z = -x * x / 4.0;
433u = 1.0;
434y = u;
435k = 1.0;
436t = 1.0;
437
438while( t > MACHEP )
439 {
440 u *= z / (k * (n+k));
441 y += u;
442 k += 1.0;
443 if( y != 0 )
444 t = fabs( u/y );
445 }
446#if DEBUG
447printf( "power series=%.5e ", y );
448#endif
449t = frexp( 0.5*x, &ex );
450ex = ex * n;
451if( (ex > -1023)
452 && (ex < 1023)
453 && (n > 0.0)
454 && (n < (MAXGAM-1.0)) )
455 {
456 t = pow( 0.5*x, n ) / gamma( n + 1.0 );
457#if DEBUG
458printf( "pow(.5*x, %.4e)/gamma(n+1)=%.5e\n", n, t );
459#endif
460 y *= t;
461 }
462else
463 {
464#if DEBUG
465 z = n * log(0.5*x);
466 k = lgam( n+1.0 );
467 t = z - k;
468 printf( "log pow=%.5e, lgam(%.4e)=%.5e\n", z, n+1.0, k );
469#else
470 t = n * log(0.5*x) - lgam(n + 1.0);
471#endif
472 if( y < 0 )
473 {
474 sgngam = -sgngam;
475 y = -y;
476 }
477 t += log(y);
478#if DEBUG
479printf( "log y=%.5e\n", log(y) );
480#endif
481 if( t < -MAXLOG )
482 {
483 return( 0.0 );
484 }
485 if( t > MAXLOG )
486 {
487 mtherr( "Jv", OVERFLOW );
488 return( MAXNUM );
489 }
490 y = sgngam * exp( t );
491 }
492return(y);
493}
494\f
495/* Hankel's asymptotic expansion
496 * for large x.
497 * AMS55 #9.2.5.
498 */
499
500static double hankel( n, x )
501double n, x;
502{
503double t, u, z, k, sign, conv;
504double p, q, j, m, pp, qq;
505int flag;
506
507m = 4.0*n*n;
508j = 1.0;
509z = 8.0 * x;
510k = 1.0;
511p = 1.0;
512u = (m - 1.0)/z;
513q = u;
514sign = 1.0;
515conv = 1.0;
516flag = 0;
517t = 1.0;
518pp = 1.0e38;
519qq = 1.0e38;
520
521while( t > MACHEP )
522 {
523 k += 2.0;
524 j += 1.0;
525 sign = -sign;
526 u *= (m - k * k)/(j * z);
527 p += sign * u;
528 k += 2.0;
529 j += 1.0;
530 u *= (m - k * k)/(j * z);
531 q += sign * u;
532 t = fabs(u/p);
533 if( t < conv )
534 {
535 conv = t;
536 qq = q;
537 pp = p;
538 flag = 1;
539 }
540/* stop if the terms start getting larger */
541 if( (flag != 0) && (t > conv) )
542 {
543#if DEBUG
544 printf( "Hankel: convergence to %.4E\n", conv );
545#endif
546 goto hank1;
547 }
548 }
549
550hank1:
551u = x - (0.5*n + 0.25) * PI;
552t = sqrt( 2.0/(PI*x) ) * ( pp * cos(u) - qq * sin(u) );
553#if DEBUG
554printf( "hank: %.6e\n", t );
555#endif
556return( t );
557}
558\f
559
560/* Asymptotic expansion for large n.
561 * AMS55 #9.3.35.
562 */
563
564static double lambda[] = {
565 1.0,
566 1.041666666666666666666667E-1,
567 8.355034722222222222222222E-2,
568 1.282265745563271604938272E-1,
569 2.918490264641404642489712E-1,
570 8.816272674437576524187671E-1,
571 3.321408281862767544702647E+0,
572 1.499576298686255465867237E+1,
573 7.892301301158651813848139E+1,
574 4.744515388682643231611949E+2,
575 3.207490090890661934704328E+3
576};
577static double mu[] = {
578 1.0,
579 -1.458333333333333333333333E-1,
580 -9.874131944444444444444444E-2,
581 -1.433120539158950617283951E-1,
582 -3.172272026784135480967078E-1,
583 -9.424291479571202491373028E-1,
584 -3.511203040826354261542798E+0,
585 -1.572726362036804512982712E+1,
586 -8.228143909718594444224656E+1,
587 -4.923553705236705240352022E+2,
588 -3.316218568547972508762102E+3
589};
590static double P1[] = {
591 -2.083333333333333333333333E-1,
592 1.250000000000000000000000E-1
593};
594static double P2[] = {
595 3.342013888888888888888889E-1,
596 -4.010416666666666666666667E-1,
597 7.031250000000000000000000E-2
598};
599static double P3[] = {
600 -1.025812596450617283950617E+0,
601 1.846462673611111111111111E+0,
602 -8.912109375000000000000000E-1,
603 7.324218750000000000000000E-2
604};
605static double P4[] = {
606 4.669584423426247427983539E+0,
607 -1.120700261622299382716049E+1,
608 8.789123535156250000000000E+0,
609 -2.364086914062500000000000E+0,
610 1.121520996093750000000000E-1
611};
612static double P5[] = {
613 -2.8212072558200244877E1,
614 8.4636217674600734632E1,
615 -9.1818241543240017361E1,
616 4.2534998745388454861E1,
617 -7.3687943594796316964E0,
618 2.27108001708984375E-1
619};
620static double P6[] = {
621 2.1257013003921712286E2,
622 -7.6525246814118164230E2,
623 1.0599904525279998779E3,
624 -6.9957962737613254123E2,
625 2.1819051174421159048E2,
626 -2.6491430486951555525E1,
627 5.7250142097473144531E-1
628};
629static double P7[] = {
630 -1.9194576623184069963E3,
631 8.0617221817373093845E3,
632 -1.3586550006434137439E4,
633 1.1655393336864533248E4,
634 -5.3056469786134031084E3,
635 1.2009029132163524628E3,
636 -1.0809091978839465550E2,
637 1.7277275025844573975E0
638};
639
640
641static double jnx( n, x )
642double n, x;
643{
644double zeta, sqz, zz, zp, np;
645double cbn, n23, t, z, sz;
646double pp, qq, z32i, zzi;
647double ak, bk, akl, bkl;
648int sign, doa, dob, nflg, k, s, tk, tkp1, m;
649static double u[8];
650static double ai, aip, bi, bip;
651
652/* Test for x very close to n.
653 * Use expansion for transition region if so.
654 */
655cbn = cbrt(n);
656z = (x - n)/cbn;
657if( fabs(z) <= 0.7 )
658 return( jnt(n,x) );
659
660z = x/n;
661zz = 1.0 - z*z;
662if( zz == 0.0 )
663 return(0.0);
664
665if( zz > 0.0 )
666 {
667 sz = sqrt( zz );
668 t = 1.5 * (log( (1.0+sz)/z ) - sz ); /* zeta ** 3/2 */
669 zeta = cbrt( t * t );
670 nflg = 1;
671 }
672else
673 {
674 sz = sqrt(-zz);
675 t = 1.5 * (sz - acos(1.0/z));
676 zeta = -cbrt( t * t );
677 nflg = -1;
678 }
679z32i = fabs(1.0/t);
680sqz = cbrt(t);
681
682/* Airy function */
683n23 = cbrt( n * n );
684t = n23 * zeta;
685
686#if DEBUG
687printf("zeta %.5E, Airy(%.5E)\n", zeta, t );
688#endif
689airy( t, &ai, &aip, &bi, &bip );
690
691/* polynomials in expansion */
692u[0] = 1.0;
693zzi = 1.0/zz;
694u[1] = polevl( zzi, P1, 1 )/sz;
695u[2] = polevl( zzi, P2, 2 )/zz;
696u[3] = polevl( zzi, P3, 3 )/(sz*zz);
697pp = zz*zz;
698u[4] = polevl( zzi, P4, 4 )/pp;
699u[5] = polevl( zzi, P5, 5 )/(pp*sz);
700pp *= zz;
701u[6] = polevl( zzi, P6, 6 )/pp;
702u[7] = polevl( zzi, P7, 7 )/(pp*sz);
703
704#if DEBUG
705for( k=0; k<=7; k++ )
706 printf( "u[%d] = %.5E\n", k, u[k] );
707#endif
708
709pp = 0.0;
710qq = 0.0;
711np = 1.0;
712/* flags to stop when terms get larger */
713doa = 1;
714dob = 1;
715akl = MAXNUM;
716bkl = MAXNUM;
717
718for( k=0; k<=3; k++ )
719 {
720 tk = 2 * k;
721 tkp1 = tk + 1;
722 zp = 1.0;
723 ak = 0.0;
724 bk = 0.0;
725 for( s=0; s<=tk; s++ )
726 {
727 if( doa )
728 {
729 if( (s & 3) > 1 )
730 sign = nflg;
731 else
732 sign = 1;
733 ak += sign * mu[s] * zp * u[tk-s];
734 }
735
736 if( dob )
737 {
738 m = tkp1 - s;
739 if( ((m+1) & 3) > 1 )
740 sign = nflg;
741 else
742 sign = 1;
743 bk += sign * lambda[s] * zp * u[m];
744 }
745 zp *= z32i;
746 }
747
748 if( doa )
749 {
750 ak *= np;
751 t = fabs(ak);
752 if( t < akl )
753 {
754 akl = t;
755 pp += ak;
756 }
757 else
758 doa = 0;
759 }
760
761 if( dob )
762 {
763 bk += lambda[tkp1] * zp * u[0];
764 bk *= -np/sqz;
765 t = fabs(bk);
766 if( t < bkl )
767 {
768 bkl = t;
769 qq += bk;
770 }
771 else
772 dob = 0;
773 }
774#if DEBUG
775 printf("a[%d] %.5E, b[%d] %.5E\n", k, ak, k, bk );
776#endif
777 if( np < MACHEP )
778 break;
779 np /= n*n;
780 }
781
782/* normalizing factor ( 4*zeta/(1 - z**2) )**1/4 */
783t = 4.0 * zeta/zz;
784t = sqrt( sqrt(t) );
785
786t *= ai*pp/cbrt(n) + aip*qq/(n23*n);
787return(t);
788}
789\f
790/* Asymptotic expansion for transition region,
791 * n large and x close to n.
792 * AMS55 #9.3.23.
793 */
794
795static double PF2[] = {
796 -9.0000000000000000000e-2,
797 8.5714285714285714286e-2
798};
799static double PF3[] = {
800 1.3671428571428571429e-1,
801 -5.4920634920634920635e-2,
802 -4.4444444444444444444e-3
803};
804static double PF4[] = {
805 1.3500000000000000000e-3,
806 -1.6036054421768707483e-1,
807 4.2590187590187590188e-2,
808 2.7330447330447330447e-3
809};
810static double PG1[] = {
811 -2.4285714285714285714e-1,
812 1.4285714285714285714e-2
813};
814static double PG2[] = {
815 -9.0000000000000000000e-3,
816 1.9396825396825396825e-1,
817 -1.1746031746031746032e-2
818};
819static double PG3[] = {
820 1.9607142857142857143e-2,
821 -1.5983694083694083694e-1,
822 6.3838383838383838384e-3
823};
824
825
826static double jnt( n, x )
827double n, x;
828{
829double z, zz, z3;
830double cbn, n23, cbtwo;
831double ai, aip, bi, bip; /* Airy functions */
832double nk, fk, gk, pp, qq;
833double F[5], G[4];
834int k;
835
836cbn = cbrt(n);
837z = (x - n)/cbn;
838cbtwo = cbrt( 2.0 );
839
840/* Airy function */
841zz = -cbtwo * z;
842airy( zz, &ai, &aip, &bi, &bip );
843
844/* polynomials in expansion */
845zz = z * z;
846z3 = zz * z;
847F[0] = 1.0;
848F[1] = -z/5.0;
849F[2] = polevl( z3, PF2, 1 ) * zz;
850F[3] = polevl( z3, PF3, 2 );
851F[4] = polevl( z3, PF4, 3 ) * z;
852G[0] = 0.3 * zz;
853G[1] = polevl( z3, PG1, 1 );
854G[2] = polevl( z3, PG2, 2 ) * z;
855G[3] = polevl( z3, PG3, 2 ) * zz;
856#if DEBUG
857for( k=0; k<=4; k++ )
858 printf( "F[%d] = %.5E\n", k, F[k] );
859for( k=0; k<=3; k++ )
860 printf( "G[%d] = %.5E\n", k, G[k] );
861#endif
862pp = 0.0;
863qq = 0.0;
864nk = 1.0;
865n23 = cbrt( n * n );
866
867for( k=0; k<=4; k++ )
868 {
869 fk = F[k]*nk;
870 pp += fk;
871 if( k != 4 )
872 {
873 gk = G[k]*nk;
874 qq += gk;
875 }
876#if DEBUG
877 printf("fk[%d] %.5E, gk[%d] %.5E\n", k, fk, k, gk );
878#endif
879 nk /= n23;
880 }
881
882fk = cbtwo * ai * pp/cbn + cbrt(4.0) * aip * qq/n;
883return(fk);
884}
This page took 0.118572 seconds and 4 git commands to generate.