78584b72c9b1f7cac495079122adb74cc7b0225e
[alexxy/gromacs.git] / src / gromacs / simd / tests / simd_floatingpoint.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #include "gmxpre.h"
36
37 #include "config.h"
38
39 #include <math.h>
40 #include "gromacs/math/utilities.h"
41
42 #include "simd.h"
43
44 namespace gmx
45 {
46 namespace test
47 {
48 namespace
49 {
50
51 /*! \cond internal */
52 /*! \addtogroup module_simd */
53 /*! \{ */
54
55 #ifdef GMX_SIMD_HAVE_REAL
56
57 /*! \brief Test fixture for floating-point tests (identical to the generic \ref SimdTest) */
58 typedef SimdTest SimdFloatingpointTest;
59
60 TEST_F(SimdFloatingpointTest, gmxSimdSetZeroR)
61 {
62     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(0.0), gmx_simd_setzero_r());
63 }
64
65 TEST_F(SimdFloatingpointTest, gmxSimdSet1R)
66 {
67     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(1.0), gmx_simd_set1_r(1.0));
68 }
69
70 TEST_F(SimdFloatingpointTest, gmxSimdLoad1R)
71 {
72     real r = 2.0;
73     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(r), gmx_simd_load1_r(&r));
74 }
75
76 TEST_F(SimdFloatingpointTest, gmxSimdAddR)
77 {
78     GMX_EXPECT_SIMD_REAL_EQ(rSimd_5_7_9,
79                             gmx_simd_add_r(rSimd_1_2_3, rSimd_4_5_6)); // 1+4=5, 2+5=7, 3+6=9
80 }
81
82 TEST_F(SimdFloatingpointTest, gmxSimdSubR)
83 {
84     GMX_EXPECT_SIMD_REAL_EQ(rSimd_4_5_6,
85                             gmx_simd_sub_r(rSimd_5_7_9, rSimd_1_2_3)); // 5-1=4, 7-2=5, 9-3=6
86 }
87
88 TEST_F(SimdFloatingpointTest, gmxSimdMulR)
89 {
90     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(4, 10, 18),
91                             gmx_simd_mul_r(rSimd_1_2_3, rSimd_4_5_6));
92 }
93
94 TEST_F(SimdFloatingpointTest, gmxSimdFmaddR)
95 {
96     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(11, 18, 27),
97                             gmx_simd_fmadd_r(rSimd_1_2_3, rSimd_4_5_6, rSimd_7_8_9)); // 1*4+7, etc.
98 }
99
100 TEST_F(SimdFloatingpointTest, gmxSimdFmsubR)
101 {
102     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-3, 2, 9),
103                             gmx_simd_fmsub_r(rSimd_1_2_3, rSimd_4_5_6, rSimd_7_8_9)); // 1*4-7, etc.
104 }
105
106 TEST_F(SimdFloatingpointTest, gmxSimdFnmaddR)
107 {
108     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(3, -2, -9),
109                             gmx_simd_fnmadd_r(rSimd_1_2_3, rSimd_4_5_6, rSimd_7_8_9)); // -1*4+7, etc.
110 }
111
112 TEST_F(SimdFloatingpointTest, gmxSimdFnmsubR)
113 {
114     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-11, -18, -27),
115                             gmx_simd_fnmsub_r(rSimd_1_2_3, rSimd_4_5_6, rSimd_7_8_9)); // -1*4-7, etc.
116 }
117
118 TEST_F(SimdFloatingpointTest, gmxSimdFabsR)
119 {
120     GMX_EXPECT_SIMD_REAL_EQ(rSimd_1_2_3, gmx_simd_fabs_r(rSimd_1_2_3));    // fabs(x)=x
121     GMX_EXPECT_SIMD_REAL_EQ(rSimd_1_2_3, gmx_simd_fabs_r(rSimd_m1_m2_m3)); // fabs(-x)=x
122 }
123
124 TEST_F(SimdFloatingpointTest, gmxSimdFnegR)
125 {
126     GMX_EXPECT_SIMD_REAL_EQ(rSimd_m1_m2_m3, gmx_simd_fneg_r(rSimd_1_2_3));    // fneg(x)=-x
127     GMX_EXPECT_SIMD_REAL_EQ(rSimd_1_2_3,    gmx_simd_fneg_r(rSimd_m1_m2_m3)); // fneg(-x)=x
128 }
129
130 #ifdef GMX_SIMD_HAVE_LOGICAL
131 /* 1.3333282470703125 has mantissa 0101010101010101 (followed by zeros)
132  * 1.79998779296875   has mantissa 1100110011001100 (followed by zeros)
133  * 1.26666259765625   has mantissa 0100010001000100 (followed by zeros)
134  * 1.8666534423828125 has mantissa 1101110111011101 (followed by zeros)
135  *
136  * Since all of them have the same exponent (2^0), the exponent will
137  * not change with AND or OR operations.
138  */
139 TEST_F(SimdFloatingpointTest, gmxSimdAndR)
140 {
141     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(1.26666259765625),
142                             gmx_simd_and_r(gmx_simd_set1_r(1.3333282470703125),
143                                            gmx_simd_set1_r(1.79998779296875)));
144 }
145
146 TEST_F(SimdFloatingpointTest, gmxSimdOrR)
147 {
148     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(1.8666534423828125),
149                             gmx_simd_or_r(gmx_simd_set1_r(1.3333282470703125),
150                                           gmx_simd_set1_r(1.79998779296875)));
151 }
152
153 TEST_F(SimdFloatingpointTest, gmxSimdXorR)
154 {
155     /* Test xor by taking xor with a number and its negative. This should result
156      * in only the sign bit being set. We then use this bit change the sign of
157      * different numbers.
158      */
159     gmx_simd_real_t signbit = gmx_simd_xor_r(gmx_simd_set1_r(1.5), gmx_simd_set1_r(-1.5));
160     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-1, 2, -3), gmx_simd_xor_r(signbit, setSimdRealFrom3R(1, -2, 3)));
161 }
162
163 TEST_F(SimdFloatingpointTest, gmxSimdAndnotR)
164 {
165     /* Use xor (which we already tested, so fix that first if both tests fail)
166      * to extract the sign bit, and then use andnot to take absolute values.
167      */
168     gmx_simd_real_t signbit = gmx_simd_xor_r(gmx_simd_set1_r(1.5), gmx_simd_set1_r(-1.5));
169     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1, 2, 3), gmx_simd_andnot_r(signbit, setSimdRealFrom3R(-1, 2, -3)));
170 }
171
172 #endif
173
174 TEST_F(SimdFloatingpointTest, gmxSimdMaxR)
175 {
176     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(3, 2, 4), gmx_simd_max_r(rSimd_1_2_3, rSimd_3_1_4));
177     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(3, 2, 4), gmx_simd_max_r(rSimd_3_1_4, rSimd_1_2_3));
178     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-1, -1, -3), gmx_simd_max_r(rSimd_m1_m2_m3, rSimd_m3_m1_m4));
179     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-1, -1, -3), gmx_simd_max_r(rSimd_m3_m1_m4, rSimd_m1_m2_m3));
180 }
181
182 TEST_F(SimdFloatingpointTest, gmxSimdMinR)
183 {
184     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1, 1, 3), gmx_simd_min_r(rSimd_1_2_3, rSimd_3_1_4));
185     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1, 1, 3), gmx_simd_min_r(rSimd_3_1_4, rSimd_1_2_3));
186     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-3, -2, -4), gmx_simd_min_r(rSimd_m1_m2_m3, rSimd_m3_m1_m4));
187     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-3, -2, -4), gmx_simd_min_r(rSimd_m3_m1_m4, rSimd_m1_m2_m3));
188 }
189
190 TEST_F(SimdFloatingpointTest, gmxSimdRoundR)
191 {
192     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(2), gmx_simd_round_r(gmx_simd_set1_r(2.25)));
193     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(4), gmx_simd_round_r(gmx_simd_set1_r(3.75)));
194     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-2), gmx_simd_round_r(gmx_simd_set1_r(-2.25)));
195     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-4), gmx_simd_round_r(gmx_simd_set1_r(-3.75)));
196 }
197
198 TEST_F(SimdFloatingpointTest, gmxSimdTruncR)
199 {
200     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(2), gmx_simd_trunc_r(rSimd_2p25));
201     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(3), gmx_simd_trunc_r(rSimd_3p75));
202     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-2), gmx_simd_trunc_r(rSimd_m2p25));
203     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-3), gmx_simd_trunc_r(rSimd_m3p75));
204 }
205
206 TEST_F(SimdFloatingpointTest, gmxSimdFractionR)
207 {
208     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(0.25), gmx_simd_fraction_r(rSimd_2p25));   // fract(2.25)=0.25
209     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(0.75), gmx_simd_fraction_r(rSimd_3p75));   // fract(3.75)=0.75
210     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-0.25), gmx_simd_fraction_r(rSimd_m2p25)); // fract(-2.25)=-0.25
211     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-0.75), gmx_simd_fraction_r(rSimd_m3p75)); // fract(-3.75)=-0.75
212 }
213
214 // We explicitly test the exponent/mantissa routines with double precision data,
215 // since these usually rely on direct manipulation and shift of the SIMD registers,
216 // where it is easy to make mistakes with single vs double precision.
217
218 TEST_F(SimdFloatingpointTest, gmxSimdGetExponentR)
219 {
220     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(60.0, -41.0, 54.0), gmx_simd_get_exponent_r(rSimd_Exp));
221 #if (defined GMX_SIMD_HAVE_DOUBLE) && (defined GMX_DOUBLE)
222     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(587.0, -462.0, 672.0), gmx_simd_get_exponent_r(rSimd_ExpDouble));
223 #endif
224 }
225
226 TEST_F(SimdFloatingpointTest, gmxSimdGetMantissaR)
227 {
228     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1.219097320577810839026256,
229                                               1.166738027848349235071623,
230                                               1.168904015004464724825084), gmx_simd_get_mantissa_r(rSimd_Exp));
231 #if (defined GMX_SIMD_HAVE_DOUBLE) && (defined GMX_DOUBLE)
232     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1.241261238952345623563251,
233                                               1.047294723759123852359232,
234                                               1.856066204750275957395734), gmx_simd_get_mantissa_r(rSimd_ExpDouble));
235 #endif
236 }
237
238 TEST_F(SimdFloatingpointTest, gmxSimdSetExponentR)
239 {
240     gmx_simd_real_t x0 = setSimdRealFrom3R(0.5, 11.5, 99.5);
241     gmx_simd_real_t x1 = setSimdRealFrom3R(-0.5, -11.5, -99.5);
242
243     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(pow(2.0, 60.0), pow(2.0, -41.0), pow(2.0, 54.0)),
244                             gmx_simd_set_exponent_r(setSimdRealFrom3R(60.0, -41.0, 54.0)));
245 #if (defined GMX_SIMD_HAVE_DOUBLE) && (defined GMX_DOUBLE)
246     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(pow(2.0, 587.0), pow(2.0, -462.0), pow(2.0, 672.0)),
247                             gmx_simd_set_exponent_r(setSimdRealFrom3R(587.0, -462.0, 672.0)));
248 #endif
249     /* Rounding mode in gmx_simd_set_exponent_r() must be consistent with gmx_simd_round_r() */
250     GMX_EXPECT_SIMD_REAL_EQ(gmx_simd_set_exponent_r(gmx_simd_round_r(x0)), gmx_simd_set_exponent_r(x0));
251     GMX_EXPECT_SIMD_REAL_EQ(gmx_simd_set_exponent_r(gmx_simd_round_r(x1)), gmx_simd_set_exponent_r(x1));
252 }
253
254 /*
255  * We do extensive 1/sqrt(x) and 1/x accuracy testing in the math module, so
256  * we just make sure the lookup instructions appear to work here
257  */
258
259 TEST_F(SimdFloatingpointTest, gmxSimdRsqrtR)
260 {
261     gmx_simd_real_t x      = setSimdRealFrom3R(4.0, M_PI, 1234567890.0);
262     gmx_simd_real_t ref    = setSimdRealFrom3R(0.5, 1.0/sqrt(M_PI), 1.0/sqrt(1234567890.0));
263
264     /* Set the allowed ulp error as 2 to the power of the number of bits in
265      * the mantissa that do not have to be correct after the table lookup.
266      */
267     setUlpTol(1LL << (std::numeric_limits<real>::digits-GMX_SIMD_RSQRT_BITS));
268
269     GMX_EXPECT_SIMD_REAL_NEAR(ref, gmx_simd_rsqrt_r(x));
270 }
271
272 TEST_F(SimdFloatingpointTest, gmxSimdRcpR)
273 {
274     gmx_simd_real_t x      = setSimdRealFrom3R(4.0, M_PI, 1234567890.0);
275     gmx_simd_real_t ref    = setSimdRealFrom3R(0.25, 1.0/M_PI, 1.0/1234567890.0);
276
277     /* Set the allowed ulp error as 2 to the power of the number of bits in
278      * the mantissa that do not have to be correct after the table lookup.
279      */
280     setUlpTol(1LL << (std::numeric_limits<real>::digits-GMX_SIMD_RCP_BITS));
281
282     GMX_EXPECT_SIMD_REAL_NEAR(ref, gmx_simd_rcp_r(x));
283 }
284
285 TEST_F(SimdFloatingpointTest, gmxSimdBoolCmpEqAndBlendZeroR)
286 {
287     gmx_simd_bool_t eq   = gmx_simd_cmpeq_r(rSimd_5_7_9, rSimd_7_8_9);
288     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0, 0, 3), gmx_simd_blendzero_r(rSimd_1_2_3, eq));
289 }
290
291 TEST_F(SimdFloatingpointTest, gmxSimdBlendNotZeroR)
292 {
293     gmx_simd_bool_t eq   = gmx_simd_cmpeq_r(rSimd_5_7_9, rSimd_7_8_9);
294     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1, 2, 0), gmx_simd_blendnotzero_r(rSimd_1_2_3, eq));
295 }
296
297 TEST_F(SimdFloatingpointTest, gmxSimdBoolCmpLER)
298 {
299     gmx_simd_bool_t le   = gmx_simd_cmple_r(rSimd_5_7_9, rSimd_7_8_9);
300     GMX_EXPECT_SIMD_REAL_EQ(rSimd_1_2_3, gmx_simd_blendzero_r(rSimd_1_2_3, le));
301 }
302
303 TEST_F(SimdFloatingpointTest, gmxSimdBoolCmpLTR)
304 {
305     gmx_simd_bool_t lt   = gmx_simd_cmplt_r(rSimd_5_7_9, rSimd_7_8_9);
306     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1, 2, 0), gmx_simd_blendzero_r(rSimd_1_2_3, lt));
307 }
308
309 TEST_F(SimdFloatingpointTest, gmxSimdBoolAndB)
310 {
311     gmx_simd_bool_t eq   = gmx_simd_cmpeq_r(rSimd_5_7_9, rSimd_7_8_9);
312     gmx_simd_bool_t le   = gmx_simd_cmple_r(rSimd_5_7_9, rSimd_7_8_9);
313     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0, 0, 3), gmx_simd_blendzero_r(rSimd_1_2_3, gmx_simd_and_b(eq, le)));
314 }
315
316 TEST_F(SimdFloatingpointTest, gmxSimdBoolOrB)
317 {
318     gmx_simd_bool_t eq   = gmx_simd_cmpeq_r(rSimd_5_7_9, rSimd_7_8_9);
319     gmx_simd_bool_t lt   = gmx_simd_cmplt_r(rSimd_5_7_9, rSimd_7_8_9);
320     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1, 2, 3), gmx_simd_blendzero_r(rSimd_1_2_3, gmx_simd_or_b(eq, lt)));
321 }
322
323 TEST_F(SimdFloatingpointTest, gmxSimdAnytrueB)
324 {
325     gmx_simd_bool_t eq;
326
327     /* this test is a bit tricky since we don't know the simd width.
328      * We cannot check for truth values for "any" element beyond the first,
329      * since that part of the data will not be used if simd width is 1.
330      */
331     eq = gmx_simd_cmpeq_r(rSimd_5_7_9, setSimdRealFrom3R(5, 0, 0));
332     EXPECT_NE(0, gmx_simd_anytrue_b(eq));
333
334     eq = gmx_simd_cmpeq_r(rSimd_1_2_3, rSimd_4_5_6);
335     EXPECT_EQ(0, gmx_simd_anytrue_b(eq));
336 }
337
338 TEST_F(SimdFloatingpointTest, gmxSimdBlendvR)
339 {
340     gmx_simd_bool_t lt   = gmx_simd_cmplt_r(rSimd_5_7_9, rSimd_7_8_9);
341     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(4, 5, 3), gmx_simd_blendv_r(rSimd_1_2_3, rSimd_4_5_6, lt));
342 }
343
344 TEST_F(SimdFloatingpointTest, gmxSimdReduceR)
345 {
346     // The horizontal sum of the SIMD variable depends on the width, so
347     // simply store it an extra time and calculate what the sum should be
348     std::vector<real> v   = simdReal2Vector(rSimd_4_5_6);
349     real              sum = 0.0;
350
351     for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
352     {
353         sum += v[i];
354     }
355
356     EXPECT_EQ(sum, gmx_simd_reduce_r(rSimd_4_5_6));
357 }
358
359 #endif      // GMX_SIMD_HAVE_REAL
360
361 /*! \} */
362 /*! \endcond */
363
364 }      // namespace
365 }      // namespace
366 }      // namespace