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_FLOAT_H
37 #define GMX_SIMD_IMPLEMENTATION_IBM_VSX_SIMD_FLOAT_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 SimdFloat(float gmx_unused f) : simdInternal_(vec_splats(f)) {}
57 // Internal utility constructor to simplify return statements
58 SimdFloat(__vector float simd) : simdInternal_(simd) {}
60 __vector float simdInternal_;
68 // gcc-4.9 does not recognize that we use the parameter
69 SimdFInt32(std::int32_t gmx_unused i) : simdInternal_(vec_splats(i)) {}
71 // Internal utility constructor to simplify return statements
72 SimdFInt32(__vector signed int simd) : simdInternal_(simd) {}
74 __vector signed int simdInternal_;
83 simdInternal_(reinterpret_cast<__vector vsxBool int>(vec_splats(b ? 0xFFFFFFFF : 0)))
87 // Internal utility constructor to simplify return statements
88 SimdFBool(__vector vsxBool int simd) : simdInternal_(simd) {}
90 __vector vsxBool int simdInternal_;
99 simdInternal_(reinterpret_cast<__vector vsxBool int>(vec_splats(b ? 0xFFFFFFFF : 0)))
103 // Internal utility constructor to simplify return statements
104 SimdFIBool(__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 SimdFloat gmx_simdcall simdLoad(const float* m, SimdFloatTag = {})
114 return { *reinterpret_cast<const __vector float*>(m) };
117 static inline void gmx_simdcall store(float* m, SimdFloat a)
119 *reinterpret_cast<__vector float*>(m) = a.simdInternal_;
122 static inline SimdFloat gmx_simdcall simdLoadU(const float* m, SimdFloatTag = {})
127 *reinterpret_cast<const __vector float*>(m)
134 static inline void gmx_simdcall storeU(float* m, SimdFloat a)
137 *reinterpret_cast<__vector float*>(m) = a.simdInternal_;
139 vec_xst(a.simdInternal_, 0, m);
143 static inline SimdFloat gmx_simdcall setZeroF()
145 return { vec_splats(0.0F) };
148 static inline SimdFInt32 gmx_simdcall simdLoad(const std::int32_t* m, SimdFInt32Tag)
150 return { *reinterpret_cast<const __vector int*>(m) };
153 static inline void gmx_simdcall store(std::int32_t* m, SimdFInt32 a)
155 *reinterpret_cast<__vector int*>(m) = a.simdInternal_;
158 static inline SimdFInt32 gmx_simdcall simdLoadU(const std::int32_t* m, SimdFInt32Tag)
163 *reinterpret_cast<const __vector int*>(m)
170 static inline void gmx_simdcall storeU(std::int32_t* m, SimdFInt32 a)
173 *reinterpret_cast<__vector int*>(m) = a.simdInternal_;
175 vec_xst(a.simdInternal_, 0, m);
179 static inline SimdFInt32 gmx_simdcall setZeroFI()
181 return { vec_splats(static_cast<int>(0)) };
184 // gcc-4.9 does not detect that vec_extract() uses its argument
186 static inline std::int32_t gmx_simdcall extract(SimdFInt32 gmx_unused a)
188 return vec_extract(a.simdInternal_, index);
191 static inline SimdFloat gmx_simdcall operator&(SimdFloat a, SimdFloat b)
193 return { vec_and(a.simdInternal_, b.simdInternal_) };
196 static inline SimdFloat gmx_simdcall andNot(SimdFloat a, SimdFloat b)
198 return { vec_andc(b.simdInternal_, a.simdInternal_) };
201 static inline SimdFloat gmx_simdcall operator|(SimdFloat a, SimdFloat b)
203 return { vec_or(a.simdInternal_, b.simdInternal_) };
206 static inline SimdFloat gmx_simdcall operator^(SimdFloat a, SimdFloat b)
208 return { vec_xor(a.simdInternal_, b.simdInternal_) };
211 static inline SimdFloat gmx_simdcall operator+(SimdFloat a, SimdFloat b)
213 return { vec_add(a.simdInternal_, b.simdInternal_) };
216 static inline SimdFloat gmx_simdcall operator-(SimdFloat a, SimdFloat b)
218 return { vec_sub(a.simdInternal_, b.simdInternal_) };
221 static inline SimdFloat gmx_simdcall operator-(SimdFloat x)
223 return { -x.simdInternal_ };
226 static inline SimdFloat gmx_simdcall operator*(SimdFloat a, SimdFloat b)
228 return { vec_mul(a.simdInternal_, b.simdInternal_) };
231 static inline SimdFloat gmx_simdcall fma(SimdFloat a, SimdFloat b, SimdFloat c)
233 return { vec_madd(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
236 static inline SimdFloat gmx_simdcall fms(SimdFloat a, SimdFloat b, SimdFloat c)
238 return { vec_msub(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
241 static inline SimdFloat gmx_simdcall fnma(SimdFloat a, SimdFloat b, SimdFloat c)
243 return { vec_nmsub(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
246 static inline SimdFloat gmx_simdcall fnms(SimdFloat a, SimdFloat b, SimdFloat c)
248 return { vec_nmadd(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
251 static inline SimdFloat gmx_simdcall rsqrt(SimdFloat x)
253 return { vec_rsqrte(x.simdInternal_) };
256 static inline SimdFloat gmx_simdcall rcp(SimdFloat x)
258 return { vec_re(x.simdInternal_) };
261 static inline SimdFloat gmx_simdcall maskAdd(SimdFloat a, SimdFloat b, SimdFBool m)
263 return { vec_add(a.simdInternal_,
264 vec_and(b.simdInternal_, reinterpret_cast<__vector float>(m.simdInternal_))) };
267 static inline SimdFloat gmx_simdcall maskzMul(SimdFloat a, SimdFloat b, SimdFBool m)
269 SimdFloat prod = a * b;
271 return { vec_and(prod.simdInternal_, reinterpret_cast<__vector float>(m.simdInternal_)) };
274 static inline SimdFloat gmx_simdcall maskzFma(SimdFloat a, SimdFloat b, SimdFloat c, SimdFBool m)
276 SimdFloat prod = fma(a, b, c);
278 return { vec_and(prod.simdInternal_, reinterpret_cast<__vector float>(m.simdInternal_)) };
281 static inline SimdFloat gmx_simdcall maskzRsqrt(SimdFloat x, SimdFBool m)
284 x.simdInternal_ = vec_sel(vec_splats(1.0F), x.simdInternal_, m.simdInternal_);
286 return { vec_and(vec_rsqrte(x.simdInternal_), reinterpret_cast<__vector float>(m.simdInternal_)) };
289 static inline SimdFloat gmx_simdcall maskzRcp(SimdFloat x, SimdFBool m)
292 x.simdInternal_ = vec_sel(vec_splats(1.0F), x.simdInternal_, m.simdInternal_);
294 return { vec_and(vec_re(x.simdInternal_), reinterpret_cast<__vector float>(m.simdInternal_)) };
297 static inline SimdFloat gmx_simdcall abs(SimdFloat x)
299 return { vec_abs(x.simdInternal_) };
302 static inline SimdFloat gmx_simdcall max(SimdFloat a, SimdFloat b)
304 return { vec_max(a.simdInternal_, b.simdInternal_) };
307 static inline SimdFloat gmx_simdcall min(SimdFloat a, SimdFloat b)
309 return { vec_min(a.simdInternal_, b.simdInternal_) };
312 static inline SimdFloat gmx_simdcall round(SimdFloat x)
314 return { vec_round(x.simdInternal_) };
317 static inline SimdFloat gmx_simdcall trunc(SimdFloat x)
319 return { vec_trunc(x.simdInternal_) };
322 static inline SimdFloat gmx_simdcall frexp(SimdFloat value, SimdFInt32* exponent)
324 const __vector float exponentMask = reinterpret_cast<__vector float>(vec_splats(0x7F800000U));
325 const __vector signed int exponentBias = vec_splats(126);
326 const __vector float half = vec_splats(0.5F);
327 __vector signed int iExponent;
329 iExponent = reinterpret_cast<__vector signed int>(vec_and(value.simdInternal_, exponentMask));
330 iExponent = vec_sub(vec_sr(iExponent, vec_splats(23U)), exponentBias);
331 exponent->simdInternal_ = iExponent;
333 return { vec_or(vec_andc(value.simdInternal_, exponentMask), half) };
336 template<MathOptimization opt = MathOptimization::Safe>
337 static inline SimdFloat gmx_simdcall ldexp(SimdFloat value, SimdFInt32 exponent)
339 const __vector signed int exponentBias = vec_splats(127);
340 __vector signed int iExponent;
342 iExponent = vec_add(exponent.simdInternal_, exponentBias);
344 if (opt == MathOptimization::Safe)
346 // Make sure biased argument is not negative
347 iExponent = vec_max(iExponent, vec_splat_s32(0));
350 iExponent = vec_sl(iExponent, vec_splats(23U));
352 return { vec_mul(value.simdInternal_, reinterpret_cast<__vector float>(iExponent)) };
355 static inline float gmx_simdcall reduce(SimdFloat x)
357 const __vector unsigned char perm1 = { 8, 9, 10, 11, 12, 13, 14, 15, 0, 1, 2, 3, 4, 5, 6, 7 };
358 const __vector unsigned char perm2 = { 4, 5, 6, 7, 0, 1, 2, 3, 4, 5, 6, 7, 0, 1, 2, 3 };
360 x.simdInternal_ = vec_add(x.simdInternal_, vec_perm(x.simdInternal_, x.simdInternal_, perm1));
361 x.simdInternal_ = vec_add(x.simdInternal_, vec_perm(x.simdInternal_, x.simdInternal_, perm2));
362 return vec_extract(x.simdInternal_, 0);
365 static inline SimdFBool gmx_simdcall operator==(SimdFloat a, SimdFloat b)
367 return { vec_cmpeq(a.simdInternal_, b.simdInternal_) };
370 static inline SimdFBool gmx_simdcall operator!=(SimdFloat a, SimdFloat b)
372 return { vec_or(vec_cmpgt(a.simdInternal_, b.simdInternal_),
373 vec_cmplt(a.simdInternal_, b.simdInternal_)) };
376 static inline SimdFBool gmx_simdcall operator<(SimdFloat a, SimdFloat b)
378 return { vec_cmplt(a.simdInternal_, b.simdInternal_) };
381 static inline SimdFBool gmx_simdcall operator<=(SimdFloat a, SimdFloat b)
383 return { vec_cmple(a.simdInternal_, b.simdInternal_) };
386 static inline SimdFBool gmx_simdcall testBits(SimdFloat a)
388 return { vec_cmpgt(reinterpret_cast<__vector unsigned int>(a.simdInternal_), vec_splats(0U)) };
391 static inline SimdFBool gmx_simdcall operator&&(SimdFBool a, SimdFBool b)
393 return { vec_and(a.simdInternal_, b.simdInternal_) };
396 static inline SimdFBool gmx_simdcall operator||(SimdFBool a, SimdFBool b)
398 return { vec_or(a.simdInternal_, b.simdInternal_) };
401 static inline bool gmx_simdcall anyTrue(SimdFBool a)
403 return vec_any_ne(a.simdInternal_, reinterpret_cast<__vector vsxBool int>(vec_splats(0)));
406 static inline SimdFloat gmx_simdcall selectByMask(SimdFloat a, SimdFBool m)
408 return { vec_and(a.simdInternal_, reinterpret_cast<__vector float>(m.simdInternal_)) };
411 static inline SimdFloat gmx_simdcall selectByNotMask(SimdFloat a, SimdFBool m)
413 return { vec_andc(a.simdInternal_, reinterpret_cast<__vector float>(m.simdInternal_)) };
416 static inline SimdFloat gmx_simdcall blend(SimdFloat a, SimdFloat b, SimdFBool sel)
418 return { vec_sel(a.simdInternal_, b.simdInternal_, sel.simdInternal_) };
421 static inline SimdFInt32 gmx_simdcall operator&(SimdFInt32 a, SimdFInt32 b)
423 return { vec_and(a.simdInternal_, b.simdInternal_) };
426 static inline SimdFInt32 gmx_simdcall andNot(SimdFInt32 a, SimdFInt32 b)
428 return { vec_andc(b.simdInternal_, a.simdInternal_) };
431 static inline SimdFInt32 gmx_simdcall operator|(SimdFInt32 a, SimdFInt32 b)
433 return { vec_or(a.simdInternal_, b.simdInternal_) };
436 static inline SimdFInt32 gmx_simdcall operator^(SimdFInt32 a, SimdFInt32 b)
438 return { vec_xor(a.simdInternal_, b.simdInternal_) };
441 static inline SimdFInt32 gmx_simdcall operator+(SimdFInt32 a, SimdFInt32 b)
443 return { vec_add(a.simdInternal_, b.simdInternal_) };
446 static inline SimdFInt32 gmx_simdcall operator-(SimdFInt32 a, SimdFInt32 b)
448 return { vec_sub(a.simdInternal_, b.simdInternal_) };
451 static inline SimdFInt32 gmx_simdcall operator*(SimdFInt32 a, SimdFInt32 b)
453 return { a.simdInternal_ * b.simdInternal_ };
456 static inline SimdFIBool gmx_simdcall operator==(SimdFInt32 a, SimdFInt32 b)
458 return { vec_cmpeq(a.simdInternal_, b.simdInternal_) };
461 static inline SimdFIBool gmx_simdcall testBits(SimdFInt32 a)
463 return { vec_cmpgt(reinterpret_cast<__vector unsigned int>(a.simdInternal_), vec_splats(0U)) };
466 static inline SimdFIBool gmx_simdcall operator<(SimdFInt32 a, SimdFInt32 b)
468 return { vec_cmplt(a.simdInternal_, b.simdInternal_) };
471 static inline SimdFIBool gmx_simdcall operator&&(SimdFIBool a, SimdFIBool b)
473 return { vec_and(a.simdInternal_, b.simdInternal_) };
476 static inline SimdFIBool gmx_simdcall operator||(SimdFIBool a, SimdFIBool b)
478 return { vec_or(a.simdInternal_, b.simdInternal_) };
481 static inline bool gmx_simdcall anyTrue(SimdFIBool a)
483 return vec_any_ne(a.simdInternal_, reinterpret_cast<__vector vsxBool int>(vec_splats(0)));
486 static inline SimdFInt32 gmx_simdcall selectByMask(SimdFInt32 a, SimdFIBool m)
488 return { vec_and(a.simdInternal_, reinterpret_cast<__vector signed int>(m.simdInternal_)) };
491 static inline SimdFInt32 gmx_simdcall selectByNotMask(SimdFInt32 a, SimdFIBool m)
493 return { vec_andc(a.simdInternal_, reinterpret_cast<__vector signed int>(m.simdInternal_)) };
496 static inline SimdFInt32 gmx_simdcall blend(SimdFInt32 a, SimdFInt32 b, SimdFIBool sel)
498 return { vec_sel(a.simdInternal_, b.simdInternal_, sel.simdInternal_) };
501 static inline SimdFInt32 gmx_simdcall cvtR2I(SimdFloat a)
503 return { vec_cts(vec_round(a.simdInternal_), 0) };
506 static inline SimdFInt32 gmx_simdcall cvttR2I(SimdFloat a)
508 return { vec_cts(a.simdInternal_, 0) };
511 static inline SimdFloat gmx_simdcall cvtI2R(SimdFInt32 a)
513 return { vec_ctf(a.simdInternal_, 0) };
516 static inline SimdFIBool gmx_simdcall cvtB2IB(SimdFBool a)
518 return { a.simdInternal_ };
521 static inline SimdFBool gmx_simdcall cvtIB2B(SimdFIBool a)
523 return { a.simdInternal_ };
526 static inline SimdFloat gmx_simdcall copysign(SimdFloat x, SimdFloat y)
528 #if defined(__GNUC__) && !defined(__ibmxl__) && !defined(__xlC__)
530 __asm__("xvcpsgnsp %x0,%x1,%x2" : "=wf"(res) : "wf"(y.simdInternal_), "wf"(x.simdInternal_));
533 return { vec_cpsgn(y.simdInternal_, x.simdInternal_) };
539 #endif // GMX_SIMD_IMPLEMENTATION_IBM_VSX_SIMD_FLOAT_H