Merge release-4-6 into release-5-0
[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 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
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 TEST_F(SimdFloatingpointTest, gmxSimdAndR)
132 {
133     GMX_EXPECT_SIMD_REAL_EQ(rSimd_Bits3, gmx_simd_and_r(rSimd_Bits1, rSimd_Bits2)); // Bits1 & Bits2 = Bits3
134 }
135
136 TEST_F(SimdFloatingpointTest, gmxSimdAndnotR)
137 {
138     GMX_EXPECT_SIMD_REAL_EQ(rSimd_Bits4, gmx_simd_andnot_r(rSimd_Bits1, rSimd_Bits2)); // (~Bits1) & Bits2 = Bits3
139 }
140
141 TEST_F(SimdFloatingpointTest, gmxSimdOrR)
142 {
143     GMX_EXPECT_SIMD_REAL_EQ(rSimd_Bits5, gmx_simd_or_r(rSimd_Bits1, rSimd_Bits2)); // Bits1 | Bits2 = Bits3
144 }
145
146 TEST_F(SimdFloatingpointTest, gmxSimdXorR)
147 {
148     GMX_EXPECT_SIMD_REAL_EQ(rSimd_Bits6, gmx_simd_xor_r(rSimd_Bits1, rSimd_Bits2)); // Bits1 ^ Bits2 = Bits3
149 }
150 #endif
151
152 TEST_F(SimdFloatingpointTest, gmxSimdMaxR)
153 {
154     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(3, 2, 4), gmx_simd_max_r(rSimd_1_2_3, rSimd_3_1_4));
155     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(3, 2, 4), gmx_simd_max_r(rSimd_3_1_4, rSimd_1_2_3));
156     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-1, -1, -3), gmx_simd_max_r(rSimd_m1_m2_m3, rSimd_m3_m1_m4));
157     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-1, -1, -3), gmx_simd_max_r(rSimd_m3_m1_m4, rSimd_m1_m2_m3));
158 }
159
160 TEST_F(SimdFloatingpointTest, gmxSimdMinR)
161 {
162     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1, 1, 3), gmx_simd_min_r(rSimd_1_2_3, rSimd_3_1_4));
163     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1, 1, 3), gmx_simd_min_r(rSimd_3_1_4, rSimd_1_2_3));
164     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-3, -2, -4), gmx_simd_min_r(rSimd_m1_m2_m3, rSimd_m3_m1_m4));
165     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-3, -2, -4), gmx_simd_min_r(rSimd_m3_m1_m4, rSimd_m1_m2_m3));
166 }
167
168 TEST_F(SimdFloatingpointTest, gmxSimdRoundR)
169 {
170     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(2), gmx_simd_round_r(gmx_simd_set1_r(2.25)));
171     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(4), gmx_simd_round_r(gmx_simd_set1_r(3.75)));
172     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-2), gmx_simd_round_r(gmx_simd_set1_r(-2.25)));
173     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-4), gmx_simd_round_r(gmx_simd_set1_r(-3.75)));
174 }
175
176 TEST_F(SimdFloatingpointTest, gmxSimdTruncR)
177 {
178     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(2), gmx_simd_trunc_r(rSimd_2p25));
179     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(3), gmx_simd_trunc_r(rSimd_3p75));
180     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-2), gmx_simd_trunc_r(rSimd_m2p25));
181     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-3), gmx_simd_trunc_r(rSimd_m3p75));
182 }
183
184 TEST_F(SimdFloatingpointTest, gmxSimdFractionR)
185 {
186     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(0.25), gmx_simd_fraction_r(rSimd_2p25));   // fract(2.25)=0.25
187     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(0.75), gmx_simd_fraction_r(rSimd_3p75));   // fract(3.75)=0.75
188     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-0.25), gmx_simd_fraction_r(rSimd_m2p25)); // fract(-2.25)=-0.25
189     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-0.75), gmx_simd_fraction_r(rSimd_m3p75)); // fract(-3.75)=-0.75
190 }
191
192 // We explicitly test the exponent/mantissa routines with double precision data,
193 // since these usually rely on direct manipulation and shift of the SIMD registers,
194 // where it is easy to make mistakes with single vs double precision.
195
196 TEST_F(SimdFloatingpointTest, gmxSimdGetExponentR)
197 {
198     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(60.0, -41.0, 54.0), gmx_simd_get_exponent_r(rSimd_Exp));
199 #if (defined GMX_SIMD_HAVE_DOUBLE) && (defined GMX_DOUBLE)
200     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(587.0, -462.0, 672.0), gmx_simd_get_exponent_r(rSimd_ExpDouble));
201 #endif
202 }
203
204 TEST_F(SimdFloatingpointTest, gmxSimdGetMantissaR)
205 {
206     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1.219097320577810839026256,
207                                               1.166738027848349235071623,
208                                               1.168904015004464724825084), gmx_simd_get_mantissa_r(rSimd_Exp));
209 #if (defined GMX_SIMD_HAVE_DOUBLE) && (defined GMX_DOUBLE)
210     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1.241261238952345623563251,
211                                               1.047294723759123852359232,
212                                               1.856066204750275957395734), gmx_simd_get_mantissa_r(rSimd_ExpDouble));
213 #endif
214 }
215
216 TEST_F(SimdFloatingpointTest, gmxSimdSetExponentR)
217 {
218     gmx_simd_real_t x0 = setSimdRealFrom3R(0.5, 11.5, 99.5);
219     gmx_simd_real_t x1 = setSimdRealFrom3R(-0.5, -11.5, -99.5);
220
221     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(pow(2.0, 60.0), pow(2.0, -41.0), pow(2.0, 54.0)),
222                             gmx_simd_set_exponent_r(setSimdRealFrom3R(60.0, -41.0, 54.0)));
223 #if (defined GMX_SIMD_HAVE_DOUBLE) && (defined GMX_DOUBLE)
224     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(pow(2.0, 587.0), pow(2.0, -462.0), pow(2.0, 672.0)),
225                             gmx_simd_set_exponent_r(setSimdRealFrom3R(587.0, -462.0, 672.0)));
226 #endif
227     /* Rounding mode in gmx_simd_set_exponent_r() must be consistent with gmx_simd_round_r() */
228     GMX_EXPECT_SIMD_REAL_EQ(gmx_simd_set_exponent_r(gmx_simd_round_r(x0)), gmx_simd_set_exponent_r(x0));
229     GMX_EXPECT_SIMD_REAL_EQ(gmx_simd_set_exponent_r(gmx_simd_round_r(x1)), gmx_simd_set_exponent_r(x1));
230 }
231
232 /*
233  * We do extensive 1/sqrt(x) and 1/x accuracy testing in the math module, so
234  * we just make sure the lookup instructions appear to work here
235  */
236
237 TEST_F(SimdFloatingpointTest, gmxSimdRsqrtR)
238 {
239     gmx_simd_real_t x      = setSimdRealFrom3R(4.0, M_PI, 1234567890.0);
240     gmx_simd_real_t ref    = setSimdRealFrom3R(0.5, 1.0/sqrt(M_PI), 1.0/sqrt(1234567890.0));
241
242     /* Set the allowed ulp error as 2 to the power of the number of bits in
243      * the mantissa that do not have to be correct after the table lookup.
244      */
245     setUlpTol(1LL << (std::numeric_limits<real>::digits-GMX_SIMD_RSQRT_BITS));
246
247     GMX_EXPECT_SIMD_REAL_NEAR(ref, gmx_simd_rsqrt_r(x));
248 }
249
250 TEST_F(SimdFloatingpointTest, gmxSimdRcpR)
251 {
252     gmx_simd_real_t x      = setSimdRealFrom3R(4.0, M_PI, 1234567890.0);
253     gmx_simd_real_t ref    = setSimdRealFrom3R(0.25, 1.0/M_PI, 1.0/1234567890.0);
254
255     /* Set the allowed ulp error as 2 to the power of the number of bits in
256      * the mantissa that do not have to be correct after the table lookup.
257      */
258     setUlpTol(1LL << (std::numeric_limits<real>::digits-GMX_SIMD_RCP_BITS));
259
260     GMX_EXPECT_SIMD_REAL_NEAR(ref, gmx_simd_rcp_r(x));
261 }
262
263 TEST_F(SimdFloatingpointTest, gmxSimdBoolCmpEqAndBlendZeroR)
264 {
265     gmx_simd_bool_t eq   = gmx_simd_cmpeq_r(rSimd_5_7_9, rSimd_7_8_9);
266     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0, 0, 3), gmx_simd_blendzero_r(rSimd_1_2_3, eq));
267 }
268
269 TEST_F(SimdFloatingpointTest, gmxSimdBlendNotZeroR)
270 {
271     gmx_simd_bool_t eq   = gmx_simd_cmpeq_r(rSimd_5_7_9, rSimd_7_8_9);
272     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1, 2, 0), gmx_simd_blendnotzero_r(rSimd_1_2_3, eq));
273 }
274
275 TEST_F(SimdFloatingpointTest, gmxSimdBoolCmpLER)
276 {
277     gmx_simd_bool_t le   = gmx_simd_cmple_r(rSimd_5_7_9, rSimd_7_8_9);
278     GMX_EXPECT_SIMD_REAL_EQ(rSimd_1_2_3, gmx_simd_blendzero_r(rSimd_1_2_3, le));
279 }
280
281 TEST_F(SimdFloatingpointTest, gmxSimdBoolCmpLTR)
282 {
283     gmx_simd_bool_t lt   = gmx_simd_cmplt_r(rSimd_5_7_9, rSimd_7_8_9);
284     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1, 2, 0), gmx_simd_blendzero_r(rSimd_1_2_3, lt));
285 }
286
287 TEST_F(SimdFloatingpointTest, gmxSimdBoolAndB)
288 {
289     gmx_simd_bool_t eq   = gmx_simd_cmpeq_r(rSimd_5_7_9, rSimd_7_8_9);
290     gmx_simd_bool_t le   = gmx_simd_cmple_r(rSimd_5_7_9, rSimd_7_8_9);
291     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0, 0, 3), gmx_simd_blendzero_r(rSimd_1_2_3, gmx_simd_and_b(eq, le)));
292 }
293
294 TEST_F(SimdFloatingpointTest, gmxSimdBoolOrB)
295 {
296     gmx_simd_bool_t eq   = gmx_simd_cmpeq_r(rSimd_5_7_9, rSimd_7_8_9);
297     gmx_simd_bool_t lt   = gmx_simd_cmplt_r(rSimd_5_7_9, rSimd_7_8_9);
298     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(1, 2, 3), gmx_simd_blendzero_r(rSimd_1_2_3, gmx_simd_or_b(eq, lt)));
299 }
300
301 TEST_F(SimdFloatingpointTest, gmxSimdAnytrueB)
302 {
303     gmx_simd_bool_t eq;
304
305     /* this test is a bit tricky since we don't know the simd width.
306      * We cannot check for truth values for "any" element beyond the first,
307      * since that part of the data will not be used if simd width is 1.
308      */
309     eq = gmx_simd_cmpeq_r(rSimd_5_7_9, setSimdRealFrom3R(5, 0, 0));
310     EXPECT_NE(0, gmx_simd_anytrue_b(eq));
311
312     eq = gmx_simd_cmpeq_r(rSimd_1_2_3, rSimd_4_5_6);
313     EXPECT_EQ(0, gmx_simd_anytrue_b(eq));
314 }
315
316 TEST_F(SimdFloatingpointTest, gmxSimdBlendvR)
317 {
318     gmx_simd_bool_t lt   = gmx_simd_cmplt_r(rSimd_5_7_9, rSimd_7_8_9);
319     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(4, 5, 3), gmx_simd_blendv_r(rSimd_1_2_3, rSimd_4_5_6, lt));
320 }
321
322 TEST_F(SimdFloatingpointTest, gmxSimdReduceR)
323 {
324     // The horizontal sum of the SIMD variable depends on the width, so
325     // simply store it an extra time and calculate what the sum should be
326     std::vector<real> v   = simdReal2Vector(rSimd_4_5_6);
327     real              sum = 0.0;
328
329     for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
330     {
331         sum += v[i];
332     }
333
334     EXPECT_EQ(sum, gmx_simd_reduce_r(rSimd_4_5_6));
335 }
336
337 #endif      // GMX_SIMD_HAVE_REAL
338
339 /*! \} */
340 /*! \endcond */
341
342 }      // namespace
343 }      // namespace
344 }      // namespace