code cleanup of named constants
[musl] / src / math / erf.c
1 /* origin: FreeBSD /usr/src/lib/msun/src/s_erf.c */
2 /*
3  * ====================================================
4  * Copyright (C) 1993 by Sun Microsystems, Inc. All rights reserved.
5  *
6  * Developed at SunPro, a Sun Microsystems, Inc. business.
7  * Permission to use, copy, modify, and distribute this
8  * software is freely granted, provided that this notice
9  * is preserved.
10  * ====================================================
11  */
12 /* double erf(double x)
13  * double erfc(double x)
14  *                           x
15  *                    2      |\
16  *     erf(x)  =  ---------  | exp(-t*t)dt
17  *                 sqrt(pi) \|
18  *                           0
19  *
20  *     erfc(x) =  1-erf(x)
21  *  Note that
22  *              erf(-x) = -erf(x)
23  *              erfc(-x) = 2 - erfc(x)
24  *
25  * Method:
26  *      1. For |x| in [0, 0.84375]
27  *          erf(x)  = x + x*R(x^2)
28  *          erfc(x) = 1 - erf(x)           if x in [-.84375,0.25]
29  *                  = 0.5 + ((0.5-x)-x*R)  if x in [0.25,0.84375]
30  *         where R = P/Q where P is an odd poly of degree 8 and
31  *         Q is an odd poly of degree 10.
32  *                                               -57.90
33  *                      | R - (erf(x)-x)/x | <= 2
34  *
35  *
36  *         Remark. The formula is derived by noting
37  *          erf(x) = (2/sqrt(pi))*(x - x^3/3 + x^5/10 - x^7/42 + ....)
38  *         and that
39  *          2/sqrt(pi) = 1.128379167095512573896158903121545171688
40  *         is close to one. The interval is chosen because the fix
41  *         point of erf(x) is near 0.6174 (i.e., erf(x)=x when x is
42  *         near 0.6174), and by some experiment, 0.84375 is chosen to
43  *         guarantee the error is less than one ulp for erf.
44  *
45  *      2. For |x| in [0.84375,1.25], let s = |x| - 1, and
46  *         c = 0.84506291151 rounded to single (24 bits)
47  *              erf(x)  = sign(x) * (c  + P1(s)/Q1(s))
48  *              erfc(x) = (1-c)  - P1(s)/Q1(s) if x > 0
49  *                        1+(c+P1(s)/Q1(s))    if x < 0
50  *              |P1/Q1 - (erf(|x|)-c)| <= 2**-59.06
51  *         Remark: here we use the taylor series expansion at x=1.
52  *              erf(1+s) = erf(1) + s*Poly(s)
53  *                       = 0.845.. + P1(s)/Q1(s)
54  *         That is, we use rational approximation to approximate
55  *                      erf(1+s) - (c = (single)0.84506291151)
56  *         Note that |P1/Q1|< 0.078 for x in [0.84375,1.25]
57  *         where
58  *              P1(s) = degree 6 poly in s
59  *              Q1(s) = degree 6 poly in s
60  *
61  *      3. For x in [1.25,1/0.35(~2.857143)],
62  *              erfc(x) = (1/x)*exp(-x*x-0.5625+R1/S1)
63  *              erf(x)  = 1 - erfc(x)
64  *         where
65  *              R1(z) = degree 7 poly in z, (z=1/x^2)
66  *              S1(z) = degree 8 poly in z
67  *
68  *      4. For x in [1/0.35,28]
69  *              erfc(x) = (1/x)*exp(-x*x-0.5625+R2/S2) if x > 0
70  *                      = 2.0 - (1/x)*exp(-x*x-0.5625+R2/S2) if -6<x<0
71  *                      = 2.0 - tiny            (if x <= -6)
72  *              erf(x)  = sign(x)*(1.0 - erfc(x)) if x < 6, else
73  *              erf(x)  = sign(x)*(1.0 - tiny)
74  *         where
75  *              R2(z) = degree 6 poly in z, (z=1/x^2)
76  *              S2(z) = degree 7 poly in z
77  *
78  *      Note1:
79  *         To compute exp(-x*x-0.5625+R/S), let s be a single
80  *         precision number and s := x; then
81  *              -x*x = -s*s + (s-x)*(s+x)
82  *              exp(-x*x-0.5626+R/S) =
83  *                      exp(-s*s-0.5625)*exp((s-x)*(s+x)+R/S);
84  *      Note2:
85  *         Here 4 and 5 make use of the asymptotic series
86  *                        exp(-x*x)
87  *              erfc(x) ~ ---------- * ( 1 + Poly(1/x^2) )
88  *                        x*sqrt(pi)
89  *         We use rational approximation to approximate
90  *              g(s)=f(1/x^2) = log(erfc(x)*x) - x*x + 0.5625
91  *         Here is the error bound for R1/S1 and R2/S2
92  *              |R1/S1 - f(x)|  < 2**(-62.57)
93  *              |R2/S2 - f(x)|  < 2**(-61.52)
94  *
95  *      5. For inf > x >= 28
96  *              erf(x)  = sign(x) *(1 - tiny)  (raise inexact)
97  *              erfc(x) = tiny*tiny (raise underflow) if x > 0
98  *                      = 2 - tiny if x<0
99  *
100  *      7. Special case:
101  *              erf(0)  = 0, erf(inf)  = 1, erf(-inf) = -1,
102  *              erfc(0) = 1, erfc(inf) = 0, erfc(-inf) = 2,
103  *              erfc/erf(NaN) is NaN
104  */
105
106 #include "libm.h"
107
108 static const double
109 tiny = 1e-300,
110 /* c = (float)0.84506291151 */
111 erx  = 8.45062911510467529297e-01, /* 0x3FEB0AC1, 0x60000000 */
112 /*
113  * Coefficients for approximation to  erf on [0,0.84375]
114  */
115 efx  =  1.28379167095512586316e-01, /* 0x3FC06EBA, 0x8214DB69 */
116 efx8 =  1.02703333676410069053e+00, /* 0x3FF06EBA, 0x8214DB69 */
117 pp0  =  1.28379167095512558561e-01, /* 0x3FC06EBA, 0x8214DB68 */
118 pp1  = -3.25042107247001499370e-01, /* 0xBFD4CD7D, 0x691CB913 */
119 pp2  = -2.84817495755985104766e-02, /* 0xBF9D2A51, 0xDBD7194F */
120 pp3  = -5.77027029648944159157e-03, /* 0xBF77A291, 0x236668E4 */
121 pp4  = -2.37630166566501626084e-05, /* 0xBEF8EAD6, 0x120016AC */
122 qq1  =  3.97917223959155352819e-01, /* 0x3FD97779, 0xCDDADC09 */
123 qq2  =  6.50222499887672944485e-02, /* 0x3FB0A54C, 0x5536CEBA */
124 qq3  =  5.08130628187576562776e-03, /* 0x3F74D022, 0xC4D36B0F */
125 qq4  =  1.32494738004321644526e-04, /* 0x3F215DC9, 0x221C1A10 */
126 qq5  = -3.96022827877536812320e-06, /* 0xBED09C43, 0x42A26120 */
127 /*
128  * Coefficients for approximation to  erf  in [0.84375,1.25]
129  */
130 pa0  = -2.36211856075265944077e-03, /* 0xBF6359B8, 0xBEF77538 */
131 pa1  =  4.14856118683748331666e-01, /* 0x3FDA8D00, 0xAD92B34D */
132 pa2  = -3.72207876035701323847e-01, /* 0xBFD7D240, 0xFBB8C3F1 */
133 pa3  =  3.18346619901161753674e-01, /* 0x3FD45FCA, 0x805120E4 */
134 pa4  = -1.10894694282396677476e-01, /* 0xBFBC6398, 0x3D3E28EC */
135 pa5  =  3.54783043256182359371e-02, /* 0x3FA22A36, 0x599795EB */
136 pa6  = -2.16637559486879084300e-03, /* 0xBF61BF38, 0x0A96073F */
137 qa1  =  1.06420880400844228286e-01, /* 0x3FBB3E66, 0x18EEE323 */
138 qa2  =  5.40397917702171048937e-01, /* 0x3FE14AF0, 0x92EB6F33 */
139 qa3  =  7.18286544141962662868e-02, /* 0x3FB2635C, 0xD99FE9A7 */
140 qa4  =  1.26171219808761642112e-01, /* 0x3FC02660, 0xE763351F */
141 qa5  =  1.36370839120290507362e-02, /* 0x3F8BEDC2, 0x6B51DD1C */
142 qa6  =  1.19844998467991074170e-02, /* 0x3F888B54, 0x5735151D */
143 /*
144  * Coefficients for approximation to  erfc in [1.25,1/0.35]
145  */
146 ra0  = -9.86494403484714822705e-03, /* 0xBF843412, 0x600D6435 */
147 ra1  = -6.93858572707181764372e-01, /* 0xBFE63416, 0xE4BA7360 */
148 ra2  = -1.05586262253232909814e+01, /* 0xC0251E04, 0x41B0E726 */
149 ra3  = -6.23753324503260060396e+01, /* 0xC04F300A, 0xE4CBA38D */
150 ra4  = -1.62396669462573470355e+02, /* 0xC0644CB1, 0x84282266 */
151 ra5  = -1.84605092906711035994e+02, /* 0xC067135C, 0xEBCCABB2 */
152 ra6  = -8.12874355063065934246e+01, /* 0xC0545265, 0x57E4D2F2 */
153 ra7  = -9.81432934416914548592e+00, /* 0xC023A0EF, 0xC69AC25C */
154 sa1  =  1.96512716674392571292e+01, /* 0x4033A6B9, 0xBD707687 */
155 sa2  =  1.37657754143519042600e+02, /* 0x4061350C, 0x526AE721 */
156 sa3  =  4.34565877475229228821e+02, /* 0x407B290D, 0xD58A1A71 */
157 sa4  =  6.45387271733267880336e+02, /* 0x40842B19, 0x21EC2868 */
158 sa5  =  4.29008140027567833386e+02, /* 0x407AD021, 0x57700314 */
159 sa6  =  1.08635005541779435134e+02, /* 0x405B28A3, 0xEE48AE2C */
160 sa7  =  6.57024977031928170135e+00, /* 0x401A47EF, 0x8E484A93 */
161 sa8  = -6.04244152148580987438e-02, /* 0xBFAEEFF2, 0xEE749A62 */
162 /*
163  * Coefficients for approximation to  erfc in [1/.35,28]
164  */
165 rb0  = -9.86494292470009928597e-03, /* 0xBF843412, 0x39E86F4A */
166 rb1  = -7.99283237680523006574e-01, /* 0xBFE993BA, 0x70C285DE */
167 rb2  = -1.77579549177547519889e+01, /* 0xC031C209, 0x555F995A */
168 rb3  = -1.60636384855821916062e+02, /* 0xC064145D, 0x43C5ED98 */
169 rb4  = -6.37566443368389627722e+02, /* 0xC083EC88, 0x1375F228 */
170 rb5  = -1.02509513161107724954e+03, /* 0xC0900461, 0x6A2E5992 */
171 rb6  = -4.83519191608651397019e+02, /* 0xC07E384E, 0x9BDC383F */
172 sb1  =  3.03380607434824582924e+01, /* 0x403E568B, 0x261D5190 */
173 sb2  =  3.25792512996573918826e+02, /* 0x40745CAE, 0x221B9F0A */
174 sb3  =  1.53672958608443695994e+03, /* 0x409802EB, 0x189D5118 */
175 sb4  =  3.19985821950859553908e+03, /* 0x40A8FFB7, 0x688C246A */
176 sb5  =  2.55305040643316442583e+03, /* 0x40A3F219, 0xCEDF3BE6 */
177 sb6  =  4.74528541206955367215e+02, /* 0x407DA874, 0xE79FE763 */
178 sb7  = -2.24409524465858183362e+01; /* 0xC03670E2, 0x42712D62 */
179
180 double erf(double x)
181 {
182         int32_t hx,ix,i;
183         double R,S,P,Q,s,y,z,r;
184
185         GET_HIGH_WORD(hx, x);
186         ix = hx & 0x7fffffff;
187         if (ix >= 0x7ff00000) {
188                 /* erf(nan)=nan, erf(+-inf)=+-1 */
189                 i = ((uint32_t)hx>>31)<<1;
190                 return (double)(1-i) + 1.0/x;
191         }
192         if (ix < 0x3feb0000) {  /* |x|<0.84375 */
193                 if (ix < 0x3e300000) {  /* |x|<2**-28 */
194                         if (ix < 0x00800000)
195                                 /* avoid underflow */
196                                 return 0.125*(8.0*x + efx8*x);
197                         return x + efx*x;
198                 }
199                 z = x*x;
200                 r = pp0+z*(pp1+z*(pp2+z*(pp3+z*pp4)));
201                 s = 1.0+z*(qq1+z*(qq2+z*(qq3+z*(qq4+z*qq5))));
202                 y = r/s;
203                 return x + x*y;
204         }
205         if (ix < 0x3ff40000) {  /* 0.84375 <= |x| < 1.25 */
206                 s = fabs(x)-1.0;
207                 P = pa0+s*(pa1+s*(pa2+s*(pa3+s*(pa4+s*(pa5+s*pa6)))));
208                 Q = 1.0+s*(qa1+s*(qa2+s*(qa3+s*(qa4+s*(qa5+s*qa6)))));
209                 if (hx >= 0)
210                         return erx + P/Q;
211                 return -erx - P/Q;
212         }
213         if (ix >= 0x40180000) {  /* inf > |x| >= 6 */
214                 if (hx >= 0)
215                         return 1.0 - tiny;
216                 return tiny - 1.0;
217         }
218         x = fabs(x);
219         s = 1.0/(x*x);
220         if (ix < 0x4006DB6E) {  /* |x| < 1/0.35 */
221                 R = ra0+s*(ra1+s*(ra2+s*(ra3+s*(ra4+s*(
222                      ra5+s*(ra6+s*ra7))))));
223                 S = 1.0+s*(sa1+s*(sa2+s*(sa3+s*(sa4+s*(
224                      sa5+s*(sa6+s*(sa7+s*sa8)))))));
225         } else {                /* |x| >= 1/0.35 */
226                 R = rb0+s*(rb1+s*(rb2+s*(rb3+s*(rb4+s*(
227                      rb5+s*rb6)))));
228                 S = 1.0+s*(sb1+s*(sb2+s*(sb3+s*(sb4+s*(
229                      sb5+s*(sb6+s*sb7))))));
230         }
231         z = x;
232         SET_LOW_WORD(z,0);
233         r = exp(-z*z-0.5625)*exp((z-x)*(z+x)+R/S);
234         if (hx >= 0)
235                 return 1.0 - r/x;
236         return r/x - 1.0;
237 }
238
239 double erfc(double x)
240 {
241         int32_t hx,ix;
242         double R,S,P,Q,s,y,z,r;
243
244         GET_HIGH_WORD(hx, x);
245         ix = hx & 0x7fffffff;
246         if (ix >= 0x7ff00000) {
247                 /* erfc(nan)=nan, erfc(+-inf)=0,2 */
248                 return (double)(((uint32_t)hx>>31)<<1) + 1.0/x;
249         }
250         if (ix < 0x3feb0000) {  /* |x| < 0.84375 */
251                 if (ix < 0x3c700000)  /* |x| < 2**-56 */
252                         return 1.0 - x;
253                 z = x*x;
254                 r = pp0+z*(pp1+z*(pp2+z*(pp3+z*pp4)));
255                 s = 1.0+z*(qq1+z*(qq2+z*(qq3+z*(qq4+z*qq5))));
256                 y = r/s;
257                 if (hx < 0x3fd00000) {  /* x < 1/4 */
258                         return 1.0 - (x+x*y);
259                 } else {
260                         r = x*y;
261                         r += x - 0.5;
262                         return 0.5 - r ;
263                 }
264         }
265         if (ix < 0x3ff40000) {  /* 0.84375 <= |x| < 1.25 */
266                 s = fabs(x)-1.0;
267                 P = pa0+s*(pa1+s*(pa2+s*(pa3+s*(pa4+s*(pa5+s*pa6)))));
268                 Q = 1.0+s*(qa1+s*(qa2+s*(qa3+s*(qa4+s*(qa5+s*qa6)))));
269                 if (hx >= 0) {
270                         z = 1.0-erx;
271                         return z - P/Q;
272                 } else {
273                         z = erx+P/Q;
274                         return 1.0 + z;
275                 }
276         }
277         if (ix < 0x403c0000) {  /* |x| < 28 */
278                 x = fabs(x);
279                 s = 1.0/(x*x);
280                 if (ix < 0x4006DB6D) {  /* |x| < 1/.35 ~ 2.857143*/
281                         R = ra0+s*(ra1+s*(ra2+s*(ra3+s*(ra4+s*(
282                              ra5+s*(ra6+s*ra7))))));
283                         S = 1.0+s*(sa1+s*(sa2+s*(sa3+s*(sa4+s*(
284                              sa5+s*(sa6+s*(sa7+s*sa8)))))));
285                 } else {                /* |x| >= 1/.35 ~ 2.857143 */
286                         if (hx < 0 && ix >= 0x40180000)  /* x < -6 */
287                                 return 2.0 - tiny;
288                         R = rb0+s*(rb1+s*(rb2+s*(rb3+s*(rb4+s*(
289                              rb5+s*rb6)))));
290                         S = 1.0+s*(sb1+s*(sb2+s*(sb3+s*(sb4+s*(
291                              sb5+s*(sb6+s*sb7))))));
292                 }
293                 z = x;
294                 SET_LOW_WORD(z, 0);
295                 r = exp(-z*z-0.5625)*exp((z-x)*(z+x)+R/S);
296                 if (hx > 0)
297                         return r/x;
298                 return 2.0 - r/x;
299         }
300         if (hx > 0)
301                 return tiny*tiny;
302         return 2.0 - tiny;
303 }