Bumped version to 4.5-beta3
[gromacs.git] / include / maths.h
blobfc95bfb66230325759d84962a6bcaec71daa61d4
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * Gromacs Runs On Most of All Computer Systems
37 #ifndef _maths_h
38 #define _maths_h
40 #ifdef HAVE_CONFIG_H
41 #include <config.h>
42 #endif
44 #include <math.h>
45 #include "types/simple.h"
46 #include "typedefs.h"
48 #ifdef __cplusplus
49 extern "C" {
50 #endif
52 #ifndef M_PI
53 #define M_PI 3.14159265358979323846
54 #endif
56 #ifndef M_PI_2
57 #define M_PI_2 1.57079632679489661923
58 #endif
60 #ifndef M_2PI
61 #define M_2PI 6.28318530717958647692
62 #endif
64 #ifndef M_SQRT2
65 #define M_SQRT2 sqrt(2.0)
66 #endif
68 /* Suzuki-Yoshida Constants, for n=3 and n=5, for symplectic integration */
69 /* for n=1, w0 = 1 */
70 /* for n=3, w0 = w2 = 1/(2-2^-(1/3)), w1 = 1-2*w0 */
71 /* for n=5, w0 = w1 = w3 = w4 = 1/(4-4^-(1/3)), w1 = 1-4*w0 */
73 #define MAX_SUZUKI_YOSHIDA_NUM 5
74 #define SUZUKI_YOSHIDA_NUM 5
76 static const double sy_const_1[] = { 1. };
77 static const double sy_const_3[] = { 0.828981543588751,-0.657963087177502,0.828981543588751 };
78 static const double sy_const_5[] = { 0.2967324292201065,0.2967324292201065,-0.186929716880426,0.2967324292201065,0.2967324292201065 };
80 static const double* sy_const[] = {
81 NULL,
82 sy_const_1,
83 NULL,
84 sy_const_3,
85 NULL,
86 sy_const_5
90 static const double sy_const[MAX_SUZUKI_YOSHIDA_NUM+1][MAX_SUZUKI_YOSHIDA_NUM+1] = {
91 {},
92 {1},
93 {},
94 {0.828981543588751,-0.657963087177502,0.828981543588751},
95 {},
96 {0.2967324292201065,0.2967324292201065,-0.186929716880426,0.2967324292201065,0.2967324292201065}
97 };*/
99 extern int gmx_nint(real a);
100 extern real sign(real x,real y);
102 extern int gmx_nint(real a);
103 extern real sign(real x,real y);
104 extern real cuberoot (real a);
105 extern real gmx_erf(real x);
106 extern real gmx_erfc(real x);
108 /*! \brief Check if two numbers are within a tolerance
110 * This routine checks if the relative difference between two numbers is
111 * approximately within the given tolerance, defined as
112 * fabs(f1-f2)<=tolerance*fabs(f1+f2).
114 * To check if two floating-point numbers are almost identical, use this routine
115 * with the tolerance GMX_REAL_EPS, or GMX_DOUBLE_EPS if the check should be
116 * done in double regardless of Gromacs precision.
118 * To check if two algorithms produce similar results you will normally need
119 * to relax the tolerance significantly since many operations (e.g. summation)
120 * accumulate floating point errors.
122 * \param f1 First number to compare
123 * \param f2 Second number to compare
124 * \param tol Tolerance to use
126 * \return 1 if the relative difference is within tolerance, 0 if not.
128 static int
129 gmx_within_tol(double f1,
130 double f2,
131 double tol)
133 /* The or-equal is important - otherwise we return false if f1==f2==0 */
134 if( fabs(f1-f2) <= tol*0.5*(fabs(f1)+fabs(f2)) )
136 return 1;
138 else
140 return 0;
146 /**
147 * Check if a number is smaller than some preset safe minimum
148 * value, currently defined as GMX_REAL_MIN/GMX_REAL_EPS.
150 * If a number is smaller than this value we risk numerical overflow
151 * if any number larger than 1.0/GMX_REAL_EPS is divided by it.
153 * \return 1 if 'almost' numerically zero, 0 otherwise.
155 static int
156 gmx_numzero(double a)
158 return gmx_within_tol(a,0.0,GMX_REAL_MIN/GMX_REAL_EPS);
162 static real
163 gmx_log2(real x)
165 const real iclog2 = 1.0/log( 2.0 );
167 return log( x ) * iclog2;
171 #ifdef __cplusplus
173 #endif
175 #endif /* _maths_h */