2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017,2018, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
39 #include <gtest/gtest.h>
41 #include "gromacs/math/utilities.h"
42 #include "gromacs/simd/simd.h"
43 #include "gromacs/utility/basedefinitions.h"
45 #include "testutils/testasserts.h"
58 /*! \addtogroup module_simd */
61 // Since these functions are mostly wrappers for similar standard library
62 // functions, we typically just test 1-2 values to make sure we call the right
65 /******************************************
66 * Default floating-point precision tests *
67 ******************************************/
69 TEST(SimdScalarMathTest
, copysign
)
71 EXPECT_EQ(real(-c1
), copysign(real(c1
), real(-c2
)));
72 EXPECT_EQ(real(c2
), copysign(real(c2
), real(c3
)));
75 TEST(SimdScalarMathTest
, invsqrtPair
)
82 invsqrtPair(x0
, x1
, &out0
, &out1
);
84 EXPECT_EQ(invsqrt(x0
), out0
);
85 EXPECT_EQ(invsqrt(x1
), out1
);
88 TEST(SimdScalarMathTest
, inv
)
92 EXPECT_EQ(real(1.0)/x0
, inv(x0
));
95 TEST(SimdScalarMathTest
, maskzInvsqrt
)
99 EXPECT_EQ(invsqrt(x0
), maskzInvsqrt(x0
, true));
100 EXPECT_EQ(real(0), maskzInvsqrt(x0
, false));
103 TEST(SimdScalarMathTest
, log
)
107 EXPECT_EQ(std::log(x0
), log(x0
));
110 TEST(SimdScalarMathTest
, exp2
)
114 EXPECT_EQ(std::exp2(x0
), exp2(x0
));
117 TEST(SimdScalarMathTest
, exp
)
121 EXPECT_EQ(std::exp(x0
), exp(x0
));
124 TEST(SimdScalarMathTest
, erf
)
128 EXPECT_EQ(std::erf(x0
), erf(x0
));
131 TEST(SimdScalarMathTest
, erfc
)
135 EXPECT_EQ(std::erfc(x0
), erfc(x0
));
138 TEST(SimdScalarMathTest
, sincos
)
145 EXPECT_EQ(std::sin(x0
), s
);
146 EXPECT_EQ(std::cos(x0
), c
);
149 TEST(SimdScalarMathTest
, sin
)
153 EXPECT_EQ(std::sin(x0
), sin(x0
));
156 TEST(SimdScalarMathTest
, cos
)
160 EXPECT_EQ(std::cos(x0
), cos(x0
));
163 TEST(SimdScalarMathTest
, tan
)
167 EXPECT_EQ(std::tan(x0
), tan(x0
));
171 TEST(SimdScalarMathTest
, asin
)
175 EXPECT_EQ(std::asin(x0
), asin(x0
));
178 TEST(SimdScalarMathTest
, acos
)
182 EXPECT_EQ(std::acos(x0
), acos(x0
));
185 TEST(SimdScalarMathTest
, atan
)
189 EXPECT_EQ(std::atan(x0
), atan(x0
));
192 TEST(SimdScalarMathTest
, atan2
)
195 real y
= std::sqrt(c0
);
198 EXPECT_EQ(std::atan2(y
, x
), atan2(y
, x
));
201 TEST(SimdScalarMathTest
, pmeForceCorrection
)
205 // Calculate reference value for z2!=0
206 real z
= std::sqrt(z2
);
207 real ref
= 2.0*std::exp(-z2
)/(std::sqrt(M_PI
)*z2
) - std::erf(z
)/(z2
*z
);
209 // Pme correction only needs to be ~1e-6 accuracy single, 1e-10 double
211 FloatingPointTolerance
tolerance(relativeToleranceAsFloatingPoint(1.0, 5e-10));
213 FloatingPointTolerance
tolerance(relativeToleranceAsFloatingPoint(1.0, 5e-6));
216 EXPECT_REAL_EQ_TOL(ref
, pmeForceCorrection(z2
), tolerance
);
219 TEST(SimdScalarMathTest
, pmePotentialCorrection
)
223 // Calculate reference value for z2!=0
224 real z
= std::sqrt(z2
);
225 real ref
= std::erf(z
)/z
;
227 // Pme correction only needs to be ~1e-6 accuracy single, 1e-10 double
229 FloatingPointTolerance
tolerance(relativeToleranceAsFloatingPoint(1.0, 5e-10));
231 FloatingPointTolerance
tolerance(relativeToleranceAsFloatingPoint(1.0, 5e-6));
234 EXPECT_REAL_EQ_TOL(ref
, pmePotentialCorrection(z2
), tolerance
);
237 /*******************************************
238 * Double precision, single accuracy tests *
239 *******************************************/
241 TEST(SimdScalarMathTest
, invsqrtPairSingleAccuracy
)
248 invsqrtPairSingleAccuracy(x0
, x1
, &out0
, &out1
);
250 EXPECT_EQ(invsqrt(static_cast<float>(x0
)), static_cast<float>(out0
));
251 EXPECT_EQ(invsqrt(static_cast<float>(x1
)), static_cast<float>(out1
));
254 TEST(SimdScalarMathTest
, invSingleAccuracy
)
258 EXPECT_EQ(1.0f
/static_cast<float>(x0
),
259 static_cast<float>(invSingleAccuracy(x0
)));
262 TEST(SimdScalarMathTest
, maskzInvsqrtSingleAccuracy
)
266 EXPECT_EQ(invsqrt(static_cast<float>(x0
)),
267 static_cast<float>(maskzInvsqrtSingleAccuracy(x0
, true)));
269 static_cast<float>(maskzInvsqrtSingleAccuracy(x0
, false)));
272 TEST(SimdScalarMathTest
, logSingleAccuracy
)
276 EXPECT_EQ(std::log(static_cast<float>(x0
)),
277 static_cast<float>(logSingleAccuracy(x0
)));
280 TEST(SimdScalarMathTest
, exp2SingleAccuracy
)
284 EXPECT_EQ(std::exp2(static_cast<float>(x0
)),
285 static_cast<float>(exp2SingleAccuracy(x0
)));
288 TEST(SimdScalarMathTest
, expSingleAccuracy
)
292 EXPECT_EQ(std::exp(static_cast<float>(x0
)),
293 static_cast<float>(expSingleAccuracy(x0
)));
296 TEST(SimdScalarMathTest
, erfSingleAccuracy
)
300 EXPECT_EQ(std::erf(static_cast<float>(x0
)), static_cast<float>(erfSingleAccuracy(x0
)));
303 TEST(SimdScalarMathTest
, erfcSingleAccuracy
)
307 EXPECT_EQ(std::erfc(static_cast<float>(x0
)), static_cast<float>(erfcSingleAccuracy(x0
)));
310 TEST(SimdScalarMathTest
, sincosSingleAccuracy
)
315 sincosSingleAccuracy(x0
, &s
, &c
);
317 EXPECT_EQ(std::sin(static_cast<float>(x0
)), static_cast<float>(s
));
318 EXPECT_EQ(std::cos(static_cast<float>(x0
)), static_cast<float>(c
));
321 TEST(SimdScalarMathTest
, sinSingleAccuracy
)
325 EXPECT_EQ(std::sin(static_cast<float>(x0
)),
326 static_cast<float>(sinSingleAccuracy(x0
)));
329 TEST(SimdScalarMathTest
, cosSingleAccuracy
)
333 EXPECT_EQ(std::cos(static_cast<float>(x0
)),
334 static_cast<float>(cosSingleAccuracy(x0
)));
337 TEST(SimdScalarMathTest
, tanSingleAccuracy
)
341 EXPECT_EQ(std::tan(static_cast<float>(x0
)),
342 static_cast<float>(tanSingleAccuracy(x0
)));
346 TEST(SimdScalarMathTest
, asinSingleAccuracy
)
350 EXPECT_EQ(std::asin(static_cast<float>(x0
)),
351 static_cast<float>(asinSingleAccuracy(x0
)));
354 TEST(SimdScalarMathTest
, acosSingleAccuracy
)
358 EXPECT_EQ(std::acos(static_cast<float>(x0
)),
359 static_cast<float>(acosSingleAccuracy(x0
)));
362 TEST(SimdScalarMathTest
, atanSingleAccuracy
)
366 EXPECT_EQ(std::atan(static_cast<float>(x0
)),
367 static_cast<float>(atanSingleAccuracy(x0
)));
370 TEST(SimdScalarMathTest
, atan2SingleAccuracy
)
373 double y
= std::sqrt(c0
);
376 EXPECT_EQ(std::atan2(static_cast<float>(y
), static_cast<float>(x
)),
377 static_cast<float>(atan2SingleAccuracy(y
, x
)));
380 TEST(SimdScalarMathTest
, pmeForceCorrectionSingleAccuracy
)
384 // Calculate reference value for z2!=0 in single precision
385 float z
= std::sqrt(static_cast<float>(z2
));
386 float ref
= 2.0*std::exp(static_cast<float>(-z2
))/(std::sqrt(static_cast<float>(M_PI
))*z2
) - std::erf(z
)/(z2
*z
);
388 // Pme correction only needs to be ~1e-6 accuracy single
389 FloatingPointTolerance
tolerance(relativeToleranceAsFloatingPoint(1.0, 5e-6));
391 EXPECT_REAL_EQ_TOL(ref
, static_cast<float>(pmeForceCorrectionSingleAccuracy(z2
)), tolerance
);
394 TEST(SimdScalarMathTest
, pmePotentialCorrectionSingleAccuracy
)
398 // Calculate reference value for z2!=0 in single precision
399 float z
= std::sqrt(static_cast<float>(z2
));
400 float ref
= std::erf(z
)/z
;
402 // Pme correction only needs to be ~1e-6 accuracy single
403 FloatingPointTolerance
tolerance(relativeToleranceAsFloatingPoint(1.0, 5e-6));
405 EXPECT_REAL_EQ_TOL(ref
, static_cast<float>(pmePotentialCorrectionSingleAccuracy(z2
)), tolerance
);
409 /*! \endcond internal */