2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017,2018,2019, 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.
36 #ifndef GMX_SIMD_IMPLEMENTATION_IBM_VSX_SIMD_DOUBLE_H
37 #define GMX_SIMD_IMPLEMENTATION_IBM_VSX_SIMD_DOUBLE_H
41 #include "gromacs/math/utilities.h"
42 #include "gromacs/utility/basedefinitions.h"
44 #include "impl_ibm_vsx_definitions.h"
54 // gcc-4.9 does not recognize that we use the parameter
55 SimdDouble(double gmx_unused d) : simdInternal_(vec_splats(d)) {}
57 // Internal utility constructor to simplify return statements
58 SimdDouble(__vector double simd) : simdInternal_(simd) {}
60 __vector double simdInternal_;
68 // gcc-4.9 does not recognize that we use the parameter
69 SimdDInt32(std::int32_t gmx_unused i) : simdInternal_(vec_splats(i)) {}
71 // Internal utility constructor to simplify return statements
72 SimdDInt32(__vector signed int simd) : simdInternal_(simd) {}
74 __vector signed int simdInternal_;
83 simdInternal_(reinterpret_cast<__vector vsxBool long long>(vec_splats(b ? 0xFFFFFFFFFFFFFFFFULL : 0)))
87 // Internal utility constructor to simplify return statements
88 SimdDBool(__vector vsxBool long long simd) : simdInternal_(simd) {}
90 __vector vsxBool long long simdInternal_;
99 simdInternal_(reinterpret_cast<__vector vsxBool int>(vec_splats(b ? 0xFFFFFFFF : 0)))
103 // Internal utility constructor to simplify return statements
104 SimdDIBool(__vector vsxBool int simd) : simdInternal_(simd) {}
106 __vector vsxBool int simdInternal_;
109 // Note that the interfaces we use here have been a mess in xlc;
110 // currently version 13.1.5 is required.
112 static inline SimdDouble gmx_simdcall simdLoad(const double* m, SimdDoubleTag = {})
116 #if defined(__ibmxl__)
120 *reinterpret_cast<const __vector double*>(m)
128 static inline void gmx_simdcall store(double* m, SimdDouble a)
130 #if defined(__ibmxl__)
131 vec_st(a.simdInternal_, 0, m);
134 *reinterpret_cast<__vector double*>(m) = a.simdInternal_;
136 vec_vsx_st(a.simdInternal_, 0, m);
141 static inline SimdDouble gmx_simdcall simdLoadU(const double* m, SimdDoubleTag = {})
145 #if defined(__ibmxl__)
149 *reinterpret_cast<const __vector double*>(m)
157 static inline void gmx_simdcall storeU(double* m, SimdDouble a)
159 #if defined(__ibmxl__)
160 vec_xst(a.simdInternal_, 0, m);
163 *reinterpret_cast<__vector double*>(m) = a.simdInternal_;
165 vec_vsx_st(a.simdInternal_, 0, m);
170 static inline SimdDouble gmx_simdcall setZeroD()
172 return { vec_splats(0.0) };
175 static inline SimdDInt32 gmx_simdcall simdLoad(const std::int32_t* m, SimdDInt32Tag)
177 __vector signed int t0, t1;
178 const __vector unsigned char perm = { 0, 1, 2, 3, 0, 1, 2, 3, 16, 17, 18, 19, 16, 17, 18, 19 };
179 t0 = vec_splats(m[0]);
180 t1 = vec_splats(m[1]);
181 return { vec_perm(t0, t1, perm) };
184 // gcc-4.9 does not understand that arguments to vec_extract() are used
185 static inline void gmx_simdcall store(std::int32_t* m, SimdDInt32 gmx_unused x)
187 m[0] = vec_extract(x.simdInternal_, 0);
188 m[1] = vec_extract(x.simdInternal_, 2);
191 static inline SimdDInt32 gmx_simdcall simdLoadU(const std::int32_t* m, SimdDInt32Tag)
193 return simdLoad(m, SimdDInt32Tag());
196 static inline void gmx_simdcall storeU(std::int32_t* m, SimdDInt32 a)
201 static inline SimdDInt32 gmx_simdcall setZeroDI()
203 return { vec_splats(static_cast<int>(0)) };
206 // gcc-4.9 does not detect that vec_extract() uses its argument
208 static inline std::int32_t gmx_simdcall extract(SimdDInt32 gmx_unused a)
210 return vec_extract(a.simdInternal_, 2 * index);
213 static inline SimdDouble gmx_simdcall operator&(SimdDouble a, SimdDouble b)
215 return { vec_and(a.simdInternal_, b.simdInternal_) };
218 static inline SimdDouble gmx_simdcall andNot(SimdDouble a, SimdDouble b)
220 return { vec_andc(b.simdInternal_, a.simdInternal_) };
223 static inline SimdDouble gmx_simdcall operator|(SimdDouble a, SimdDouble b)
225 return { vec_or(a.simdInternal_, b.simdInternal_) };
228 static inline SimdDouble gmx_simdcall operator^(SimdDouble a, SimdDouble b)
230 return { vec_xor(a.simdInternal_, b.simdInternal_) };
233 static inline SimdDouble gmx_simdcall operator+(SimdDouble a, SimdDouble b)
235 return { vec_add(a.simdInternal_, b.simdInternal_) };
238 static inline SimdDouble gmx_simdcall operator-(SimdDouble a, SimdDouble b)
240 return { vec_sub(a.simdInternal_, b.simdInternal_) };
243 static inline SimdDouble gmx_simdcall operator-(SimdDouble x)
245 return { -x.simdInternal_ };
248 static inline SimdDouble gmx_simdcall operator*(SimdDouble a, SimdDouble b)
250 return { vec_mul(a.simdInternal_, b.simdInternal_) };
253 static inline SimdDouble gmx_simdcall fma(SimdDouble a, SimdDouble b, SimdDouble c)
255 return { vec_madd(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
258 static inline SimdDouble gmx_simdcall fms(SimdDouble a, SimdDouble b, SimdDouble c)
260 return { vec_msub(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
263 static inline SimdDouble gmx_simdcall fnma(SimdDouble a, SimdDouble b, SimdDouble c)
265 return { vec_nmsub(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
268 static inline SimdDouble gmx_simdcall fnms(SimdDouble a, SimdDouble b, SimdDouble c)
270 return { vec_nmadd(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
273 static inline SimdDouble gmx_simdcall rsqrt(SimdDouble x)
275 return { vec_rsqrte(x.simdInternal_) };
278 static inline SimdDouble gmx_simdcall rcp(SimdDouble x)
280 return { vec_re(x.simdInternal_) };
283 static inline SimdDouble gmx_simdcall maskAdd(SimdDouble a, SimdDouble b, SimdDBool m)
285 return { vec_add(a.simdInternal_,
286 vec_and(b.simdInternal_, reinterpret_cast<__vector double>(m.simdInternal_))) };
289 static inline SimdDouble gmx_simdcall maskzMul(SimdDouble a, SimdDouble b, SimdDBool m)
291 SimdDouble prod = a * b;
293 return { vec_and(prod.simdInternal_, reinterpret_cast<__vector double>(m.simdInternal_)) };
296 static inline SimdDouble gmx_simdcall maskzFma(SimdDouble a, SimdDouble b, SimdDouble c, SimdDBool m)
298 SimdDouble prod = fma(a, b, c);
300 return { vec_and(prod.simdInternal_, reinterpret_cast<__vector double>(m.simdInternal_)) };
303 static inline SimdDouble gmx_simdcall maskzRsqrt(SimdDouble x, SimdDBool m)
306 x.simdInternal_ = vec_sel(vec_splats(1.0), x.simdInternal_, m.simdInternal_);
308 return { vec_and(vec_rsqrte(x.simdInternal_), reinterpret_cast<__vector double>(m.simdInternal_)) };
311 static inline SimdDouble gmx_simdcall maskzRcp(SimdDouble x, SimdDBool m)
314 x.simdInternal_ = vec_sel(vec_splats(1.0), x.simdInternal_, m.simdInternal_);
316 return { vec_and(vec_re(x.simdInternal_), reinterpret_cast<__vector double>(m.simdInternal_)) };
319 static inline SimdDouble gmx_simdcall abs(SimdDouble x)
321 return { vec_abs(x.simdInternal_) };
324 static inline SimdDouble gmx_simdcall max(SimdDouble a, SimdDouble b)
326 return { vec_max(a.simdInternal_, b.simdInternal_) };
329 static inline SimdDouble gmx_simdcall min(SimdDouble a, SimdDouble b)
331 return { vec_min(a.simdInternal_, b.simdInternal_) };
334 static inline SimdDouble gmx_simdcall round(SimdDouble x)
336 #if defined(__GNUC__) && !defined(__ibmxl__) && !defined(__xlC__)
337 // gcc up to at least version 4.9 does not have vec_round() in double precision - use inline asm
339 __asm__("xvrdpi %x0,%x1" : "=wd"(res) : "wd"(x.simdInternal_));
342 return { vec_round(x.simdInternal_) };
346 static inline SimdDouble gmx_simdcall trunc(SimdDouble x)
348 return { vec_trunc(x.simdInternal_) };
351 static inline SimdDouble frexp(SimdDouble value, SimdDInt32* exponent)
353 const __vector double exponentMask =
354 reinterpret_cast<__vector double>(vec_splats(0x7FF0000000000000ULL));
355 const __vector signed int exponentBias = vec_splats(1022);
356 const __vector double half = vec_splats(0.5);
357 __vector signed int iExponent;
359 iExponent = reinterpret_cast<__vector signed int>(vec_and(value.simdInternal_, exponentMask));
360 // The data is in the upper half of each double (corresponding to elements 1 and 3).
361 // First shift 52-32=20bits, and then permute to swap element 0 with 1 and element 2 with 3
362 // For big endian they are in opposite order, so then we simply skip the swap.
363 iExponent = vec_sr(iExponent, vec_splats(20U));
364 #ifndef __BIG_ENDIAN__
365 const __vector unsigned char perm = { 4, 5, 6, 7, 0, 1, 2, 3, 12, 13, 14, 15, 8, 9, 10, 11 };
366 iExponent = vec_perm(iExponent, iExponent, perm);
368 iExponent = vec_sub(iExponent, exponentBias);
369 exponent->simdInternal_ = iExponent;
371 return { vec_or(vec_andc(value.simdInternal_, exponentMask), half) };
374 template<MathOptimization opt = MathOptimization::Safe>
375 static inline SimdDouble ldexp(SimdDouble value, SimdDInt32 exponent)
377 const __vector signed int exponentBias = vec_splats(1023);
378 __vector signed int iExponent;
379 #ifdef __BIG_ENDIAN__
380 const __vector unsigned char perm = { 0, 1, 2, 3, 16, 17, 18, 19, 8, 9, 10, 11, 16, 17, 18, 19 };
382 const __vector unsigned char perm = { 16, 17, 18, 19, 0, 1, 2, 3, 16, 17, 18, 19, 8, 9, 10, 11 };
385 iExponent = vec_add(exponent.simdInternal_, exponentBias);
387 if (opt == MathOptimization::Safe)
389 // Make sure biased argument is not negative
390 iExponent = vec_max(iExponent, vec_splat_s32(0));
393 // exponent is now present in pairs of integers; 0011.
394 // Elements 0/2 already correspond to the upper half of each double,
395 // so we only need to shift by another 52-32=20 bits.
396 // The remaining elements are set to zero.
397 iExponent = vec_sl(iExponent, vec_splats(20U));
398 iExponent = vec_perm(iExponent, vec_splats(0), perm);
400 return { vec_mul(value.simdInternal_, reinterpret_cast<__vector double>(iExponent)) };
403 static inline double gmx_simdcall reduce(SimdDouble x)
405 const __vector unsigned char perm = { 8, 9, 10, 11, 12, 13, 14, 15, 0, 1, 2, 3, 4, 5, 6, 7 };
407 /* old xlc version 12 does not understand vec_perm() with double arguments */
408 x.simdInternal_ = vec_add(
409 x.simdInternal_, reinterpret_cast<__vector double>(vec_perm(
410 reinterpret_cast<__vector signed int>(x.simdInternal_),
411 reinterpret_cast<__vector signed int>(x.simdInternal_), perm)));
413 x.simdInternal_ = vec_add(x.simdInternal_, vec_perm(x.simdInternal_, x.simdInternal_, perm));
415 return vec_extract(x.simdInternal_, 0);
418 static inline SimdDBool gmx_simdcall operator==(SimdDouble a, SimdDouble b)
420 return { vec_cmpeq(a.simdInternal_, b.simdInternal_) };
423 static inline SimdDBool gmx_simdcall operator!=(SimdDouble a, SimdDouble b)
425 return { reinterpret_cast<__vector vsxBool long long>(vec_or(
426 reinterpret_cast<__vector signed int>(vec_cmpgt(a.simdInternal_, b.simdInternal_)),
427 reinterpret_cast<__vector signed int>(vec_cmplt(a.simdInternal_, b.simdInternal_)))) };
430 static inline SimdDBool gmx_simdcall operator<(SimdDouble a, SimdDouble b)
432 return { vec_cmplt(a.simdInternal_, b.simdInternal_) };
435 static inline SimdDBool gmx_simdcall operator<=(SimdDouble a, SimdDouble b)
437 return { vec_cmple(a.simdInternal_, b.simdInternal_) };
440 static inline SimdDBool gmx_simdcall testBits(SimdDouble a)
442 #ifdef __POWER8_VECTOR__
443 // Power8 VSX has proper support for operations on long long integers
444 return { vec_cmpgt(reinterpret_cast<__vector unsigned long long>(a.simdInternal_), vec_splats(0ULL)) };
446 // No support for long long operations.
447 // Start with comparing 32-bit subfields bitwise by casting to integers
448 __vector vsxBool int tmp =
449 vec_cmpgt(reinterpret_cast<__vector unsigned int>(a.simdInternal_), vec_splats(0U));
451 // Shuffle low/high 32-bit fields of tmp into tmp2
452 const __vector unsigned char perm = { 4, 5, 6, 7, 0, 1, 2, 3, 12, 13, 14, 15, 8, 9, 10, 11 };
453 __vector vsxBool int tmp2 = vec_perm(tmp, tmp, perm);
455 // Return the or:d parts of tmp & tmp2
456 return { reinterpret_cast<__vector vsxBool long long>(vec_or(tmp, tmp2)) };
460 static inline SimdDBool gmx_simdcall operator&&(SimdDBool a, SimdDBool b)
462 return { reinterpret_cast<__vector vsxBool long long>(
463 vec_and(reinterpret_cast<__vector signed int>(a.simdInternal_),
464 reinterpret_cast<__vector signed int>(b.simdInternal_))) };
467 static inline SimdDBool gmx_simdcall operator||(SimdDBool a, SimdDBool b)
469 return { reinterpret_cast<__vector vsxBool long long>(
470 vec_or(reinterpret_cast<__vector signed int>(a.simdInternal_),
471 reinterpret_cast<__vector signed int>(b.simdInternal_))) };
474 static inline bool gmx_simdcall anyTrue(SimdDBool a)
476 return vec_any_ne(reinterpret_cast<__vector vsxBool int>(a.simdInternal_),
477 reinterpret_cast<__vector vsxBool int>(vec_splats(0)));
480 static inline SimdDouble gmx_simdcall selectByMask(SimdDouble a, SimdDBool m)
482 return { vec_and(a.simdInternal_, reinterpret_cast<__vector double>(m.simdInternal_)) };
485 static inline SimdDouble gmx_simdcall selectByNotMask(SimdDouble a, SimdDBool m)
487 return { vec_andc(a.simdInternal_, reinterpret_cast<__vector double>(m.simdInternal_)) };
490 static inline SimdDouble gmx_simdcall blend(SimdDouble a, SimdDouble b, SimdDBool sel)
492 return { vec_sel(a.simdInternal_, b.simdInternal_, sel.simdInternal_) };
495 static inline SimdDInt32 gmx_simdcall operator&(SimdDInt32 a, SimdDInt32 b)
497 return { vec_and(a.simdInternal_, b.simdInternal_) };
500 static inline SimdDInt32 gmx_simdcall andNot(SimdDInt32 a, SimdDInt32 b)
502 return { vec_andc(b.simdInternal_, a.simdInternal_) };
505 static inline SimdDInt32 gmx_simdcall operator|(SimdDInt32 a, SimdDInt32 b)
507 return { vec_or(a.simdInternal_, b.simdInternal_) };
510 static inline SimdDInt32 gmx_simdcall operator^(SimdDInt32 a, SimdDInt32 b)
512 return { vec_xor(a.simdInternal_, b.simdInternal_) };
515 static inline SimdDInt32 gmx_simdcall operator+(SimdDInt32 a, SimdDInt32 b)
517 return { vec_add(a.simdInternal_, b.simdInternal_) };
520 static inline SimdDInt32 gmx_simdcall operator-(SimdDInt32 a, SimdDInt32 b)
522 return { vec_sub(a.simdInternal_, b.simdInternal_) };
525 static inline SimdDInt32 gmx_simdcall operator*(SimdDInt32 a, SimdDInt32 b)
527 return { a.simdInternal_ * b.simdInternal_ };
530 static inline SimdDIBool gmx_simdcall operator==(SimdDInt32 a, SimdDInt32 b)
532 return { vec_cmpeq(a.simdInternal_, b.simdInternal_) };
535 static inline SimdDIBool gmx_simdcall testBits(SimdDInt32 a)
537 return { vec_cmpgt(reinterpret_cast<__vector unsigned int>(a.simdInternal_), vec_splats(0U)) };
540 static inline SimdDIBool gmx_simdcall operator<(SimdDInt32 a, SimdDInt32 b)
542 return { vec_cmplt(a.simdInternal_, b.simdInternal_) };
545 static inline SimdDIBool gmx_simdcall operator&&(SimdDIBool a, SimdDIBool b)
547 return { vec_and(a.simdInternal_, b.simdInternal_) };
550 static inline SimdDIBool gmx_simdcall operator||(SimdDIBool a, SimdDIBool b)
552 return { vec_or(a.simdInternal_, b.simdInternal_) };
555 static inline bool gmx_simdcall anyTrue(SimdDIBool a)
557 return vec_any_ne(a.simdInternal_, reinterpret_cast<__vector vsxBool int>(vec_splats(0)));
560 static inline SimdDInt32 gmx_simdcall selectByMask(SimdDInt32 a, SimdDIBool m)
562 return { vec_and(a.simdInternal_, reinterpret_cast<__vector signed int>(m.simdInternal_)) };
565 static inline SimdDInt32 gmx_simdcall selectByNotMask(SimdDInt32 a, SimdDIBool m)
567 return { vec_andc(a.simdInternal_, reinterpret_cast<__vector signed int>(m.simdInternal_)) };
570 static inline SimdDInt32 gmx_simdcall blend(SimdDInt32 a, SimdDInt32 b, SimdDIBool sel)
572 return { vec_sel(a.simdInternal_, b.simdInternal_, sel.simdInternal_) };
575 static inline SimdDInt32 gmx_simdcall cvttR2I(SimdDouble a)
577 #if defined(__GNUC__) && !defined(__ibmxl__) && !defined(__xlC__)
578 // gcc up to at least version 6.1 is missing intrinsics for converting double to/from int - use inline asm
579 const __vector unsigned char perm = { 4, 5, 6, 7, 0, 1, 2, 3, 12, 13, 14, 15, 8, 9, 10, 11 };
582 __asm__("xvcvdpsxws %x0,%x1" : "=wa"(ix) : "wd"(a.simdInternal_));
584 return { reinterpret_cast<__vector signed int>(vec_perm(ix, ix, perm)) };
586 return { vec_cts(a.simdInternal_, 0) };
590 static inline SimdDInt32 gmx_simdcall cvtR2I(SimdDouble a)
592 return cvttR2I(round(a));
595 static inline SimdDouble gmx_simdcall cvtI2R(SimdDInt32 a)
597 #if defined(__GNUC__) && !defined(__ibmxl__) && !defined(__xlC__)
598 // gcc up to at least version 4.9 is missing intrinsics for converting double to/from int - use inline asm
600 # ifndef __BIG_ENDIAN__
601 const __vector unsigned char perm = { 4, 5, 6, 7, 0, 1, 2, 3, 12, 13, 14, 15, 8, 9, 10, 11 };
602 a.simdInternal_ = vec_perm(a.simdInternal_, a.simdInternal_, perm);
605 __asm__("xvcvsxwdp %x0,%x1" : "=wd"(x) : "wa"(a.simdInternal_));
609 return { vec_ctd(a.simdInternal_, 0) };
613 static inline SimdDIBool gmx_simdcall cvtB2IB(SimdDBool a)
615 return { reinterpret_cast<__vector vsxBool int>(a.simdInternal_) };
618 static inline SimdDBool gmx_simdcall cvtIB2B(SimdDIBool a)
620 return { reinterpret_cast<__vector vsxBool long long>(a.simdInternal_) };
623 static inline void gmx_simdcall cvtF2DD(SimdFloat f, SimdDouble* d0, SimdDouble* d1)
625 __vector float fA, fB;
626 fA = vec_mergeh(f.simdInternal_, f.simdInternal_); /* 0011 */
627 fB = vec_mergel(f.simdInternal_, f.simdInternal_); /* 2233 */
628 #if defined(__GNUC__) && !defined(__ibmxl__) && !defined(__xlC__)
629 // gcc-4.9 is missing double-to-float/float-to-double conversions.
630 __asm__("xvcvspdp %x0,%x1" : "=wd"(d0->simdInternal_) : "wf"(fA));
631 __asm__("xvcvspdp %x0,%x1" : "=wd"(d1->simdInternal_) : "wf"(fB));
633 d0->simdInternal_ = vec_cvf(fA); /* 01 */
634 d1->simdInternal_ = vec_cvf(fB); /* 23 */
638 static inline SimdFloat gmx_simdcall cvtDD2F(SimdDouble d0, SimdDouble d1)
640 __vector float fA, fB, fC, fD, fE;
641 #if defined(__GNUC__) && !defined(__ibmxl__) && !defined(__xlC__)
642 // gcc-4.9 is missing double-to-float/float-to-double conversions.
643 __asm__("xvcvdpsp %x0,%x1" : "=wf"(fA) : "wd"(d0.simdInternal_));
644 __asm__("xvcvdpsp %x0,%x1" : "=wf"(fB) : "wd"(d1.simdInternal_));
646 fA = vec_cvf(d0.simdInternal_); /* 0x1x */
647 fB = vec_cvf(d1.simdInternal_); /* 2x3x */
649 fC = vec_mergeh(fA, fB); /* 02xx */
650 fD = vec_mergel(fA, fB); /* 13xx */
651 fE = vec_mergeh(fC, fD); /* 0123 */
655 static inline SimdDouble gmx_simdcall copysign(SimdDouble x, SimdDouble y)
657 #if defined(__GNUC__) && !defined(__ibmxl__) && !defined(__xlC__)
659 __asm__("xvcpsgndp %x0,%x1,%x2" : "=wd"(res) : "wd"(y.simdInternal_), "wd"(x.simdInternal_));
662 return { vec_cpsgn(y.simdInternal_, x.simdInternal_) };
668 #endif // GMX_SIMD_IMPLEMENTATION_IBM_VSX_SIMD_DOUBLE_H