Disable instruction fusion on Power8
[gromacs.git] / src / gromacs / math / utilities.h
blobde262687267b7c536083ff7bfc3a87e4ddb69482
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017, 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.
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.
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.
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.
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.
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.
37 #ifndef GMX_MATH_UTILITIES_H
38 #define GMX_MATH_UTILITIES_H
40 #include <limits.h>
42 #include <cmath>
44 #include "gromacs/utility/basedefinitions.h"
45 #include "gromacs/utility/real.h"
47 #ifdef __cplusplus
48 extern "C" {
49 #endif
51 #ifndef M_PI
52 #define M_PI 3.14159265358979323846
53 #endif
55 #ifndef M_PI_2
56 #define M_PI_2 1.57079632679489661923
57 #endif
59 #ifndef M_2PI
60 #define M_2PI 6.28318530717958647692
61 #endif
63 #ifndef M_SQRT2
64 #define M_SQRT2 sqrt(2.0)
65 #endif
67 #ifndef M_1_PI
68 #define M_1_PI 0.31830988618379067154
69 #endif
71 #ifndef M_FLOAT_1_SQRTPI /* used in GPU kernels */
72 /* 1.0 / sqrt(M_PI) */
73 #define M_FLOAT_1_SQRTPI 0.564189583547756f
74 #endif
76 #ifndef M_1_SQRTPI
77 /* 1.0 / sqrt(M_PI) */
78 #define M_1_SQRTPI 0.564189583547756
79 #endif
81 #ifndef M_2_SQRTPI
82 /* 2.0 / sqrt(M_PI) */
83 #define M_2_SQRTPI 1.128379167095513
84 #endif
86 /*! \brief Enum to select safe or highly unsafe (faster) math functions.
88 * Normally all the Gromacs math functions should apply reasonable care with
89 * input arguments. While we do not necessarily adhere strictly to IEEE
90 * (in particular not for arguments that might result in NaN, inf, etc.), the
91 * functions should return reasonable values or e.g. clamp results to zero.
93 * However, in a few cases where we are extremely performance-sensitive it
94 * makes sense to forego these checks too in cases where we know the exact
95 * properties if the input data, and we really need to save every cycle we can.
97 * This class is typically used as a template parameter to such calls to enable
98 * the caller to select the level of aggressiveness. We should always use the
99 * safe alternative as the default value, and document carefully what might
100 * happen with the unsafe alternative.
102 enum class MathOptimization
104 Safe, //!< Don't do unsafe optimizations. This should always be default.
105 Unsafe //!< Allow optimizations that can be VERY dangerous for general code.
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.
129 gmx_within_tol(double f1,
130 double f2,
131 double tol);
134 * \brief Check if a number is smaller than some preset safe minimum
135 * value, currently defined as GMX_REAL_MIN/GMX_REAL_EPS.
137 * If a number is smaller than this value we risk numerical overflow
138 * if any number larger than 1.0/GMX_REAL_EPS is divided by it.
140 * \return 1 if 'almost' numerically zero, 0 otherwise.
143 gmx_numzero(double a);
145 /*! \brief Multiply two large ints
147 * \return False iff overflow occurred
149 gmx_bool
150 check_int_multiply_for_overflow(gmx_int64_t a,
151 gmx_int64_t b,
152 gmx_int64_t *result);
154 /*! \brief Find greatest common divisor of two numbers
156 * \return GCD of the two inputs
159 gmx_greatest_common_divisor(int p, int q);
162 /*! \brief Enable floating-point exceptions if supported on OS
164 * Enables division-by-zero, invalid, and overflow.
166 int gmx_feenableexcept();
168 /*! \brief Return cut-off to use
170 * Takes the max of two cut-offs. However a cut-off of 0
171 * signifies that the cut-off in fact is infinite, and
172 * this requires this special routine.
173 * \param[in] cutoff1 The first cutoff (e.g. coulomb)
174 * \param[in] cutoff2 The second cutoff (e.g. vdw)
175 * \return 0 if either is 0, the normal max of the two otherwise.
177 real max_cutoff(real cutoff1, real cutoff2);
179 #ifdef __cplusplus
181 #endif
183 #endif