Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / math / utilities.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gromacs/math/utilities.h"
38
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42
43 #include <assert.h>
44 #include <math.h>
45 #include <limits.h>
46 #ifdef HAVE__FINITE
47 #include <float.h>
48 #endif
49
50 int gmx_nint(real a)
51 {
52     const real half = .5;
53     int        result;
54
55     result = (a < 0.) ? ((int)(a - half)) : ((int)(a + half));
56     return result;
57 }
58
59 real cuberoot(real x)
60 {
61     if (x < 0)
62     {
63         return (-pow(-x, 1.0/3.0));
64     }
65     else
66     {
67         return (pow(x, 1.0/3.0));
68     }
69 }
70
71 real sign(real x, real y)
72 {
73     if (y < 0)
74     {
75         return -fabs(x);
76     }
77     else
78     {
79         return +fabs(x);
80     }
81 }
82
83 /* Double and single precision erf() and erfc() from
84  * the Sun Freely Distributable Math Library FDLIBM.
85  * See http://www.netlib.org/fdlibm
86  * Specific file used: s_erf.c, version 1.3 95/01/18
87  */
88 /*
89  * ====================================================
90  * Copyright (C) 1993 by Sun Microsystems, Inc. All rights reserved.
91  *
92  * Developed at SunSoft, a Sun Microsystems, Inc. business.
93  * Permission to use, copy, modify, and distribute this
94  * software is freely granted, provided that this notice
95  * is preserved.
96  * ====================================================
97  */
98
99 static const double
100     tiny        = 1e-300,
101     half        =  5.00000000000000000000e-01, /* 0x3FE00000, 0x00000000 */
102     one         =  1.00000000000000000000e+00, /* 0x3FF00000, 0x00000000 */
103     two         =  2.00000000000000000000e+00, /* 0x40000000, 0x00000000 */
104 /* c = (float)0.84506291151 */
105     erx =  8.45062911510467529297e-01,         /* 0x3FEB0AC1, 0x60000000 */
106 /*
107  * Coefficients for approximation to  erf on [0,0.84375]
108  */
109     efx  =  1.28379167095512586316e-01,        /* 0x3FC06EBA, 0x8214DB69 */
110     efx8 =  1.02703333676410069053e+00,        /* 0x3FF06EBA, 0x8214DB69 */
111     pp0  =  1.28379167095512558561e-01,        /* 0x3FC06EBA, 0x8214DB68 */
112     pp1  = -3.25042107247001499370e-01,        /* 0xBFD4CD7D, 0x691CB913 */
113     pp2  = -2.84817495755985104766e-02,        /* 0xBF9D2A51, 0xDBD7194F */
114     pp3  = -5.77027029648944159157e-03,        /* 0xBF77A291, 0x236668E4 */
115     pp4  = -2.37630166566501626084e-05,        /* 0xBEF8EAD6, 0x120016AC */
116     qq1  =  3.97917223959155352819e-01,        /* 0x3FD97779, 0xCDDADC09 */
117     qq2  =  6.50222499887672944485e-02,        /* 0x3FB0A54C, 0x5536CEBA */
118     qq3  =  5.08130628187576562776e-03,        /* 0x3F74D022, 0xC4D36B0F */
119     qq4  =  1.32494738004321644526e-04,        /* 0x3F215DC9, 0x221C1A10 */
120     qq5  = -3.96022827877536812320e-06,        /* 0xBED09C43, 0x42A26120 */
121 /*
122  * Coefficients for approximation to  erf  in [0.84375,1.25]
123  */
124     pa0  = -2.36211856075265944077e-03,        /* 0xBF6359B8, 0xBEF77538 */
125     pa1  =  4.14856118683748331666e-01,        /* 0x3FDA8D00, 0xAD92B34D */
126     pa2  = -3.72207876035701323847e-01,        /* 0xBFD7D240, 0xFBB8C3F1 */
127     pa3  =  3.18346619901161753674e-01,        /* 0x3FD45FCA, 0x805120E4 */
128     pa4  = -1.10894694282396677476e-01,        /* 0xBFBC6398, 0x3D3E28EC */
129     pa5  =  3.54783043256182359371e-02,        /* 0x3FA22A36, 0x599795EB */
130     pa6  = -2.16637559486879084300e-03,        /* 0xBF61BF38, 0x0A96073F */
131     qa1  =  1.06420880400844228286e-01,        /* 0x3FBB3E66, 0x18EEE323 */
132     qa2  =  5.40397917702171048937e-01,        /* 0x3FE14AF0, 0x92EB6F33 */
133     qa3  =  7.18286544141962662868e-02,        /* 0x3FB2635C, 0xD99FE9A7 */
134     qa4  =  1.26171219808761642112e-01,        /* 0x3FC02660, 0xE763351F */
135     qa5  =  1.36370839120290507362e-02,        /* 0x3F8BEDC2, 0x6B51DD1C */
136     qa6  =  1.19844998467991074170e-02,        /* 0x3F888B54, 0x5735151D */
137 /*
138  * Coefficients for approximation to  erfc in [1.25,1/0.35]
139  */
140     ra0  = -9.86494403484714822705e-03,        /* 0xBF843412, 0x600D6435 */
141     ra1  = -6.93858572707181764372e-01,        /* 0xBFE63416, 0xE4BA7360 */
142     ra2  = -1.05586262253232909814e+01,        /* 0xC0251E04, 0x41B0E726 */
143     ra3  = -6.23753324503260060396e+01,        /* 0xC04F300A, 0xE4CBA38D */
144     ra4  = -1.62396669462573470355e+02,        /* 0xC0644CB1, 0x84282266 */
145     ra5  = -1.84605092906711035994e+02,        /* 0xC067135C, 0xEBCCABB2 */
146     ra6  = -8.12874355063065934246e+01,        /* 0xC0545265, 0x57E4D2F2 */
147     ra7  = -9.81432934416914548592e+00,        /* 0xC023A0EF, 0xC69AC25C */
148     sa1  =  1.96512716674392571292e+01,        /* 0x4033A6B9, 0xBD707687 */
149     sa2  =  1.37657754143519042600e+02,        /* 0x4061350C, 0x526AE721 */
150     sa3  =  4.34565877475229228821e+02,        /* 0x407B290D, 0xD58A1A71 */
151     sa4  =  6.45387271733267880336e+02,        /* 0x40842B19, 0x21EC2868 */
152     sa5  =  4.29008140027567833386e+02,        /* 0x407AD021, 0x57700314 */
153     sa6  =  1.08635005541779435134e+02,        /* 0x405B28A3, 0xEE48AE2C */
154     sa7  =  6.57024977031928170135e+00,        /* 0x401A47EF, 0x8E484A93 */
155     sa8  = -6.04244152148580987438e-02,        /* 0xBFAEEFF2, 0xEE749A62 */
156 /*
157  * Coefficients for approximation to  erfc in [1/.35,28]
158  */
159     rb0  = -9.86494292470009928597e-03,        /* 0xBF843412, 0x39E86F4A */
160     rb1  = -7.99283237680523006574e-01,        /* 0xBFE993BA, 0x70C285DE */
161     rb2  = -1.77579549177547519889e+01,        /* 0xC031C209, 0x555F995A */
162     rb3  = -1.60636384855821916062e+02,        /* 0xC064145D, 0x43C5ED98 */
163     rb4  = -6.37566443368389627722e+02,        /* 0xC083EC88, 0x1375F228 */
164     rb5  = -1.02509513161107724954e+03,        /* 0xC0900461, 0x6A2E5992 */
165     rb6  = -4.83519191608651397019e+02,        /* 0xC07E384E, 0x9BDC383F */
166     sb1  =  3.03380607434824582924e+01,        /* 0x403E568B, 0x261D5190 */
167     sb2  =  3.25792512996573918826e+02,        /* 0x40745CAE, 0x221B9F0A */
168     sb3  =  1.53672958608443695994e+03,        /* 0x409802EB, 0x189D5118 */
169     sb4  =  3.19985821950859553908e+03,        /* 0x40A8FFB7, 0x688C246A */
170     sb5  =  2.55305040643316442583e+03,        /* 0x40A3F219, 0xCEDF3BE6 */
171     sb6  =  4.74528541206955367215e+02,        /* 0x407DA874, 0xE79FE763 */
172     sb7  = -2.24409524465858183362e+01;        /* 0xC03670E2, 0x42712D62 */
173
174 double gmx_erfd(double x)
175 {
176     gmx_int32_t hx, ix, i;
177     double      R, S, P, Q, s, y, z, r;
178
179     union
180     {
181         double d;
182         int    i[2];
183     }
184     conv;
185
186     conv.d = x;
187
188 #ifdef GMX_IEEE754_BIG_ENDIAN_WORD_ORDER
189     hx = conv.i[0];
190 #else
191     hx = conv.i[1];
192 #endif
193
194     ix = hx&0x7fffffff;
195     if (ix >= 0x7ff00000)
196     {
197         /* erf(nan)=nan */
198         i = ((gmx_uint32_t)hx>>31)<<1;
199         return (double)(1-i)+one/x; /* erf(+-inf)=+-1 */
200     }
201
202     if (ix < 0x3feb0000)
203     {
204         /* |x|<0.84375 */
205         if (ix < 0x3e300000)
206         {
207             /* |x|<2**-28 */
208             if (ix < 0x00800000)
209             {
210                 return 0.125*(8.0*x+efx8*x);  /*avoid underflow */
211             }
212             return x + efx*x;
213         }
214         z = x*x;
215         r = pp0+z*(pp1+z*(pp2+z*(pp3+z*pp4)));
216         s = one+z*(qq1+z*(qq2+z*(qq3+z*(qq4+z*qq5))));
217         y = r/s;
218         return x + x*y;
219     }
220     if (ix < 0x3ff40000)
221     {
222         /* 0.84375 <= |x| < 1.25 */
223         s = fabs(x)-one;
224         P = pa0+s*(pa1+s*(pa2+s*(pa3+s*(pa4+s*(pa5+s*pa6)))));
225         Q = one+s*(qa1+s*(qa2+s*(qa3+s*(qa4+s*(qa5+s*qa6)))));
226         if (hx >= 0)
227         {
228             return erx + P/Q;
229         }
230         else
231         {
232             return -erx - P/Q;
233         }
234     }
235     if (ix >= 0x40180000)
236     {
237         /* inf>|x|>=6 */
238         if (hx >= 0)
239         {
240             return one-tiny;
241         }
242         else
243         {
244             return tiny-one;
245         }
246     }
247     x = fabs(x);
248     s = one/(x*x);
249     if (ix < 0x4006DB6E)
250     {
251         /* |x| < 1/0.35 */
252         R = ra0+s*(ra1+s*(ra2+s*(ra3+s*(ra4+s*(ra5+s*(ra6+s*ra7))))));
253         S = one+s*(sa1+s*(sa2+s*(sa3+s*(sa4+s*(sa5+s*(sa6+s*(sa7+s*sa8)))))));
254     }
255     else
256     {
257         /* |x| >= 1/0.35 */
258         R = rb0+s*(rb1+s*(rb2+s*(rb3+s*(rb4+s*(rb5+s*rb6)))));
259         S = one+s*(sb1+s*(sb2+s*(sb3+s*(sb4+s*(sb5+s*(sb6+s*sb7))))));
260     }
261
262     conv.d = x;
263
264 #ifdef GMX_IEEE754_BIG_ENDIAN_WORD_ORDER
265     conv.i[1] = 0;
266 #else
267     conv.i[0] = 0;
268 #endif
269
270     z = conv.d;
271
272     r  =  exp(-z*z-0.5625)*exp((z-x)*(z+x)+R/S);
273     if (hx >= 0)
274     {
275         return one-r/x;
276     }
277     else
278     {
279         return r/x-one;
280     }
281 }
282
283
284 double gmx_erfcd(double x)
285 {
286     gmx_int32_t hx, ix;
287     double      R, S, P, Q, s, y, z, r;
288
289     union
290     {
291         double d;
292         int    i[2];
293     }
294     conv;
295
296     conv.d = x;
297
298 #ifdef GMX_IEEE754_BIG_ENDIAN_WORD_ORDER
299     hx = conv.i[0];
300 #else
301     hx = conv.i[1];
302 #endif
303
304     ix = hx&0x7fffffff;
305     if (ix >= 0x7ff00000)
306     {
307         /* erfc(nan)=nan */
308         /* erfc(+-inf)=0,2 */
309         return (double)(((gmx_uint32_t)hx>>31)<<1)+one/x;
310     }
311
312     if (ix < 0x3feb0000)
313     {
314         /* |x|<0.84375 */
315         double r1, r2, s1, s2, s3, z2, z4;
316         if (ix < 0x3c700000)     /* |x|<2**-56 */
317         {
318             return one-x;
319         }
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)
325         {
326             /* x<1/4 */
327             return one-(x+x*y);
328         }
329         else
330         {
331             r  = x*y;
332             r += (x-half);
333             return half - r;
334         }
335     }
336
337     if (ix < 0x3ff40000)
338     {
339         /* 0.84375 <= |x| < 1.25 */
340         s = fabs(x)-one;
341         P = pa0+s*(pa1+s*(pa2+s*(pa3+s*(pa4+s*(pa5+s*pa6)))));
342         Q = one+s*(qa1+s*(qa2+s*(qa3+s*(qa4+s*(qa5+s*qa6)))));
343         if (hx >= 0)
344         {
345             z  = one-erx; return z - P/Q;
346         }
347         else
348         {
349             z = erx+P/Q; return one+z;
350         }
351     }
352     if (ix < 0x403c0000)
353     {
354         /* |x|<28 */
355         x = fabs(x);
356         s = one/(x*x);
357         if (ix < 0x4006DB6D)
358         {
359             /* |x| < 1/.35 ~ 2.857143*/
360             R = ra0+s*(ra1+s*(ra2+s*(ra3+s*(ra4+s*(ra5+s*(ra6+s*ra7))))));
361             S = one+s*(sa1+s*(sa2+s*(sa3+s*(sa4+s*(sa5+s*(sa6+s*(sa7+s*sa8)))))));
362         }
363         else
364         {
365             /* |x| >= 1/.35 ~ 2.857143 */
366             if (hx < 0 && ix >= 0x40180000)
367             {
368                 return two-tiny; /* x < -6 */
369             }
370             R = rb0+s*(rb1+s*(rb2+s*(rb3+s*(rb4+s*(rb5+s*rb6)))));
371             S = one+s*(sb1+s*(sb2+s*(sb3+s*(sb4+s*(sb5+s*(sb6+s*sb7))))));
372         }
373
374         conv.d = x;
375
376 #ifdef GMX_IEEE754_BIG_ENDIAN_WORD_ORDER
377         conv.i[1] = 0;
378 #else
379         conv.i[0] = 0;
380 #endif
381
382         z = conv.d;
383
384         r  =  exp(-z*z-0.5625)*exp((z-x)*(z+x)+R/S);
385
386         if (hx > 0)
387         {
388             return r/x;
389         }
390         else
391         {
392             return two-r/x;
393         }
394     }
395     else
396     {
397         if (hx > 0)
398         {
399             return tiny*tiny;
400         }
401         else
402         {
403             return two-tiny;
404         }
405     }
406 }
407
408
409 static const float
410     tinyf =  1e-30,
411     halff =  5.0000000000e-01, /* 0x3F000000 */
412     onef  =  1.0000000000e+00, /* 0x3F800000 */
413     twof  =  2.0000000000e+00, /* 0x40000000 */
414 /* c = (subfloat)0.84506291151 */
415     erxf =  8.4506291151e-01,  /* 0x3f58560b */
416 /*
417  * Coefficients for approximation to  erf on [0,0.84375]
418  */
419     efxf  =  1.2837916613e-01, /* 0x3e0375d4 */
420     efx8f =  1.0270333290e+00, /* 0x3f8375d4 */
421     pp0f  =  1.2837916613e-01, /* 0x3e0375d4 */
422     pp1f  = -3.2504209876e-01, /* 0xbea66beb */
423     pp2f  = -2.8481749818e-02, /* 0xbce9528f */
424     pp3f  = -5.7702702470e-03, /* 0xbbbd1489 */
425     pp4f  = -2.3763017452e-05, /* 0xb7c756b1 */
426     qq1f  =  3.9791721106e-01, /* 0x3ecbbbce */
427     qq2f  =  6.5022252500e-02, /* 0x3d852a63 */
428     qq3f  =  5.0813062117e-03, /* 0x3ba68116 */
429     qq4f  =  1.3249473704e-04, /* 0x390aee49 */
430     qq5f  = -3.9602282413e-06, /* 0xb684e21a */
431 /*
432  * Coefficients for approximation to  erf  in [0.84375,1.25]
433  */
434     pa0f = -2.3621185683e-03,  /* 0xbb1acdc6 */
435     pa1f =  4.1485610604e-01,  /* 0x3ed46805 */
436     pa2f = -3.7220788002e-01,  /* 0xbebe9208 */
437     pa3f =  3.1834661961e-01,  /* 0x3ea2fe54 */
438     pa4f = -1.1089469492e-01,  /* 0xbde31cc2 */
439     pa5f =  3.5478305072e-02,  /* 0x3d1151b3 */
440     pa6f = -2.1663755178e-03,  /* 0xbb0df9c0 */
441     qa1f =  1.0642088205e-01,  /* 0x3dd9f331 */
442     qa2f =  5.4039794207e-01,  /* 0x3f0a5785 */
443     qa3f =  7.1828655899e-02,  /* 0x3d931ae7 */
444     qa4f =  1.2617121637e-01,  /* 0x3e013307 */
445     qa5f =  1.3637083583e-02,  /* 0x3c5f6e13 */
446     qa6f =  1.1984500103e-02,  /* 0x3c445aa3 */
447 /*
448  * Coefficients for approximation to  erfc in [1.25,1/0.35]
449  */
450     ra0f = -9.8649440333e-03,  /* 0xbc21a093 */
451     ra1f = -6.9385856390e-01,  /* 0xbf31a0b7 */
452     ra2f = -1.0558626175e+01,  /* 0xc128f022 */
453     ra3f = -6.2375331879e+01,  /* 0xc2798057 */
454     ra4f = -1.6239666748e+02,  /* 0xc322658c */
455     ra5f = -1.8460508728e+02,  /* 0xc3389ae7 */
456     ra6f = -8.1287437439e+01,  /* 0xc2a2932b */
457     ra7f = -9.8143291473e+00,  /* 0xc11d077e */
458     sa1f =  1.9651271820e+01,  /* 0x419d35ce */
459     sa2f =  1.3765776062e+02,  /* 0x4309a863 */
460     sa3f =  4.3456588745e+02,  /* 0x43d9486f */
461     sa4f =  6.4538726807e+02,  /* 0x442158c9 */
462     sa5f =  4.2900814819e+02,  /* 0x43d6810b */
463     sa6f =  1.0863500214e+02,  /* 0x42d9451f */
464     sa7f =  6.5702495575e+00,  /* 0x40d23f7c */
465     sa8f = -6.0424413532e-02,  /* 0xbd777f97 */
466 /*
467  * Coefficients for approximation to  erfc in [1/.35,28]
468  */
469     rb0f = -9.8649431020e-03,  /* 0xbc21a092 */
470     rb1f = -7.9928326607e-01,  /* 0xbf4c9dd4 */
471     rb2f = -1.7757955551e+01,  /* 0xc18e104b */
472     rb3f = -1.6063638306e+02,  /* 0xc320a2ea */
473     rb4f = -6.3756646729e+02,  /* 0xc41f6441 */
474     rb5f = -1.0250950928e+03,  /* 0xc480230b */
475     rb6f = -4.8351919556e+02,  /* 0xc3f1c275 */
476     sb1f =  3.0338060379e+01,  /* 0x41f2b459 */
477     sb2f =  3.2579251099e+02,  /* 0x43a2e571 */
478     sb3f =  1.5367296143e+03,  /* 0x44c01759 */
479     sb4f =  3.1998581543e+03,  /* 0x4547fdbb */
480     sb5f =  2.5530502930e+03,  /* 0x451f90ce */
481     sb6f =  4.7452853394e+02,  /* 0x43ed43a7 */
482     sb7f = -2.2440952301e+01;  /* 0xc1b38712 */
483
484
485 typedef union
486 {
487     float         value;
488     gmx_uint32_t  word;
489 } ieee_float_shape_type;
490
491 #define GET_FLOAT_WORD(i, d)                 \
492     do {                                \
493         ieee_float_shape_type gf_u;                   \
494         gf_u.value = (d);                     \
495         (i)        = gf_u.word;                      \
496     } while (0)
497
498
499 #define SET_FLOAT_WORD(d, i)                 \
500     do {                                \
501         ieee_float_shape_type sf_u;                   \
502         sf_u.word = (i);                      \
503         (d)       = sf_u.value;                     \
504     } while (0)
505
506
507 float gmx_erff(float x)
508 {
509     gmx_int32_t hx, ix, i;
510     float       R, S, P, Q, s, y, z, r;
511
512     union
513     {
514         float  f;
515         int    i;
516     }
517     conv;
518
519     conv.f = x;
520     hx     = conv.i;
521
522     ix = hx&0x7fffffff;
523     if (ix >= 0x7f800000)
524     {
525         /* erf(nan)=nan */
526         i = ((gmx_uint32_t)hx>>31)<<1;
527         return (float)(1-i)+onef/x; /* erf(+-inf)=+-1 */
528     }
529
530     if (ix < 0x3f580000)
531     {
532         /* |x|<0.84375 */
533         if (ix < 0x31800000)
534         {
535             /* |x|<2**-28 */
536             if (ix < 0x04000000)
537             {
538                 return (float)0.125*((float)8.0*x+efx8f*x);             /*avoid underflow */
539             }
540             return x + efxf*x;
541         }
542         z = x*x;
543         r = pp0f+z*(pp1f+z*(pp2f+z*(pp3f+z*pp4f)));
544         s = onef+z*(qq1f+z*(qq2f+z*(qq3f+z*(qq4f+z*qq5f))));
545         y = r/s;
546         return x + x*y;
547     }
548     if (ix < 0x3fa00000)
549     {
550         /* 0.84375 <= |x| < 1.25 */
551         s = fabs(x)-onef;
552         P = pa0f+s*(pa1f+s*(pa2f+s*(pa3f+s*(pa4f+s*(pa5f+s*pa6f)))));
553         Q = onef+s*(qa1f+s*(qa2f+s*(qa3f+s*(qa4f+s*(qa5f+s*qa6f)))));
554         if (hx >= 0)
555         {
556             return erxf + P/Q;
557         }
558         else
559         {
560             return -erxf - P/Q;
561         }
562     }
563     if (ix >= 0x40c00000)
564     {
565         /* inf>|x|>=6 */
566         if (hx >= 0)
567         {
568             return onef-tinyf;
569         }
570         else
571         {
572             return tinyf-onef;
573         }
574     }
575     x = fabs(x);
576     s = onef/(x*x);
577     if (ix < 0x4036DB6E)
578     {
579         /* |x| < 1/0.35 */
580         R = ra0f+s*(ra1f+s*(ra2f+s*(ra3f+s*(ra4f+s*(ra5f+s*(ra6f+s*ra7f))))));
581         S = onef+s*(sa1f+s*(sa2f+s*(sa3f+s*(sa4f+s*(sa5f+s*(sa6f+s*(sa7f+s*sa8f)))))));
582     }
583     else
584     {
585         /* |x| >= 1/0.35 */
586         R = rb0f+s*(rb1f+s*(rb2f+s*(rb3f+s*(rb4f+s*(rb5f+s*rb6f)))));
587         S = onef+s*(sb1f+s*(sb2f+s*(sb3f+s*(sb4f+s*(sb5f+s*(sb6f+s*sb7f))))));
588     }
589
590     conv.f = x;
591     conv.i = conv.i & 0xfffff000;
592     z      = conv.f;
593
594     r  =  exp(-z*z-(float)0.5625)*exp((z-x)*(z+x)+R/S);
595     if (hx >= 0)
596     {
597         return onef-r/x;
598     }
599     else
600     {
601         return r/x-onef;
602     }
603 }
604
605 float gmx_erfcf(float x)
606 {
607     gmx_int32_t hx, ix;
608     float       R, S, P, Q, s, y, z, r;
609
610     union
611     {
612         float  f;
613         int    i;
614     }
615     conv;
616
617     conv.f = x;
618     hx     = conv.i;
619
620     ix = hx&0x7fffffff;
621     if (ix >= 0x7f800000)
622     {
623         /* erfc(nan)=nan */
624         /* erfc(+-inf)=0,2 */
625         return (float)(((gmx_uint32_t)hx>>31)<<1)+onef/x;
626     }
627
628     if (ix < 0x3f580000)
629     {
630         /* |x|<0.84375 */
631         if (ix < 0x23800000)
632         {
633             return onef-x;  /* |x|<2**-56 */
634         }
635         z = x*x;
636         r = pp0f+z*(pp1f+z*(pp2f+z*(pp3f+z*pp4f)));
637         s = onef+z*(qq1f+z*(qq2f+z*(qq3f+z*(qq4f+z*qq5f))));
638         y = r/s;
639         if (hx < 0x3e800000)
640         {
641             /* x<1/4 */
642             return onef-(x+x*y);
643         }
644         else
645         {
646             r  = x*y;
647             r += (x-halff);
648             return halff - r;
649         }
650     }
651     if (ix < 0x3fa00000)
652     {
653         /* 0.84375 <= |x| < 1.25 */
654         s = fabs(x)-onef;
655         P = pa0f+s*(pa1f+s*(pa2f+s*(pa3f+s*(pa4f+s*(pa5f+s*pa6f)))));
656         Q = onef+s*(qa1f+s*(qa2f+s*(qa3f+s*(qa4f+s*(qa5f+s*qa6f)))));
657         if (hx >= 0)
658         {
659             z  = onef-erxf; return z - P/Q;
660         }
661         else
662         {
663             z = erxf+P/Q; return onef+z;
664         }
665     }
666     if (ix < 0x41e00000)
667     {
668         /* |x|<28 */
669         x = fabs(x);
670         s = onef/(x*x);
671         if (ix < 0x4036DB6D)
672         {
673             /* |x| < 1/.35 ~ 2.857143*/
674             R = ra0f+s*(ra1f+s*(ra2f+s*(ra3f+s*(ra4f+s*(ra5f+s*(ra6f+s*ra7f))))));
675             S = onef+s*(sa1f+s*(sa2f+s*(sa3f+s*(sa4f+s*(sa5f+s*(sa6f+s*(sa7f+s*sa8f)))))));
676         }
677         else
678         {
679             /* |x| >= 1/.35 ~ 2.857143 */
680             if (hx < 0 && ix >= 0x40c00000)
681             {
682                 return twof-tinyf;                     /* x < -6 */
683             }
684             R = rb0f+s*(rb1f+s*(rb2f+s*(rb3f+s*(rb4f+s*(rb5f+s*rb6f)))));
685             S = onef+s*(sb1f+s*(sb2f+s*(sb3f+s*(sb4f+s*(sb5f+s*(sb6f+s*sb7f))))));
686         }
687
688         conv.f = x;
689         conv.i = conv.i & 0xfffff000;
690         z      = conv.f;
691
692         r  =  exp(-z*z-(float)0.5625)*exp((z-x)*(z+x)+R/S);
693         if (hx > 0)
694         {
695             return r/x;
696         }
697         else
698         {
699             return twof-r/x;
700         }
701     }
702     else
703     {
704         if (hx > 0)
705         {
706             return tinyf*tinyf;
707         }
708         else
709         {
710             return twof-tinyf;
711         }
712     }
713 }
714
715
716 gmx_bool gmx_isfinite(real gmx_unused x)
717 {
718     gmx_bool returnval = TRUE;
719     /* If no suitable function was found, assume the value is
720      * finite. */
721
722 #ifdef HAVE__FINITE
723     returnval = _finite(x);
724 #elif defined HAVE_ISFINITE
725     returnval = isfinite(x);
726 #elif defined HAVE__ISFINITE
727     returnval = _isfinite(x);
728 #endif
729     return returnval;
730 }
731
732 gmx_bool gmx_isnan(real x)
733 {
734     return x != x;
735 }
736
737 int
738 gmx_within_tol(double   f1,
739                double   f2,
740                double   tol)
741 {
742     /* The or-equal is important - otherwise we return false if f1==f2==0 */
743     if (fabs(f1-f2) <= tol*0.5*(fabs(f1)+fabs(f2)) )
744     {
745         return 1;
746     }
747     else
748     {
749         return 0;
750     }
751 }
752
753 int
754 gmx_numzero(double a)
755 {
756     return gmx_within_tol(a, 0.0, GMX_REAL_MIN/GMX_REAL_EPS);
757 }
758
759 unsigned int
760 gmx_log2i(unsigned int n)
761 {
762     assert(n != 0); /* behavior differs for 0 */
763 #if defined(__INTEL_COMPILER)
764     return _bit_scan_reverse(n);
765 #elif defined(__GNUC__) && UINT_MAX == 4294967295U /*also for clang*/
766     return __builtin_clz(n) ^ 31U;                 /* xor gets optimized out */
767 #elif defined(_MSC_VER) && _MSC_VER >= 1400
768     {
769         unsigned long i;
770         _BitScanReverse(&i, n);
771         return i;
772     }
773 #elif defined(__xlC__)
774     return 31 - __cntlz4(n);
775 #else
776     /* http://graphics.stanford.edu/~seander/bithacks.html#IntegerLogLookup */
777 #define LT(n) n, n, n, n, n, n, n, n, n, n, n, n, n, n, n, n
778     static const char     LogTable256[256] = {
779         -1, 0, 1, 1, 2, 2, 2, 2, 3, 3, 3, 3, 3, 3, 3, 3,
780         LT(4), LT(5), LT(5), LT(6), LT(6), LT(6), LT(6),
781         LT(7), LT(7), LT(7), LT(7), LT(7), LT(7), LT(7), LT(7)
782     };
783 #undef LT
784
785     unsigned int r;     /* r will be lg(n) */
786     unsigned int t, tt; /* temporaries */
787
788     if ((tt = n >> 16) != 0)
789     {
790         r = ((t = tt >> 8) != 0) ? 24 + LogTable256[t] : 16 + LogTable256[tt];
791     }
792     else
793     {
794         r = ((t = n >> 8) != 0) ? 8 + LogTable256[t] : LogTable256[n];
795     }
796     return r;
797 #endif
798 }
799
800 gmx_bool
801 check_int_multiply_for_overflow(gmx_int64_t  a,
802                                 gmx_int64_t  b,
803                                 gmx_int64_t *result)
804 {
805     gmx_int64_t sign = 1;
806     if ((0 == a) || (0 == b))
807     {
808         *result = 0;
809         return TRUE;
810     }
811     if (a < 0)
812     {
813         a    = -a;
814         sign = -sign;
815     }
816     if (b < 0)
817     {
818         b    = -b;
819         sign = -sign;
820     }
821     if (GMX_INT64_MAX / b < a)
822     {
823         *result = (sign > 0) ? GMX_INT64_MAX : GMX_INT64_MIN;
824         return FALSE;
825     }
826     *result = sign * a * b;
827     return TRUE;
828 }
829
830 int gmx_greatest_common_divisor(int p, int q)
831 {
832     int tmp;
833     while (q != 0)
834     {
835         tmp = q;
836         q   = p % q;
837         p   = tmp;
838     }
839     return p;
840 }