jv.c 15 KB

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