1
2 /* @(#)s_erf.c 5.1 93/09/24 */
3 /*
4 * ====================================================
5 * Copyright (C) 1993 by Sun Microsystems, Inc. All rights reserved.
6 *
7 * Developed at SunPro, a Sun Microsystems, Inc. business.
8 * Permission to use, copy, modify, and distribute this
9 * software is freely granted, provided that this notice
10 * is preserved.
11 * ====================================================
12 */
13
14 /*
15 FUNCTION
16 <<erf>>, <<erff>>, <<erfc>>, <<erfcf>>---error function
17 INDEX
18 erf
19 INDEX
20 erff
21 INDEX
22 erfc
23 INDEX
24 erfcf
25
26 SYNOPSIS
27 #include <math.h>
28 double erf(double <[x]>);
29 float erff(float <[x]>);
30 double erfc(double <[x]>);
31 float erfcf(float <[x]>);
32
33 DESCRIPTION
34 <<erf>> calculates an approximation to the ``error function'',
35 which estimates the probability that an observation will fall within
36 <[x]> standard deviations of the mean (assuming a normal
37 distribution).
38 @tex
39 The error function is defined as
40 $${2\over\sqrt\pi}\times\int_0^x e^{-t^2}dt$$
41 @end tex
42
43 <<erfc>> calculates the complementary probability; that is,
44 <<erfc(<[x]>)>> is <<1 - erf(<[x]>)>>. <<erfc>> is computed directly,
45 so that you can use it to avoid the loss of precision that would
46 result from subtracting large probabilities (on large <[x]>) from 1.
47
48 <<erff>> and <<erfcf>> differ from <<erf>> and <<erfc>> only in the
49 argument and result types.
50
51 RETURNS
52 For positive arguments, <<erf>> and all its variants return a
53 probability---a number between 0 and 1.
54
55 PORTABILITY
56 None of the variants of <<erf>> are ANSI C.
57 */
58
59 /* double erf(double x)
60 * double erfc(double x)
61 * x
62 * 2 |\
63 * erf(x) = --------- | exp(-t*t)dt
64 * sqrt(pi) \|
65 * 0
66 *
67 * erfc(x) = 1-erf(x)
68 * Note that
69 * erf(-x) = -erf(x)
70 * erfc(-x) = 2 - erfc(x)
71 *
72 * Method:
73 * 1. For |x| in [0, 0.84375]
74 * erf(x) = x + x*R(x^2)
75 * erfc(x) = 1 - erf(x) if x in [-.84375,0.25]
76 * = 0.5 + ((0.5-x)-x*R) if x in [0.25,0.84375]
77 * where R = P/Q where P is an odd poly of degree 8 and
78 * Q is an odd poly of degree 10.
79 * -57.90
80 * | R - (erf(x)-x)/x | <= 2
81 *
82 *
83 * Remark. The formula is derived by noting
84 * erf(x) = (2/sqrt(pi))*(x - x^3/3 + x^5/10 - x^7/42 + ....)
85 * and that
86 * 2/sqrt(pi) = 1.128379167095512573896158903121545171688
87 * is close to one. The interval is chosen because the fix
88 * point of erf(x) is near 0.6174 (i.e., erf(x)=x when x is
89 * near 0.6174), and by some experiment, 0.84375 is chosen to
90 * guarantee the error is less than one ulp for erf.
91 *
92 * 2. For |x| in [0.84375,1.25], let s = |x| - 1, and
93 * c = 0.84506291151 rounded to single (24 bits)
94 * erf(x) = sign(x) * (c + P1(s)/Q1(s))
95 * erfc(x) = (1-c) - P1(s)/Q1(s) if x > 0
96 * 1+(c+P1(s)/Q1(s)) if x < 0
97 * |P1/Q1 - (erf(|x|)-c)| <= 2**-59.06
98 * Remark: here we use the taylor series expansion at x=1.
99 * erf(1+s) = erf(1) + s*Poly(s)
100 * = 0.845.. + P1(s)/Q1(s)
101 * That is, we use rational approximation to approximate
102 * erf(1+s) - (c = (single)0.84506291151)
103 * Note that |P1/Q1|< 0.078 for x in [0.84375,1.25]
104 * where
105 * P1(s) = degree 6 poly in s
106 * Q1(s) = degree 6 poly in s
107 *
108 * 3. For x in [1.25,1/0.35(~2.857143)],
109 * erfc(x) = (1/x)*exp(-x*x-0.5625+R1/S1)
110 * erf(x) = 1 - erfc(x)
111 * where
112 * R1(z) = degree 7 poly in z, (z=1/x^2)
113 * S1(z) = degree 8 poly in z
114 *
115 * 4. For x in [1/0.35,28]
116 * erfc(x) = (1/x)*exp(-x*x-0.5625+R2/S2) if x > 0
117 * = 2.0 - (1/x)*exp(-x*x-0.5625+R2/S2) if -6<x<0
118 * = 2.0 - tiny (if x <= -6)
119 * erf(x) = sign(x)*(1.0 - erfc(x)) if x < 6, else
120 * erf(x) = sign(x)*(1.0 - tiny)
121 * where
122 * R2(z) = degree 6 poly in z, (z=1/x^2)
123 * S2(z) = degree 7 poly in z
124 *
125 * Note1:
126 * To compute exp(-x*x-0.5625+R/S), let s be a single
127 * precision number and s := x; then
128 * -x*x = -s*s + (s-x)*(s+x)
129 * exp(-x*x-0.5626+R/S) =
130 * exp(-s*s-0.5625)*exp((s-x)*(s+x)+R/S);
131 * Note2:
132 * Here 4 and 5 make use of the asymptotic series
133 * exp(-x*x)
134 * erfc(x) ~ ---------- * ( 1 + Poly(1/x^2) )
135 * x*sqrt(pi)
136 * We use rational approximation to approximate
137 * g(s)=f(1/x^2) = log(erfc(x)*x) - x*x + 0.5625
138 * Here is the error bound for R1/S1 and R2/S2
139 * |R1/S1 - f(x)| < 2**(-62.57)
140 * |R2/S2 - f(x)| < 2**(-61.52)
141 *
142 * 5. For inf > x >= 28
143 * erf(x) = sign(x) *(1 - tiny) (raise inexact)
144 * erfc(x) = tiny*tiny (raise underflow) if x > 0
145 * = 2 - tiny if x<0
146 *
147 * 7. Special case:
148 * erf(0) = 0, erf(inf) = 1, erf(-inf) = -1,
149 * erfc(0) = 1, erfc(inf) = 0, erfc(-inf) = 2,
150 * erfc/erf(NaN) is NaN
151 */
152
153 #include "fdlibm.h"
154
155 #ifdef _NEED_FLOAT64
156
157 static const __float64
158 tiny = _F_64(1e-300),
159 half = _F_64(5.00000000000000000000e-01), /* 0x3FE00000, 0x00000000 */
160 one = _F_64(1.00000000000000000000e+00), /* 0x3FF00000, 0x00000000 */
161 two = _F_64(2.00000000000000000000e+00), /* 0x40000000, 0x00000000 */
162 /* c = (float)0.84506291151 */
163 erx = _F_64(8.45062911510467529297e-01), /* 0x3FEB0AC1, 0x60000000 */
164 /*
165 * Coefficients for approximation to erf on [0,0.84375]
166 */
167 efx = _F_64(1.28379167095512586316e-01), /* 0x3FC06EBA, 0x8214DB69 */
168 efx8 = _F_64(1.02703333676410069053e+00), /* 0x3FF06EBA, 0x8214DB69 */
169 pp0 = _F_64(1.28379167095512558561e-01), /* 0x3FC06EBA, 0x8214DB68 */
170 pp1 = _F_64(-3.25042107247001499370e-01), /* 0xBFD4CD7D, 0x691CB913 */
171 pp2 = _F_64(-2.84817495755985104766e-02), /* 0xBF9D2A51, 0xDBD7194F */
172 pp3 = _F_64(-5.77027029648944159157e-03), /* 0xBF77A291, 0x236668E4 */
173 pp4 = _F_64(-2.37630166566501626084e-05), /* 0xBEF8EAD6, 0x120016AC */
174 qq1 = _F_64(3.97917223959155352819e-01), /* 0x3FD97779, 0xCDDADC09 */
175 qq2 = _F_64(6.50222499887672944485e-02), /* 0x3FB0A54C, 0x5536CEBA */
176 qq3 = _F_64(5.08130628187576562776e-03), /* 0x3F74D022, 0xC4D36B0F */
177 qq4 = _F_64(1.32494738004321644526e-04), /* 0x3F215DC9, 0x221C1A10 */
178 qq5 = _F_64(-3.96022827877536812320e-06), /* 0xBED09C43, 0x42A26120 */
179 /*
180 * Coefficients for approximation to erf in [0.84375,1.25]
181 */
182 pa0 = _F_64(-2.36211856075265944077e-03), /* 0xBF6359B8, 0xBEF77538 */
183 pa1 = _F_64(4.14856118683748331666e-01), /* 0x3FDA8D00, 0xAD92B34D */
184 pa2 = _F_64(-3.72207876035701323847e-01), /* 0xBFD7D240, 0xFBB8C3F1 */
185 pa3 = _F_64(3.18346619901161753674e-01), /* 0x3FD45FCA, 0x805120E4 */
186 pa4 = _F_64(-1.10894694282396677476e-01), /* 0xBFBC6398, 0x3D3E28EC */
187 pa5 = _F_64(3.54783043256182359371e-02), /* 0x3FA22A36, 0x599795EB */
188 pa6 = _F_64(-2.16637559486879084300e-03), /* 0xBF61BF38, 0x0A96073F */
189 qa1 = _F_64(1.06420880400844228286e-01), /* 0x3FBB3E66, 0x18EEE323 */
190 qa2 = _F_64(5.40397917702171048937e-01), /* 0x3FE14AF0, 0x92EB6F33 */
191 qa3 = _F_64(7.18286544141962662868e-02), /* 0x3FB2635C, 0xD99FE9A7 */
192 qa4 = _F_64(1.26171219808761642112e-01), /* 0x3FC02660, 0xE763351F */
193 qa5 = _F_64(1.36370839120290507362e-02), /* 0x3F8BEDC2, 0x6B51DD1C */
194 qa6 = _F_64(1.19844998467991074170e-02), /* 0x3F888B54, 0x5735151D */
195 /*
196 * Coefficients for approximation to erfc in [1.25,1/0.35]
197 */
198 ra0 = _F_64(-9.86494403484714822705e-03), /* 0xBF843412, 0x600D6435 */
199 ra1 = _F_64(-6.93858572707181764372e-01), /* 0xBFE63416, 0xE4BA7360 */
200 ra2 = _F_64(-1.05586262253232909814e+01), /* 0xC0251E04, 0x41B0E726 */
201 ra3 = _F_64(-6.23753324503260060396e+01), /* 0xC04F300A, 0xE4CBA38D */
202 ra4 = _F_64(-1.62396669462573470355e+02), /* 0xC0644CB1, 0x84282266 */
203 ra5 = _F_64(-1.84605092906711035994e+02), /* 0xC067135C, 0xEBCCABB2 */
204 ra6 = _F_64(-8.12874355063065934246e+01), /* 0xC0545265, 0x57E4D2F2 */
205 ra7 = _F_64(-9.81432934416914548592e+00), /* 0xC023A0EF, 0xC69AC25C */
206 sa1 = _F_64(1.96512716674392571292e+01), /* 0x4033A6B9, 0xBD707687 */
207 sa2 = _F_64(1.37657754143519042600e+02), /* 0x4061350C, 0x526AE721 */
208 sa3 = _F_64(4.34565877475229228821e+02), /* 0x407B290D, 0xD58A1A71 */
209 sa4 = _F_64(6.45387271733267880336e+02), /* 0x40842B19, 0x21EC2868 */
210 sa5 = _F_64(4.29008140027567833386e+02), /* 0x407AD021, 0x57700314 */
211 sa6 = _F_64(1.08635005541779435134e+02), /* 0x405B28A3, 0xEE48AE2C */
212 sa7 = _F_64(6.57024977031928170135e+00), /* 0x401A47EF, 0x8E484A93 */
213 sa8 = _F_64(-6.04244152148580987438e-02), /* 0xBFAEEFF2, 0xEE749A62 */
214 /*
215 * Coefficients for approximation to erfc in [1/.35,28]
216 */
217 rb0 = _F_64(-9.86494292470009928597e-03), /* 0xBF843412, 0x39E86F4A */
218 rb1 = _F_64(-7.99283237680523006574e-01), /* 0xBFE993BA, 0x70C285DE */
219 rb2 = _F_64(-1.77579549177547519889e+01), /* 0xC031C209, 0x555F995A */
220 rb3 = _F_64(-1.60636384855821916062e+02), /* 0xC064145D, 0x43C5ED98 */
221 rb4 = _F_64(-6.37566443368389627722e+02), /* 0xC083EC88, 0x1375F228 */
222 rb5 = _F_64(-1.02509513161107724954e+03), /* 0xC0900461, 0x6A2E5992 */
223 rb6 = _F_64(-4.83519191608651397019e+02), /* 0xC07E384E, 0x9BDC383F */
224 sb1 = _F_64(3.03380607434824582924e+01), /* 0x403E568B, 0x261D5190 */
225 sb2 = _F_64(3.25792512996573918826e+02), /* 0x40745CAE, 0x221B9F0A */
226 sb3 = _F_64(1.53672958608443695994e+03), /* 0x409802EB, 0x189D5118 */
227 sb4 = _F_64(3.19985821950859553908e+03), /* 0x40A8FFB7, 0x688C246A */
228 sb5 = _F_64(2.55305040643316442583e+03), /* 0x40A3F219, 0xCEDF3BE6 */
229 sb6 = _F_64(4.74528541206955367215e+02), /* 0x407DA874, 0xE79FE763 */
230 sb7 = _F_64(-2.24409524465858183362e+01); /* 0xC03670E2, 0x42712D62 */
231
232 __float64
erf64(__float64 x)233 erf64(__float64 x)
234 {
235 __int32_t hx, ix, i;
236 __float64 R, S, P, Q, s, y, z, r;
237 GET_HIGH_WORD(hx, x);
238 ix = hx & 0x7fffffff;
239 if (ix >= 0x7ff00000) { /* erf(nan)=nan */
240 i = ((__uint32_t)hx >> 31) << 1;
241 return (__float64)(1 - i) + one / x; /* erf(+-inf)=+-1 */
242 }
243
244 if (ix < 0x3feb0000) { /* |x|<0.84375 */
245 if (ix < 0x3e300000) { /* |x|<2**-28 */
246 if (ix < 0x00800000)
247 return _F_64(0.125) * (_F_64(8.0) * x + efx8 * x); /*avoid underflow */
248 return x + efx * x;
249 }
250 z = x * x;
251 r = pp0 + z * (pp1 + z * (pp2 + z * (pp3 + z * pp4)));
252 s = one + z * (qq1 + z * (qq2 + z * (qq3 + z * (qq4 + z * qq5))));
253 y = r / s;
254 return x + x * y;
255 }
256 if (ix < 0x3ff40000) { /* 0.84375 <= |x| < 1.25 */
257 s = fabs64(x) - one;
258 P = pa0 +
259 s * (pa1 + s * (pa2 + s * (pa3 + s * (pa4 + s * (pa5 + s * pa6)))));
260 Q = one +
261 s * (qa1 + s * (qa2 + s * (qa3 + s * (qa4 + s * (qa5 + s * qa6)))));
262 if (hx >= 0)
263 return erx + P / Q;
264 else
265 return -erx - P / Q;
266 }
267 if (ix >= 0x40180000) { /* inf>|x|>=6 */
268 if (hx >= 0)
269 return one - tiny;
270 else
271 return tiny - one;
272 }
273 x = fabs64(x);
274 s = one / (x * x);
275 if (ix < 0x4006DB6E) { /* |x| < 1/0.35 */
276 R = ra0 +
277 s * (ra1 +
278 s * (ra2 +
279 s * (ra3 + s * (ra4 + s * (ra5 + s * (ra6 + s * ra7))))));
280 S = one +
281 s * (sa1 +
282 s * (sa2 +
283 s * (sa3 +
284 s * (sa4 +
285 s * (sa5 + s * (sa6 + s * (sa7 + s * sa8)))))));
286 } else { /* |x| >= 1/0.35 */
287 R = rb0 +
288 s * (rb1 + s * (rb2 + s * (rb3 + s * (rb4 + s * (rb5 + s * rb6)))));
289 S = one +
290 s * (sb1 +
291 s * (sb2 +
292 s * (sb3 + s * (sb4 + s * (sb5 + s * (sb6 + s * sb7))))));
293 }
294 z = x;
295 SET_LOW_WORD(z, 0);
296 r = exp(-z * z - _F_64(0.5625)) * exp((z - x) * (z + x) + R / S);
297 if (hx >= 0)
298 return one - r / x;
299 else
300 return r / x - one;
301 }
302
_MATH_ALIAS_d_d(erf)303 _MATH_ALIAS_d_d(erf)
304
305 __float64
306 erfc64(__float64 x)
307 {
308 __int32_t hx, ix;
309 __float64 R, S, P, Q, s, y, z, r;
310 GET_HIGH_WORD(hx, x);
311 ix = hx & 0x7fffffff;
312 if (ix >= 0x7ff00000) { /* erfc(nan)=nan */
313 /* erfc(+-inf)=0,2 */
314 return (__float64)(((__uint32_t)hx >> 31) << 1) + one / x;
315 }
316
317 if (ix < 0x3feb0000) { /* |x|<0.84375 */
318 if (ix < 0x3c700000) /* |x|<2**-56 */
319 return one - x;
320 z = x * x;
321 r = pp0 + z * (pp1 + z * (pp2 + z * (pp3 + z * pp4)));
322 s = one + z * (qq1 + z * (qq2 + z * (qq3 + z * (qq4 + z * qq5))));
323 y = r / s;
324 if (hx < 0x3fd00000) { /* x<1/4 */
325 return one - (x + x * y);
326 } else {
327 r = x * y;
328 r += (x - half);
329 return half - r;
330 }
331 }
332 if (ix < 0x3ff40000) { /* 0.84375 <= |x| < 1.25 */
333 s = fabs64(x) - one;
334 P = pa0 +
335 s * (pa1 + s * (pa2 + s * (pa3 + s * (pa4 + s * (pa5 + s * pa6)))));
336 Q = one +
337 s * (qa1 + s * (qa2 + s * (qa3 + s * (qa4 + s * (qa5 + s * qa6)))));
338 if (hx >= 0) {
339 z = one - erx;
340 return z - P / Q;
341 } else {
342 z = erx + P / Q;
343 return one + z;
344 }
345 }
346 if (ix < 0x403c0000) { /* |x|<28 */
347 x = fabs64(x);
348 s = one / (x * x);
349 if (ix < 0x4006DB6D) { /* |x| < 1/.35 ~ 2.857143*/
350 R = ra0 +
351 s * (ra1 +
352 s * (ra2 +
353 s * (ra3 +
354 s * (ra4 + s * (ra5 + s * (ra6 + s * ra7))))));
355 S = one +
356 s * (sa1 +
357 s * (sa2 +
358 s * (sa3 +
359 s * (sa4 +
360 s * (sa5 +
361 s * (sa6 + s * (sa7 + s * sa8)))))));
362 } else { /* |x| >= 1/.35 ~ 2.857143 */
363 if (hx < 0 && ix >= 0x40180000)
364 return two - tiny; /* x < -6 */
365 R = rb0 +
366 s * (rb1 +
367 s * (rb2 + s * (rb3 + s * (rb4 + s * (rb5 + s * rb6)))));
368 S = one +
369 s * (sb1 +
370 s * (sb2 +
371 s * (sb3 +
372 s * (sb4 + s * (sb5 + s * (sb6 + s * sb7))))));
373 }
374 z = x;
375 SET_LOW_WORD(z, 0);
376 r = exp(-z * z - _F_64(0.5625)) * exp((z - x) * (z + x) + R / S);
377 if (hx > 0)
378 return r / x;
379 else
380 return two - r / x;
381 } else {
382 if (hx > 0)
383 return __math_uflow(0);
384 else
385 return two - tiny;
386 }
387 }
388
389 _MATH_ALIAS_d_d(erfc)
390
391 #endif /* _NEED_FLOAT64 */
392