Improved SIMD test data to use all bits
[gromacs.git] / src / gromacs / simd / tests / simd4_floatingpoint.cpp
blobb2300ec8be733abe5c03a77103b17332224c47ad
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017, 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.
35 #include "gmxpre.h"
37 #include <cmath>
39 #include "gromacs/math/utilities.h"
40 #include "gromacs/simd/simd.h"
41 #include "gromacs/utility/basedefinitions.h"
43 #include "testutils/testasserts.h"
45 #include "data.h"
46 #include "simd4.h"
48 #if GMX_SIMD
50 namespace gmx
52 namespace test
54 namespace
57 /*! \cond internal */
58 /*! \addtogroup module_simd */
59 /*! \{ */
61 #if GMX_SIMD4_HAVE_REAL
63 /*! \brief Test fixture for SIMD4 floating-point operations (identical to the SIMD4 \ref Simd4Test) */
64 typedef Simd4Test Simd4FloatingpointTest;
66 TEST_F(Simd4FloatingpointTest, setZero)
68 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom1R(0.0), setZero());
71 TEST_F(Simd4FloatingpointTest, set)
73 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom1R(c1), Simd4Real(c1));
76 TEST_F(Simd4FloatingpointTest, add)
78 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(c0 + c3, c1 + c4, c2 + c5 ),
79 rSimd4_c0c1c2 + rSimd4_c3c4c5);
83 TEST_F(Simd4FloatingpointTest, sub)
85 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(c0 - c3, c1 - c4, c2 - c5 ),
86 rSimd4_c0c1c2 - rSimd4_c3c4c5);
89 TEST_F(Simd4FloatingpointTest, mul)
91 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(c0 * c3, c1 * c4, c2 * c5 ),
92 rSimd4_c0c1c2 * rSimd4_c3c4c5);
95 TEST_F(Simd4FloatingpointTest, fma)
97 // The last bit of FMA operations depends on hardware, so we don't require exact match
98 GMX_EXPECT_SIMD4_REAL_NEAR(setSimd4RealFrom3R(c0 * c3 + c6, c1 * c4 + c7, c2 * c5 + c8),
99 fma(rSimd4_c0c1c2, rSimd4_c3c4c5, rSimd4_c6c7c8));
102 TEST_F(Simd4FloatingpointTest, fms)
104 // The last bit of FMA operations depends on hardware, so we don't require exact match
105 GMX_EXPECT_SIMD4_REAL_NEAR(setSimd4RealFrom3R(c0 * c3 - c6, c1 * c4 - c7, c2 * c5 - c8),
106 fms(rSimd4_c0c1c2, rSimd4_c3c4c5, rSimd4_c6c7c8));
109 TEST_F(Simd4FloatingpointTest, fnma)
111 // The last bit of FMA operations depends on hardware, so we don't require exact match
112 GMX_EXPECT_SIMD4_REAL_NEAR(setSimd4RealFrom3R(c6 - c0 * c3, c7 - c1 * c4, c8 - c2 * c5),
113 fnma(rSimd4_c0c1c2, rSimd4_c3c4c5, rSimd4_c6c7c8));
116 TEST_F(Simd4FloatingpointTest, fnms)
118 // The last bit of FMA operations depends on hardware, so we don't require exact match
119 GMX_EXPECT_SIMD4_REAL_NEAR(setSimd4RealFrom3R(-c0 * c3 - c6, -c1 * c4 - c7, -c2 * c5 - c8),
120 fnms(rSimd4_c0c1c2, rSimd4_c3c4c5, rSimd4_c6c7c8));
123 TEST_F(Simd4FloatingpointTest, abs)
125 GMX_EXPECT_SIMD4_REAL_EQ(rSimd4_c0c1c2, abs(rSimd4_c0c1c2)); // fabs(x)=x
126 GMX_EXPECT_SIMD4_REAL_EQ(rSimd4_c0c1c2, abs(rSimd4_m0m1m2)); // fabs(-x)=x
129 TEST_F(Simd4FloatingpointTest, neg)
131 GMX_EXPECT_SIMD4_REAL_EQ(rSimd4_m0m1m2, -(rSimd4_c0c1c2)); // fneg(x)=-x
132 GMX_EXPECT_SIMD4_REAL_EQ(rSimd4_c0c1c2, -(rSimd4_m0m1m2)); // fneg(-x)=x
135 #if GMX_SIMD_HAVE_LOGICAL
136 TEST_F(Simd4FloatingpointTest, and)
138 GMX_EXPECT_SIMD4_REAL_EQ(rSimd4_logicalResultAnd,
139 (rSimd4_logicalA & rSimd4_logicalB));
142 TEST_F(Simd4FloatingpointTest, or)
144 GMX_EXPECT_SIMD4_REAL_EQ(rSimd4_logicalResultOr,
145 (rSimd4_logicalA | rSimd4_logicalB));
148 TEST_F(Simd4FloatingpointTest, xor)
150 /* Test xor by taking xor with a number and its negative. This should result
151 * in only the sign bit being set. We then use this bit change the sign of
152 * different numbers.
154 Simd4Real signbit = Simd4Real(c1) ^ Simd4Real(-c1);
155 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(-c2, c3, -c4), signbit ^ setSimd4RealFrom3R(c2, -c3, c4));
158 TEST_F(Simd4FloatingpointTest, andNot)
160 /* Use xor (which we already tested, so fix that first if both tests fail)
161 * to extract the sign bit, and then use andnot to take absolute values.
163 Simd4Real signbit = Simd4Real(c1) ^ Simd4Real(-c1);
164 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(c2, c3, c4), andNot(signbit, setSimd4RealFrom3R(-c2, c3, -c4)));
167 #endif
169 TEST_F(Simd4FloatingpointTest, max)
171 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R( c3, c1, c4), max(rSimd4_c0c1c2, rSimd4_c3c0c4));
172 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R( c3, c1, c4), max(rSimd4_c3c0c4, rSimd4_c0c1c2));
173 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(-c0, -c0, -c2), max(rSimd4_m0m1m2, rSimd4_m3m0m4));
174 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(-c0, -c0, -c2), max(rSimd4_m3m0m4, rSimd4_m0m1m2));
177 TEST_F(Simd4FloatingpointTest, min)
179 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R( c0, c0, c2), min(rSimd4_c0c1c2, rSimd4_c3c0c4));
180 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R( c0, c0, c2), min(rSimd4_c3c0c4, rSimd4_c0c1c2));
181 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(-c3, -c1, -c4), min(rSimd4_m0m1m2, rSimd4_m3m0m4));
182 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(-c3, -c1, -c4), min(rSimd4_m3m0m4, rSimd4_m0m1m2));
185 TEST_F(Simd4FloatingpointTest, round)
187 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom1R(2), round(Simd4Real(2.25)));
188 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom1R(4), round(Simd4Real(3.75)));
189 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom1R(-2), round(Simd4Real(-2.25)));
190 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom1R(-4), round(Simd4Real(-3.75)));
193 TEST_F(Simd4FloatingpointTest, trunc)
195 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom1R(2), trunc(rSimd4_2p25));
196 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom1R(3), trunc(rSimd4_3p75));
197 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom1R(-2), trunc(rSimd4_m2p25));
198 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom1R(-3), trunc(rSimd4_m3p75));
201 /* We do extensive 1/sqrt(x) and 1/x accuracy testing in the tests for
202 * the SIMD math functions, so we just make sure the lookup instructions
203 * appear to work for a few values here.
205 TEST_F(Simd4FloatingpointTest, gmxSimd4RsqrtR)
207 Simd4Real x = setSimd4RealFrom3R(4.0, M_PI, 1234567890.0);
208 Simd4Real ref = setSimd4RealFrom3R(0.5, 1.0/std::sqrt(M_PI), 1.0/std::sqrt(1234567890.0));
209 int shiftbits = std::numeric_limits<real>::digits-GMX_SIMD_RSQRT_BITS;
211 if (shiftbits < 0)
213 shiftbits = 0;
216 // The allowed Ulp deviation is 2 to the power of the number of mantissa
217 // digits, minus the number of bits provided by the table lookup
218 setUlpTol(1LL << shiftbits);
219 GMX_EXPECT_SIMD4_REAL_NEAR(ref, rsqrt(x));
222 TEST_F(Simd4FloatingpointTest, cmpEqAndSelectByMask)
224 Simd4Bool eq = rSimd4_c4c6c8 == rSimd4_c6c7c8;
225 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(0, 0, c2), selectByMask(rSimd4_c0c1c2, eq));
228 TEST_F(Simd4FloatingpointTest, selectByNotMask)
230 Simd4Bool eq = rSimd4_c4c6c8 == rSimd4_c6c7c8;
231 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(c0, c1, 0), selectByNotMask(rSimd4_c0c1c2, eq));
234 TEST_F(Simd4FloatingpointTest, cmpNe)
236 Simd4Bool eq = rSimd4_c4c6c8 != rSimd4_c6c7c8;
237 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(c0, c1, 0), selectByMask(rSimd4_c0c1c2, eq));
240 TEST_F(Simd4FloatingpointTest, cmpLe)
242 Simd4Bool le = rSimd4_c4c6c8 <= rSimd4_c6c7c8;
243 GMX_EXPECT_SIMD4_REAL_EQ(rSimd4_c0c1c2, selectByMask(rSimd4_c0c1c2, le));
246 TEST_F(Simd4FloatingpointTest, cmpLt)
248 Simd4Bool lt = rSimd4_c4c6c8 < rSimd4_c6c7c8;
249 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(c0, c1, 0), selectByMask(rSimd4_c0c1c2, lt));
252 TEST_F(Simd4FloatingpointTest, andB)
254 Simd4Bool eq = rSimd4_c4c6c8 == rSimd4_c6c7c8;
255 Simd4Bool le = rSimd4_c4c6c8 <= rSimd4_c6c7c8;
256 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(0, 0, c2), selectByMask(rSimd4_c0c1c2, (eq && le)));
259 TEST_F(Simd4FloatingpointTest, orB)
261 Simd4Bool eq = rSimd4_c4c6c8 == rSimd4_c6c7c8;
262 Simd4Bool lt = rSimd4_c4c6c8 < rSimd4_c6c7c8;
263 GMX_EXPECT_SIMD4_REAL_EQ(rSimd4_c0c1c2, selectByMask(rSimd4_c0c1c2, (eq || lt)));
266 TEST_F(Simd4FloatingpointTest, anyTrue)
268 Simd4Bool eq;
270 /* this test is a bit tricky since we don't know the simd width.
271 * We cannot check for truth values for "any" element beyond the first,
272 * since that part of the data will not be used if simd width is 1.
274 eq = (rSimd4_c4c6c8 == setSimd4RealFrom3R(c4, 0, 0));
275 EXPECT_TRUE(anyTrue(eq));
277 eq = (rSimd4_c0c1c2 == rSimd4_c3c4c5);
278 EXPECT_FALSE(anyTrue(eq));
281 TEST_F(Simd4FloatingpointTest, blend)
283 Simd4Bool lt = rSimd4_c4c6c8 < rSimd4_c6c7c8;
284 GMX_EXPECT_SIMD4_REAL_EQ(setSimd4RealFrom3R(c3, c4, c2), blend(rSimd4_c0c1c2, rSimd4_c3c4c5, lt));
287 TEST_F(Simd4FloatingpointTest, reduce)
289 // The horizontal sum of the SIMD variable depends on the width, so
290 // simply store it an extra time and calculate what the sum should be
291 std::vector<real> v = simd4Real2Vector(rSimd4_c3c4c5);
292 real sum = 0.0;
294 for (int i = 0; i < GMX_SIMD4_WIDTH; i++)
296 sum += v[i];
299 EXPECT_REAL_EQ_TOL(sum, reduce(rSimd4_c3c4c5), defaultRealTolerance() );
303 TEST_F(Simd4FloatingpointTest, dotProduct)
305 real res = c0*c3 + c1*c4 + c2*c5;
307 EXPECT_REAL_EQ_TOL(res, dotProduct(rSimd4_c0c1c2, rSimd4_c3c4c5), defaultRealTolerance());
310 TEST_F(Simd4FloatingpointTest, transpose)
312 Simd4Real v0, v1, v2, v3;
313 int i;
314 // aligned pointers
315 GMX_ALIGNED(real, GMX_SIMD4_WIDTH) p0[4*GMX_SIMD4_WIDTH];
316 real * p1 = p0 + GMX_SIMD4_WIDTH;
317 real * p2 = p0 + 2*GMX_SIMD4_WIDTH;
318 real * p3 = p0 + 3*GMX_SIMD4_WIDTH;
320 // Assign data with tens as row, single-digit as column
321 for (i = 0; i < 4; i++)
323 // Scale by 1+100*eps to use low bits tii
324 p0[i] = (0*10 + i*1) * (1.0 + 100*GMX_REAL_EPS);
325 p1[i] = (1*10 + i*1) * (1.0 + 100*GMX_REAL_EPS);
326 p2[i] = (2*10 + i*1) * (1.0 + 100*GMX_REAL_EPS);
327 p3[i] = (3*10 + i*1) * (1.0 + 100*GMX_REAL_EPS);
330 v0 = load4(p0);
331 v1 = load4(p1);
332 v2 = load4(p2);
333 v3 = load4(p3);
335 transpose(&v0, &v1, &v2, &v3);
337 store4(p0, v0);
338 store4(p1, v1);
339 store4(p2, v2);
340 store4(p3, v3);
342 for (i = 0; i < 4; i++)
344 EXPECT_REAL_EQ_TOL( (i*10+0) * (1.0 + 100*GMX_REAL_EPS), p0[i], defaultRealTolerance());
345 EXPECT_REAL_EQ_TOL( (i*10+1) * (1.0 + 100*GMX_REAL_EPS), p1[i], defaultRealTolerance());
346 EXPECT_REAL_EQ_TOL( (i*10+2) * (1.0 + 100*GMX_REAL_EPS), p2[i], defaultRealTolerance());
347 EXPECT_REAL_EQ_TOL( (i*10+3) * (1.0 + 100*GMX_REAL_EPS), p3[i], defaultRealTolerance());
351 #endif // GMX_SIMD4_HAVE_REAL
353 /*! \} */
354 /*! \endcond */
356 } // namespace
357 } // namespace
358 } // namespace
360 #endif // GMX_SIMD