Fix bugs in double Simd4 implementations
[gromacs.git] / src / gromacs / simd / impl_x86_mic / impl_x86_mic_simd4_double.h
blob29320adf412280aafbac9a5769caf139e2d27a48
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2017, 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.
36 #ifndef GMX_SIMD_IMPL_X86_MIC_SIMD4_DOUBLE_H
37 #define GMX_SIMD_IMPL_X86_MIC_SIMD4_DOUBLE_H
39 #include "config.h"
41 #include <cassert>
43 #include <immintrin.h>
45 #include "gromacs/utility/basedefinitions.h"
47 #include "impl_x86_mic_simd_double.h"
49 namespace gmx
52 class Simd4Double
54 public:
55 Simd4Double() {}
57 Simd4Double(double d) : simdInternal_(_mm512_set1_pd(d)) {}
59 // Internal utility constructor to simplify return statements
60 Simd4Double(__m512d simd) : simdInternal_(simd) {}
62 __m512d simdInternal_;
65 class Simd4DBool
67 public:
68 Simd4DBool() {}
70 // Internal utility constructor to simplify return statements
71 Simd4DBool(__mmask16 simd) : simdInternal_(simd) {}
73 __mmask16 simdInternal_;
76 static inline Simd4Double gmx_simdcall
77 load4(const double *m)
79 assert(size_t(m) % 32 == 0);
80 return {
81 _mm512_mask_extload_pd(_mm512_undefined_pd(), _mm512_int2mask(0xF), m, _MM_UPCONV_PD_NONE, _MM_BROADCAST_4X8, _MM_HINT_NONE)
85 static inline void gmx_simdcall
86 store4(double *m, Simd4Double a)
88 assert(size_t(m) % 32 == 0);
89 _mm512_mask_packstorelo_pd(m, _mm512_int2mask(0xF), a.simdInternal_);
92 static inline Simd4Double gmx_simdcall
93 load4U(const double *m)
95 return {
96 _mm512_mask_loadunpackhi_pd(_mm512_mask_loadunpacklo_pd(_mm512_undefined_pd(), _mm512_int2mask(0xF), m), _mm512_int2mask(0xF), m+8)
100 static inline void gmx_simdcall
101 store4U(double *m, Simd4Double a)
103 _mm512_mask_packstorelo_pd(m, _mm512_int2mask(0xF), a.simdInternal_);
104 _mm512_mask_packstorehi_pd(m+8, _mm512_int2mask(0xF), a.simdInternal_);
107 static inline Simd4Double gmx_simdcall
108 simd4SetZeroD()
110 return {
111 _mm512_setzero_pd()
115 static inline Simd4Double gmx_simdcall
116 operator&(Simd4Double a, Simd4Double b)
118 return {
119 _mm512_castsi512_pd(_mm512_mask_and_epi32(_mm512_undefined_epi32(), _mm512_int2mask(0x00FF), _mm512_castpd_si512(a.simdInternal_),
120 _mm512_castpd_si512(b.simdInternal_)))
124 static inline Simd4Double gmx_simdcall
125 andNot(Simd4Double a, Simd4Double b)
127 return {
128 _mm512_castsi512_pd(_mm512_mask_andnot_epi32(_mm512_undefined_epi32(), _mm512_int2mask(0x00FF), _mm512_castpd_si512(a.simdInternal_),
129 _mm512_castpd_si512(b.simdInternal_)))
133 static inline Simd4Double gmx_simdcall
134 operator|(Simd4Double a, Simd4Double b)
136 return {
137 _mm512_castsi512_pd(_mm512_mask_or_epi32(_mm512_undefined_epi32(), _mm512_int2mask(0x00FF), _mm512_castpd_si512(a.simdInternal_),
138 _mm512_castpd_si512(b.simdInternal_)))
142 static inline Simd4Double gmx_simdcall
143 operator^(Simd4Double a, Simd4Double b)
145 return {
146 _mm512_castsi512_pd(_mm512_mask_xor_epi32(_mm512_undefined_epi32(), _mm512_int2mask(0x00FF), _mm512_castpd_si512(a.simdInternal_),
147 _mm512_castpd_si512(b.simdInternal_)))
151 static inline Simd4Double gmx_simdcall
152 operator+(Simd4Double a, Simd4Double b)
154 return {
155 _mm512_mask_add_pd(_mm512_undefined_pd(), _mm512_int2mask(0xF), a.simdInternal_, b.simdInternal_)
159 static inline Simd4Double gmx_simdcall
160 operator-(Simd4Double a, Simd4Double b)
162 return {
163 _mm512_mask_sub_pd(_mm512_undefined_pd(), _mm512_int2mask(0xF), a.simdInternal_, b.simdInternal_)
167 static inline Simd4Double gmx_simdcall
168 operator-(Simd4Double x)
170 return {
171 _mm512_mask_addn_pd(_mm512_undefined_pd(), _mm512_int2mask(0xF), x.simdInternal_, _mm512_setzero_pd())
175 static inline Simd4Double gmx_simdcall
176 operator*(Simd4Double a, Simd4Double b)
178 return {
179 _mm512_mask_mul_pd(_mm512_undefined_pd(), _mm512_int2mask(0xF), a.simdInternal_, b.simdInternal_)
183 static inline Simd4Double gmx_simdcall
184 fma(Simd4Double a, Simd4Double b, Simd4Double c)
186 return {
187 _mm512_mask_fmadd_pd(a.simdInternal_, _mm512_int2mask(0xF), b.simdInternal_, c.simdInternal_)
191 static inline Simd4Double gmx_simdcall
192 fms(Simd4Double a, Simd4Double b, Simd4Double c)
194 return {
195 _mm512_mask_fmsub_pd(a.simdInternal_, _mm512_int2mask(0xF), b.simdInternal_, c.simdInternal_)
199 static inline Simd4Double gmx_simdcall
200 fnma(Simd4Double a, Simd4Double b, Simd4Double c)
202 return {
203 _mm512_mask_fnmadd_pd(a.simdInternal_, _mm512_int2mask(0xF), b.simdInternal_, c.simdInternal_)
207 static inline Simd4Double gmx_simdcall
208 fnms(Simd4Double a, Simd4Double b, Simd4Double c)
210 return {
211 _mm512_mask_fnmsub_pd(a.simdInternal_, _mm512_int2mask(0xF), b.simdInternal_, c.simdInternal_)
215 static inline Simd4Double gmx_simdcall
216 rsqrt(Simd4Double x)
218 return {
219 _mm512_mask_cvtpslo_pd(_mm512_undefined_pd(),
220 _mm512_int2mask(0xF),
221 _mm512_mask_rsqrt23_ps(_mm512_undefined_ps(),
222 _mm512_int2mask(0xF),
223 _mm512_mask_cvtpd_pslo(_mm512_undefined_ps(),
224 _mm512_int2mask(0xF), x.simdInternal_)))
228 static inline Simd4Double gmx_simdcall
229 abs(Simd4Double x)
231 return {
232 _mm512_castsi512_pd(_mm512_mask_andnot_epi32(_mm512_undefined_epi32(), _mm512_int2mask(0x00FF),
233 _mm512_castpd_si512(_mm512_set1_pd(GMX_DOUBLE_NEGZERO)),
234 _mm512_castpd_si512(x.simdInternal_)))
239 static inline Simd4Double gmx_simdcall
240 max(Simd4Double a, Simd4Double b)
242 return {
243 _mm512_mask_gmax_pd(_mm512_undefined_pd(), _mm512_int2mask(0xF), a.simdInternal_, b.simdInternal_)
247 static inline Simd4Double gmx_simdcall
248 min(Simd4Double a, Simd4Double b)
250 return {
251 _mm512_mask_gmin_pd(_mm512_undefined_pd(), _mm512_int2mask(0xF), a.simdInternal_, b.simdInternal_)
255 static inline Simd4Double gmx_simdcall
256 round(Simd4Double x)
258 return {
259 _mm512_mask_roundfxpnt_adjust_pd(_mm512_undefined_pd(), _mm512_int2mask(0xF), x.simdInternal_, _MM_FROUND_TO_NEAREST_INT, _MM_EXPADJ_NONE)
263 static inline Simd4Double gmx_simdcall
264 trunc(Simd4Double x)
266 return {
267 _mm512_mask_roundfxpnt_adjust_pd(_mm512_undefined_pd(), _mm512_int2mask(0xF), x.simdInternal_, _MM_FROUND_TO_ZERO, _MM_EXPADJ_NONE)
271 static inline double gmx_simdcall
272 dotProduct(Simd4Double a, Simd4Double b)
274 return _mm512_mask_reduce_add_pd(_mm512_int2mask(7),
275 _mm512_mask_mul_pd(_mm512_undefined_pd(), _mm512_int2mask(7),
276 a.simdInternal_, b.simdInternal_));
279 static inline void gmx_simdcall
280 transpose(Simd4Double * v0, Simd4Double * v1,
281 Simd4Double * v2, Simd4Double * v3)
283 __m512i t0 = _mm512_mask_permute4f128_epi32(_mm512_castpd_si512(v0->simdInternal_), 0xFF00,
284 _mm512_castpd_si512(v1->simdInternal_), _MM_PERM_BABA);
285 __m512i t1 = _mm512_mask_permute4f128_epi32(_mm512_castpd_si512(v2->simdInternal_), 0xFF00,
286 _mm512_castpd_si512(v3->simdInternal_), _MM_PERM_BABA);
288 t0 = _mm512_permutevar_epi32(_mm512_set_epi32(15, 14, 7, 6, 13, 12, 5, 4, 11, 10, 3, 2, 9, 8, 1, 0), t0);
289 t1 = _mm512_permutevar_epi32(_mm512_set_epi32(15, 14, 7, 6, 13, 12, 5, 4, 11, 10, 3, 2, 9, 8, 1, 0), t1);
291 v0->simdInternal_ = _mm512_mask_swizzle_pd(_mm512_castsi512_pd(t0), _mm512_int2mask(0xCC),
292 _mm512_castsi512_pd(t1), _MM_SWIZ_REG_BADC);
293 v1->simdInternal_ = _mm512_mask_swizzle_pd(_mm512_castsi512_pd(t1), _mm512_int2mask(0x33),
294 _mm512_castsi512_pd(t0), _MM_SWIZ_REG_BADC);
296 v2->simdInternal_ = _mm512_castps_pd(_mm512_permute4f128_ps(_mm512_castpd_ps(v0->simdInternal_), _MM_PERM_DCDC));
297 v3->simdInternal_ = _mm512_castps_pd(_mm512_permute4f128_ps(_mm512_castpd_ps(v1->simdInternal_), _MM_PERM_DCDC));
300 // Picky, picky, picky:
301 // icc-16 complains about "Illegal value of immediate argument to intrinsic"
302 // unless we use
303 // 1) Ordered-quiet for ==
304 // 2) Unordered-quiet for !=
305 // 3) Ordered-signaling for < and <=
307 static inline Simd4DBool gmx_simdcall
308 operator==(Simd4Double a, Simd4Double b)
310 return {
311 _mm512_mask_cmp_pd_mask(_mm512_int2mask(0xF), a.simdInternal_, b.simdInternal_, _CMP_EQ_OQ)
315 static inline Simd4DBool gmx_simdcall
316 operator!=(Simd4Double a, Simd4Double b)
318 return {
319 _mm512_mask_cmp_pd_mask(_mm512_int2mask(0xF), a.simdInternal_, b.simdInternal_, _CMP_NEQ_UQ)
323 static inline Simd4DBool gmx_simdcall
324 operator<(Simd4Double a, Simd4Double b)
326 return {
327 _mm512_mask_cmp_pd_mask(_mm512_int2mask(0xF), a.simdInternal_, b.simdInternal_, _CMP_LT_OS)
331 static inline Simd4DBool gmx_simdcall
332 operator<=(Simd4Double a, Simd4Double b)
334 return {
335 _mm512_mask_cmp_pd_mask(_mm512_int2mask(0xF), a.simdInternal_, b.simdInternal_, _CMP_LE_OS)
339 static inline Simd4DBool gmx_simdcall
340 operator&&(Simd4DBool a, Simd4DBool b)
342 return {
343 _mm512_kand(a.simdInternal_, b.simdInternal_)
347 static inline Simd4DBool gmx_simdcall
348 operator||(Simd4DBool a, Simd4DBool b)
350 return {
351 _mm512_kor(a.simdInternal_, b.simdInternal_)
355 static inline bool gmx_simdcall
356 anyTrue(Simd4DBool a)
358 return (_mm512_mask2int(a.simdInternal_) & 0xF) != 0;
361 static inline Simd4Double gmx_simdcall
362 selectByMask(Simd4Double a, Simd4DBool m)
364 return {
365 _mm512_mask_mov_pd(_mm512_setzero_pd(), m.simdInternal_, a.simdInternal_)
369 static inline Simd4Double gmx_simdcall
370 selectByNotMask(Simd4Double a, Simd4DBool m)
372 return {
373 _mm512_mask_mov_pd(_mm512_setzero_pd(), _mm512_knot(m.simdInternal_), a.simdInternal_)
377 static inline Simd4Double gmx_simdcall
378 blend(Simd4Double a, Simd4Double b, Simd4DBool sel)
380 return {
381 _mm512_mask_blend_pd(sel.simdInternal_, a.simdInternal_, b.simdInternal_)
385 static inline double gmx_simdcall
386 reduce(Simd4Double a)
388 return _mm512_mask_reduce_add_pd(_mm512_int2mask(0xF), a.simdInternal_);
391 } // namespace gmx
393 #endif // GMX_SIMD_IMPL_X86_MIC_SIMD4_DOUBLE_H