825df83b6fdd07565eb5d5126a287b0681f6ae92
[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, 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 <cmath>
38
39 #include "gromacs/math/utilities.h"
40 #include "gromacs/simd/simd.h"
41 #include "gromacs/utility/basedefinitions.h"
42
43 #include "testutils/testasserts.h"
44
45 #include "data.h"
46 #include "simd.h"
47
48 #if GMX_SIMD
49
50 namespace gmx
51 {
52 namespace test
53 {
54
55 namespace
56 {
57
58 /*! \cond internal */
59 /*! \addtogroup module_simd */
60 /*! \{ */
61
62 #if GMX_SIMD_HAVE_REAL
63
64 /*! \brief Test fixture for floating-point tests (identical to the generic \ref SimdTest) */
65 typedef SimdTest SimdFloatingpointTest;
66
67 TEST_F(SimdFloatingpointTest, setZero)
68 {
69     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(0.0), setZero());
70 }
71
72 TEST_F(SimdFloatingpointTest, set)
73 {
74     const real *p  = &c0;
75     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(c1), SimdReal(c1));
76     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(c0), SimdReal(*p));
77 }
78
79 TEST_F(SimdFloatingpointTest, add)
80 {
81     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0 + c3, c1 + c4, c2 + c5 ),
82                               rSimd_c0c1c2 + rSimd_c3c4c5);
83 }
84
85 TEST_F(SimdFloatingpointTest, maskAdd)
86 {
87     SimdBool m = setSimdRealFrom3R(c6, 0, c7) != setZero();
88     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0 + c3, c1 + 0.0, c2 + c5 ),
89                               maskAdd(rSimd_c0c1c2, rSimd_c3c4c5, m));
90 }
91
92 TEST_F(SimdFloatingpointTest, sub)
93 {
94     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0 - c3, c1 - c4, c2 - c5 ),
95                               rSimd_c0c1c2 - rSimd_c3c4c5);
96 }
97
98 TEST_F(SimdFloatingpointTest, mul)
99 {
100     GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0 * c3, c1 * c4, c2 * c5 ),
101                               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,
164                             (rSimd_logicalA & rSimd_logicalB));
165 }
166
167 TEST_F(SimdFloatingpointTest, or)
168 {
169     GMX_EXPECT_SIMD_REAL_EQ(rSimd_logicalResultOr,
170                             (rSimd_logicalA | rSimd_logicalB));
171 }
172
173 TEST_F(SimdFloatingpointTest, xor)
174 {
175     /* Test xor by taking xor with a number and its negative. This should result
176      * in only the sign bit being set. We then use this bit change the sign of
177      * different numbers.
178      */
179     SimdReal signbit = SimdReal(c1) ^ SimdReal(-c1);
180     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-c2, c3, -c4), (signbit ^ setSimdRealFrom3R(c2, -c3, c4)));
181 }
182
183 TEST_F(SimdFloatingpointTest, andNot)
184 {
185     /* Use xor (which we already tested, so fix that first if both tests fail)
186      * to extract the sign bit, and then use andnot to take absolute values.
187      */
188     SimdReal signbit = SimdReal(c1) ^ SimdReal(-c1);
189     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c2, c3, c4), andNot(signbit, setSimdRealFrom3R(-c2, c3, -c4)));
190 }
191
192 #endif
193
194 TEST_F(SimdFloatingpointTest, max)
195 {
196     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R( c3,  c1,  c4), max(rSimd_c0c1c2, rSimd_c3c0c4));
197     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R( c3,  c1,  c4), max(rSimd_c3c0c4, rSimd_c0c1c2));
198     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-c0, -c0, -c2), max(rSimd_m0m1m2, rSimd_m3m0m4));
199     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-c0, -c0, -c2), max(rSimd_m3m0m4, rSimd_m0m1m2));
200 }
201
202 TEST_F(SimdFloatingpointTest, min)
203 {
204     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R( c0,  c0,  c2), min(rSimd_c0c1c2, rSimd_c3c0c4));
205     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R( c0,  c0,  c2), min(rSimd_c3c0c4, rSimd_c0c1c2));
206     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-c3, -c1, -c4), min(rSimd_m0m1m2, rSimd_m3m0m4));
207     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-c3, -c1, -c4), min(rSimd_m3m0m4, rSimd_m0m1m2));
208 }
209
210 TEST_F(SimdFloatingpointTest, round)
211 {
212     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(2), round(rSimd_2p25));
213     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(4), round(rSimd_3p75));
214     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-2), round(rSimd_m2p25));
215     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-4), round(rSimd_m3p75));
216 }
217
218 TEST_F(SimdFloatingpointTest, roundMode)
219 {
220     /* Rounding mode needs to be consistent between round and cvtR2I */
221     SimdReal x0  = setSimdRealFrom3R(0.5, 11.5, 99.5);
222     SimdReal x1  = setSimdRealFrom3R(-0.5, -11.5, -99.5);
223
224     GMX_EXPECT_SIMD_REAL_EQ(round(x0), cvtI2R(cvtR2I(x0)));
225     GMX_EXPECT_SIMD_REAL_EQ(round(x1), cvtI2R(cvtR2I(x1)));
226 }
227
228 TEST_F(SimdFloatingpointTest, trunc)
229 {
230     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(2), trunc(rSimd_2p25));
231     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(3), trunc(rSimd_3p75));
232     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-2), trunc(rSimd_m2p25));
233     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-3), trunc(rSimd_m3p75));
234 }
235
236 // We explicitly test the exponent/mantissa routines with double precision data,
237 // since these usually rely on direct manipulation and shift of the SIMD registers,
238 // where it is easy to make mistakes with single vs double precision.
239
240 TEST_F(SimdFloatingpointTest, frexp)
241 {
242     SimdReal  fraction;
243     SimdInt32 exponent;
244
245     fraction = frexp(rSimd_Exp, &exponent);
246
247     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0.609548660288905419513128,
248                                               0.5833690139241746175358116,
249                                               -0.584452007502232362412542),
250                             fraction);
251     GMX_EXPECT_SIMD_INT_EQ(setSimdIntFrom3I(61, -40, 55), exponent);
252
253
254 #if GMX_SIMD_HAVE_DOUBLE && GMX_DOUBLE
255     fraction = frexp(rSimd_ExpDouble, &exponent);
256
257     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0.6206306194761728178832527,
258                                               0.5236473618795619566768096,
259                                               -0.9280331023751380303821179),
260                             fraction);
261     GMX_EXPECT_SIMD_INT_EQ(setSimdIntFrom3I(588, -461, 673), exponent);
262 #endif
263 }
264
265 TEST_F(SimdFloatingpointTest, ldexp)
266 {
267     SimdReal one = setSimdRealFrom1R(1.0);
268
269     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(pow(2.0, 60.0), pow(2.0, -41.0), pow(2.0, 54.0)),
270                             ldexp<MathOptimization::Unsafe>(one, setSimdIntFrom3I(60, -41, 54)));
271 #if GMX_SIMD_HAVE_DOUBLE && GMX_DOUBLE
272     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(pow(2.0, 587.0), pow(2.0, -462.0), pow(2.0, 672.0)),
273                             ldexp<MathOptimization::Unsafe>(one, setSimdIntFrom3I(587, -462, 672)));
274 #endif
275     // The default safe version must be able to handle very negative arguments too
276     GMX_EXPECT_SIMD_REAL_EQ(setZero(), ldexp(one, setSimdIntFrom3I(-2000, -1000000, -1000000000)));
277 }
278
279 /*
280  * We do extensive 1/sqrt(x) and 1/x accuracy testing in the math module, so
281  * we just make sure the lookup instructions appear to work here
282  */
283
284 TEST_F(SimdFloatingpointTest, rsqrt)
285 {
286     SimdReal        x                  = setSimdRealFrom3R(4.0, M_PI, 1234567890.0);
287     SimdReal        ref                = setSimdRealFrom3R(0.5, 1.0/std::sqrt(M_PI), 1.0/std::sqrt(1234567890.0));
288     int             shiftbits          = std::numeric_limits<real>::digits-GMX_SIMD_RSQRT_BITS;
289
290     if (shiftbits < 0)
291     {
292         shiftbits = 0;
293     }
294
295     /* Set the allowed ulp error as 2 to the power of the number of bits in
296      * the mantissa that do not have to be correct after the table lookup.
297      */
298     setUlpTol(1LL << shiftbits);
299     GMX_EXPECT_SIMD_REAL_NEAR(ref, rsqrt(x));
300 }
301
302 TEST_F(SimdFloatingpointTest, maskzRsqrt)
303 {
304     SimdReal        x                  = setSimdRealFrom3R(M_PI, -4.0, 0.0);
305     // simdCmpLe is tested separately further down
306     SimdBool        m                  = setZero() < x;
307     SimdReal        ref                = setSimdRealFrom3R(1.0/std::sqrt(M_PI), 0.0, 0.0);
308     int             shiftbits          = std::numeric_limits<real>::digits-GMX_SIMD_RSQRT_BITS;
309
310     if (shiftbits < 0)
311     {
312         shiftbits = 0;
313     }
314
315     /* Set the allowed ulp error as 2 to the power of the number of bits in
316      * the mantissa that do not have to be correct after the table lookup.
317      */
318     setUlpTol(1LL << shiftbits);
319     GMX_EXPECT_SIMD_REAL_NEAR(ref, maskzRsqrt(x, m));
320 }
321
322 TEST_F(SimdFloatingpointTest, rcp)
323 {
324     SimdReal        x                  = setSimdRealFrom3R(4.0, M_PI, 1234567890.0);
325     SimdReal        ref                = setSimdRealFrom3R(0.25, 1.0/M_PI, 1.0/1234567890.0);
326     int             shiftbits          = std::numeric_limits<real>::digits-GMX_SIMD_RCP_BITS;
327
328     if (shiftbits < 0)
329     {
330         shiftbits = 0;
331     }
332
333     /* Set the allowed ulp error as 2 to the power of the number of bits in
334      * the mantissa that do not have to be correct after the table lookup.
335      */
336     setUlpTol(1LL << shiftbits);
337     GMX_EXPECT_SIMD_REAL_NEAR(ref, rcp(x));
338 }
339
340 TEST_F(SimdFloatingpointTest, maskzRcp)
341 {
342     SimdReal        x                  = setSimdRealFrom3R(M_PI, 0.0, -1234567890.0);
343     SimdBool        m                  = (x != setZero());
344     SimdReal        ref                = setSimdRealFrom3R(1.0/M_PI, 0.0, -1.0/1234567890.0);
345     int             shiftbits          = std::numeric_limits<real>::digits-GMX_SIMD_RCP_BITS;
346
347     if (shiftbits < 0)
348     {
349         shiftbits = 0;
350     }
351
352     /* Set the allowed ulp error as 2 to the power of the number of bits in
353      * the mantissa that do not have to be correct after the table lookup.
354      */
355     setUlpTol(1LL << shiftbits);
356     GMX_EXPECT_SIMD_REAL_NEAR(ref, maskzRcp(x, m));
357 }
358
359 TEST_F(SimdFloatingpointTest, cmpEqAndSelectByMask)
360 {
361     SimdBool eq   = rSimd_c4c6c8 == rSimd_c6c7c8;
362     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0, 0, c2), selectByMask(rSimd_c0c1c2, eq));
363 }
364
365 TEST_F(SimdFloatingpointTest, selectByNotMask)
366 {
367     SimdBool eq   = rSimd_c4c6c8 == rSimd_c6c7c8;
368     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0, c1, 0), selectByNotMask(rSimd_c0c1c2, eq));
369 }
370
371 TEST_F(SimdFloatingpointTest, cmpNe)
372 {
373     SimdBool eq   = rSimd_c4c6c8 != rSimd_c6c7c8;
374     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0, c1, 0), selectByMask(rSimd_c0c1c2, eq));
375 }
376
377 TEST_F(SimdFloatingpointTest, cmpLe)
378 {
379     SimdBool le   = rSimd_c4c6c8 <= rSimd_c6c7c8;
380     GMX_EXPECT_SIMD_REAL_EQ(rSimd_c0c1c2, selectByMask(rSimd_c0c1c2, le));
381 }
382
383 TEST_F(SimdFloatingpointTest, cmpLt)
384 {
385     SimdBool lt   = rSimd_c4c6c8 < rSimd_c6c7c8;
386     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0, c1, 0), selectByMask(rSimd_c0c1c2, lt));
387 }
388
389 #if GMX_SIMD_HAVE_INT32_LOGICAL || GMX_SIMD_HAVE_LOGICAL
390 TEST_F(SimdFloatingpointTest, testBits)
391 {
392     SimdBool eq   = testBits(setSimdRealFrom3R(c1, 0, c1));
393     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0, 0, c2), selectByMask(rSimd_c0c1c2, eq));
394
395     // Test if we detect only the sign bit being set
396     eq            = testBits(setSimdRealFrom1R(GMX_REAL_NEGZERO));
397     GMX_EXPECT_SIMD_REAL_EQ(rSimd_c0c1c2, selectByMask(rSimd_c0c1c2, eq));
398 }
399 #endif
400
401 TEST_F(SimdFloatingpointTest, andB)
402 {
403     SimdBool eq   = rSimd_c4c6c8 == rSimd_c6c7c8;
404     SimdBool le   = rSimd_c4c6c8 <= rSimd_c6c7c8;
405     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0, 0, c2), selectByMask(rSimd_c0c1c2, (eq && le)));
406 }
407
408 TEST_F(SimdFloatingpointTest, orB)
409 {
410     SimdBool eq   = rSimd_c4c6c8 == rSimd_c6c7c8;
411     SimdBool lt   = rSimd_c4c6c8  < rSimd_c6c7c8;
412     GMX_EXPECT_SIMD_REAL_EQ(rSimd_c0c1c2, selectByMask(rSimd_c0c1c2, (eq || lt)));
413 }
414
415 TEST_F(SimdFloatingpointTest, anyTrueB)
416 {
417     SimdBool eq;
418
419     /* this test is a bit tricky since we don't know the simd width.
420      * We cannot check for truth values for "any" element beyond the first,
421      * since that part of the data will not be used if simd width is 1.
422      */
423     eq = rSimd_c4c6c8 == setSimdRealFrom3R(c4, 0, 0);
424     EXPECT_TRUE(anyTrue(eq));
425
426     eq = rSimd_c0c1c2 == rSimd_c3c4c5;
427     EXPECT_FALSE(anyTrue(eq));
428 }
429
430 TEST_F(SimdFloatingpointTest, blend)
431 {
432     SimdBool lt   = rSimd_c4c6c8 < rSimd_c6c7c8;
433     GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c3, c4, c2), blend(rSimd_c0c1c2, rSimd_c3c4c5, lt));
434 }
435
436 TEST_F(SimdFloatingpointTest, reduce)
437 {
438     // The horizontal sum of the SIMD variable depends on the width, so
439     // simply store it an extra time and calculate what the sum should be
440     std::vector<real> v   = simdReal2Vector(rSimd_c3c4c5);
441     real              sum = 0.0;
442
443     for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
444     {
445         sum += v[i];
446     }
447
448     EXPECT_REAL_EQ_TOL(sum, reduce(rSimd_c3c4c5), defaultRealTolerance() );
449 }
450
451 #endif      // GMX_SIMD_HAVE_REAL
452
453 #if GMX_SIMD_HAVE_FLOAT && GMX_SIMD_HAVE_DOUBLE
454 TEST_F(SimdFloatingpointTest, cvtFloat2Double)
455 {
456     alignas(GMX_SIMD_ALIGNMENT) float   f[GMX_SIMD_FLOAT_WIDTH];
457     alignas(GMX_SIMD_ALIGNMENT) double  d[GMX_SIMD_FLOAT_WIDTH];  // Yes, double array length should be same as float
458
459     int                               i;
460     SimdFloat                         vf;
461     SimdDouble                        vd0;
462     FloatingPointTolerance            tolerance(defaultRealTolerance());
463
464     for (i = 0; i < GMX_SIMD_FLOAT_WIDTH; i++)
465     {
466         // Scale by 1+100*eps to use low bits too.
467         // Due to the conversions we want to avoid being too sensitive to fluctuations in last bit
468         f[i] = i * (1.0 + 100*GMX_FLOAT_EPS);
469     }
470
471     vf = load<SimdFloat>(f);
472 #if (GMX_SIMD_FLOAT_WIDTH == 2*GMX_SIMD_DOUBLE_WIDTH)
473     SimdDouble vd1;
474     cvtF2DD(vf, &vd0, &vd1);
475     store(d + GMX_SIMD_DOUBLE_WIDTH, vd1); // Store upper part halfway through array
476 #elif (GMX_SIMD_FLOAT_WIDTH == GMX_SIMD_DOUBLE_WIDTH)
477     vd0 = cvtF2D(vf);
478 #else
479 #    error Width of float SIMD must either be identical to double, or twice the width.
480 #endif
481     store(d, vd0); // store lower (or whole) part from start of vector
482
483     for (i = 0; i < GMX_SIMD_FLOAT_WIDTH; i++)
484     {
485         EXPECT_REAL_EQ_TOL(f[i], d[i], tolerance);
486     }
487 }
488
489 TEST_F(SimdFloatingpointTest, cvtDouble2Float)
490 {
491     alignas(GMX_SIMD_ALIGNMENT) float     f[GMX_SIMD_FLOAT_WIDTH];
492     alignas(GMX_SIMD_ALIGNMENT) double    d[GMX_SIMD_FLOAT_WIDTH];  // Yes, double array length should be same as float
493     int                               i;
494     SimdFloat                         vf;
495     SimdDouble                        vd0;
496     FloatingPointTolerance            tolerance(defaultRealTolerance());
497
498     // This fills elements for pd1 too when double width is 2*single width
499     for (i = 0; i < GMX_SIMD_FLOAT_WIDTH; i++)
500     {
501         // Scale by 1+eps to use low bits too.
502         // Due to the conversions we want to avoid being too sensitive to fluctuations in last bit
503         d[i] = i * (1.0 + 100*GMX_FLOAT_EPS);
504     }
505
506     vd0 = load<SimdDouble>(d);
507 #if (GMX_SIMD_FLOAT_WIDTH == 2*GMX_SIMD_DOUBLE_WIDTH)
508     SimdDouble vd1 = load<SimdDouble>(d + GMX_SIMD_DOUBLE_WIDTH); // load upper half of data
509     vf = cvtDD2F(vd0, vd1);
510 #elif (GMX_SIMD_FLOAT_WIDTH == GMX_SIMD_DOUBLE_WIDTH)
511     vf = cvtD2F(vd0);
512 #else
513 #    error Width of float SIMD must either be identical to double, or twice the width.
514 #endif
515     store(f, vf);
516
517     // This will check elements in pd1 too when double width is 2*single width
518     for (i = 0; i < GMX_SIMD_FLOAT_WIDTH; i++)
519     {
520         EXPECT_FLOAT_EQ_TOL(d[i], f[i], tolerance);
521     }
522 }
523 #endif      // GMX_SIMD_HAVE_FLOAT && GMX_SIMD_HAVE_DOUBLE
524
525 /*! \} */
526 /*! \endcond */
527
528 }      // namespace
529 }      // namespace
530 }      // namespace
531
532 #endif // GMX_SIMD