2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017, 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.
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.
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.
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.
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.
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.
39 #include "gromacs/math/utilities.h"
40 #include "gromacs/simd/simd.h"
41 #include "gromacs/utility/basedefinitions.h"
43 #include "testutils/testasserts.h"
59 /*! \addtogroup module_simd */
62 #if GMX_SIMD_HAVE_REAL
64 /*! \brief Test fixture for floating-point tests (identical to the generic \ref SimdTest) */
65 typedef SimdTest SimdFloatingpointTest
;
67 TEST_F(SimdFloatingpointTest
, setZero
)
69 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(0.0), setZero());
72 TEST_F(SimdFloatingpointTest
, set
)
75 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(c1
), SimdReal(c1
));
76 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(c0
), SimdReal(*p
));
79 TEST_F(SimdFloatingpointTest
, add
)
81 GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0
+ c3
, c1
+ c4
, c2
+ c5
),
82 rSimd_c0c1c2
+ rSimd_c3c4c5
);
85 TEST_F(SimdFloatingpointTest
, maskAdd
)
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
));
92 TEST_F(SimdFloatingpointTest
, sub
)
94 GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0
- c3
, c1
- c4
, c2
- c5
),
95 rSimd_c0c1c2
- rSimd_c3c4c5
);
98 TEST_F(SimdFloatingpointTest
, mul
)
100 GMX_EXPECT_SIMD_REAL_NEAR(setSimdRealFrom3R(c0
* c3
, c1
* c4
, c2
* c5
),
101 rSimd_c0c1c2
* rSimd_c3c4c5
);
104 TEST_F(SimdFloatingpointTest
, maskzMul
)
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
));
111 TEST_F(SimdFloatingpointTest
, fma
)
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
));
119 TEST_F(SimdFloatingpointTest
, maskzFma
)
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
));
127 TEST_F(SimdFloatingpointTest
, fms
)
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
));
134 TEST_F(SimdFloatingpointTest
, fnma
)
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
));
141 TEST_F(SimdFloatingpointTest
, fnms
)
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
));
148 TEST_F(SimdFloatingpointTest
, abs
)
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
154 TEST_F(SimdFloatingpointTest
, neg
)
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
160 #if GMX_SIMD_HAVE_LOGICAL
161 TEST_F(SimdFloatingpointTest
, and)
163 GMX_EXPECT_SIMD_REAL_EQ(rSimd_logicalResultAnd
,
164 (rSimd_logicalA
& rSimd_logicalB
));
167 TEST_F(SimdFloatingpointTest
, or)
169 GMX_EXPECT_SIMD_REAL_EQ(rSimd_logicalResultOr
,
170 (rSimd_logicalA
| rSimd_logicalB
));
173 TEST_F(SimdFloatingpointTest
, xor)
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
179 SimdReal signbit
= SimdReal(c1
) ^ SimdReal(-c1
);
180 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(-c2
, c3
, -c4
), (signbit
^ setSimdRealFrom3R(c2
, -c3
, c4
)));
183 TEST_F(SimdFloatingpointTest
, andNot
)
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.
188 SimdReal signbit
= SimdReal(c1
) ^ SimdReal(-c1
);
189 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c2
, c3
, c4
), andNot(signbit
, setSimdRealFrom3R(-c2
, c3
, -c4
)));
194 TEST_F(SimdFloatingpointTest
, max
)
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
));
202 TEST_F(SimdFloatingpointTest
, min
)
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
));
210 TEST_F(SimdFloatingpointTest
, round
)
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
));
218 TEST_F(SimdFloatingpointTest
, trunc
)
220 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(2), trunc(rSimd_2p25
));
221 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(3), trunc(rSimd_3p75
));
222 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-2), trunc(rSimd_m2p25
));
223 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom1R(-3), trunc(rSimd_m3p75
));
226 // We explicitly test the exponent/mantissa routines with double precision data,
227 // since these usually rely on direct manipulation and shift of the SIMD registers,
228 // where it is easy to make mistakes with single vs double precision.
230 TEST_F(SimdFloatingpointTest
, frexp
)
235 fraction
= frexp(rSimd_Exp
, &exponent
);
237 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0.609548660288905419513128,
238 0.5833690139241746175358116,
239 -0.584452007502232362412542),
241 GMX_EXPECT_SIMD_INT_EQ(setSimdIntFrom3I(61, -40, 55), exponent
);
244 #if GMX_SIMD_HAVE_DOUBLE && GMX_DOUBLE
245 fraction
= frexp(rSimd_ExpDouble
, &exponent
);
247 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0.6206306194761728178832527,
248 0.5236473618795619566768096,
249 -0.9280331023751380303821179),
251 GMX_EXPECT_SIMD_INT_EQ(setSimdIntFrom3I(588, -461, 673), exponent
);
255 TEST_F(SimdFloatingpointTest
, ldexp
)
257 SimdReal x0
= setSimdRealFrom3R(0.5, 11.5, 99.5);
258 SimdReal x1
= setSimdRealFrom3R(-0.5, -11.5, -99.5);
259 SimdReal one
= setSimdRealFrom1R(1.0);
261 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(pow(2.0, 60.0), pow(2.0, -41.0), pow(2.0, 54.0)),
262 ldexp(one
, setSimdIntFrom3I(60, -41, 54)));
263 #if GMX_SIMD_HAVE_DOUBLE && GMX_DOUBLE
264 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(pow(2.0, 587.0), pow(2.0, -462.0), pow(2.0, 672.0)),
265 ldexp(one
, setSimdIntFrom3I(587, -462, 672)));
267 /* Rounding mode in conversions must be consistent with simdRound() for SetExponent() to work */
268 GMX_EXPECT_SIMD_REAL_EQ(ldexp(one
, cvtR2I(round(x0
))), ldexp(one
, cvtR2I(x0
)));
269 GMX_EXPECT_SIMD_REAL_EQ(ldexp(one
, cvtR2I(round(x1
))), ldexp(one
, cvtR2I(x1
)));
273 * We do extensive 1/sqrt(x) and 1/x accuracy testing in the math module, so
274 * we just make sure the lookup instructions appear to work here
277 TEST_F(SimdFloatingpointTest
, rsqrt
)
279 SimdReal x
= setSimdRealFrom3R(4.0, M_PI
, 1234567890.0);
280 SimdReal ref
= setSimdRealFrom3R(0.5, 1.0/std::sqrt(M_PI
), 1.0/std::sqrt(1234567890.0));
281 int shiftbits
= std::numeric_limits
<real
>::digits
-GMX_SIMD_RSQRT_BITS
;
288 /* Set the allowed ulp error as 2 to the power of the number of bits in
289 * the mantissa that do not have to be correct after the table lookup.
291 setUlpTol(1LL << shiftbits
);
292 GMX_EXPECT_SIMD_REAL_NEAR(ref
, rsqrt(x
));
295 TEST_F(SimdFloatingpointTest
, maskzRsqrt
)
297 SimdReal x
= setSimdRealFrom3R(M_PI
, -4.0, 0.0);
298 // simdCmpLe is tested separately further down
299 SimdBool m
= setZero() < x
;
300 SimdReal ref
= setSimdRealFrom3R(1.0/std::sqrt(M_PI
), 0.0, 0.0);
301 int shiftbits
= std::numeric_limits
<real
>::digits
-GMX_SIMD_RSQRT_BITS
;
308 /* Set the allowed ulp error as 2 to the power of the number of bits in
309 * the mantissa that do not have to be correct after the table lookup.
311 setUlpTol(1LL << shiftbits
);
312 GMX_EXPECT_SIMD_REAL_NEAR(ref
, maskzRsqrt(x
, m
));
315 TEST_F(SimdFloatingpointTest
, rcp
)
317 SimdReal x
= setSimdRealFrom3R(4.0, M_PI
, 1234567890.0);
318 SimdReal ref
= setSimdRealFrom3R(0.25, 1.0/M_PI
, 1.0/1234567890.0);
319 int shiftbits
= std::numeric_limits
<real
>::digits
-GMX_SIMD_RCP_BITS
;
326 /* Set the allowed ulp error as 2 to the power of the number of bits in
327 * the mantissa that do not have to be correct after the table lookup.
329 setUlpTol(1LL << shiftbits
);
330 GMX_EXPECT_SIMD_REAL_NEAR(ref
, rcp(x
));
333 TEST_F(SimdFloatingpointTest
, maskzRcp
)
335 SimdReal x
= setSimdRealFrom3R(M_PI
, 0.0, -1234567890.0);
336 SimdBool m
= (x
!= setZero());
337 SimdReal ref
= setSimdRealFrom3R(1.0/M_PI
, 0.0, -1.0/1234567890.0);
338 int shiftbits
= std::numeric_limits
<real
>::digits
-GMX_SIMD_RCP_BITS
;
345 /* Set the allowed ulp error as 2 to the power of the number of bits in
346 * the mantissa that do not have to be correct after the table lookup.
348 setUlpTol(1LL << shiftbits
);
349 GMX_EXPECT_SIMD_REAL_NEAR(ref
, maskzRcp(x
, m
));
352 TEST_F(SimdFloatingpointTest
, cmpEqAndSelectByMask
)
354 SimdBool eq
= rSimd_c4c6c8
== rSimd_c6c7c8
;
355 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0, 0, c2
), selectByMask(rSimd_c0c1c2
, eq
));
358 TEST_F(SimdFloatingpointTest
, selectByNotMask
)
360 SimdBool eq
= rSimd_c4c6c8
== rSimd_c6c7c8
;
361 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0
, c1
, 0), selectByNotMask(rSimd_c0c1c2
, eq
));
364 TEST_F(SimdFloatingpointTest
, cmpNe
)
366 SimdBool eq
= rSimd_c4c6c8
!= rSimd_c6c7c8
;
367 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0
, c1
, 0), selectByMask(rSimd_c0c1c2
, eq
));
370 TEST_F(SimdFloatingpointTest
, cmpLe
)
372 SimdBool le
= rSimd_c4c6c8
<= rSimd_c6c7c8
;
373 GMX_EXPECT_SIMD_REAL_EQ(rSimd_c0c1c2
, selectByMask(rSimd_c0c1c2
, le
));
376 TEST_F(SimdFloatingpointTest
, cmpLt
)
378 SimdBool lt
= rSimd_c4c6c8
< rSimd_c6c7c8
;
379 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0
, c1
, 0), selectByMask(rSimd_c0c1c2
, lt
));
382 #if GMX_SIMD_HAVE_INT32_LOGICAL || GMX_SIMD_HAVE_LOGICAL
383 TEST_F(SimdFloatingpointTest
, testBits
)
385 SimdBool eq
= testBits(setSimdRealFrom3R(c1
, 0, c1
));
386 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c0
, 0, c2
), selectByMask(rSimd_c0c1c2
, eq
));
388 // Test if we detect only the sign bit being set
389 eq
= testBits(setSimdRealFrom1R(GMX_REAL_NEGZERO
));
390 GMX_EXPECT_SIMD_REAL_EQ(rSimd_c0c1c2
, selectByMask(rSimd_c0c1c2
, eq
));
394 TEST_F(SimdFloatingpointTest
, andB
)
396 SimdBool eq
= rSimd_c4c6c8
== rSimd_c6c7c8
;
397 SimdBool le
= rSimd_c4c6c8
<= rSimd_c6c7c8
;
398 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(0, 0, c2
), selectByMask(rSimd_c0c1c2
, (eq
&& le
)));
401 TEST_F(SimdFloatingpointTest
, orB
)
403 SimdBool eq
= rSimd_c4c6c8
== rSimd_c6c7c8
;
404 SimdBool lt
= rSimd_c4c6c8
< rSimd_c6c7c8
;
405 GMX_EXPECT_SIMD_REAL_EQ(rSimd_c0c1c2
, selectByMask(rSimd_c0c1c2
, (eq
|| lt
)));
408 TEST_F(SimdFloatingpointTest
, anyTrueB
)
412 /* this test is a bit tricky since we don't know the simd width.
413 * We cannot check for truth values for "any" element beyond the first,
414 * since that part of the data will not be used if simd width is 1.
416 eq
= rSimd_c4c6c8
== setSimdRealFrom3R(c4
, 0, 0);
417 EXPECT_TRUE(anyTrue(eq
));
419 eq
= rSimd_c0c1c2
== rSimd_c3c4c5
;
420 EXPECT_FALSE(anyTrue(eq
));
423 TEST_F(SimdFloatingpointTest
, blend
)
425 SimdBool lt
= rSimd_c4c6c8
< rSimd_c6c7c8
;
426 GMX_EXPECT_SIMD_REAL_EQ(setSimdRealFrom3R(c3
, c4
, c2
), blend(rSimd_c0c1c2
, rSimd_c3c4c5
, lt
));
429 TEST_F(SimdFloatingpointTest
, reduce
)
431 // The horizontal sum of the SIMD variable depends on the width, so
432 // simply store it an extra time and calculate what the sum should be
433 std::vector
<real
> v
= simdReal2Vector(rSimd_c3c4c5
);
436 for (int i
= 0; i
< GMX_SIMD_REAL_WIDTH
; i
++)
441 EXPECT_REAL_EQ_TOL(sum
, reduce(rSimd_c3c4c5
), defaultRealTolerance() );
444 #endif // GMX_SIMD_HAVE_REAL
446 #if GMX_SIMD_HAVE_FLOAT && GMX_SIMD_HAVE_DOUBLE
447 TEST_F(SimdFloatingpointTest
, cvtFloat2Double
)
449 GMX_ALIGNED(float, GMX_SIMD_FLOAT_WIDTH
) f
[GMX_SIMD_FLOAT_WIDTH
];
450 GMX_ALIGNED(double, GMX_SIMD_DOUBLE_WIDTH
) d
[GMX_SIMD_FLOAT_WIDTH
]; // Yes, double array length should be same as float
455 FloatingPointTolerance
tolerance(defaultRealTolerance());
457 for (i
= 0; i
< GMX_SIMD_FLOAT_WIDTH
; i
++)
459 // Scale by 1+100*eps to use low bits too.
460 // Due to the conversions we want to avoid being too sensitive to fluctuations in last bit
461 f
[i
] = i
* (1.0 + 100*GMX_FLOAT_EPS
);
465 #if (GMX_SIMD_FLOAT_WIDTH == 2*GMX_SIMD_DOUBLE_WIDTH)
467 cvtF2DD(vf
, &vd0
, &vd1
);
468 store(d
+ GMX_SIMD_DOUBLE_WIDTH
, vd1
); // Store upper part halfway through array
469 #elif (GMX_SIMD_FLOAT_WIDTH == GMX_SIMD_DOUBLE_WIDTH)
472 # error Width of float SIMD must either be identical to double, or twice the width.
474 store(d
, vd0
); // store lower (or whole) part from start of vector
476 for (i
= 0; i
< GMX_SIMD_FLOAT_WIDTH
; i
++)
478 EXPECT_REAL_EQ_TOL(f
[i
], d
[i
], tolerance
);
482 TEST_F(SimdFloatingpointTest
, cvtDouble2Float
)
484 GMX_ALIGNED(float, GMX_SIMD_FLOAT_WIDTH
) f
[GMX_SIMD_FLOAT_WIDTH
];
485 GMX_ALIGNED(double, GMX_SIMD_DOUBLE_WIDTH
) d
[GMX_SIMD_FLOAT_WIDTH
]; // Yes, double array length should be same as float
489 FloatingPointTolerance
tolerance(defaultRealTolerance());
491 // This fills elements for pd1 too when double width is 2*single width
492 for (i
= 0; i
< GMX_SIMD_FLOAT_WIDTH
; i
++)
494 // Scale by 1+eps to use low bits too.
495 // Due to the conversions we want to avoid being too sensitive to fluctuations in last bit
496 d
[i
] = i
* (1.0 + 100*GMX_FLOAT_EPS
);
500 #if (GMX_SIMD_FLOAT_WIDTH == 2*GMX_SIMD_DOUBLE_WIDTH)
501 SimdDouble vd1
= load(d
+ GMX_SIMD_DOUBLE_WIDTH
); // load upper half of data
502 vf
= cvtDD2F(vd0
, vd1
);
503 #elif (GMX_SIMD_FLOAT_WIDTH == GMX_SIMD_DOUBLE_WIDTH)
506 # error Width of float SIMD must either be identical to double, or twice the width.
510 // This will check elements in pd1 too when double width is 2*single width
511 for (i
= 0; i
< GMX_SIMD_FLOAT_WIDTH
; i
++)
513 EXPECT_REAL_EQ_TOL(d
[i
], f
[i
], tolerance
);
516 #endif // GMX_SIMD_HAVE_FLOAT && GMX_SIMD_HAVE_DOUBLE