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