28d07cb7713e90a33504890ae4f0bbac0f158d30
[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,2015,2016,2017,2018 by the GROMACS development team.
5  * Copyright (c) 2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 #include "gmxpre.h"
37
38 #include <cmath>
39
40 #include <array>
41
42 #include "gromacs/math/utilities.h"
43 #include "gromacs/simd/simd.h"
44 #include "gromacs/utility/basedefinitions.h"
45
46 #include "testutils/testasserts.h"
47
48 #include "data.h"
49 #include "simd.h"
50
51 #if GMX_SIMD
52
53 namespace gmx
54 {
55 namespace test
56 {
57
58 namespace
59 {
60
61 /*! \cond internal */
62 /*! \addtogroup module_simd */
63 /*! \{ */
64
65 #    if GMX_SIMD_HAVE_REAL
66
67 /*! \brief Test fixture for floating-point tests (identical to the generic \ref SimdTest) */
68 typedef SimdTest SimdFloatingpointTest;
69
70 TEST_F(SimdFloatingpointTest, setZero)
71 {
72     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(0.0), setZero());
73 }
74
75 TEST_F(SimdFloatingpointTest, set)
76 {
77     const real* p = &c0;
78     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(c1), SimdReal(c1));
79     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(c0), SimdReal(*p));
80 }
81
82 TEST_F(SimdFloatingpointTest, add)
83 {
84     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0 + c3, c1 + c4, c2 + c5), rSimd_c0c1c2 + rSimd_c3c4c5);
85 }
86
87 TEST_F(SimdFloatingpointTest, maskAdd)
88 {
89     SimdBool m = setSimdRealFrom3R(c6, 0, c7) != setZero();
90     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0 + c3, c1 + 0.0, c2 + c5),
91                               maskAdd(rSimd_c0c1c2, rSimd_c3c4c5, m));
92 }
93
94 TEST_F(SimdFloatingpointTest, sub)
95 {
96     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0 - c3, c1 - c4, c2 - c5), rSimd_c0c1c2 - rSimd_c3c4c5);
97 }
98
99 TEST_F(SimdFloatingpointTest, mul)
100 {
101     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0 * c3, c1 * c4, c2 * c5), rSimd_c0c1c2 * rSimd_c3c4c5);
102 }
103
104 TEST_F(SimdFloatingpointTest, maskzMul)
105 {
106     SimdBool m = setSimdRealFrom3R(c1, 0, c1) != setZero();
107     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0 * c3, 0.0, c2 * c5),
108                               maskzMul(rSimd_c0c1c2, rSimd_c3c4c5, m));
109 }
110
111 TEST_F(SimdFloatingpointTest, fma)
112 {
113     // The last bit of FMA operations depends on hardware, so we don't require exact match
114     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0 * c3 + c6, c1 * c4 + c7, c2 * c5 + c8),
115                               fma(rSimd_c0c1c2, rSimd_c3c4c5, rSimd_c6c7c8));
116 }
117
118
119 TEST_F(SimdFloatingpointTest, maskzFma)
120 {
121     SimdBool m = setSimdRealFrom3R(c2, 0, c3) != setZero();
122     // The last bit of FMA operations depends on hardware, so we don't require exact match
123     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0 * c3 + c6, 0.0, c2 * c5 + c8),
124                               maskzFma(rSimd_c0c1c2, rSimd_c3c4c5, rSimd_c6c7c8, m));
125 }
126
127 TEST_F(SimdFloatingpointTest, fms)
128 {
129     // The last bit of FMA operations depends on hardware, so we don't require exact match
130     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0 * c3 - c6, c1 * c4 - c7, c2 * c5 - c8),
131                               fms(rSimd_c0c1c2, rSimd_c3c4c5, rSimd_c6c7c8));
132 }
133
134 TEST_F(SimdFloatingpointTest, fnma)
135 {
136     // The last bit of FMA operations depends on hardware, so we don't require exact match
137     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c6 - c0 * c3, c7 - c1 * c4, c8 - c2 * c5),
138                               fnma(rSimd_c0c1c2, rSimd_c3c4c5, rSimd_c6c7c8));
139 }
140
141 TEST_F(SimdFloatingpointTest, fnms)
142 {
143     // The last bit of FMA operations depends on hardware, so we don't require exact match
144     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(-c0 * c3 - c6, -c1 * c4 - c7, -c2 * c5 - c8),
145                               fnms(rSimd_c0c1c2, rSimd_c3c4c5, rSimd_c6c7c8));
146 }
147
148 TEST_F(SimdFloatingpointTest, abs)
149 {
150     GMX_EXPECT_SIMD_REAL_EQ(rSimd_c0c1c2, abs(rSimd_c0c1c2)); // fabs(x)=x
151     GMX_EXPECT_SIMD_REAL_EQ(rSimd_c0c1c2, abs(rSimd_m0m1m2)); // fabs(-x)=x
152 }
153
154 TEST_F(SimdFloatingpointTest, neg)
155 {
156     GMX_EXPECT_SIMD_REAL_EQ(rSimd_m0m1m2, -(rSimd_c0c1c2)); // fneg(x)=-x
157     GMX_EXPECT_SIMD_REAL_EQ(rSimd_c0c1c2, -(rSimd_m0m1m2)); // fneg(-x)=x
158 }
159
160 #        if GMX_SIMD_HAVE_LOGICAL
161 TEST_F(SimdFloatingpointTest, and)
162 {
163     GMX_EXPECT_SIMD_REAL_EQ(rSimd_logicalResultAnd, (rSimd_logicalA & rSimd_logicalB));
164 }
165
166 TEST_F(SimdFloatingpointTest, or)
167 {
168     GMX_EXPECT_SIMD_REAL_EQ(rSimd_logicalResultOr, (rSimd_logicalA | rSimd_logicalB));
169 }
170
171 TEST_F(SimdFloatingpointTest, xor)
172 {
173     /* Test xor by taking xor with a number and its negative. This should result
174      * in only the sign bit being set. We then use this bit change the sign of
175      * different numbers.
176      */
177     SimdReal signbit = SimdReal(c1) ^ SimdReal(-c1);
178     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-c2, c3, -c4), (signbit ^ setSimdRealFrom3R(c2, -c3, c4)));
179 }
180
181 TEST_F(SimdFloatingpointTest, andNot)
182 {
183     /* Use xor (which we already tested, so fix that first if both tests fail)
184      * to extract the sign bit, and then use andnot to take absolute values.
185      */
186     SimdReal signbit = SimdReal(c1) ^ SimdReal(-c1);
187     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c2, c3, c4),
188                             andNot(signbit, setSimdRealFrom3R(-c2, c3, -c4)));
189 }
190
191 #        endif
192
193 TEST_F(SimdFloatingpointTest, max)
194 {
195     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c3, c1, c4), max(rSimd_c0c1c2, rSimd_c3c0c4));
196     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c3, c1, c4), max(rSimd_c3c0c4, rSimd_c0c1c2));
197     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-c0, -c0, -c2), max(rSimd_m0m1m2, rSimd_m3m0m4));
198     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-c0, -c0, -c2), max(rSimd_m3m0m4, rSimd_m0m1m2));
199 }
200
201 TEST_F(SimdFloatingpointTest, min)
202 {
203     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0, c0, c2), min(rSimd_c0c1c2, rSimd_c3c0c4));
204     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0, c0, c2), min(rSimd_c3c0c4, rSimd_c0c1c2));
205     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-c3, -c1, -c4), min(rSimd_m0m1m2, rSimd_m3m0m4));
206     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-c3, -c1, -c4), min(rSimd_m3m0m4, rSimd_m0m1m2));
207 }
208
209 TEST_F(SimdFloatingpointTest, round)
210 {
211     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(2), round(rSimd_2p25));
212     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(4), round(rSimd_3p75));
213     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-2), round(rSimd_m2p25));
214     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-4), round(rSimd_m3p75));
215 }
216
217 TEST_F(SimdFloatingpointTest, roundMode)
218 {
219     /* Rounding mode needs to be consistent between round and cvtR2I */
220     SimdReal x0 = setSimdRealFrom3R(0.5, 11.5, 99.5);
221     SimdReal x1 = setSimdRealFrom3R(-0.5, -11.5, -99.5);
222
223     GMX_EXPECT_SIMD_REAL_EQ(round(x0), cvtI2R(cvtR2I(x0)));
224     GMX_EXPECT_SIMD_REAL_EQ(round(x1), cvtI2R(cvtR2I(x1)));
225 }
226
227 TEST_F(SimdFloatingpointTest, trunc)
228 {
229     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(2), trunc(rSimd_2p25));
230     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(3), trunc(rSimd_3p75));
231     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-2), trunc(rSimd_m2p25));
232     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-3), trunc(rSimd_m3p75));
233 }
234
235 // We explicitly test the exponent/mantissa routines with double precision data,
236 // since these usually rely on direct manipulation and shift of the SIMD registers,
237 // where it is easy to make mistakes with single vs double precision.
238
239 TEST_F(SimdFloatingpointTest, frexp)
240 {
241     SimdReal  fraction;
242     SimdInt32 exponent;
243
244
245     fraction = frexp(rSimd_Exp, &exponent);
246     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0.609548660288905419513128, 0.5833690139241746175358116,
247                                               -0.584452007502232362412542),
248                             fraction);
249     GMX_EXPECT_SIMD_INT_EQ(setSimdIntFrom3I(61, -40, 55), exponent);
250
251     // Test the unsafe flavor too, in case they use different branches
252     fraction = frexp<MathOptimization::Unsafe>(rSimd_Exp, &exponent);
253     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0.609548660288905419513128, 0.5833690139241746175358116,
254                                               -0.584452007502232362412542),
255                             fraction);
256     GMX_EXPECT_SIMD_INT_EQ(setSimdIntFrom3I(61, -40, 55), exponent);
257
258     // Use ulp testing with 0 bit ulp tolerance for testing to separate 0.0 and -0.0
259     setUlpTol(0);
260     fraction = frexp(setSimdRealFrom1R(0.0), &exponent);
261     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom1R(0.0), fraction);
262     GMX_EXPECT_SIMD_INT_EQ(setSimdIntFrom1I(0), exponent);
263
264     // Second -0.0.
265     fraction = frexp(setSimdRealFrom1R(-0.0), &exponent);
266     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom1R(-0.0), fraction);
267     GMX_EXPECT_SIMD_INT_EQ(setSimdIntFrom1I(0), exponent);
268
269     // Reset to default ulp tolerance
270     setUlpTol(defaultRealUlpTol());
271
272 #        if GMX_SIMD_HAVE_DOUBLE && GMX_DOUBLE
273     // Test exponents larger than what fit in single precision, as well as mixtures of 0 and non-zero values, to
274     // make sure the shuffling operations in the double-precision implementations don't do anything bad.
275     fraction = frexp(rSimd_ExpDouble1, &exponent);
276     GMX_EXPECT_SIMD_REAL_EQ(
277             setSimdRealFrom3R(0.0, 0.5236473618795619566768096, -0.9280331023751380303821179), fraction);
278     GMX_EXPECT_SIMD_INT_EQ(setSimdIntFrom3I(0, -461, 673), exponent);
279
280     fraction = frexp(rSimd_ExpDouble2, &exponent);
281     GMX_EXPECT_SIMD_REAL_EQ(
282             setSimdRealFrom3R(0.6206306194761728178832527, 0.0, -0.9280331023751380303821179), fraction);
283     GMX_EXPECT_SIMD_INT_EQ(setSimdIntFrom3I(588, 0, 673), exponent);
284 #        endif
285 }
286
287 TEST_F(SimdFloatingpointTest, ldexp)
288 {
289     SimdReal one = setSimdRealFrom1R(1.0);
290
291     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(pow(2.0, 60.0), pow(2.0, -41.0), pow(2.0, 54.0)),
292                             ldexp<MathOptimization::Unsafe>(one, setSimdIntFrom3I(60, -41, 54)));
293 #        if GMX_SIMD_HAVE_DOUBLE && GMX_DOUBLE
294     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(pow(2.0, 587.0), pow(2.0, -462.0), pow(2.0, 672.0)),
295                             ldexp<MathOptimization::Unsafe>(one, setSimdIntFrom3I(587, -462, 672)));
296 #        endif
297     // The default safe version must be able to handle very negative arguments too
298     GMX_EXPECT_SIMD_REAL_EQ(setZero(), ldexp(one, setSimdIntFrom3I(-2000, -1000000, -1000000000)));
299 }
300
301 /*
302  * We do extensive 1/sqrt(x) and 1/x accuracy testing in the math module, so
303  * we just make sure the lookup instructions appear to work here
304  */
305
306 TEST_F(SimdFloatingpointTest, rsqrt)
307 {
308     SimdReal x   = setSimdRealFrom3R(4.0, M_PI, 1234567890.0);
309     SimdReal ref = setSimdRealFrom3R(0.5, 1.0 / std::sqrt(M_PI), 1.0 / std::sqrt(1234567890.0));
310     int      shiftbits = std::numeric_limits<real>::digits - GMX_SIMD_RSQRT_BITS;
311
312     if (shiftbits < 0)
313     {
314         shiftbits = 0;
315     }
316
317     /* Set the allowed ulp error as 2 to the power of the number of bits in
318      * the mantissa that do not have to be correct after the table lookup.
319      */
320     setUlpTol(1LL << shiftbits);
321     GMX_EXPECT_SIMD_REAL_NEAR(ref, rsqrt(x));
322 }
323
324 TEST_F(SimdFloatingpointTest, maskzRsqrt)
325 {
326     SimdReal x = setSimdRealFrom3R(M_PI, -4.0, 0.0);
327     // simdCmpLe is tested separately further down
328     SimdBool m         = setZero() < x;
329     SimdReal ref       = setSimdRealFrom3R(1.0 / std::sqrt(M_PI), 0.0, 0.0);
330     int      shiftbits = std::numeric_limits<real>::digits - GMX_SIMD_RSQRT_BITS;
331
332     if (shiftbits < 0)
333     {
334         shiftbits = 0;
335     }
336
337     /* Set the allowed ulp error as 2 to the power of the number of bits in
338      * the mantissa that do not have to be correct after the table lookup.
339      */
340     setUlpTol(1LL << shiftbits);
341     GMX_EXPECT_SIMD_REAL_NEAR(ref, maskzRsqrt(x, m));
342 }
343
344 TEST_F(SimdFloatingpointTest, rcp)
345 {
346     SimdReal x         = setSimdRealFrom3R(4.0, M_PI, 1234567890.0);
347     SimdReal ref       = setSimdRealFrom3R(0.25, 1.0 / M_PI, 1.0 / 1234567890.0);
348     int      shiftbits = std::numeric_limits<real>::digits - GMX_SIMD_RCP_BITS;
349
350     if (shiftbits < 0)
351     {
352         shiftbits = 0;
353     }
354
355     /* Set the allowed ulp error as 2 to the power of the number of bits in
356      * the mantissa that do not have to be correct after the table lookup.
357      */
358     setUlpTol(1LL << shiftbits);
359     GMX_EXPECT_SIMD_REAL_NEAR(ref, rcp(x));
360 }
361
362 TEST_F(SimdFloatingpointTest, maskzRcp)
363 {
364     SimdReal x         = setSimdRealFrom3R(M_PI, 0.0, -1234567890.0);
365     SimdBool m         = (x != setZero());
366     SimdReal ref       = setSimdRealFrom3R(1.0 / M_PI, 0.0, -1.0 / 1234567890.0);
367     int      shiftbits = std::numeric_limits<real>::digits - GMX_SIMD_RCP_BITS;
368
369     if (shiftbits < 0)
370     {
371         shiftbits = 0;
372     }
373
374     /* Set the allowed ulp error as 2 to the power of the number of bits in
375      * the mantissa that do not have to be correct after the table lookup.
376      */
377     setUlpTol(1LL << shiftbits);
378     GMX_EXPECT_SIMD_REAL_NEAR(ref, maskzRcp(x, m));
379 }
380
381 TEST_F(SimdFloatingpointTest, cmpEqAndSelectByMask)
382 {
383     SimdBool eq = rSimd_c4c6c8 == rSimd_c6c7c8;
384     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0, 0, c2), selectByMask(rSimd_c0c1c2, eq));
385 }
386
387 TEST_F(SimdFloatingpointTest, selectByNotMask)
388 {
389     SimdBool eq = rSimd_c4c6c8 == rSimd_c6c7c8;
390     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0, c1, 0), selectByNotMask(rSimd_c0c1c2, eq));
391 }
392
393 TEST_F(SimdFloatingpointTest, cmpNe)
394 {
395     SimdBool eq = rSimd_c4c6c8 != rSimd_c6c7c8;
396     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0, c1, 0), selectByMask(rSimd_c0c1c2, eq));
397 }
398
399 TEST_F(SimdFloatingpointTest, cmpLe)
400 {
401     SimdBool le = rSimd_c4c6c8 <= rSimd_c6c7c8;
402     GMX_EXPECT_SIMD_REAL_EQ(rSimd_c0c1c2, selectByMask(rSimd_c0c1c2, le));
403 }
404
405 TEST_F(SimdFloatingpointTest, cmpLt)
406 {
407     SimdBool lt = rSimd_c4c6c8 < rSimd_c6c7c8;
408     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0, c1, 0), selectByMask(rSimd_c0c1c2, lt));
409 }
410
411 #        if GMX_SIMD_HAVE_INT32_LOGICAL || GMX_SIMD_HAVE_LOGICAL
412 TEST_F(SimdFloatingpointTest, testBits)
413 {
414     SimdBool eq = testBits(setSimdRealFrom3R(c1, 0, c1));
415     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0, 0, c2), selectByMask(rSimd_c0c1c2, eq));
416
417     // Test if we detect only the sign bit being set
418     eq = testBits(setSimdRealFrom1R(GMX_REAL_NEGZERO));
419     GMX_EXPECT_SIMD_REAL_EQ(rSimd_c0c1c2, selectByMask(rSimd_c0c1c2, eq));
420 }
421 #        endif
422
423 TEST_F(SimdFloatingpointTest, andB)
424 {
425     SimdBool eq = rSimd_c4c6c8 == rSimd_c6c7c8;
426     SimdBool le = rSimd_c4c6c8 <= rSimd_c6c7c8;
427     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0, 0, c2), selectByMask(rSimd_c0c1c2, (eq && le)));
428 }
429
430 TEST_F(SimdFloatingpointTest, orB)
431 {
432     SimdBool eq = rSimd_c4c6c8 == rSimd_c6c7c8;
433     SimdBool lt = rSimd_c4c6c8 < rSimd_c6c7c8;
434     GMX_EXPECT_SIMD_REAL_EQ(rSimd_c0c1c2, selectByMask(rSimd_c0c1c2, (eq || lt)));
435 }
436
437 TEST_F(SimdFloatingpointTest, anyTrueB)
438 {
439     alignas(GMX_SIMD_ALIGNMENT) std::array<real, GMX_SIMD_REAL_WIDTH> mem{};
440
441     // Test the false case
442     EXPECT_FALSE(anyTrue(setZero() < load<SimdReal>(mem.data())));
443
444     // Test each bit (these should all be true)
445     for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
446     {
447         mem.fill(0.0);
448         mem[i] = 1.0;
449         EXPECT_TRUE(anyTrue(setZero() < load<SimdReal>(mem.data())))
450                 << "Not detecting true in element " << i;
451     }
452 }
453
454 TEST_F(SimdFloatingpointTest, blend)
455 {
456     SimdBool lt = rSimd_c4c6c8 < rSimd_c6c7c8;
457     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c3, c4, c2), blend(rSimd_c0c1c2, rSimd_c3c4c5, lt));
458 }
459
460 TEST_F(SimdFloatingpointTest, reduce)
461 {
462     // The horizontal sum of the SIMD variable depends on the width, so
463     // simply store it an extra time and calculate what the sum should be
464     std::vector<real> v   = simdReal2Vector(rSimd_c3c4c5);
465     real              sum = 0.0;
466
467     for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
468     {
469         sum += v[i];
470     }
471
472     EXPECT_REAL_EQ_TOL(sum, reduce(rSimd_c3c4c5), defaultRealTolerance());
473 }
474
475 #    endif // GMX_SIMD_HAVE_REAL
476
477 #    if GMX_SIMD_HAVE_FLOAT && GMX_SIMD_HAVE_DOUBLE
478 TEST_F(SimdFloatingpointTest, cvtFloat2Double)
479 {
480     alignas(GMX_SIMD_ALIGNMENT) float f[GMX_SIMD_FLOAT_WIDTH];
481     alignas(GMX_SIMD_ALIGNMENT) double d[GMX_SIMD_FLOAT_WIDTH]; // Yes, double array length should be same as float
482
483     int                    i;
484     SimdFloat              vf;
485     SimdDouble             vd0;
486     FloatingPointTolerance tolerance(defaultRealTolerance());
487
488     for (i = 0; i < GMX_SIMD_FLOAT_WIDTH; i++)
489     {
490         // Scale by 1+100*eps to use low bits too.
491         // Due to the conversions we want to avoid being too sensitive to fluctuations in last bit
492         f[i] = i * (1.0 + 100 * GMX_FLOAT_EPS);
493     }
494
495     vf = load<SimdFloat>(f);
496 #        if (GMX_SIMD_FLOAT_WIDTH == 2 * GMX_SIMD_DOUBLE_WIDTH)
497     SimdDouble vd1;
498     cvtF2DD(vf, &vd0, &vd1);
499     store(d + GMX_SIMD_DOUBLE_WIDTH, vd1); // Store upper part halfway through array
500 #        elif (GMX_SIMD_FLOAT_WIDTH == GMX_SIMD_DOUBLE_WIDTH)
501     vd0 = cvtF2D(vf);
502 #        else
503 #            error Width of float SIMD must either be identical to double, or twice the width.
504 #        endif
505     store(d, vd0); // store lower (or whole) part from start of vector
506
507     for (i = 0; i < GMX_SIMD_FLOAT_WIDTH; i++)
508     {
509         EXPECT_REAL_EQ_TOL(f[i], d[i], tolerance);
510     }
511 }
512
513 TEST_F(SimdFloatingpointTest, cvtDouble2Float)
514 {
515     alignas(GMX_SIMD_ALIGNMENT) float f[GMX_SIMD_FLOAT_WIDTH];
516     alignas(GMX_SIMD_ALIGNMENT) double d[GMX_SIMD_FLOAT_WIDTH]; // Yes, double array length should be same as float
517     int                    i;
518     SimdFloat              vf;
519     SimdDouble             vd0;
520     FloatingPointTolerance tolerance(defaultRealTolerance());
521
522     // This fills elements for pd1 too when double width is 2*single width
523     for (i = 0; i < GMX_SIMD_FLOAT_WIDTH; i++)
524     {
525         // Scale by 1+eps to use low bits too.
526         // Due to the conversions we want to avoid being too sensitive to fluctuations in last bit
527         d[i] = i * (1.0 + 100 * GMX_FLOAT_EPS);
528     }
529
530     vd0 = load<SimdDouble>(d);
531 #        if (GMX_SIMD_FLOAT_WIDTH == 2 * GMX_SIMD_DOUBLE_WIDTH)
532     SimdDouble vd1 = load<SimdDouble>(d + GMX_SIMD_DOUBLE_WIDTH); // load upper half of data
533     vf             = cvtDD2F(vd0, vd1);
534 #        elif (GMX_SIMD_FLOAT_WIDTH == GMX_SIMD_DOUBLE_WIDTH)
535     vf  = cvtD2F(vd0);
536 #        else
537 #            error Width of float SIMD must either be identical to double, or twice the width.
538 #        endif
539     store(f, vf);
540
541     // This will check elements in pd1 too when double width is 2*single width
542     for (i = 0; i < GMX_SIMD_FLOAT_WIDTH; i++)
543     {
544         EXPECT_FLOAT_EQ_TOL(d[i], f[i], tolerance);
545     }
546 }
547 #    endif // GMX_SIMD_HAVE_FLOAT && GMX_SIMD_HAVE_DOUBLE
548
549 /*! \} */
550 /*! \endcond */
551
552 } // namespace
553 } // namespace test
554 } // namespace gmx
555
556 #endif // GMX_SIMD