Update instructions in containers.rst
[gromacs.git] / src / gromacs / simd / tests / scalar_math.cpp
blob048386b148d573350e587da53145acecbcd011d3
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017,2018 by the GROMACS development team.
5 * Copyright (c) 2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
36 #include "gmxpre.h"
38 #include <cmath>
40 #include <gtest/gtest.h>
42 #include "gromacs/math/utilities.h"
43 #include "gromacs/simd/simd.h"
44 #include "gromacs/utility/basedefinitions.h"
46 #include "testutils/testasserts.h"
48 #include "data.h"
50 namespace gmx
52 namespace test
55 namespace
58 /*! \cond internal */
59 /*! \addtogroup module_simd */
60 /*! \{ */
62 // Since these functions are mostly wrappers for similar standard library
63 // functions, we typically just test 1-2 values to make sure we call the right
64 // function.
66 /******************************************
67 * Default floating-point precision tests *
68 ******************************************/
70 TEST(SimdScalarMathTest, copysign)
72 EXPECT_EQ(real(-c1), copysign(real(c1), real(-c2)));
73 EXPECT_EQ(real(c2), copysign(real(c2), real(c3)));
76 TEST(SimdScalarMathTest, invsqrtPair)
78 real x0 = c1;
79 real x1 = c2;
81 real out0, out1;
83 invsqrtPair(x0, x1, &out0, &out1);
85 EXPECT_EQ(invsqrt(x0), out0);
86 EXPECT_EQ(invsqrt(x1), out1);
89 TEST(SimdScalarMathTest, inv)
91 real x0 = c0;
93 EXPECT_EQ(real(1.0) / x0, inv(x0));
96 TEST(SimdScalarMathTest, maskzInvsqrt)
98 real x0 = c0;
100 EXPECT_EQ(invsqrt(x0), maskzInvsqrt(x0, true));
101 EXPECT_EQ(real(0), maskzInvsqrt(x0, false));
104 TEST(SimdScalarMathTest, log)
106 real x0 = c0;
108 EXPECT_EQ(std::log(x0), log(x0));
111 TEST(SimdScalarMathTest, exp2)
113 real x0 = c0;
115 EXPECT_EQ(std::exp2(x0), exp2(x0));
118 TEST(SimdScalarMathTest, exp)
120 real x0 = c0;
122 EXPECT_EQ(std::exp(x0), exp(x0));
125 TEST(SimdScalarMathTest, erf)
127 real x0 = c0;
129 EXPECT_EQ(std::erf(x0), erf(x0));
132 TEST(SimdScalarMathTest, erfc)
134 real x0 = c0;
136 EXPECT_EQ(std::erfc(x0), erfc(x0));
139 TEST(SimdScalarMathTest, sincos)
141 real x0 = c0;
142 real s, c;
144 sincos(x0, &s, &c);
146 EXPECT_EQ(std::sin(x0), s);
147 EXPECT_EQ(std::cos(x0), c);
150 TEST(SimdScalarMathTest, sin)
152 real x0 = c0;
154 EXPECT_EQ(std::sin(x0), sin(x0));
157 TEST(SimdScalarMathTest, cos)
159 real x0 = c0;
161 EXPECT_EQ(std::cos(x0), cos(x0));
164 TEST(SimdScalarMathTest, tan)
166 real x0 = c0;
168 EXPECT_EQ(std::tan(x0), tan(x0));
172 TEST(SimdScalarMathTest, asin)
174 real x0 = c0;
176 EXPECT_EQ(std::asin(x0), asin(x0));
179 TEST(SimdScalarMathTest, acos)
181 real x0 = c0;
183 EXPECT_EQ(std::acos(x0), acos(x0));
186 TEST(SimdScalarMathTest, atan)
188 real x0 = c0;
190 EXPECT_EQ(std::atan(x0), atan(x0));
193 TEST(SimdScalarMathTest, atan2)
195 real x = c0;
196 real y = std::sqrt(c0);
199 EXPECT_EQ(std::atan2(y, x), atan2(y, x));
202 TEST(SimdScalarMathTest, pmeForceCorrection)
204 real z2 = c0;
206 // Calculate reference value for z2!=0
207 real z = std::sqrt(z2);
208 real ref = 2.0 * std::exp(-z2) / (std::sqrt(M_PI) * z2) - std::erf(z) / (z2 * z);
210 // Pme correction only needs to be ~1e-6 accuracy single, 1e-10 double
211 #if GMX_DOUBLE
212 FloatingPointTolerance tolerance(relativeToleranceAsFloatingPoint(1.0, 5e-10));
213 #else
214 FloatingPointTolerance tolerance(relativeToleranceAsFloatingPoint(1.0, 5e-6));
215 #endif
217 EXPECT_REAL_EQ_TOL(ref, pmeForceCorrection(z2), tolerance);
220 TEST(SimdScalarMathTest, pmePotentialCorrection)
222 real z2 = c0;
224 // Calculate reference value for z2!=0
225 real z = std::sqrt(z2);
226 real ref = std::erf(z) / z;
228 // Pme correction only needs to be ~1e-6 accuracy single, 1e-10 double
229 #if GMX_DOUBLE
230 FloatingPointTolerance tolerance(relativeToleranceAsFloatingPoint(1.0, 5e-10));
231 #else
232 FloatingPointTolerance tolerance(relativeToleranceAsFloatingPoint(1.0, 5e-6));
233 #endif
235 EXPECT_REAL_EQ_TOL(ref, pmePotentialCorrection(z2), tolerance);
238 /*******************************************
239 * Double precision, single accuracy tests *
240 *******************************************/
242 TEST(SimdScalarMathTest, invsqrtPairSingleAccuracy)
244 double x0 = c1;
245 double x1 = c2;
247 double out0, out1;
249 invsqrtPairSingleAccuracy(x0, x1, &out0, &out1);
251 EXPECT_EQ(invsqrt(static_cast<float>(x0)), static_cast<float>(out0));
252 EXPECT_EQ(invsqrt(static_cast<float>(x1)), static_cast<float>(out1));
255 TEST(SimdScalarMathTest, invSingleAccuracy)
257 double x0 = c1;
259 EXPECT_EQ(1.0F / static_cast<float>(x0), static_cast<float>(invSingleAccuracy(x0)));
262 TEST(SimdScalarMathTest, maskzInvsqrtSingleAccuracy)
264 double x0 = c1;
266 EXPECT_EQ(invsqrt(static_cast<float>(x0)), static_cast<float>(maskzInvsqrtSingleAccuracy(x0, true)));
267 EXPECT_EQ(0.0F, static_cast<float>(maskzInvsqrtSingleAccuracy(x0, false)));
270 TEST(SimdScalarMathTest, logSingleAccuracy)
272 double x0 = c1;
274 EXPECT_EQ(std::log(static_cast<float>(x0)), static_cast<float>(logSingleAccuracy(x0)));
277 TEST(SimdScalarMathTest, exp2SingleAccuracy)
279 double x0 = c1;
281 EXPECT_EQ(std::exp2(static_cast<float>(x0)), static_cast<float>(exp2SingleAccuracy(x0)));
284 TEST(SimdScalarMathTest, expSingleAccuracy)
286 double x0 = c1;
288 EXPECT_EQ(std::exp(static_cast<float>(x0)), static_cast<float>(expSingleAccuracy(x0)));
291 TEST(SimdScalarMathTest, erfSingleAccuracy)
293 double x0 = c0;
295 EXPECT_EQ(std::erf(static_cast<float>(x0)), static_cast<float>(erfSingleAccuracy(x0)));
298 TEST(SimdScalarMathTest, erfcSingleAccuracy)
300 double x0 = c0;
302 EXPECT_EQ(std::erfc(static_cast<float>(x0)), static_cast<float>(erfcSingleAccuracy(x0)));
305 TEST(SimdScalarMathTest, sincosSingleAccuracy)
307 double x0 = c0;
308 double s, c;
310 sincosSingleAccuracy(x0, &s, &c);
312 EXPECT_EQ(std::sin(static_cast<float>(x0)), static_cast<float>(s));
313 EXPECT_EQ(std::cos(static_cast<float>(x0)), static_cast<float>(c));
316 TEST(SimdScalarMathTest, sinSingleAccuracy)
318 double x0 = c0;
320 EXPECT_EQ(std::sin(static_cast<float>(x0)), static_cast<float>(sinSingleAccuracy(x0)));
323 TEST(SimdScalarMathTest, cosSingleAccuracy)
325 double x0 = c0;
327 EXPECT_EQ(std::cos(static_cast<float>(x0)), static_cast<float>(cosSingleAccuracy(x0)));
330 TEST(SimdScalarMathTest, tanSingleAccuracy)
332 double x0 = c0;
334 EXPECT_EQ(std::tan(static_cast<float>(x0)), static_cast<float>(tanSingleAccuracy(x0)));
338 TEST(SimdScalarMathTest, asinSingleAccuracy)
340 double x0 = c0;
342 EXPECT_EQ(std::asin(static_cast<float>(x0)), static_cast<float>(asinSingleAccuracy(x0)));
345 TEST(SimdScalarMathTest, acosSingleAccuracy)
347 double x0 = c0;
349 EXPECT_EQ(std::acos(static_cast<float>(x0)), static_cast<float>(acosSingleAccuracy(x0)));
352 TEST(SimdScalarMathTest, atanSingleAccuracy)
354 double x0 = c0;
356 EXPECT_EQ(std::atan(static_cast<float>(x0)), static_cast<float>(atanSingleAccuracy(x0)));
359 TEST(SimdScalarMathTest, atan2SingleAccuracy)
361 double x = c0;
362 double y = std::sqrt(c0);
365 EXPECT_EQ(std::atan2(static_cast<float>(y), static_cast<float>(x)),
366 static_cast<float>(atan2SingleAccuracy(y, x)));
369 TEST(SimdScalarMathTest, pmeForceCorrectionSingleAccuracy)
371 double z2 = c0;
373 // Calculate reference value for z2!=0 in single precision
374 float z = std::sqrt(static_cast<float>(z2));
375 float ref = 2.0 * std::exp(static_cast<float>(-z2)) / (std::sqrt(static_cast<float>(M_PI)) * z2)
376 - std::erf(z) / (z2 * z);
378 // Pme correction only needs to be ~1e-6 accuracy single
379 FloatingPointTolerance tolerance(relativeToleranceAsFloatingPoint(1.0, 5e-6));
381 EXPECT_REAL_EQ_TOL(ref, static_cast<float>(pmeForceCorrectionSingleAccuracy(z2)), tolerance);
384 TEST(SimdScalarMathTest, pmePotentialCorrectionSingleAccuracy)
386 double z2 = c0;
388 // Calculate reference value for z2!=0 in single precision
389 float z = std::sqrt(static_cast<float>(z2));
390 float ref = std::erf(z) / z;
392 // Pme correction only needs to be ~1e-6 accuracy single
393 FloatingPointTolerance tolerance(relativeToleranceAsFloatingPoint(1.0, 5e-6));
395 EXPECT_REAL_EQ_TOL(ref, static_cast<float>(pmePotentialCorrectionSingleAccuracy(z2)), tolerance);
398 /*! \} */
399 /*! \endcond internal */
401 } // namespace
402 } // namespace test
403 } // namespace gmx