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