2 * This file is part of the GROMACS molecular simulation package.
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.
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_UTIL_FLOAT_H
37 #define GMX_SIMD_IMPL_X86_MIC_UTIL_FLOAT_H
44 #include <immintrin.h>
46 #include "gromacs/utility/basedefinitions.h"
48 #include "impl_x86_mic_simd_float.h"
53 // On MIC it is better to use scatter operations, so we define the load routines
54 // that use a SIMD offset variable first.
57 static inline void gmx_simdcall gatherLoadBySimdIntTranspose(const float* base,
58 SimdFInt32 simdoffset,
64 assert(std::size_t(base) % 16 == 0);
65 assert(align % 4 == 0);
67 // All instructions might be latency ~4 on MIC, so we use shifts where we
68 // only need a single instruction (since the shift parameter is an immediate),
69 // but multiplication otherwise.
72 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 2);
76 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 3);
80 simdoffset = simdoffset * SimdFInt32(align);
83 v0->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base, sizeof(float));
84 v1->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base + 1, sizeof(float));
85 v2->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base + 2, sizeof(float));
86 v3->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base + 3, sizeof(float));
90 static inline void gmx_simdcall
91 gatherLoadUBySimdIntTranspose(const float* base, SimdFInt32 simdoffset, SimdFloat* v0, SimdFloat* v1)
93 // All instructions might be latency ~4 on MIC, so we use shifts where we
94 // only need a single instruction (since the shift parameter is an immediate),
95 // but multiplication otherwise.
96 // For align == 2 we can merge the constant into the scale parameter,
97 // which can take constants up to 8 in total.
100 v0->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base, align * sizeof(float));
102 _mm512_i32gather_ps(simdoffset.simdInternal_, base + 1, align * sizeof(float));
108 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 2);
112 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 3);
116 simdoffset = simdoffset * SimdFInt32(align);
118 v0->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base, sizeof(float));
119 v1->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base + 1, sizeof(float));
124 static inline void gmx_simdcall
125 gatherLoadBySimdIntTranspose(const float* base, SimdFInt32 simdoffset, SimdFloat* v0, SimdFloat* v1)
127 assert(std::size_t(base) % 8 == 0);
128 assert(align % 2 == 0);
129 gatherLoadUBySimdIntTranspose<align>(base, simdoffset, v0, v1);
133 static inline void gmx_simdcall gatherLoadTranspose(const float* base,
134 const std::int32_t offset[],
140 gatherLoadBySimdIntTranspose<align>(base, simdLoad(offset, SimdFInt32Tag()), v0, v1, v2, v3);
144 static inline void gmx_simdcall
145 gatherLoadTranspose(const float* base, const std::int32_t offset[], SimdFloat* v0, SimdFloat* v1)
147 gatherLoadBySimdIntTranspose<align>(base, simdLoad(offset, SimdFInt32Tag()), v0, v1);
150 static const int c_simdBestPairAlignmentFloat = 2;
153 static inline void gmx_simdcall gatherLoadUTranspose(const float* base,
154 const std::int32_t offset[],
159 SimdFInt32 simdoffset;
161 assert(std::size_t(offset) % 64 == 0);
163 simdoffset = simdLoad(offset, SimdFInt32Tag());
165 // All instructions might be latency ~4 on MIC, so we use shifts where we
166 // only need a single instruction (since the shift parameter is an immediate),
167 // but multiplication otherwise.
170 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 2);
174 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 3);
178 simdoffset = simdoffset * SimdFInt32(align);
181 v0->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base, sizeof(float));
182 v1->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base + 1, sizeof(float));
183 v2->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base + 2, sizeof(float));
188 static inline void gmx_simdcall
189 transposeScatterStoreU(float* base, const std::int32_t offset[], SimdFloat v0, SimdFloat v1, SimdFloat v2)
191 SimdFInt32 simdoffset;
193 assert(std::size_t(offset) % 64 == 0);
195 simdoffset = simdLoad(offset, SimdFInt32Tag());
197 // All instructions might be latency ~4 on MIC, so we use shifts where we
198 // only need a single instruction (since the shift parameter is an immediate),
199 // but multiplication otherwise.
202 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 2);
206 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 3);
210 simdoffset = simdoffset * SimdFInt32(align);
213 _mm512_i32scatter_ps(base, simdoffset.simdInternal_, v0.simdInternal_, sizeof(float));
214 _mm512_i32scatter_ps(base + 1, simdoffset.simdInternal_, v1.simdInternal_, sizeof(float));
215 _mm512_i32scatter_ps(base + 2, simdoffset.simdInternal_, v2.simdInternal_, sizeof(float));
220 static inline void gmx_simdcall
221 transposeScatterIncrU(float* base, const std::int32_t offset[], SimdFloat v0, SimdFloat v1, SimdFloat v2)
223 alignas(GMX_SIMD_ALIGNMENT) float rdata0[GMX_SIMD_FLOAT_WIDTH];
224 alignas(GMX_SIMD_ALIGNMENT) float rdata1[GMX_SIMD_FLOAT_WIDTH];
225 alignas(GMX_SIMD_ALIGNMENT) float rdata2[GMX_SIMD_FLOAT_WIDTH];
231 for (int i = 0; i < GMX_SIMD_FLOAT_WIDTH; i++)
233 base[align * offset[i] + 0] += rdata0[i];
234 base[align * offset[i] + 1] += rdata1[i];
235 base[align * offset[i] + 2] += rdata2[i];
240 static inline void gmx_simdcall
241 transposeScatterDecrU(float* base, const std::int32_t offset[], SimdFloat v0, SimdFloat v1, SimdFloat v2)
243 alignas(GMX_SIMD_ALIGNMENT) float rdata0[GMX_SIMD_FLOAT_WIDTH];
244 alignas(GMX_SIMD_ALIGNMENT) float rdata1[GMX_SIMD_FLOAT_WIDTH];
245 alignas(GMX_SIMD_ALIGNMENT) float rdata2[GMX_SIMD_FLOAT_WIDTH];
251 for (int i = 0; i < GMX_SIMD_FLOAT_WIDTH; i++)
253 base[align * offset[i] + 0] -= rdata0[i];
254 base[align * offset[i] + 1] -= rdata1[i];
255 base[align * offset[i] + 2] -= rdata2[i];
259 static inline void gmx_simdcall expandScalarsToTriplets(SimdFloat scalar,
260 SimdFloat* triplets0,
261 SimdFloat* triplets1,
262 SimdFloat* triplets2)
264 triplets0->simdInternal_ = _mm512_castsi512_ps(
265 _mm512_permutevar_epi32(_mm512_set_epi32(5, 4, 4, 4, 3, 3, 3, 2, 2, 2, 1, 1, 1, 0, 0, 0),
266 _mm512_castps_si512(scalar.simdInternal_)));
267 triplets1->simdInternal_ = _mm512_castsi512_ps(_mm512_permutevar_epi32(
268 _mm512_set_epi32(10, 10, 9, 9, 9, 8, 8, 8, 7, 7, 7, 6, 6, 6, 5, 5),
269 _mm512_castps_si512(scalar.simdInternal_)));
270 triplets2->simdInternal_ = _mm512_castsi512_ps(_mm512_permutevar_epi32(
271 _mm512_set_epi32(15, 15, 15, 14, 14, 14, 13, 13, 13, 12, 12, 12, 11, 11, 11, 10),
272 _mm512_castps_si512(scalar.simdInternal_)));
276 static inline float gmx_simdcall reduceIncr4ReturnSum(float* m, SimdFloat v0, SimdFloat v1, SimdFloat v2, SimdFloat v3)
279 __m512 t0, t1, t2, t3;
281 assert(std::size_t(m) % 16 == 0);
283 t0 = _mm512_add_ps(v0.simdInternal_, _mm512_swizzle_ps(v0.simdInternal_, _MM_SWIZ_REG_BADC));
284 t0 = _mm512_mask_add_ps(t0, _mm512_int2mask(0xCCCC), v2.simdInternal_,
285 _mm512_swizzle_ps(v2.simdInternal_, _MM_SWIZ_REG_BADC));
286 t1 = _mm512_add_ps(v1.simdInternal_, _mm512_swizzle_ps(v1.simdInternal_, _MM_SWIZ_REG_BADC));
287 t1 = _mm512_mask_add_ps(t1, _mm512_int2mask(0xCCCC), v3.simdInternal_,
288 _mm512_swizzle_ps(v3.simdInternal_, _MM_SWIZ_REG_BADC));
289 t2 = _mm512_add_ps(t0, _mm512_swizzle_ps(t0, _MM_SWIZ_REG_CDAB));
290 t2 = _mm512_mask_add_ps(t2, _mm512_int2mask(0xAAAA), t1, _mm512_swizzle_ps(t1, _MM_SWIZ_REG_CDAB));
292 t2 = _mm512_add_ps(t2, _mm512_permute4f128_ps(t2, _MM_PERM_BADC));
293 t2 = _mm512_add_ps(t2, _mm512_permute4f128_ps(t2, _MM_PERM_CDAB));
295 t0 = _mm512_mask_extload_ps(_mm512_undefined_ps(), _mm512_int2mask(0xF), m, _MM_UPCONV_PS_NONE,
296 _MM_BROADCAST_4X16, _MM_HINT_NONE);
297 t0 = _mm512_add_ps(t0, t2);
298 _mm512_mask_packstorelo_ps(m, _mm512_int2mask(0xF), t0);
300 t2 = _mm512_add_ps(t2, _mm512_swizzle_ps(t2, _MM_SWIZ_REG_BADC));
301 t2 = _mm512_add_ps(t2, _mm512_swizzle_ps(t2, _MM_SWIZ_REG_CDAB));
303 _mm512_mask_packstorelo_ps(&f, _mm512_mask2int(0x1), t2);
307 static inline SimdFloat gmx_simdcall loadDualHsimd(const float* m0, const float* m1)
309 assert(std::size_t(m0) % 32 == 0);
310 assert(std::size_t(m1) % 32 == 0);
312 return _mm512_castpd_ps(_mm512_mask_extload_pd(
313 _mm512_extload_pd(reinterpret_cast<const double*>(m0), _MM_UPCONV_PD_NONE,
314 _MM_BROADCAST_4X8, _MM_HINT_NONE),
315 _mm512_int2mask(0xF0), reinterpret_cast<const double*>(m1), _MM_UPCONV_PD_NONE,
316 _MM_BROADCAST_4X8, _MM_HINT_NONE));
319 static inline SimdFloat gmx_simdcall loadDuplicateHsimd(const float* m)
321 assert(std::size_t(m) % 32 == 0);
323 return _mm512_castpd_ps(_mm512_extload_pd(reinterpret_cast<const double*>(m),
324 _MM_UPCONV_PD_NONE, _MM_BROADCAST_4X8, _MM_HINT_NONE));
327 static inline SimdFloat gmx_simdcall loadU1DualHsimd(const float* m)
329 return _mm512_mask_extload_ps(
330 _mm512_extload_ps(m, _MM_UPCONV_PS_NONE, _MM_BROADCAST_1X16, _MM_HINT_NONE),
331 _mm512_int2mask(0xFF00), m + 1, _MM_UPCONV_PS_NONE, _MM_BROADCAST_1X16, _MM_HINT_NONE);
335 static inline void gmx_simdcall storeDualHsimd(float* m0, float* m1, SimdFloat a)
339 assert(std::size_t(m0) % 32 == 0);
340 assert(std::size_t(m1) % 32 == 0);
342 _mm512_mask_packstorelo_ps(m0, _mm512_int2mask(0x00FF), a.simdInternal_);
343 _mm512_mask_packstorelo_ps(m1, _mm512_int2mask(0xFF00), a.simdInternal_);
346 static inline void gmx_simdcall incrDualHsimd(float* m0, float* m1, SimdFloat a)
348 assert(std::size_t(m0) % 32 == 0);
349 assert(std::size_t(m1) % 32 == 0);
354 x = _mm512_castpd_ps(_mm512_extload_pd(reinterpret_cast<const double*>(m0), _MM_UPCONV_PD_NONE,
355 _MM_BROADCAST_4X8, _MM_HINT_NONE));
356 x = _mm512_add_ps(x, a.simdInternal_);
357 _mm512_mask_packstorelo_ps(m0, _mm512_int2mask(0x00FF), x);
360 x = _mm512_castpd_ps(_mm512_extload_pd(reinterpret_cast<const double*>(m1), _MM_UPCONV_PD_NONE,
361 _MM_BROADCAST_4X8, _MM_HINT_NONE));
362 x = _mm512_add_ps(x, a.simdInternal_);
363 _mm512_mask_packstorelo_ps(m1, _mm512_int2mask(0xFF00), x);
366 static inline void gmx_simdcall decrHsimd(float* m, SimdFloat a)
370 assert(std::size_t(m) % 32 == 0);
372 t = _mm512_castpd_ps(_mm512_extload_pd(reinterpret_cast<const double*>(m), _MM_UPCONV_PD_NONE,
373 _MM_BROADCAST_4X8, _MM_HINT_NONE));
374 a = _mm512_add_ps(a.simdInternal_, _mm512_permute4f128_ps(a.simdInternal_, _MM_PERM_BADC));
375 t = _mm512_sub_ps(t, a.simdInternal_);
376 _mm512_mask_packstorelo_ps(m, _mm512_int2mask(0x00FF), t);
381 static inline void gmx_simdcall gatherLoadTransposeHsimd(const float* base0,
383 const std::int32_t offset[],
387 __m512i idx0, idx1, idx;
390 assert(std::size_t(offset) % 32 == 0);
391 assert(std::size_t(base0) % 8 == 0);
392 assert(std::size_t(base1) % 8 == 0);
393 assert(std::size_t(align) % 2 == 0);
395 idx0 = _mm512_loadunpacklo_epi32(_mm512_undefined_epi32(), offset);
397 idx0 = _mm512_mullo_epi32(idx0, _mm512_set1_epi32(align));
398 idx1 = _mm512_add_epi32(idx0, _mm512_set1_epi32(1));
400 idx = _mm512_mask_permute4f128_epi32(idx0, _mm512_int2mask(0xFF00), idx1, _MM_PERM_BABA);
402 tmp1 = _mm512_i32gather_ps(idx, base0, sizeof(float));
403 tmp2 = _mm512_i32gather_ps(idx, base1, sizeof(float));
405 v0->simdInternal_ = _mm512_mask_permute4f128_ps(tmp1, _mm512_int2mask(0xFF00), tmp2, _MM_PERM_BABA);
406 v1->simdInternal_ = _mm512_mask_permute4f128_ps(tmp2, _mm512_int2mask(0x00FF), tmp1, _MM_PERM_DCDC);
409 static inline float gmx_simdcall reduceIncr4ReturnSumHsimd(float* m, SimdFloat v0, SimdFloat v1)
414 assert(std::size_t(m) % 32 == 0);
416 t0 = _mm512_add_ps(v0.simdInternal_, _mm512_swizzle_ps(v0.simdInternal_, _MM_SWIZ_REG_BADC));
417 t0 = _mm512_mask_add_ps(t0, _mm512_int2mask(0xCCCC), v1.simdInternal_,
418 _mm512_swizzle_ps(v1.simdInternal_, _MM_SWIZ_REG_BADC));
419 t0 = _mm512_add_ps(t0, _mm512_swizzle_ps(t0, _MM_SWIZ_REG_CDAB));
420 t0 = _mm512_add_ps(t0, _mm512_castpd_ps(_mm512_swizzle_pd(_mm512_castps_pd(t0), _MM_SWIZ_REG_BADC)));
421 t0 = _mm512_mask_permute4f128_ps(t0, _mm512_int2mask(0xAAAA), t0, _MM_PERM_BADC);
422 t1 = _mm512_mask_extload_ps(_mm512_undefined_ps(), _mm512_int2mask(0xF), m, _MM_UPCONV_PS_NONE,
423 _MM_BROADCAST_4X16, _MM_HINT_NONE);
424 t1 = _mm512_add_ps(t1, t0);
425 _mm512_mask_packstorelo_ps(m, _mm512_int2mask(0xF), t1);
427 t0 = _mm512_add_ps(t0, _mm512_swizzle_ps(t0, _MM_SWIZ_REG_BADC));
428 t0 = _mm512_add_ps(t0, _mm512_swizzle_ps(t0, _MM_SWIZ_REG_CDAB));
430 _mm512_mask_packstorelo_ps(&f, _mm512_mask2int(0x1), t0);
436 #endif // GMX_SIMD_IMPL_X86_MIC_UTIL_FLOAT_H