Split lines with many copyright years
[alexxy/gromacs.git] / src / gromacs / simd / impl_x86_avx_512 / impl_x86_avx_512_simd_float.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2016,2017,2018 by the GROMACS development team.
5  * Copyright (c) 2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 #ifndef GMX_SIMD_IMPL_X86_AVX_512_SIMD_FLOAT_H
38 #define GMX_SIMD_IMPL_X86_AVX_512_SIMD_FLOAT_H
39
40 #include "config.h"
41
42 #include <cassert>
43 #include <cstdint>
44
45 #include <immintrin.h>
46
47 #include "gromacs/math/utilities.h"
48 #include "gromacs/utility/real.h"
49
50 #include "impl_x86_avx_512_general.h"
51
52 namespace gmx
53 {
54
55 class SimdFloat
56 {
57 public:
58     SimdFloat() {}
59
60     SimdFloat(float f) : simdInternal_(_mm512_set1_ps(f)) {}
61
62     // Internal utility constructor to simplify return statements
63     SimdFloat(__m512 simd) : simdInternal_(simd) {}
64
65     __m512 simdInternal_;
66 };
67
68 class SimdFInt32
69 {
70 public:
71     SimdFInt32() {}
72
73     SimdFInt32(std::int32_t i) : simdInternal_(_mm512_set1_epi32(i)) {}
74
75     // Internal utility constructor to simplify return statements
76     SimdFInt32(__m512i simd) : simdInternal_(simd) {}
77
78     __m512i simdInternal_;
79 };
80
81 class SimdFBool
82 {
83 public:
84     SimdFBool() {}
85
86     // Internal utility constructor to simplify return statements
87     SimdFBool(__mmask16 simd) : simdInternal_(simd) {}
88
89     __mmask16 simdInternal_;
90 };
91
92 class SimdFIBool
93 {
94 public:
95     SimdFIBool() {}
96
97     // Internal utility constructor to simplify return statements
98     SimdFIBool(__mmask16 simd) : simdInternal_(simd) {}
99
100     __mmask16 simdInternal_;
101 };
102
103 static inline SimdFloat gmx_simdcall simdLoad(const float* m, SimdFloatTag = {})
104 {
105     assert(std::size_t(m) % 64 == 0);
106     return { _mm512_load_ps(m) };
107 }
108
109 static inline void gmx_simdcall store(float* m, SimdFloat a)
110 {
111     assert(std::size_t(m) % 64 == 0);
112     _mm512_store_ps(m, a.simdInternal_);
113 }
114
115 static inline SimdFloat gmx_simdcall simdLoadU(const float* m, SimdFloatTag = {})
116 {
117     return { _mm512_loadu_ps(m) };
118 }
119
120 static inline void gmx_simdcall storeU(float* m, SimdFloat a)
121 {
122     _mm512_storeu_ps(m, a.simdInternal_);
123 }
124
125 static inline SimdFloat gmx_simdcall setZeroF()
126 {
127     return { _mm512_setzero_ps() };
128 }
129
130 static inline SimdFInt32 gmx_simdcall simdLoad(const std::int32_t* m, SimdFInt32Tag)
131 {
132     assert(std::size_t(m) % 64 == 0);
133     return { _mm512_load_si512(m) };
134 }
135
136 static inline void gmx_simdcall store(std::int32_t* m, SimdFInt32 a)
137 {
138     assert(std::size_t(m) % 64 == 0);
139     _mm512_store_si512(m, a.simdInternal_);
140 }
141
142 static inline SimdFInt32 gmx_simdcall simdLoadU(const std::int32_t* m, SimdFInt32Tag)
143 {
144     return { _mm512_loadu_si512(m) };
145 }
146
147 static inline void gmx_simdcall storeU(std::int32_t* m, SimdFInt32 a)
148 {
149     _mm512_storeu_si512(m, a.simdInternal_);
150 }
151
152 static inline SimdFInt32 gmx_simdcall setZeroFI()
153 {
154     return { _mm512_setzero_si512() };
155 }
156
157
158 static inline SimdFloat gmx_simdcall operator&(SimdFloat a, SimdFloat b)
159 {
160     return { _mm512_castsi512_ps(_mm512_and_epi32(_mm512_castps_si512(a.simdInternal_),
161                                                   _mm512_castps_si512(b.simdInternal_))) };
162 }
163
164 static inline SimdFloat gmx_simdcall andNot(SimdFloat a, SimdFloat b)
165 {
166     return { _mm512_castsi512_ps(_mm512_andnot_epi32(_mm512_castps_si512(a.simdInternal_),
167                                                      _mm512_castps_si512(b.simdInternal_))) };
168 }
169
170 static inline SimdFloat gmx_simdcall operator|(SimdFloat a, SimdFloat b)
171 {
172     return { _mm512_castsi512_ps(_mm512_or_epi32(_mm512_castps_si512(a.simdInternal_),
173                                                  _mm512_castps_si512(b.simdInternal_))) };
174 }
175
176 static inline SimdFloat gmx_simdcall operator^(SimdFloat a, SimdFloat b)
177 {
178     return { _mm512_castsi512_ps(_mm512_xor_epi32(_mm512_castps_si512(a.simdInternal_),
179                                                   _mm512_castps_si512(b.simdInternal_))) };
180 }
181
182 static inline SimdFloat gmx_simdcall operator+(SimdFloat a, SimdFloat b)
183 {
184     return { _mm512_add_ps(a.simdInternal_, b.simdInternal_) };
185 }
186
187 static inline SimdFloat gmx_simdcall operator-(SimdFloat a, SimdFloat b)
188 {
189     return { _mm512_sub_ps(a.simdInternal_, b.simdInternal_) };
190 }
191
192 static inline SimdFloat gmx_simdcall operator-(SimdFloat x)
193 {
194     return { _mm512_castsi512_ps(_mm512_xor_epi32(_mm512_castps_si512(x.simdInternal_),
195                                                   _mm512_castps_si512(_mm512_set1_ps(GMX_FLOAT_NEGZERO)))) };
196 }
197
198 static inline SimdFloat gmx_simdcall operator*(SimdFloat a, SimdFloat b)
199 {
200     return { _mm512_mul_ps(a.simdInternal_, b.simdInternal_) };
201 }
202
203 static inline SimdFloat gmx_simdcall fma(SimdFloat a, SimdFloat b, SimdFloat c)
204 {
205     return { _mm512_fmadd_ps(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
206 }
207
208 static inline SimdFloat gmx_simdcall fms(SimdFloat a, SimdFloat b, SimdFloat c)
209 {
210     return { _mm512_fmsub_ps(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
211 }
212
213 static inline SimdFloat gmx_simdcall fnma(SimdFloat a, SimdFloat b, SimdFloat c)
214 {
215     return { _mm512_fnmadd_ps(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
216 }
217
218 static inline SimdFloat gmx_simdcall fnms(SimdFloat a, SimdFloat b, SimdFloat c)
219 {
220     return { _mm512_fnmsub_ps(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
221 }
222
223 // Override for AVX-512-KNL
224 #if GMX_SIMD_X86_AVX_512
225 static inline SimdFloat gmx_simdcall rsqrt(SimdFloat x)
226 {
227     return { _mm512_rsqrt14_ps(x.simdInternal_) };
228 }
229
230 static inline SimdFloat gmx_simdcall rcp(SimdFloat x)
231 {
232     return { _mm512_rcp14_ps(x.simdInternal_) };
233 }
234 #endif
235
236 static inline SimdFloat gmx_simdcall maskAdd(SimdFloat a, SimdFloat b, SimdFBool m)
237 {
238     return { _mm512_mask_add_ps(a.simdInternal_, m.simdInternal_, a.simdInternal_, b.simdInternal_) };
239 }
240
241 static inline SimdFloat gmx_simdcall maskzMul(SimdFloat a, SimdFloat b, SimdFBool m)
242 {
243     return { _mm512_maskz_mul_ps(m.simdInternal_, a.simdInternal_, b.simdInternal_) };
244 }
245
246 static inline SimdFloat gmx_simdcall maskzFma(SimdFloat a, SimdFloat b, SimdFloat c, SimdFBool m)
247 {
248     return { _mm512_maskz_fmadd_ps(m.simdInternal_, a.simdInternal_, b.simdInternal_, c.simdInternal_) };
249 }
250
251 // Override for AVX-512-KNL
252 #if GMX_SIMD_X86_AVX_512
253 static inline SimdFloat gmx_simdcall maskzRsqrt(SimdFloat x, SimdFBool m)
254 {
255     return { _mm512_maskz_rsqrt14_ps(m.simdInternal_, x.simdInternal_) };
256 }
257
258 static inline SimdFloat gmx_simdcall maskzRcp(SimdFloat x, SimdFBool m)
259 {
260     return { _mm512_maskz_rcp14_ps(m.simdInternal_, x.simdInternal_) };
261 }
262 #endif
263
264 static inline SimdFloat gmx_simdcall abs(SimdFloat x)
265 {
266     return { _mm512_castsi512_ps(_mm512_andnot_epi32(_mm512_castps_si512(_mm512_set1_ps(GMX_FLOAT_NEGZERO)),
267                                                      _mm512_castps_si512(x.simdInternal_))) };
268 }
269
270 static inline SimdFloat gmx_simdcall max(SimdFloat a, SimdFloat b)
271 {
272     return { _mm512_max_ps(a.simdInternal_, b.simdInternal_) };
273 }
274
275 static inline SimdFloat gmx_simdcall min(SimdFloat a, SimdFloat b)
276 {
277     return { _mm512_min_ps(a.simdInternal_, b.simdInternal_) };
278 }
279
280 static inline SimdFloat gmx_simdcall round(SimdFloat x)
281 {
282     return { _mm512_roundscale_ps(x.simdInternal_, 0) };
283 }
284
285 static inline SimdFloat gmx_simdcall trunc(SimdFloat x)
286 {
287 #if defined(__INTEL_COMPILER) || defined(__ECC)
288     return { _mm512_trunc_ps(x.simdInternal_) };
289 #else
290     return { _mm512_cvtepi32_ps(_mm512_cvttps_epi32(x.simdInternal_)) };
291 #endif
292 }
293
294 static inline SimdFloat gmx_simdcall frexp(SimdFloat value, SimdFInt32* exponent)
295 {
296     __m512  rExponent = _mm512_getexp_ps(value.simdInternal_);
297     __m512i iExponent = _mm512_cvtps_epi32(rExponent);
298
299     exponent->simdInternal_ = _mm512_add_epi32(iExponent, _mm512_set1_epi32(1));
300
301     return { _mm512_getmant_ps(value.simdInternal_, _MM_MANT_NORM_p5_1, _MM_MANT_SIGN_src) };
302 }
303
304 template<MathOptimization opt = MathOptimization::Safe>
305 static inline SimdFloat gmx_simdcall ldexp(SimdFloat value, SimdFInt32 exponent)
306 {
307     const __m512i exponentBias = _mm512_set1_epi32(127);
308     __m512i       iExponent    = _mm512_add_epi32(exponent.simdInternal_, exponentBias);
309
310     if (opt == MathOptimization::Safe)
311     {
312         // Make sure biased argument is not negative
313         iExponent = _mm512_max_epi32(iExponent, _mm512_setzero_epi32());
314     }
315
316     iExponent = _mm512_slli_epi32(iExponent, 23);
317
318     return { _mm512_mul_ps(value.simdInternal_, _mm512_castsi512_ps(iExponent)) };
319 }
320
321 static inline float gmx_simdcall reduce(SimdFloat a)
322 {
323     __m512 x = a.simdInternal_;
324     x        = _mm512_add_ps(x, _mm512_shuffle_f32x4(x, x, 0xEE));
325     x        = _mm512_add_ps(x, _mm512_shuffle_f32x4(x, x, 0x11));
326     x        = _mm512_add_ps(x, _mm512_permute_ps(x, 0xEE));
327     x        = _mm512_add_ps(x, _mm512_permute_ps(x, 0x11));
328     return *reinterpret_cast<float*>(&x);
329 }
330
331 static inline SimdFBool gmx_simdcall operator==(SimdFloat a, SimdFloat b)
332 {
333     return { _mm512_cmp_ps_mask(a.simdInternal_, b.simdInternal_, _CMP_EQ_OQ) };
334 }
335
336 static inline SimdFBool gmx_simdcall operator!=(SimdFloat a, SimdFloat b)
337 {
338     return { _mm512_cmp_ps_mask(a.simdInternal_, b.simdInternal_, _CMP_NEQ_OQ) };
339 }
340
341 static inline SimdFBool gmx_simdcall operator<(SimdFloat a, SimdFloat b)
342 {
343     return { _mm512_cmp_ps_mask(a.simdInternal_, b.simdInternal_, _CMP_LT_OQ) };
344 }
345
346 static inline SimdFBool gmx_simdcall operator<=(SimdFloat a, SimdFloat b)
347 {
348     return { _mm512_cmp_ps_mask(a.simdInternal_, b.simdInternal_, _CMP_LE_OQ) };
349 }
350
351 static inline SimdFBool gmx_simdcall testBits(SimdFloat a)
352 {
353     return { _mm512_test_epi32_mask(_mm512_castps_si512(a.simdInternal_),
354                                     _mm512_castps_si512(a.simdInternal_)) };
355 }
356
357 static inline SimdFBool gmx_simdcall operator&&(SimdFBool a, SimdFBool b)
358 {
359     return { _mm512_kand(a.simdInternal_, b.simdInternal_) };
360 }
361
362 static inline SimdFBool gmx_simdcall operator||(SimdFBool a, SimdFBool b)
363 {
364     return { _mm512_kor(a.simdInternal_, b.simdInternal_) };
365 }
366
367 static inline bool gmx_simdcall anyTrue(SimdFBool a)
368 {
369     return (avx512Mask2Int(a.simdInternal_) != 0);
370 }
371
372 static inline SimdFloat gmx_simdcall selectByMask(SimdFloat a, SimdFBool m)
373 {
374     return { _mm512_mask_mov_ps(_mm512_setzero_ps(), m.simdInternal_, a.simdInternal_) };
375 }
376
377 static inline SimdFloat gmx_simdcall selectByNotMask(SimdFloat a, SimdFBool m)
378 {
379     return { _mm512_mask_mov_ps(a.simdInternal_, m.simdInternal_, _mm512_setzero_ps()) };
380 }
381
382 static inline SimdFloat gmx_simdcall blend(SimdFloat a, SimdFloat b, SimdFBool sel)
383 {
384     return { _mm512_mask_blend_ps(sel.simdInternal_, a.simdInternal_, b.simdInternal_) };
385 }
386
387 static inline SimdFloat gmx_simdcall copysign(SimdFloat a, SimdFloat b)
388 {
389     return { _mm512_castsi512_ps(_mm512_ternarylogic_epi32(_mm512_castps_si512(a.simdInternal_),
390                                                            _mm512_castps_si512(b.simdInternal_),
391                                                            _mm512_set1_epi32(INT32_MIN), 0xD8)) };
392 }
393
394 static inline SimdFInt32 gmx_simdcall operator&(SimdFInt32 a, SimdFInt32 b)
395 {
396     return { _mm512_and_epi32(a.simdInternal_, b.simdInternal_) };
397 }
398
399 static inline SimdFInt32 gmx_simdcall andNot(SimdFInt32 a, SimdFInt32 b)
400 {
401     return { _mm512_andnot_epi32(a.simdInternal_, b.simdInternal_) };
402 }
403
404 static inline SimdFInt32 gmx_simdcall operator|(SimdFInt32 a, SimdFInt32 b)
405 {
406     return { _mm512_or_epi32(a.simdInternal_, b.simdInternal_) };
407 }
408
409 static inline SimdFInt32 gmx_simdcall operator^(SimdFInt32 a, SimdFInt32 b)
410 {
411     return { _mm512_xor_epi32(a.simdInternal_, b.simdInternal_) };
412 }
413
414 static inline SimdFInt32 gmx_simdcall operator+(SimdFInt32 a, SimdFInt32 b)
415 {
416     return { _mm512_add_epi32(a.simdInternal_, b.simdInternal_) };
417 }
418
419 static inline SimdFInt32 gmx_simdcall operator-(SimdFInt32 a, SimdFInt32 b)
420 {
421     return { _mm512_sub_epi32(a.simdInternal_, b.simdInternal_) };
422 }
423
424 static inline SimdFInt32 gmx_simdcall operator*(SimdFInt32 a, SimdFInt32 b)
425 {
426     return { _mm512_mullo_epi32(a.simdInternal_, b.simdInternal_) };
427 }
428
429 static inline SimdFIBool gmx_simdcall operator==(SimdFInt32 a, SimdFInt32 b)
430 {
431     return { _mm512_cmp_epi32_mask(a.simdInternal_, b.simdInternal_, _MM_CMPINT_EQ) };
432 }
433
434 static inline SimdFIBool gmx_simdcall testBits(SimdFInt32 a)
435 {
436     return { _mm512_test_epi32_mask(a.simdInternal_, a.simdInternal_) };
437 }
438
439 static inline SimdFIBool gmx_simdcall operator<(SimdFInt32 a, SimdFInt32 b)
440 {
441     return { _mm512_cmp_epi32_mask(a.simdInternal_, b.simdInternal_, _MM_CMPINT_LT) };
442 }
443
444 static inline SimdFIBool gmx_simdcall operator&&(SimdFIBool a, SimdFIBool b)
445 {
446     return { _mm512_kand(a.simdInternal_, b.simdInternal_) };
447 }
448
449 static inline SimdFIBool gmx_simdcall operator||(SimdFIBool a, SimdFIBool b)
450 {
451     return { _mm512_kor(a.simdInternal_, b.simdInternal_) };
452 }
453
454 static inline bool gmx_simdcall anyTrue(SimdFIBool a)
455 {
456     return (avx512Mask2Int(a.simdInternal_) != 0);
457 }
458
459 static inline SimdFInt32 gmx_simdcall selectByMask(SimdFInt32 a, SimdFIBool m)
460 {
461     return { _mm512_mask_mov_epi32(_mm512_setzero_epi32(), m.simdInternal_, a.simdInternal_) };
462 }
463
464 static inline SimdFInt32 gmx_simdcall selectByNotMask(SimdFInt32 a, SimdFIBool m)
465 {
466     return { _mm512_mask_mov_epi32(a.simdInternal_, m.simdInternal_, _mm512_setzero_epi32()) };
467 }
468
469 static inline SimdFInt32 gmx_simdcall blend(SimdFInt32 a, SimdFInt32 b, SimdFIBool sel)
470 {
471     return { _mm512_mask_blend_epi32(sel.simdInternal_, a.simdInternal_, b.simdInternal_) };
472 }
473
474 static inline SimdFInt32 gmx_simdcall cvtR2I(SimdFloat a)
475 {
476     return { _mm512_cvtps_epi32(a.simdInternal_) };
477 }
478
479 static inline SimdFInt32 gmx_simdcall cvttR2I(SimdFloat a)
480 {
481     return { _mm512_cvttps_epi32(a.simdInternal_) };
482 }
483
484 static inline SimdFloat gmx_simdcall cvtI2R(SimdFInt32 a)
485 {
486     return { _mm512_cvtepi32_ps(a.simdInternal_) };
487 }
488
489 static inline SimdFIBool gmx_simdcall cvtB2IB(SimdFBool a)
490 {
491     return { a.simdInternal_ };
492 }
493
494 static inline SimdFBool gmx_simdcall cvtIB2B(SimdFIBool a)
495 {
496     return { a.simdInternal_ };
497 }
498
499 } // namespace gmx
500
501 #endif // GMX_SIMD_IMPL_X86_AVX_512_SIMD_FLOAT_H