2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017,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.
35 #ifndef GMX_SIMD_IMPL_X86_SSE2_SIMD_FLOAT_H
36 #define GMX_SIMD_IMPL_X86_SSE2_SIMD_FLOAT_H
44 #include <emmintrin.h>
46 #include "gromacs/math/utilities.h"
56 SimdFloat(float f) : simdInternal_(_mm_set1_ps(f)) {}
58 // Internal utility constructor to simplify return statements
59 SimdFloat(__m128 simd) : simdInternal_(simd) {}
69 SimdFInt32(std::int32_t i) : simdInternal_(_mm_set1_epi32(i)) {}
71 // Internal utility constructor to simplify return statements
72 SimdFInt32(__m128i simd) : simdInternal_(simd) {}
74 __m128i simdInternal_;
82 SimdFBool(bool b) : simdInternal_(_mm_castsi128_ps(_mm_set1_epi32(b ? 0xFFFFFFFF : 0))) {}
84 // Internal utility constructor to simplify return statements
85 SimdFBool(__m128 simd) : simdInternal_(simd) {}
95 SimdFIBool(bool b) : simdInternal_(_mm_set1_epi32(b ? 0xFFFFFFFF : 0)) {}
97 // Internal utility constructor to simplify return statements
98 SimdFIBool(__m128i simd) : simdInternal_(simd) {}
100 __m128i simdInternal_;
103 static inline SimdFloat gmx_simdcall simdLoad(const float* m, SimdFloatTag = {})
105 assert(std::size_t(m) % 16 == 0);
106 return { _mm_load_ps(m) };
109 static inline void gmx_simdcall store(float* m, SimdFloat a)
111 assert(std::size_t(m) % 16 == 0);
112 _mm_store_ps(m, a.simdInternal_);
115 static inline SimdFloat gmx_simdcall simdLoadU(const float* m, SimdFloatTag = {})
117 return { _mm_loadu_ps(m) };
120 static inline void gmx_simdcall storeU(float* m, SimdFloat a)
122 _mm_storeu_ps(m, a.simdInternal_);
125 static inline SimdFloat gmx_simdcall setZeroF()
127 return { _mm_setzero_ps() };
130 static inline SimdFInt32 gmx_simdcall simdLoad(const std::int32_t* m, SimdFInt32Tag)
132 assert(std::size_t(m) % 16 == 0);
133 return { _mm_load_si128(reinterpret_cast<const __m128i*>(m)) };
136 static inline void gmx_simdcall store(std::int32_t* m, SimdFInt32 a)
138 assert(std::size_t(m) % 16 == 0);
139 _mm_store_si128(reinterpret_cast<__m128i*>(m), a.simdInternal_);
142 static inline SimdFInt32 gmx_simdcall simdLoadU(const std::int32_t* m, SimdFInt32Tag)
144 return { _mm_loadu_si128(reinterpret_cast<const __m128i*>(m)) };
147 static inline void gmx_simdcall storeU(std::int32_t* m, SimdFInt32 a)
149 _mm_storeu_si128(reinterpret_cast<__m128i*>(m), a.simdInternal_);
152 static inline SimdFInt32 gmx_simdcall setZeroFI()
154 return { _mm_setzero_si128() };
158 // Override for SSE4.1 and higher
159 #if GMX_SIMD_X86_SSE2
161 static inline std::int32_t gmx_simdcall extract(SimdFInt32 a)
163 return _mm_cvtsi128_si32(_mm_srli_si128(a.simdInternal_, 4 * index));
167 static inline SimdFloat gmx_simdcall operator&(SimdFloat a, SimdFloat b)
169 return { _mm_and_ps(a.simdInternal_, b.simdInternal_) };
172 static inline SimdFloat gmx_simdcall andNot(SimdFloat a, SimdFloat b)
174 return { _mm_andnot_ps(a.simdInternal_, b.simdInternal_) };
177 static inline SimdFloat gmx_simdcall operator|(SimdFloat a, SimdFloat b)
179 return { _mm_or_ps(a.simdInternal_, b.simdInternal_) };
182 static inline SimdFloat gmx_simdcall operator^(SimdFloat a, SimdFloat b)
184 return { _mm_xor_ps(a.simdInternal_, b.simdInternal_) };
187 static inline SimdFloat gmx_simdcall operator+(SimdFloat a, SimdFloat b)
189 return { _mm_add_ps(a.simdInternal_, b.simdInternal_) };
192 static inline SimdFloat gmx_simdcall operator-(SimdFloat a, SimdFloat b)
194 return { _mm_sub_ps(a.simdInternal_, b.simdInternal_) };
197 static inline SimdFloat gmx_simdcall operator-(SimdFloat x)
199 return { _mm_xor_ps(x.simdInternal_, _mm_set1_ps(GMX_FLOAT_NEGZERO)) };
202 static inline SimdFloat gmx_simdcall operator*(SimdFloat a, SimdFloat b)
204 return { _mm_mul_ps(a.simdInternal_, b.simdInternal_) };
207 // Override for AVX-128-FMA and higher
208 #if GMX_SIMD_X86_SSE2 || GMX_SIMD_X86_SSE4_1
209 static inline SimdFloat gmx_simdcall fma(SimdFloat a, SimdFloat b, SimdFloat c)
211 return { _mm_add_ps(_mm_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_) };
214 static inline SimdFloat gmx_simdcall fms(SimdFloat a, SimdFloat b, SimdFloat c)
216 return { _mm_sub_ps(_mm_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_) };
219 static inline SimdFloat gmx_simdcall fnma(SimdFloat a, SimdFloat b, SimdFloat c)
221 return { _mm_sub_ps(c.simdInternal_, _mm_mul_ps(a.simdInternal_, b.simdInternal_)) };
224 static inline SimdFloat gmx_simdcall fnms(SimdFloat a, SimdFloat b, SimdFloat c)
226 return { _mm_sub_ps(_mm_setzero_ps(),
227 _mm_add_ps(_mm_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_)) };
231 static inline SimdFloat gmx_simdcall rsqrt(SimdFloat x)
233 return { _mm_rsqrt_ps(x.simdInternal_) };
236 static inline SimdFloat gmx_simdcall rcp(SimdFloat x)
238 return { _mm_rcp_ps(x.simdInternal_) };
241 static inline SimdFloat gmx_simdcall maskAdd(SimdFloat a, SimdFloat b, SimdFBool m)
243 return { _mm_add_ps(a.simdInternal_, _mm_and_ps(b.simdInternal_, m.simdInternal_)) };
246 static inline SimdFloat gmx_simdcall maskzMul(SimdFloat a, SimdFloat b, SimdFBool m)
248 return { _mm_and_ps(_mm_mul_ps(a.simdInternal_, b.simdInternal_), m.simdInternal_) };
251 static inline SimdFloat gmx_simdcall maskzFma(SimdFloat a, SimdFloat b, SimdFloat c, SimdFBool m)
253 return { _mm_and_ps(_mm_add_ps(_mm_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_),
257 // Override for SSE4.1 and higher
258 #if GMX_SIMD_X86_SSE2
259 static inline SimdFloat gmx_simdcall maskzRsqrt(SimdFloat x, SimdFBool m)
262 x.simdInternal_ = _mm_or_ps(_mm_andnot_ps(m.simdInternal_, _mm_set1_ps(1.0F)),
263 _mm_and_ps(m.simdInternal_, x.simdInternal_));
265 return { _mm_and_ps(_mm_rsqrt_ps(x.simdInternal_), m.simdInternal_) };
268 static inline SimdFloat gmx_simdcall maskzRcp(SimdFloat x, SimdFBool m)
271 x.simdInternal_ = _mm_or_ps(_mm_andnot_ps(m.simdInternal_, _mm_set1_ps(1.0F)),
272 _mm_and_ps(m.simdInternal_, x.simdInternal_));
274 return { _mm_and_ps(_mm_rcp_ps(x.simdInternal_), m.simdInternal_) };
278 static inline SimdFloat gmx_simdcall abs(SimdFloat x)
280 return { _mm_andnot_ps(_mm_set1_ps(GMX_FLOAT_NEGZERO), x.simdInternal_) };
283 static inline SimdFloat gmx_simdcall max(SimdFloat a, SimdFloat b)
285 return { _mm_max_ps(a.simdInternal_, b.simdInternal_) };
288 static inline SimdFloat gmx_simdcall min(SimdFloat a, SimdFloat b)
290 return { _mm_min_ps(a.simdInternal_, b.simdInternal_) };
293 // Override for SSE4.1 and higher
294 #if GMX_SIMD_X86_SSE2
295 static inline SimdFloat gmx_simdcall round(SimdFloat x)
297 return { _mm_cvtepi32_ps(_mm_cvtps_epi32(x.simdInternal_)) };
300 static inline SimdFloat gmx_simdcall trunc(SimdFloat x)
302 return { _mm_cvtepi32_ps(_mm_cvttps_epi32(x.simdInternal_)) };
307 static inline SimdFloat gmx_simdcall frexp(SimdFloat value, SimdFInt32* exponent)
309 const __m128 exponentMask = _mm_castsi128_ps(_mm_set1_epi32(0x7F800000));
310 const __m128 mantissaMask = _mm_castsi128_ps(_mm_set1_epi32(0x807FFFFF));
311 const __m128i exponentBias = _mm_set1_epi32(126); // add 1 to make our definition identical to frexp()
312 const __m128 half = _mm_set1_ps(0.5F);
315 iExponent = _mm_castps_si128(_mm_and_ps(value.simdInternal_, exponentMask));
316 iExponent = _mm_sub_epi32(_mm_srli_epi32(iExponent, 23), exponentBias);
317 exponent->simdInternal_ = iExponent;
319 return { _mm_or_ps(_mm_and_ps(value.simdInternal_, mantissaMask), half) };
322 // Override for SSE4.1
323 #if GMX_SIMD_X86_SSE2
324 template<MathOptimization opt = MathOptimization::Safe>
325 static inline SimdFloat gmx_simdcall ldexp(SimdFloat value, SimdFInt32 exponent)
327 const __m128i exponentBias = _mm_set1_epi32(127);
330 iExponent = _mm_add_epi32(exponent.simdInternal_, exponentBias);
332 if (opt == MathOptimization::Safe)
334 // Make sure biased argument is not negative
335 iExponent = _mm_and_si128(iExponent, _mm_cmpgt_epi32(iExponent, _mm_setzero_si128()));
338 iExponent = _mm_slli_epi32(iExponent, 23);
340 return { _mm_mul_ps(value.simdInternal_, _mm_castsi128_ps(iExponent)) };
344 // Override for AVX-128-FMA and higher
345 #if GMX_SIMD_X86_SSE2 || GMX_SIMD_X86_SSE4_1
346 static inline float gmx_simdcall reduce(SimdFloat a)
348 // Shuffle has latency 1/throughput 1, followed by add with latency 3, t-put 1.
349 // This is likely faster than using _mm_hadd_ps, which has latency 5, t-put 2.
350 a.simdInternal_ = _mm_add_ps(
351 a.simdInternal_, _mm_shuffle_ps(a.simdInternal_, a.simdInternal_, _MM_SHUFFLE(1, 0, 3, 2)));
352 a.simdInternal_ = _mm_add_ss(
353 a.simdInternal_, _mm_shuffle_ps(a.simdInternal_, a.simdInternal_, _MM_SHUFFLE(0, 3, 2, 1)));
354 return *reinterpret_cast<float*>(&a);
358 static inline SimdFBool gmx_simdcall operator==(SimdFloat a, SimdFloat b)
360 return { _mm_cmpeq_ps(a.simdInternal_, b.simdInternal_) };
363 static inline SimdFBool gmx_simdcall operator!=(SimdFloat a, SimdFloat b)
365 return { _mm_cmpneq_ps(a.simdInternal_, b.simdInternal_) };
368 static inline SimdFBool gmx_simdcall operator<(SimdFloat a, SimdFloat b)
370 return { _mm_cmplt_ps(a.simdInternal_, b.simdInternal_) };
373 static inline SimdFBool gmx_simdcall operator<=(SimdFloat a, SimdFloat b)
375 return { _mm_cmple_ps(a.simdInternal_, b.simdInternal_) };
378 static inline SimdFBool gmx_simdcall testBits(SimdFloat a)
380 __m128i ia = _mm_castps_si128(a.simdInternal_);
381 __m128i res = _mm_andnot_si128(_mm_cmpeq_epi32(ia, _mm_setzero_si128()), _mm_cmpeq_epi32(ia, ia));
383 return { _mm_castsi128_ps(res) };
386 static inline SimdFBool gmx_simdcall operator&&(SimdFBool a, SimdFBool b)
388 return { _mm_and_ps(a.simdInternal_, b.simdInternal_) };
391 static inline SimdFBool gmx_simdcall operator||(SimdFBool a, SimdFBool b)
393 return { _mm_or_ps(a.simdInternal_, b.simdInternal_) };
396 static inline bool gmx_simdcall anyTrue(SimdFBool a)
398 return _mm_movemask_ps(a.simdInternal_) != 0;
401 static inline SimdFloat gmx_simdcall selectByMask(SimdFloat a, SimdFBool mask)
403 return { _mm_and_ps(a.simdInternal_, mask.simdInternal_) };
406 static inline SimdFloat gmx_simdcall selectByNotMask(SimdFloat a, SimdFBool mask)
408 return { _mm_andnot_ps(mask.simdInternal_, a.simdInternal_) };
411 // Override for SSE4.1 and higher
412 #if GMX_SIMD_X86_SSE2
413 static inline SimdFloat gmx_simdcall blend(SimdFloat a, SimdFloat b, SimdFBool sel)
415 return { _mm_or_ps(_mm_andnot_ps(sel.simdInternal_, a.simdInternal_),
416 _mm_and_ps(sel.simdInternal_, b.simdInternal_)) };
420 static inline SimdFInt32 gmx_simdcall operator&(SimdFInt32 a, SimdFInt32 b)
422 return { _mm_and_si128(a.simdInternal_, b.simdInternal_) };
425 static inline SimdFInt32 gmx_simdcall andNot(SimdFInt32 a, SimdFInt32 b)
427 return { _mm_andnot_si128(a.simdInternal_, b.simdInternal_) };
430 static inline SimdFInt32 gmx_simdcall operator|(SimdFInt32 a, SimdFInt32 b)
432 return { _mm_or_si128(a.simdInternal_, b.simdInternal_) };
435 static inline SimdFInt32 gmx_simdcall operator^(SimdFInt32 a, SimdFInt32 b)
437 return { _mm_xor_si128(a.simdInternal_, b.simdInternal_) };
440 static inline SimdFInt32 gmx_simdcall operator+(SimdFInt32 a, SimdFInt32 b)
442 return { _mm_add_epi32(a.simdInternal_, b.simdInternal_) };
445 static inline SimdFInt32 gmx_simdcall operator-(SimdFInt32 a, SimdFInt32 b)
447 return { _mm_sub_epi32(a.simdInternal_, b.simdInternal_) };
450 // Override for SSE4.1 and higher
451 #if GMX_SIMD_X86_SSE2
452 static inline SimdFInt32 gmx_simdcall operator*(SimdFInt32 a, SimdFInt32 b)
454 __m128i a1 = _mm_srli_si128(a.simdInternal_, 4); // - a[3] a[2] a[1]
455 __m128i b1 = _mm_srli_si128(b.simdInternal_, 4); // - b[3] b[2] b[1]
456 __m128i c = _mm_mul_epu32(a.simdInternal_, b.simdInternal_);
457 __m128i c1 = _mm_mul_epu32(a1, b1);
459 c = _mm_shuffle_epi32(c, _MM_SHUFFLE(3, 1, 2, 0)); // - - a[2]*b[2] a[0]*b[0]
460 c1 = _mm_shuffle_epi32(c1, _MM_SHUFFLE(3, 1, 2, 0)); // - - a[3]*b[3] a[1]*b[1]
462 return { _mm_unpacklo_epi32(c, c1) };
466 static inline SimdFIBool gmx_simdcall operator==(SimdFInt32 a, SimdFInt32 b)
468 return { _mm_cmpeq_epi32(a.simdInternal_, b.simdInternal_) };
471 static inline SimdFIBool gmx_simdcall testBits(SimdFInt32 a)
473 __m128i x = a.simdInternal_;
474 __m128i res = _mm_andnot_si128(_mm_cmpeq_epi32(x, _mm_setzero_si128()), _mm_cmpeq_epi32(x, x));
479 static inline SimdFIBool gmx_simdcall operator<(SimdFInt32 a, SimdFInt32 b)
481 return { _mm_cmplt_epi32(a.simdInternal_, b.simdInternal_) };
484 static inline SimdFIBool gmx_simdcall operator&&(SimdFIBool a, SimdFIBool b)
486 return { _mm_and_si128(a.simdInternal_, b.simdInternal_) };
489 static inline SimdFIBool gmx_simdcall operator||(SimdFIBool a, SimdFIBool b)
491 return { _mm_or_si128(a.simdInternal_, b.simdInternal_) };
494 static inline bool gmx_simdcall anyTrue(SimdFIBool a)
496 return _mm_movemask_epi8(a.simdInternal_) != 0;
499 static inline SimdFInt32 gmx_simdcall selectByMask(SimdFInt32 a, SimdFIBool mask)
501 return { _mm_and_si128(a.simdInternal_, mask.simdInternal_) };
504 static inline SimdFInt32 gmx_simdcall selectByNotMask(SimdFInt32 a, SimdFIBool mask)
506 return { _mm_andnot_si128(mask.simdInternal_, a.simdInternal_) };
509 // Override for SSE4.1 and higher
510 #if GMX_SIMD_X86_SSE2
511 static inline SimdFInt32 gmx_simdcall blend(SimdFInt32 a, SimdFInt32 b, SimdFIBool sel)
513 return { _mm_or_si128(_mm_andnot_si128(sel.simdInternal_, a.simdInternal_),
514 _mm_and_si128(sel.simdInternal_, b.simdInternal_)) };
518 static inline SimdFInt32 gmx_simdcall cvtR2I(SimdFloat a)
520 return { _mm_cvtps_epi32(a.simdInternal_) };
523 static inline SimdFInt32 gmx_simdcall cvttR2I(SimdFloat a)
525 return { _mm_cvttps_epi32(a.simdInternal_) };
528 static inline SimdFloat gmx_simdcall cvtI2R(SimdFInt32 a)
530 return { _mm_cvtepi32_ps(a.simdInternal_) };
533 static inline SimdFIBool gmx_simdcall cvtB2IB(SimdFBool a)
535 return { _mm_castps_si128(a.simdInternal_) };
538 static inline SimdFBool gmx_simdcall cvtIB2B(SimdFIBool a)
540 return { _mm_castsi128_ps(a.simdInternal_) };
545 #endif // GMX_SIMD_IMPL_X86_SSE2_SIMD_FLOAT_H