xref: /csrg-svn/lib/libm/common_source/j0.c (revision 57151)
148402Sbostic /*-
2*57151Sbostic  * Copyright (c) 1992 The Regents of the University of California.
348402Sbostic  * All rights reserved.
448402Sbostic  *
5*57151Sbostic  * %sccs.include.redist.c%
634117Sbostic  */
734117Sbostic 
824597Szliu #ifndef lint
9*57151Sbostic static char sccsid[] = "@(#)j0.c	5.5 (Berkeley) 12/16/92";
1034117Sbostic #endif /* not lint */
1124597Szliu 
1224597Szliu /*
13*57151Sbostic  * 16 December 1992
14*57151Sbostic  * Minor modifications by Peter McIlroy to adapt non-IEEE architecture.
15*57151Sbostic  */
1624597Szliu 
17*57151Sbostic /*
18*57151Sbostic  * ====================================================
19*57151Sbostic  * Copyright (C) 1992 by Sun Microsystems, Inc.
20*57151Sbostic  *
21*57151Sbostic  * Developed at SunPro, a Sun Microsystems, Inc. business.
22*57151Sbostic  * Permission to use, copy, modify, and distribute this
23*57151Sbostic  * software is freely granted, provided that this notice
24*57151Sbostic  * is preserved.
25*57151Sbostic  * ====================================================
26*57151Sbostic  *
27*57151Sbostic  * ******************* WARNING ********************
28*57151Sbostic  * This is an alpha version of SunPro's FDLIBM (Freely
29*57151Sbostic  * Distributable Math Library) for IEEE double precision
30*57151Sbostic  * arithmetic. FDLIBM is a basic math library written
31*57151Sbostic  * in C that runs on machines that conform to IEEE
32*57151Sbostic  * Standard 754/854. This alpha version is distributed
33*57151Sbostic  * for testing purpose. Those who use this software
34*57151Sbostic  * should report any bugs to
35*57151Sbostic  *
36*57151Sbostic  *		fdlibm-comments@sunpro.eng.sun.com
37*57151Sbostic  *
38*57151Sbostic  * -- K.C. Ng, Oct 12, 1992
39*57151Sbostic  * ************************************************
40*57151Sbostic  */
4124597Szliu 
42*57151Sbostic /* double j0(double x), y0(double x)
43*57151Sbostic  * Bessel function of the first and second kinds of order zero.
44*57151Sbostic  * Method -- j0(x):
45*57151Sbostic  *	1. For tiny x, we use j0(x) = 1 - x^2/4 + x^4/64 - ...
46*57151Sbostic  *	2. Reduce x to |x| since j0(x)=j0(-x),  and
47*57151Sbostic  *	   for x in (0,2)
48*57151Sbostic  *		j0(x) = 1-z/4+ z^2*R0/S0,  where z = x*x;
49*57151Sbostic  *	   (precision:  |j0-1+z/4-z^2R0/S0 |<2**-63.67 )
50*57151Sbostic  *	   for x in (2,inf)
51*57151Sbostic  * 		j0(x) = sqrt(2/(pi*x))*(p0(x)*cos(x0)-q0(x)*sin(x0))
52*57151Sbostic  * 	   where x0 = x-pi/4. It is better to compute sin(x0),cos(x0)
53*57151Sbostic  *	   as follow:
54*57151Sbostic  *		cos(x0) = cos(x)cos(pi/4)+sin(x)sin(pi/4)
55*57151Sbostic  *			= 1/sqrt(2) * (cos(x) + sin(x))
56*57151Sbostic  *		sin(x0) = sin(x)cos(pi/4)-cos(x)sin(pi/4)
57*57151Sbostic  *			= 1/sqrt(2) * (sin(x) - cos(x))
58*57151Sbostic  * 	   (To avoid cancellation, use
59*57151Sbostic  *		sin(x) +- cos(x) = -cos(2x)/(sin(x) -+ cos(x))
60*57151Sbostic  * 	    to compute the worse one.)
61*57151Sbostic  *
62*57151Sbostic  *	3 Special cases
63*57151Sbostic  *		j0(nan)= nan
64*57151Sbostic  *		j0(0) = 1
65*57151Sbostic  *		j0(inf) = 0
66*57151Sbostic  *
67*57151Sbostic  * Method -- y0(x):
68*57151Sbostic  *	1. For x<2.
69*57151Sbostic  *	   Since
70*57151Sbostic  *		y0(x) = 2/pi*(j0(x)*(ln(x/2)+Euler) + x^2/4 - ...)
71*57151Sbostic  *	   therefore y0(x)-2/pi*j0(x)*ln(x) is an even function.
72*57151Sbostic  *	   We use the following function to approximate y0,
73*57151Sbostic  *		y0(x) = U(z)/V(z) + (2/pi)*(j0(x)*ln(x)), z= x^2
74*57151Sbostic  *	   where
75*57151Sbostic  *		U(z) = u0 + u1*z + ... + u6*z^6
76*57151Sbostic  *		V(z) = 1  + v1*z + ... + v4*z^4
77*57151Sbostic  *	   with absolute approximation error bounded by 2**-72.
78*57151Sbostic  *	   Note: For tiny x, U/V = u0 and j0(x)~1, hence
79*57151Sbostic  *		y0(tiny) = u0 + (2/pi)*ln(tiny), (choose tiny<2**-27)
80*57151Sbostic  *	2. For x>=2.
81*57151Sbostic  * 		y0(x) = sqrt(2/(pi*x))*(p0(x)*cos(x0)+q0(x)*sin(x0))
82*57151Sbostic  * 	   where x0 = x-pi/4. It is better to compute sin(x0),cos(x0)
83*57151Sbostic  *	   by the method mentioned above.
84*57151Sbostic  *	3. Special cases: y0(0)=-inf, y0(x<0)=NaN, y0(inf)=0.
85*57151Sbostic  */
8624597Szliu 
87*57151Sbostic #include <math.h>
88*57151Sbostic #include <float.h>
89*57151Sbostic #if defined(vax) || defined(tahoe)
90*57151Sbostic #define _IEEE	0
91*57151Sbostic #else
92*57151Sbostic #define _IEEE	1
93*57151Sbostic #define infnan(x) (0.0)
94*57151Sbostic #endif
9524597Szliu 
96*57151Sbostic static double pzero __P((double)), qzero __P((double));
9724597Szliu 
98*57151Sbostic static double
99*57151Sbostic huge 	= 1e300,
100*57151Sbostic zero    = 0.0,
101*57151Sbostic one	= 1.0,
102*57151Sbostic invsqrtpi= 5.641895835477562869480794515607725858441e-0001,
103*57151Sbostic tpi	= 0.636619772367581343075535053490057448,
104*57151Sbostic  		/* R0/S0 on [0, 2.00] */
105*57151Sbostic r02 =   1.562499999999999408594634421055018003102e-0002,
106*57151Sbostic r03 =  -1.899792942388547334476601771991800712355e-0004,
107*57151Sbostic r04 =   1.829540495327006565964161150603950916854e-0006,
108*57151Sbostic r05 =  -4.618326885321032060803075217804816988758e-0009,
109*57151Sbostic s01 =   1.561910294648900170180789369288114642057e-0002,
110*57151Sbostic s02 =   1.169267846633374484918570613449245536323e-0004,
111*57151Sbostic s03 =   5.135465502073181376284426245689510134134e-0007,
112*57151Sbostic s04 =   1.166140033337900097836930825478674320464e-0009;
11324597Szliu 
114*57151Sbostic double
115*57151Sbostic j0(x)
116*57151Sbostic 	double x;
117*57151Sbostic {
118*57151Sbostic 	double z, s,c,ss,cc,r,u,v;
11924597Szliu 
120*57151Sbostic 	if (!finite(x))
121*57151Sbostic 		if (_IEEE) return one/(x*x);
122*57151Sbostic 		else return (0);
123*57151Sbostic 	x = fabs(x);
124*57151Sbostic 	if (x >= 2.0) {	/* |x| >= 2.0 */
125*57151Sbostic 		s = sin(x);
126*57151Sbostic 		c = cos(x);
127*57151Sbostic 		ss = s-c;
128*57151Sbostic 		cc = s+c;
129*57151Sbostic 		if (x < .5 * DBL_MAX) {  /* make sure x+x not overflow */
130*57151Sbostic 		    z = -cos(x+x);
131*57151Sbostic 		    if ((s*c)<zero) cc = z/ss;
132*57151Sbostic 		    else 	    ss = z/cc;
133*57151Sbostic 		}
134*57151Sbostic 	/*
135*57151Sbostic 	 * j0(x) = 1/sqrt(pi) * (P(0,x)*cc - Q(0,x)*ss) / sqrt(x)
136*57151Sbostic 	 * y0(x) = 1/sqrt(pi) * (P(0,x)*ss + Q(0,x)*cc) / sqrt(x)
137*57151Sbostic 	 */
138*57151Sbostic 		if (_IEEE && x> 6.80564733841876927e+38) /* 2^129 */
139*57151Sbostic 			z = (invsqrtpi*cc)/sqrt(x);
140*57151Sbostic 		else {
141*57151Sbostic 		    u = pzero(x); v = qzero(x);
142*57151Sbostic 		    z = invsqrtpi*(u*cc-v*ss)/sqrt(x);
143*57151Sbostic 		}
144*57151Sbostic 		return z;
145*57151Sbostic 	}
146*57151Sbostic 	if (x < 1.220703125e-004) {		   /* |x| < 2**-13 */
147*57151Sbostic 	    if (huge+x > one) {			   /* raise inexact if x != 0 */
148*57151Sbostic 	        if (x < 7.450580596923828125e-009) /* |x|<2**-27 */
149*57151Sbostic 			return one;
150*57151Sbostic 	        else return (one - 0.25*x*x);
151*57151Sbostic 	    }
152*57151Sbostic 	}
153*57151Sbostic 	z = x*x;
154*57151Sbostic 	r =  z*(r02+z*(r03+z*(r04+z*r05)));
155*57151Sbostic 	s =  one+z*(s01+z*(s02+z*(s03+z*s04)));
156*57151Sbostic 	if (x < one) {			/* |x| < 1.00 */
157*57151Sbostic 	    return (one + z*(-0.25+(r/s)));
158*57151Sbostic 	} else {
159*57151Sbostic 	    u = 0.5*x;
160*57151Sbostic 	    return ((one+u)*(one-u)+z*(r/s));
161*57151Sbostic 	}
162*57151Sbostic }
16324597Szliu 
164*57151Sbostic static double
165*57151Sbostic u00 =  -7.380429510868722527422411862872999615628e-0002,
166*57151Sbostic u01 =   1.766664525091811069896442906220827182707e-0001,
167*57151Sbostic u02 =  -1.381856719455968955440002438182885835344e-0002,
168*57151Sbostic u03 =   3.474534320936836562092566861515617053954e-0004,
169*57151Sbostic u04 =  -3.814070537243641752631729276103284491172e-0006,
170*57151Sbostic u05 =   1.955901370350229170025509706510038090009e-0008,
171*57151Sbostic u06 =  -3.982051941321034108350630097330144576337e-0011,
172*57151Sbostic v01 =   1.273048348341237002944554656529224780561e-0002,
173*57151Sbostic v02 =   7.600686273503532807462101309675806839635e-0005,
174*57151Sbostic v03 =   2.591508518404578033173189144579208685163e-0007,
175*57151Sbostic v04 =   4.411103113326754838596529339004302243157e-0010;
17624597Szliu 
177*57151Sbostic double
178*57151Sbostic y0(x)
179*57151Sbostic 	double x;
180*57151Sbostic {
181*57151Sbostic 	double z, s,c,ss,cc,u,v,j0();
182*57151Sbostic     /* Y0(NaN) is NaN, y0(-inf) is Nan, y0(inf) is 0  */
183*57151Sbostic 	if (!finite(x))
184*57151Sbostic 		if (_IEEE)
185*57151Sbostic 			return (one/(x+x*x));
186*57151Sbostic 		else
187*57151Sbostic 			return (0);
188*57151Sbostic         if (x == 0)
189*57151Sbostic 		if (_IEEE)	return (-one/zero);
190*57151Sbostic 		else		return(infnan(-ERANGE));
191*57151Sbostic         if (x<0)
192*57151Sbostic 		if (_IEEE)	return (zero/zero);
193*57151Sbostic 		else		return (infnan(EDOM));
194*57151Sbostic         if (x >= 2.00) {	/* |x| >= 2.0 */
195*57151Sbostic         /* y0(x) = sqrt(2/(pi*x))*(p0(x)*sin(x0)+q0(x)*cos(x0))
196*57151Sbostic          * where x0 = x-pi/4
197*57151Sbostic          *      Better formula:
198*57151Sbostic          *              cos(x0) = cos(x)cos(pi/4)+sin(x)sin(pi/4)
199*57151Sbostic          *                      =  1/sqrt(2) * (sin(x) + cos(x))
200*57151Sbostic          *              sin(x0) = sin(x)cos(3pi/4)-cos(x)sin(3pi/4)
201*57151Sbostic          *                      =  1/sqrt(2) * (sin(x) - cos(x))
202*57151Sbostic          * To avoid cancellation, use
203*57151Sbostic          *              sin(x) +- cos(x) = -cos(2x)/(sin(x) -+ cos(x))
204*57151Sbostic          * to compute the worse one.
205*57151Sbostic          */
206*57151Sbostic                 s = sin(x);
207*57151Sbostic                 c = cos(x);
208*57151Sbostic                 ss = s-c;
209*57151Sbostic                 cc = s+c;
210*57151Sbostic 	/*
211*57151Sbostic 	 * j0(x) = 1/sqrt(pi) * (P(0,x)*cc - Q(0,x)*ss) / sqrt(x)
212*57151Sbostic 	 * y0(x) = 1/sqrt(pi) * (P(0,x)*ss + Q(0,x)*cc) / sqrt(x)
213*57151Sbostic 	 */
214*57151Sbostic                 if (x < .5 * DBL_MAX) {  /* make sure x+x not overflow */
215*57151Sbostic                     z = -cos(x+x);
216*57151Sbostic                     if ((s*c)<zero) cc = z/ss;
217*57151Sbostic                     else            ss = z/cc;
218*57151Sbostic                 }
219*57151Sbostic                 if (_IEEE && x > 6.80564733841876927e+38) /* > 2^129 */
220*57151Sbostic 			z = (invsqrtpi*ss)/sqrt(x);
221*57151Sbostic                 else {
222*57151Sbostic                     u = pzero(x); v = qzero(x);
223*57151Sbostic                     z = invsqrtpi*(u*ss+v*cc)/sqrt(x);
224*57151Sbostic                 }
225*57151Sbostic                 return z;
226*57151Sbostic 	}
227*57151Sbostic 	if (x <= 7.450580596923828125e-009) {		/* x < 2**-27 */
228*57151Sbostic 	    return (u00 + tpi*log(x));
229*57151Sbostic 	}
230*57151Sbostic 	z = x*x;
231*57151Sbostic 	u = u00+z*(u01+z*(u02+z*(u03+z*(u04+z*(u05+z*u06)))));
232*57151Sbostic 	v = one+z*(v01+z*(v02+z*(v03+z*v04)));
233*57151Sbostic 	return (u/v + tpi*(j0(x)*log(x)));
234*57151Sbostic }
23535679Sbostic 
236*57151Sbostic /* The asymptotic expansions of pzero is
237*57151Sbostic  *	1 - 9/128 s^2 + 11025/98304 s^4 - ...,	where s = 1/x.
238*57151Sbostic  * For x >= 2, We approximate pzero by
239*57151Sbostic  * 	pzero(x) = 1 + (R/S)
240*57151Sbostic  * where  R = pr0 + pr1*s^2 + pr2*s^4 + ... + pr5*s^10
241*57151Sbostic  * 	  S = 1 + ps0*s^2 + ... + ps4*s^10
242*57151Sbostic  * and
243*57151Sbostic  *	| pzero(x)-1-R/S | <= 2  ** ( -60.26)
244*57151Sbostic  */
245*57151Sbostic static double pr8[6] = { /* for x in [inf, 8]=1/[0,0.125] */
246*57151Sbostic    0.0,
247*57151Sbostic   -7.031249999999003994151563066182798210142e-0002,
248*57151Sbostic   -8.081670412753498508883963849859423939871e+0000,
249*57151Sbostic   -2.570631056797048755890526455854482662510e+0002,
250*57151Sbostic   -2.485216410094288379417154382189125598962e+0003,
251*57151Sbostic   -5.253043804907295692946647153614119665649e+0003,
25224597Szliu };
253*57151Sbostic static double ps8[5] = {
254*57151Sbostic    1.165343646196681758075176077627332052048e+0002,
255*57151Sbostic    3.833744753641218451213253490882686307027e+0003,
256*57151Sbostic    4.059785726484725470626341023967186966531e+0004,
257*57151Sbostic    1.167529725643759169416844015694440325519e+0005,
258*57151Sbostic    4.762772841467309430100106254805711722972e+0004,
25924597Szliu };
260*57151Sbostic 
261*57151Sbostic static double pr5[6] = { /* for x in [8,4.5454]=1/[0.125,0.22001] */
262*57151Sbostic   -1.141254646918944974922813501362824060117e-0011,
263*57151Sbostic   -7.031249408735992804117367183001996028304e-0002,
264*57151Sbostic   -4.159610644705877925119684455252125760478e+0000,
265*57151Sbostic   -6.767476522651671942610538094335912346253e+0001,
266*57151Sbostic   -3.312312996491729755731871867397057689078e+0002,
267*57151Sbostic   -3.464333883656048910814187305901796723256e+0002,
26824597Szliu };
269*57151Sbostic static double ps5[5] = {
270*57151Sbostic    6.075393826923003305967637195319271932944e+0001,
271*57151Sbostic    1.051252305957045869801410979087427910437e+0003,
272*57151Sbostic    5.978970943338558182743915287887408780344e+0003,
273*57151Sbostic    9.625445143577745335793221135208591603029e+0003,
274*57151Sbostic    2.406058159229391070820491174867406875471e+0003,
27524597Szliu };
276*57151Sbostic 
277*57151Sbostic static double pr3[6] = {/* for x in [4.547,2.8571]=1/[0.2199,0.35001] */
278*57151Sbostic   -2.547046017719519317420607587742992297519e-0009,
279*57151Sbostic   -7.031196163814817199050629727406231152464e-0002,
280*57151Sbostic   -2.409032215495295917537157371488126555072e+0000,
281*57151Sbostic   -2.196597747348830936268718293366935843223e+0001,
282*57151Sbostic   -5.807917047017375458527187341817239891940e+0001,
283*57151Sbostic   -3.144794705948885090518775074177485744176e+0001,
28424597Szliu };
285*57151Sbostic static double ps3[5] = {
286*57151Sbostic    3.585603380552097167919946472266854507059e+0001,
287*57151Sbostic    3.615139830503038919981567245265266294189e+0002,
288*57151Sbostic    1.193607837921115243628631691509851364715e+0003,
289*57151Sbostic    1.127996798569074250675414186814529958010e+0003,
290*57151Sbostic    1.735809308133357510239737333055228118910e+0002,
29124597Szliu };
292*57151Sbostic 
293*57151Sbostic static double pr2[6] = {/* for x in [2.8570,2]=1/[0.3499,0.5] */
294*57151Sbostic   -8.875343330325263874525704514800809730145e-0008,
295*57151Sbostic   -7.030309954836247756556445443331044338352e-0002,
296*57151Sbostic   -1.450738467809529910662233622603401167409e+0000,
297*57151Sbostic   -7.635696138235277739186371273434739292491e+0000,
298*57151Sbostic   -1.119316688603567398846655082201614524650e+0001,
299*57151Sbostic   -3.233645793513353260006821113608134669030e+0000,
30024597Szliu };
301*57151Sbostic static double ps2[5] = {
302*57151Sbostic    2.222029975320888079364901247548798910952e+0001,
303*57151Sbostic    1.362067942182152109590340823043813120940e+0002,
304*57151Sbostic    2.704702786580835044524562897256790293238e+0002,
305*57151Sbostic    1.538753942083203315263554770476850028583e+0002,
306*57151Sbostic    1.465761769482561965099880599279699314477e+0001,
30724597Szliu };
30824597Szliu 
309*57151Sbostic static double pzero(x)
310*57151Sbostic 	double x;
311*57151Sbostic {
312*57151Sbostic 	double *p,*q,z,r,s;
313*57151Sbostic 	if (x >= 8.00)			   {p = pr8; q= ps8;}
314*57151Sbostic 	else if (x >= 4.54545211791992188) {p = pr5; q= ps5;}
315*57151Sbostic 	else if (x >= 2.85714149475097656) {p = pr3; q= ps3;}
316*57151Sbostic 	else if (x >= 2.00)		   {p = pr2; q= ps2;}
317*57151Sbostic 	z = one/(x*x);
318*57151Sbostic 	r = p[0]+z*(p[1]+z*(p[2]+z*(p[3]+z*(p[4]+z*p[5]))));
319*57151Sbostic 	s = one+z*(q[0]+z*(q[1]+z*(q[2]+z*(q[3]+z*q[4]))));
320*57151Sbostic 	return one+ r/s;
321*57151Sbostic }
322*57151Sbostic 
32335679Sbostic 
324*57151Sbostic /* For x >= 8, the asymptotic expansions of qzero is
325*57151Sbostic  *	-1/8 s + 75/1024 s^3 - ..., where s = 1/x.
326*57151Sbostic  * We approximate pzero by
327*57151Sbostic  * 	qzero(x) = s*(-1.25 + (R/S))
328*57151Sbostic  * where  R = qr0 + qr1*s^2 + qr2*s^4 + ... + qr5*s^10
329*57151Sbostic  * 	  S = 1 + qs0*s^2 + ... + qs5*s^12
330*57151Sbostic  * and
331*57151Sbostic  *	| qzero(x)/s +1.25-R/S | <= 2  ** ( -61.22)
332*57151Sbostic  */
333*57151Sbostic static double qr8[6] = { /* for x in [inf, 8]=1/[0,0.125] */
334*57151Sbostic    0.0,
335*57151Sbostic    7.324218749999350414479738504551775297096e-0002,
336*57151Sbostic    1.176820646822526933903301695932765232456e+0001,
337*57151Sbostic    5.576733802564018422407734683549251364365e+0002,
338*57151Sbostic    8.859197207564685717547076568608235802317e+0003,
339*57151Sbostic    3.701462677768878501173055581933725704809e+0004,
340*57151Sbostic };
341*57151Sbostic static double qs8[6] = {
342*57151Sbostic    1.637760268956898345680262381842235272369e+0002,
343*57151Sbostic    8.098344946564498460163123708054674227492e+0003,
344*57151Sbostic    1.425382914191204905277585267143216379136e+0005,
345*57151Sbostic    8.033092571195144136565231198526081387047e+0005,
346*57151Sbostic    8.405015798190605130722042369969184811488e+0005,
347*57151Sbostic   -3.438992935378666373204500729736454421006e+0005,
348*57151Sbostic };
34924597Szliu 
350*57151Sbostic static double qr5[6] = { /* for x in [8,4.5454]=1/[0.125,0.22001] */
351*57151Sbostic    1.840859635945155400568380711372759921179e-0011,
352*57151Sbostic    7.324217666126847411304688081129741939255e-0002,
353*57151Sbostic    5.835635089620569401157245917610984757296e+0000,
354*57151Sbostic    1.351115772864498375785526599119895942361e+0002,
355*57151Sbostic    1.027243765961641042977177679021711341529e+0003,
356*57151Sbostic    1.989977858646053872589042328678602481924e+0003,
357*57151Sbostic };
358*57151Sbostic static double qs5[6] = {
359*57151Sbostic    8.277661022365377058749454444343415524509e+0001,
360*57151Sbostic    2.077814164213929827140178285401017305309e+0003,
361*57151Sbostic    1.884728877857180787101956800212453218179e+0004,
362*57151Sbostic    5.675111228949473657576693406600265778689e+0004,
363*57151Sbostic    3.597675384251145011342454247417399490174e+0004,
364*57151Sbostic   -5.354342756019447546671440667961399442388e+0003,
365*57151Sbostic };
36624597Szliu 
367*57151Sbostic static double qr3[6] = {/* for x in [4.547,2.8571]=1/[0.2199,0.35001] */
368*57151Sbostic    4.377410140897386263955149197672576223054e-0009,
369*57151Sbostic    7.324111800429115152536250525131924283018e-0002,
370*57151Sbostic    3.344231375161707158666412987337679317358e+0000,
371*57151Sbostic    4.262184407454126175974453269277100206290e+0001,
372*57151Sbostic    1.708080913405656078640701512007621675724e+0002,
373*57151Sbostic    1.667339486966511691019925923456050558293e+0002,
374*57151Sbostic };
375*57151Sbostic static double qs3[6] = {
376*57151Sbostic    4.875887297245871932865584382810260676713e+0001,
377*57151Sbostic    7.096892210566060535416958362640184894280e+0002,
378*57151Sbostic    3.704148226201113687434290319905207398682e+0003,
379*57151Sbostic    6.460425167525689088321109036469797462086e+0003,
380*57151Sbostic    2.516333689203689683999196167394889715078e+0003,
381*57151Sbostic   -1.492474518361563818275130131510339371048e+0002,
382*57151Sbostic };
38324597Szliu 
384*57151Sbostic static double qr2[6] = {/* for x in [2.8570,2]=1/[0.3499,0.5] */
385*57151Sbostic    1.504444448869832780257436041633206366087e-0007,
386*57151Sbostic    7.322342659630792930894554535717104926902e-0002,
387*57151Sbostic    1.998191740938159956838594407540292600331e+0000,
388*57151Sbostic    1.449560293478857407645853071687125850962e+0001,
389*57151Sbostic    3.166623175047815297062638132537957315395e+0001,
390*57151Sbostic    1.625270757109292688799540258329430963726e+0001,
391*57151Sbostic };
392*57151Sbostic static double qs2[6] = {
393*57151Sbostic    3.036558483552191922522729838478169383969e+0001,
394*57151Sbostic    2.693481186080498724211751445725708524507e+0002,
395*57151Sbostic    8.447837575953201460013136756723746023736e+0002,
396*57151Sbostic    8.829358451124885811233995083187666981299e+0002,
397*57151Sbostic    2.126663885117988324180482985363624996652e+0002,
398*57151Sbostic   -5.310954938826669402431816125780738924463e+0000,
399*57151Sbostic };
40024597Szliu 
401*57151Sbostic static double qzero(x)
402*57151Sbostic 	double x;
403*57151Sbostic {
404*57151Sbostic 	double *p,*q, s,r,z;
405*57151Sbostic 	if (x >= 8.00)			   {p = qr8; q= qs8;}
406*57151Sbostic 	else if (x >= 4.54545211791992188) {p = qr5; q= qs5;}
407*57151Sbostic 	else if (x >= 2.85714149475097656) {p = qr3; q= qs3;}
408*57151Sbostic 	else if (x >= 2.00)		   {p = qr2; q= qs2;}
409*57151Sbostic 	z = one/(x*x);
410*57151Sbostic 	r = p[0]+z*(p[1]+z*(p[2]+z*(p[3]+z*(p[4]+z*p[5]))));
411*57151Sbostic 	s = one+z*(q[0]+z*(q[1]+z*(q[2]+z*(q[3]+z*(q[4]+z*q[5])))));
412*57151Sbostic 	return (-.125 + r/s)/x;
41324597Szliu }
414