Split lines with many copyright years
[alexxy/gromacs.git] / src / gromacs / simd / impl_x86_avx_512 / impl_x86_avx_512_simd_double.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_DOUBLE_H
38 #define GMX_SIMD_IMPL_X86_AVX_512_SIMD_DOUBLE_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/basedefinitions.h"
49
50 #include "impl_x86_avx_512_general.h"
51 #include "impl_x86_avx_512_simd_float.h"
52
53 namespace gmx
54 {
55
56 class SimdDouble
57 {
58 public:
59     SimdDouble() {}
60
61     SimdDouble(double d) : simdInternal_(_mm512_set1_pd(d)) {}
62
63     // Internal utility constructor to simplify return statements
64     SimdDouble(__m512d simd) : simdInternal_(simd) {}
65
66     __m512d simdInternal_;
67 };
68
69 class SimdDInt32
70 {
71 public:
72     SimdDInt32() {}
73
74     SimdDInt32(std::int32_t i) : simdInternal_(_mm256_set1_epi32(i)) {}
75
76     // Internal utility constructor to simplify return statements
77     SimdDInt32(__m256i simd) : simdInternal_(simd) {}
78
79     __m256i simdInternal_;
80 };
81
82 class SimdDBool
83 {
84 public:
85     SimdDBool() {}
86
87     // Internal utility constructor to simplify return statements
88     SimdDBool(__mmask8 simd) : simdInternal_(simd) {}
89
90     __mmask8 simdInternal_;
91 };
92
93 class SimdDIBool
94 {
95 public:
96     SimdDIBool() {}
97
98     // Internal utility constructor to simplify return statements
99     SimdDIBool(__mmask16 simd) : simdInternal_(simd) {}
100
101     __mmask16 simdInternal_;
102 };
103
104 static inline SimdDouble gmx_simdcall simdLoad(const double* m, SimdDoubleTag = {})
105 {
106     assert(std::size_t(m) % 64 == 0);
107     return { _mm512_load_pd(m) };
108 }
109
110 static inline void gmx_simdcall store(double* m, SimdDouble a)
111 {
112     assert(std::size_t(m) % 64 == 0);
113     _mm512_store_pd(m, a.simdInternal_);
114 }
115
116 static inline SimdDouble gmx_simdcall simdLoadU(const double* m, SimdDoubleTag = {})
117 {
118     return { _mm512_loadu_pd(m) };
119 }
120
121 static inline void gmx_simdcall storeU(double* m, SimdDouble a)
122 {
123     _mm512_storeu_pd(m, a.simdInternal_);
124 }
125
126 static inline SimdDouble gmx_simdcall setZeroD()
127 {
128     return { _mm512_setzero_pd() };
129 }
130
131 static inline SimdDInt32 gmx_simdcall simdLoad(const std::int32_t* m, SimdDInt32Tag)
132 {
133     assert(std::size_t(m) % 32 == 0);
134     return { _mm256_load_si256(reinterpret_cast<const __m256i*>(m)) };
135 }
136
137 static inline void gmx_simdcall store(std::int32_t* m, SimdDInt32 a)
138 {
139     assert(std::size_t(m) % 32 == 0);
140     _mm256_store_si256(reinterpret_cast<__m256i*>(m), a.simdInternal_);
141 }
142
143 static inline SimdDInt32 gmx_simdcall simdLoadU(const std::int32_t* m, SimdDInt32Tag)
144 {
145     return { _mm256_loadu_si256(reinterpret_cast<const __m256i*>(m)) };
146 }
147
148 static inline void gmx_simdcall storeU(std::int32_t* m, SimdDInt32 a)
149 {
150     _mm256_storeu_si256(reinterpret_cast<__m256i*>(m), a.simdInternal_);
151 }
152
153 static inline SimdDInt32 gmx_simdcall setZeroDI()
154 {
155     return { _mm256_setzero_si256() };
156 }
157
158 static inline SimdDouble gmx_simdcall operator&(SimdDouble a, SimdDouble b)
159 {
160     return { _mm512_castsi512_pd(_mm512_and_epi32(_mm512_castpd_si512(a.simdInternal_),
161                                                   _mm512_castpd_si512(b.simdInternal_))) };
162 }
163
164 static inline SimdDouble gmx_simdcall andNot(SimdDouble a, SimdDouble b)
165 {
166     return { _mm512_castsi512_pd(_mm512_andnot_epi32(_mm512_castpd_si512(a.simdInternal_),
167                                                      _mm512_castpd_si512(b.simdInternal_))) };
168 }
169
170 static inline SimdDouble gmx_simdcall operator|(SimdDouble a, SimdDouble b)
171 {
172     return { _mm512_castsi512_pd(_mm512_or_epi32(_mm512_castpd_si512(a.simdInternal_),
173                                                  _mm512_castpd_si512(b.simdInternal_))) };
174 }
175
176 static inline SimdDouble gmx_simdcall operator^(SimdDouble a, SimdDouble b)
177 {
178     return { _mm512_castsi512_pd(_mm512_xor_epi32(_mm512_castpd_si512(a.simdInternal_),
179                                                   _mm512_castpd_si512(b.simdInternal_))) };
180 }
181
182 static inline SimdDouble gmx_simdcall operator+(SimdDouble a, SimdDouble b)
183 {
184     return { _mm512_add_pd(a.simdInternal_, b.simdInternal_) };
185 }
186
187 static inline SimdDouble gmx_simdcall operator-(SimdDouble a, SimdDouble b)
188 {
189     return { _mm512_sub_pd(a.simdInternal_, b.simdInternal_) };
190 }
191
192 static inline SimdDouble gmx_simdcall operator-(SimdDouble x)
193 {
194     return { _mm512_castsi512_pd(_mm512_xor_epi32(_mm512_castpd_si512(x.simdInternal_),
195                                                   _mm512_castpd_si512(_mm512_set1_pd(GMX_DOUBLE_NEGZERO)))) };
196 }
197
198 static inline SimdDouble gmx_simdcall operator*(SimdDouble a, SimdDouble b)
199 {
200     return { _mm512_mul_pd(a.simdInternal_, b.simdInternal_) };
201 }
202
203 static inline SimdDouble gmx_simdcall fma(SimdDouble a, SimdDouble b, SimdDouble c)
204 {
205     return { _mm512_fmadd_pd(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
206 }
207
208 static inline SimdDouble gmx_simdcall fms(SimdDouble a, SimdDouble b, SimdDouble c)
209 {
210     return { _mm512_fmsub_pd(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
211 }
212
213 static inline SimdDouble gmx_simdcall fnma(SimdDouble a, SimdDouble b, SimdDouble c)
214 {
215     return { _mm512_fnmadd_pd(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
216 }
217
218 static inline SimdDouble gmx_simdcall fnms(SimdDouble a, SimdDouble b, SimdDouble c)
219 {
220     return { _mm512_fnmsub_pd(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
221 }
222
223 // Override for AVX-512-KNL
224 #if GMX_SIMD_X86_AVX_512
225 static inline SimdDouble gmx_simdcall rsqrt(SimdDouble x)
226 {
227     return { _mm512_rsqrt14_pd(x.simdInternal_) };
228 }
229
230 static inline SimdDouble gmx_simdcall rcp(SimdDouble x)
231 {
232     return { _mm512_rcp14_pd(x.simdInternal_) };
233 }
234 #endif
235
236 static inline SimdDouble gmx_simdcall maskAdd(SimdDouble a, SimdDouble b, SimdDBool m)
237 {
238     return { _mm512_mask_add_pd(a.simdInternal_, m.simdInternal_, a.simdInternal_, b.simdInternal_) };
239 }
240
241 static inline SimdDouble gmx_simdcall maskzMul(SimdDouble a, SimdDouble b, SimdDBool m)
242 {
243     return { _mm512_maskz_mul_pd(m.simdInternal_, a.simdInternal_, b.simdInternal_) };
244 }
245
246 static inline SimdDouble gmx_simdcall maskzFma(SimdDouble a, SimdDouble b, SimdDouble c, SimdDBool m)
247 {
248     return { _mm512_maskz_fmadd_pd(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 SimdDouble gmx_simdcall maskzRsqrt(SimdDouble x, SimdDBool m)
254 {
255     return { _mm512_maskz_rsqrt14_pd(m.simdInternal_, x.simdInternal_) };
256 }
257
258 static inline SimdDouble gmx_simdcall maskzRcp(SimdDouble x, SimdDBool m)
259 {
260     return { _mm512_maskz_rcp14_pd(m.simdInternal_, x.simdInternal_) };
261 }
262 #endif
263
264 static inline SimdDouble gmx_simdcall abs(SimdDouble x)
265 {
266     return { _mm512_castsi512_pd(_mm512_andnot_epi32(_mm512_castpd_si512(_mm512_set1_pd(GMX_DOUBLE_NEGZERO)),
267                                                      _mm512_castpd_si512(x.simdInternal_))) };
268 }
269
270 static inline SimdDouble gmx_simdcall max(SimdDouble a, SimdDouble b)
271 {
272     return { _mm512_max_pd(a.simdInternal_, b.simdInternal_) };
273 }
274
275 static inline SimdDouble gmx_simdcall min(SimdDouble a, SimdDouble b)
276 {
277     return { _mm512_min_pd(a.simdInternal_, b.simdInternal_) };
278 }
279
280 static inline SimdDouble gmx_simdcall round(SimdDouble x)
281 {
282     return { _mm512_roundscale_pd(x.simdInternal_, 0) };
283 }
284
285 static inline SimdDouble gmx_simdcall trunc(SimdDouble x)
286 {
287 #if defined(__INTEL_COMPILER) || defined(__ECC)
288     return { _mm512_trunc_pd(x.simdInternal_) };
289 #else
290     return { _mm512_cvtepi32_pd(_mm512_cvttpd_epi32(x.simdInternal_)) };
291 #endif
292 }
293
294 static inline SimdDouble frexp(SimdDouble value, SimdDInt32* exponent)
295 {
296     __m512d rExponent = _mm512_getexp_pd(value.simdInternal_);
297     __m256i iExponent = _mm512_cvtpd_epi32(rExponent);
298
299     exponent->simdInternal_ = _mm256_add_epi32(iExponent, _mm256_set1_epi32(1));
300
301     return { _mm512_getmant_pd(value.simdInternal_, _MM_MANT_NORM_p5_1, _MM_MANT_SIGN_src) };
302 }
303
304 template<MathOptimization opt = MathOptimization::Safe>
305 static inline SimdDouble ldexp(SimdDouble value, SimdDInt32 exponent)
306 {
307     const __m256i exponentBias = _mm256_set1_epi32(1023);
308     __m256i       iExponent    = _mm256_add_epi32(exponent.simdInternal_, exponentBias);
309     __m512i       iExponent512;
310
311     if (opt == MathOptimization::Safe)
312     {
313         // Make sure biased argument is not negative
314         iExponent = _mm256_max_epi32(iExponent, _mm256_setzero_si256());
315     }
316
317     iExponent512 =
318             _mm512_permutexvar_epi32(_mm512_set_epi32(7, 7, 6, 6, 5, 5, 4, 4, 3, 3, 2, 2, 1, 1, 0, 0),
319                                      _mm512_castsi256_si512(iExponent));
320     iExponent512 =
321             _mm512_mask_slli_epi32(_mm512_setzero_epi32(), avx512Int2Mask(0xAAAA), iExponent512, 20);
322     return _mm512_mul_pd(_mm512_castsi512_pd(iExponent512), value.simdInternal_);
323 }
324
325 static inline double gmx_simdcall reduce(SimdDouble a)
326 {
327     __m512d x = a.simdInternal_;
328     x         = _mm512_add_pd(x, _mm512_shuffle_f64x2(x, x, 0xEE));
329     x         = _mm512_add_pd(x, _mm512_shuffle_f64x2(x, x, 0x11));
330     x         = _mm512_add_pd(x, _mm512_permute_pd(x, 0x01));
331     return *reinterpret_cast<double*>(&x);
332 }
333
334 static inline SimdDBool gmx_simdcall operator==(SimdDouble a, SimdDouble b)
335 {
336     return { _mm512_cmp_pd_mask(a.simdInternal_, b.simdInternal_, _CMP_EQ_OQ) };
337 }
338
339 static inline SimdDBool gmx_simdcall operator!=(SimdDouble a, SimdDouble b)
340 {
341     return { _mm512_cmp_pd_mask(a.simdInternal_, b.simdInternal_, _CMP_NEQ_OQ) };
342 }
343
344 static inline SimdDBool gmx_simdcall operator<(SimdDouble a, SimdDouble b)
345 {
346     return { _mm512_cmp_pd_mask(a.simdInternal_, b.simdInternal_, _CMP_LT_OQ) };
347 }
348
349 static inline SimdDBool gmx_simdcall operator<=(SimdDouble a, SimdDouble b)
350 {
351     return { _mm512_cmp_pd_mask(a.simdInternal_, b.simdInternal_, _CMP_LE_OQ) };
352 }
353
354 static inline SimdDBool gmx_simdcall testBits(SimdDouble a)
355 {
356     return { _mm512_test_epi64_mask(_mm512_castpd_si512(a.simdInternal_),
357                                     _mm512_castpd_si512(a.simdInternal_)) };
358 }
359
360 static inline SimdDBool gmx_simdcall operator&&(SimdDBool a, SimdDBool b)
361 {
362     return { static_cast<__mmask8>(_mm512_kand(a.simdInternal_, b.simdInternal_)) };
363 }
364
365 static inline SimdDBool gmx_simdcall operator||(SimdDBool a, SimdDBool b)
366 {
367     return { static_cast<__mmask8>(_mm512_kor(a.simdInternal_, b.simdInternal_)) };
368 }
369
370 static inline bool gmx_simdcall anyTrue(SimdDBool a)
371 {
372     return (avx512Mask2Int(a.simdInternal_) != 0);
373 }
374
375 static inline SimdDouble gmx_simdcall selectByMask(SimdDouble a, SimdDBool m)
376 {
377     return { _mm512_mask_mov_pd(_mm512_setzero_pd(), m.simdInternal_, a.simdInternal_) };
378 }
379
380 static inline SimdDouble gmx_simdcall selectByNotMask(SimdDouble a, SimdDBool m)
381 {
382     return { _mm512_mask_mov_pd(a.simdInternal_, m.simdInternal_, _mm512_setzero_pd()) };
383 }
384
385 static inline SimdDouble gmx_simdcall blend(SimdDouble a, SimdDouble b, SimdDBool sel)
386 {
387     return { _mm512_mask_blend_pd(sel.simdInternal_, a.simdInternal_, b.simdInternal_) };
388 }
389
390 static inline SimdDouble gmx_simdcall copysign(SimdDouble a, SimdDouble b)
391 {
392     return { _mm512_castsi512_pd(_mm512_ternarylogic_epi64(_mm512_castpd_si512(a.simdInternal_),
393                                                            _mm512_castpd_si512(b.simdInternal_),
394                                                            _mm512_set1_epi64(INT64_MIN), 0xD8)) };
395 }
396
397 static inline SimdDInt32 gmx_simdcall operator&(SimdDInt32 a, SimdDInt32 b)
398 {
399     return { _mm256_and_si256(a.simdInternal_, b.simdInternal_) };
400 }
401
402 static inline SimdDInt32 gmx_simdcall andNot(SimdDInt32 a, SimdDInt32 b)
403 {
404     return { _mm256_andnot_si256(a.simdInternal_, b.simdInternal_) };
405 }
406
407 static inline SimdDInt32 gmx_simdcall operator|(SimdDInt32 a, SimdDInt32 b)
408 {
409     return { _mm256_or_si256(a.simdInternal_, b.simdInternal_) };
410 }
411
412 static inline SimdDInt32 gmx_simdcall operator^(SimdDInt32 a, SimdDInt32 b)
413 {
414     return { _mm256_xor_si256(a.simdInternal_, b.simdInternal_) };
415 }
416
417 static inline SimdDInt32 gmx_simdcall operator+(SimdDInt32 a, SimdDInt32 b)
418 {
419     return { _mm256_add_epi32(a.simdInternal_, b.simdInternal_) };
420 }
421
422 static inline SimdDInt32 gmx_simdcall operator-(SimdDInt32 a, SimdDInt32 b)
423 {
424     return { _mm256_sub_epi32(a.simdInternal_, b.simdInternal_) };
425 }
426
427 static inline SimdDInt32 gmx_simdcall operator*(SimdDInt32 a, SimdDInt32 b)
428 {
429     return { _mm256_mullo_epi32(a.simdInternal_, b.simdInternal_) };
430 }
431
432 static inline SimdDIBool gmx_simdcall operator==(SimdDInt32 a, SimdDInt32 b)
433 {
434     return { _mm512_mask_cmp_epi32_mask(avx512Int2Mask(0xFF), _mm512_castsi256_si512(a.simdInternal_),
435                                         _mm512_castsi256_si512(b.simdInternal_), _MM_CMPINT_EQ) };
436 }
437
438 static inline SimdDIBool gmx_simdcall testBits(SimdDInt32 a)
439 {
440     return { _mm512_mask_test_epi32_mask(avx512Int2Mask(0xFF), _mm512_castsi256_si512(a.simdInternal_),
441                                          _mm512_castsi256_si512(a.simdInternal_)) };
442 }
443
444 static inline SimdDIBool gmx_simdcall operator<(SimdDInt32 a, SimdDInt32 b)
445 {
446     return { _mm512_mask_cmp_epi32_mask(avx512Int2Mask(0xFF), _mm512_castsi256_si512(a.simdInternal_),
447                                         _mm512_castsi256_si512(b.simdInternal_), _MM_CMPINT_LT) };
448 }
449
450 static inline SimdDIBool gmx_simdcall operator&&(SimdDIBool a, SimdDIBool b)
451 {
452     return { _mm512_kand(a.simdInternal_, b.simdInternal_) };
453 }
454
455 static inline SimdDIBool gmx_simdcall operator||(SimdDIBool a, SimdDIBool b)
456 {
457     return { _mm512_kor(a.simdInternal_, b.simdInternal_) };
458 }
459
460 static inline bool gmx_simdcall anyTrue(SimdDIBool a)
461 {
462     return (avx512Mask2Int(a.simdInternal_) & 0xFF) != 0;
463 }
464
465 static inline SimdDInt32 gmx_simdcall selectByMask(SimdDInt32 a, SimdDIBool m)
466 {
467     return { _mm512_castsi512_si256(_mm512_mask_mov_epi32(
468             _mm512_setzero_si512(), m.simdInternal_, _mm512_castsi256_si512(a.simdInternal_))) };
469 }
470
471 static inline SimdDInt32 gmx_simdcall selectByNotMask(SimdDInt32 a, SimdDIBool m)
472 {
473     return { _mm512_castsi512_si256(_mm512_mask_mov_epi32(
474             _mm512_castsi256_si512(a.simdInternal_), m.simdInternal_, _mm512_setzero_si512())) };
475 }
476
477 static inline SimdDInt32 gmx_simdcall blend(SimdDInt32 a, SimdDInt32 b, SimdDIBool sel)
478 {
479     return { _mm512_castsi512_si256(
480             _mm512_mask_blend_epi32(sel.simdInternal_, _mm512_castsi256_si512(a.simdInternal_),
481                                     _mm512_castsi256_si512(b.simdInternal_))) };
482 }
483
484 static inline SimdDInt32 gmx_simdcall cvtR2I(SimdDouble a)
485 {
486     return { _mm512_cvtpd_epi32(a.simdInternal_) };
487 }
488
489 static inline SimdDInt32 gmx_simdcall cvttR2I(SimdDouble a)
490 {
491     return { _mm512_cvttpd_epi32(a.simdInternal_) };
492 }
493
494 static inline SimdDouble gmx_simdcall cvtI2R(SimdDInt32 a)
495 {
496     return { _mm512_cvtepi32_pd(a.simdInternal_) };
497 }
498
499 static inline SimdDIBool gmx_simdcall cvtB2IB(SimdDBool a)
500 {
501     return { a.simdInternal_ };
502 }
503
504 static inline SimdDBool gmx_simdcall cvtIB2B(SimdDIBool a)
505 {
506     return { static_cast<__mmask8>(a.simdInternal_) };
507 }
508
509 static inline void gmx_simdcall cvtF2DD(SimdFloat f, SimdDouble* d0, SimdDouble* d1)
510 {
511     d0->simdInternal_ = _mm512_cvtps_pd(_mm512_castps512_ps256(f.simdInternal_));
512     d1->simdInternal_ = _mm512_cvtps_pd(
513             _mm512_castps512_ps256(_mm512_shuffle_f32x4(f.simdInternal_, f.simdInternal_, 0xEE)));
514 }
515
516 static inline SimdFloat gmx_simdcall cvtDD2F(SimdDouble d0, SimdDouble d1)
517 {
518     __m512 f0 = _mm512_castps256_ps512(_mm512_cvtpd_ps(d0.simdInternal_));
519     __m512 f1 = _mm512_castps256_ps512(_mm512_cvtpd_ps(d1.simdInternal_));
520     return { _mm512_shuffle_f32x4(f0, f1, 0x44) };
521 }
522
523 } // namespace gmx
524
525 #endif // GMX_SIMD_IMPL_X86_AVX_512_SIMD_DOUBLE_H