Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / simd / impl_x86_sse2 / impl_x86_sse2_simd_float.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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 #ifndef GMX_SIMD_IMPL_X86_SSE2_SIMD_FLOAT_H
36 #define GMX_SIMD_IMPL_X86_SSE2_SIMD_FLOAT_H
37
38 #include "config.h"
39
40 #include <cassert>
41 #include <cstddef>
42 #include <cstdint>
43
44 #include <emmintrin.h>
45
46 #include "gromacs/math/utilities.h"
47
48 namespace gmx
49 {
50
51 class SimdFloat
52 {
53 public:
54     SimdFloat() {}
55
56     SimdFloat(float f) : simdInternal_(_mm_set1_ps(f)) {}
57
58     // Internal utility constructor to simplify return statements
59     SimdFloat(__m128 simd) : simdInternal_(simd) {}
60
61     __m128 simdInternal_;
62 };
63
64 class SimdFInt32
65 {
66 public:
67     SimdFInt32() {}
68
69     SimdFInt32(std::int32_t i) : simdInternal_(_mm_set1_epi32(i)) {}
70
71     // Internal utility constructor to simplify return statements
72     SimdFInt32(__m128i simd) : simdInternal_(simd) {}
73
74     __m128i simdInternal_;
75 };
76
77 class SimdFBool
78 {
79 public:
80     SimdFBool() {}
81
82     SimdFBool(bool b) : simdInternal_(_mm_castsi128_ps(_mm_set1_epi32(b ? 0xFFFFFFFF : 0))) {}
83
84     // Internal utility constructor to simplify return statements
85     SimdFBool(__m128 simd) : simdInternal_(simd) {}
86
87     __m128 simdInternal_;
88 };
89
90 class SimdFIBool
91 {
92 public:
93     SimdFIBool() {}
94
95     SimdFIBool(bool b) : simdInternal_(_mm_set1_epi32(b ? 0xFFFFFFFF : 0)) {}
96
97     // Internal utility constructor to simplify return statements
98     SimdFIBool(__m128i simd) : simdInternal_(simd) {}
99
100     __m128i simdInternal_;
101 };
102
103 static inline SimdFloat gmx_simdcall simdLoad(const float* m, SimdFloatTag = {})
104 {
105     assert(std::size_t(m) % 16 == 0);
106     return { _mm_load_ps(m) };
107 }
108
109 static inline void gmx_simdcall store(float* m, SimdFloat a)
110 {
111     assert(std::size_t(m) % 16 == 0);
112     _mm_store_ps(m, a.simdInternal_);
113 }
114
115 static inline SimdFloat gmx_simdcall simdLoadU(const float* m, SimdFloatTag = {})
116 {
117     return { _mm_loadu_ps(m) };
118 }
119
120 static inline void gmx_simdcall storeU(float* m, SimdFloat a)
121 {
122     _mm_storeu_ps(m, a.simdInternal_);
123 }
124
125 static inline SimdFloat gmx_simdcall setZeroF()
126 {
127     return { _mm_setzero_ps() };
128 }
129
130 static inline SimdFInt32 gmx_simdcall simdLoad(const std::int32_t* m, SimdFInt32Tag)
131 {
132     assert(std::size_t(m) % 16 == 0);
133     return { _mm_load_si128(reinterpret_cast<const __m128i*>(m)) };
134 }
135
136 static inline void gmx_simdcall store(std::int32_t* m, SimdFInt32 a)
137 {
138     assert(std::size_t(m) % 16 == 0);
139     _mm_store_si128(reinterpret_cast<__m128i*>(m), a.simdInternal_);
140 }
141
142 static inline SimdFInt32 gmx_simdcall simdLoadU(const std::int32_t* m, SimdFInt32Tag)
143 {
144     return { _mm_loadu_si128(reinterpret_cast<const __m128i*>(m)) };
145 }
146
147 static inline void gmx_simdcall storeU(std::int32_t* m, SimdFInt32 a)
148 {
149     _mm_storeu_si128(reinterpret_cast<__m128i*>(m), a.simdInternal_);
150 }
151
152 static inline SimdFInt32 gmx_simdcall setZeroFI()
153 {
154     return { _mm_setzero_si128() };
155 }
156
157
158 // Override for SSE4.1 and higher
159 #if GMX_SIMD_X86_SSE2
160 template<int index>
161 static inline std::int32_t gmx_simdcall extract(SimdFInt32 a)
162 {
163     return _mm_cvtsi128_si32(_mm_srli_si128(a.simdInternal_, 4 * index));
164 }
165 #endif
166
167 static inline SimdFloat gmx_simdcall operator&(SimdFloat a, SimdFloat b)
168 {
169     return { _mm_and_ps(a.simdInternal_, b.simdInternal_) };
170 }
171
172 static inline SimdFloat gmx_simdcall andNot(SimdFloat a, SimdFloat b)
173 {
174     return { _mm_andnot_ps(a.simdInternal_, b.simdInternal_) };
175 }
176
177 static inline SimdFloat gmx_simdcall operator|(SimdFloat a, SimdFloat b)
178 {
179     return { _mm_or_ps(a.simdInternal_, b.simdInternal_) };
180 }
181
182 static inline SimdFloat gmx_simdcall operator^(SimdFloat a, SimdFloat b)
183 {
184     return { _mm_xor_ps(a.simdInternal_, b.simdInternal_) };
185 }
186
187 static inline SimdFloat gmx_simdcall operator+(SimdFloat a, SimdFloat b)
188 {
189     return { _mm_add_ps(a.simdInternal_, b.simdInternal_) };
190 }
191
192 static inline SimdFloat gmx_simdcall operator-(SimdFloat a, SimdFloat b)
193 {
194     return { _mm_sub_ps(a.simdInternal_, b.simdInternal_) };
195 }
196
197 static inline SimdFloat gmx_simdcall operator-(SimdFloat x)
198 {
199     return { _mm_xor_ps(x.simdInternal_, _mm_set1_ps(GMX_FLOAT_NEGZERO)) };
200 }
201
202 static inline SimdFloat gmx_simdcall operator*(SimdFloat a, SimdFloat b)
203 {
204     return { _mm_mul_ps(a.simdInternal_, b.simdInternal_) };
205 }
206
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)
210 {
211     return { _mm_add_ps(_mm_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_) };
212 }
213
214 static inline SimdFloat gmx_simdcall fms(SimdFloat a, SimdFloat b, SimdFloat c)
215 {
216     return { _mm_sub_ps(_mm_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_) };
217 }
218
219 static inline SimdFloat gmx_simdcall fnma(SimdFloat a, SimdFloat b, SimdFloat c)
220 {
221     return { _mm_sub_ps(c.simdInternal_, _mm_mul_ps(a.simdInternal_, b.simdInternal_)) };
222 }
223
224 static inline SimdFloat gmx_simdcall fnms(SimdFloat a, SimdFloat b, SimdFloat c)
225 {
226     return { _mm_sub_ps(_mm_setzero_ps(),
227                         _mm_add_ps(_mm_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_)) };
228 }
229 #endif
230
231 static inline SimdFloat gmx_simdcall rsqrt(SimdFloat x)
232 {
233     return { _mm_rsqrt_ps(x.simdInternal_) };
234 }
235
236 static inline SimdFloat gmx_simdcall rcp(SimdFloat x)
237 {
238     return { _mm_rcp_ps(x.simdInternal_) };
239 }
240
241 static inline SimdFloat gmx_simdcall maskAdd(SimdFloat a, SimdFloat b, SimdFBool m)
242 {
243     return { _mm_add_ps(a.simdInternal_, _mm_and_ps(b.simdInternal_, m.simdInternal_)) };
244 }
245
246 static inline SimdFloat gmx_simdcall maskzMul(SimdFloat a, SimdFloat b, SimdFBool m)
247 {
248     return { _mm_and_ps(_mm_mul_ps(a.simdInternal_, b.simdInternal_), m.simdInternal_) };
249 }
250
251 static inline SimdFloat gmx_simdcall maskzFma(SimdFloat a, SimdFloat b, SimdFloat c, SimdFBool m)
252 {
253     return { _mm_and_ps(_mm_add_ps(_mm_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_),
254                         m.simdInternal_) };
255 }
256
257 // Override for SSE4.1 and higher
258 #if GMX_SIMD_X86_SSE2
259 static inline SimdFloat gmx_simdcall maskzRsqrt(SimdFloat x, SimdFBool m)
260 {
261 #    ifndef NDEBUG
262     x.simdInternal_ = _mm_or_ps(_mm_andnot_ps(m.simdInternal_, _mm_set1_ps(1.0F)),
263                                 _mm_and_ps(m.simdInternal_, x.simdInternal_));
264 #    endif
265     return { _mm_and_ps(_mm_rsqrt_ps(x.simdInternal_), m.simdInternal_) };
266 }
267
268 static inline SimdFloat gmx_simdcall maskzRcp(SimdFloat x, SimdFBool m)
269 {
270 #    ifndef NDEBUG
271     x.simdInternal_ = _mm_or_ps(_mm_andnot_ps(m.simdInternal_, _mm_set1_ps(1.0F)),
272                                 _mm_and_ps(m.simdInternal_, x.simdInternal_));
273 #    endif
274     return { _mm_and_ps(_mm_rcp_ps(x.simdInternal_), m.simdInternal_) };
275 }
276 #endif
277
278 static inline SimdFloat gmx_simdcall abs(SimdFloat x)
279 {
280     return { _mm_andnot_ps(_mm_set1_ps(GMX_FLOAT_NEGZERO), x.simdInternal_) };
281 }
282
283 static inline SimdFloat gmx_simdcall max(SimdFloat a, SimdFloat b)
284 {
285     return { _mm_max_ps(a.simdInternal_, b.simdInternal_) };
286 }
287
288 static inline SimdFloat gmx_simdcall min(SimdFloat a, SimdFloat b)
289 {
290     return { _mm_min_ps(a.simdInternal_, b.simdInternal_) };
291 }
292
293 // Override for SSE4.1 and higher
294 #if GMX_SIMD_X86_SSE2
295 static inline SimdFloat gmx_simdcall round(SimdFloat x)
296 {
297     return { _mm_cvtepi32_ps(_mm_cvtps_epi32(x.simdInternal_)) };
298 }
299
300 static inline SimdFloat gmx_simdcall trunc(SimdFloat x)
301 {
302     return { _mm_cvtepi32_ps(_mm_cvttps_epi32(x.simdInternal_)) };
303 }
304
305 #endif
306
307 static inline SimdFloat gmx_simdcall frexp(SimdFloat value, SimdFInt32* exponent)
308 {
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);
313     __m128i      iExponent;
314
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;
318
319     return { _mm_or_ps(_mm_and_ps(value.simdInternal_, mantissaMask), half) };
320 }
321
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)
326 {
327     const __m128i exponentBias = _mm_set1_epi32(127);
328     __m128i       iExponent;
329
330     iExponent = _mm_add_epi32(exponent.simdInternal_, exponentBias);
331
332     if (opt == MathOptimization::Safe)
333     {
334         // Make sure biased argument is not negative
335         iExponent = _mm_and_si128(iExponent, _mm_cmpgt_epi32(iExponent, _mm_setzero_si128()));
336     }
337
338     iExponent = _mm_slli_epi32(iExponent, 23);
339
340     return { _mm_mul_ps(value.simdInternal_, _mm_castsi128_ps(iExponent)) };
341 }
342 #endif
343
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)
347 {
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);
355 }
356 #endif
357
358 static inline SimdFBool gmx_simdcall operator==(SimdFloat a, SimdFloat b)
359 {
360     return { _mm_cmpeq_ps(a.simdInternal_, b.simdInternal_) };
361 }
362
363 static inline SimdFBool gmx_simdcall operator!=(SimdFloat a, SimdFloat b)
364 {
365     return { _mm_cmpneq_ps(a.simdInternal_, b.simdInternal_) };
366 }
367
368 static inline SimdFBool gmx_simdcall operator<(SimdFloat a, SimdFloat b)
369 {
370     return { _mm_cmplt_ps(a.simdInternal_, b.simdInternal_) };
371 }
372
373 static inline SimdFBool gmx_simdcall operator<=(SimdFloat a, SimdFloat b)
374 {
375     return { _mm_cmple_ps(a.simdInternal_, b.simdInternal_) };
376 }
377
378 static inline SimdFBool gmx_simdcall testBits(SimdFloat a)
379 {
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));
382
383     return { _mm_castsi128_ps(res) };
384 }
385
386 static inline SimdFBool gmx_simdcall operator&&(SimdFBool a, SimdFBool b)
387 {
388     return { _mm_and_ps(a.simdInternal_, b.simdInternal_) };
389 }
390
391 static inline SimdFBool gmx_simdcall operator||(SimdFBool a, SimdFBool b)
392 {
393     return { _mm_or_ps(a.simdInternal_, b.simdInternal_) };
394 }
395
396 static inline bool gmx_simdcall anyTrue(SimdFBool a)
397 {
398     return _mm_movemask_ps(a.simdInternal_) != 0;
399 }
400
401 static inline SimdFloat gmx_simdcall selectByMask(SimdFloat a, SimdFBool mask)
402 {
403     return { _mm_and_ps(a.simdInternal_, mask.simdInternal_) };
404 }
405
406 static inline SimdFloat gmx_simdcall selectByNotMask(SimdFloat a, SimdFBool mask)
407 {
408     return { _mm_andnot_ps(mask.simdInternal_, a.simdInternal_) };
409 }
410
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)
414 {
415     return { _mm_or_ps(_mm_andnot_ps(sel.simdInternal_, a.simdInternal_),
416                        _mm_and_ps(sel.simdInternal_, b.simdInternal_)) };
417 }
418 #endif
419
420 static inline SimdFInt32 gmx_simdcall operator&(SimdFInt32 a, SimdFInt32 b)
421 {
422     return { _mm_and_si128(a.simdInternal_, b.simdInternal_) };
423 }
424
425 static inline SimdFInt32 gmx_simdcall andNot(SimdFInt32 a, SimdFInt32 b)
426 {
427     return { _mm_andnot_si128(a.simdInternal_, b.simdInternal_) };
428 }
429
430 static inline SimdFInt32 gmx_simdcall operator|(SimdFInt32 a, SimdFInt32 b)
431 {
432     return { _mm_or_si128(a.simdInternal_, b.simdInternal_) };
433 }
434
435 static inline SimdFInt32 gmx_simdcall operator^(SimdFInt32 a, SimdFInt32 b)
436 {
437     return { _mm_xor_si128(a.simdInternal_, b.simdInternal_) };
438 }
439
440 static inline SimdFInt32 gmx_simdcall operator+(SimdFInt32 a, SimdFInt32 b)
441 {
442     return { _mm_add_epi32(a.simdInternal_, b.simdInternal_) };
443 }
444
445 static inline SimdFInt32 gmx_simdcall operator-(SimdFInt32 a, SimdFInt32 b)
446 {
447     return { _mm_sub_epi32(a.simdInternal_, b.simdInternal_) };
448 }
449
450 // Override for SSE4.1 and higher
451 #if GMX_SIMD_X86_SSE2
452 static inline SimdFInt32 gmx_simdcall operator*(SimdFInt32 a, SimdFInt32 b)
453 {
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);
458
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]
461
462     return { _mm_unpacklo_epi32(c, c1) };
463 }
464 #endif
465
466 static inline SimdFIBool gmx_simdcall operator==(SimdFInt32 a, SimdFInt32 b)
467 {
468     return { _mm_cmpeq_epi32(a.simdInternal_, b.simdInternal_) };
469 }
470
471 static inline SimdFIBool gmx_simdcall testBits(SimdFInt32 a)
472 {
473     __m128i x   = a.simdInternal_;
474     __m128i res = _mm_andnot_si128(_mm_cmpeq_epi32(x, _mm_setzero_si128()), _mm_cmpeq_epi32(x, x));
475
476     return { res };
477 }
478
479 static inline SimdFIBool gmx_simdcall operator<(SimdFInt32 a, SimdFInt32 b)
480 {
481     return { _mm_cmplt_epi32(a.simdInternal_, b.simdInternal_) };
482 }
483
484 static inline SimdFIBool gmx_simdcall operator&&(SimdFIBool a, SimdFIBool b)
485 {
486     return { _mm_and_si128(a.simdInternal_, b.simdInternal_) };
487 }
488
489 static inline SimdFIBool gmx_simdcall operator||(SimdFIBool a, SimdFIBool b)
490 {
491     return { _mm_or_si128(a.simdInternal_, b.simdInternal_) };
492 }
493
494 static inline bool gmx_simdcall anyTrue(SimdFIBool a)
495 {
496     return _mm_movemask_epi8(a.simdInternal_) != 0;
497 }
498
499 static inline SimdFInt32 gmx_simdcall selectByMask(SimdFInt32 a, SimdFIBool mask)
500 {
501     return { _mm_and_si128(a.simdInternal_, mask.simdInternal_) };
502 }
503
504 static inline SimdFInt32 gmx_simdcall selectByNotMask(SimdFInt32 a, SimdFIBool mask)
505 {
506     return { _mm_andnot_si128(mask.simdInternal_, a.simdInternal_) };
507 }
508
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)
512 {
513     return { _mm_or_si128(_mm_andnot_si128(sel.simdInternal_, a.simdInternal_),
514                           _mm_and_si128(sel.simdInternal_, b.simdInternal_)) };
515 }
516 #endif
517
518 static inline SimdFInt32 gmx_simdcall cvtR2I(SimdFloat a)
519 {
520     return { _mm_cvtps_epi32(a.simdInternal_) };
521 }
522
523 static inline SimdFInt32 gmx_simdcall cvttR2I(SimdFloat a)
524 {
525     return { _mm_cvttps_epi32(a.simdInternal_) };
526 }
527
528 static inline SimdFloat gmx_simdcall cvtI2R(SimdFInt32 a)
529 {
530     return { _mm_cvtepi32_ps(a.simdInternal_) };
531 }
532
533 static inline SimdFIBool gmx_simdcall cvtB2IB(SimdFBool a)
534 {
535     return { _mm_castps_si128(a.simdInternal_) };
536 }
537
538 static inline SimdFBool gmx_simdcall cvtIB2B(SimdFIBool a)
539 {
540     return { _mm_castsi128_ps(a.simdInternal_) };
541 }
542
543 } // namespace gmx
544
545 #endif // GMX_SIMD_IMPL_X86_SSE2_SIMD_FLOAT_H