2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 #ifndef GMX_SIMD_IMPL_X86_MIC_UTIL_FLOAT_H
38 #define GMX_SIMD_IMPL_X86_MIC_UTIL_FLOAT_H
45 #include <immintrin.h>
47 #include "gromacs/utility/basedefinitions.h"
49 #include "impl_x86_mic_simd_float.h"
54 // On MIC it is better to use scatter operations, so we define the load routines
55 // that use a SIMD offset variable first.
58 static inline void gmx_simdcall gatherLoadBySimdIntTranspose(const float* base,
59 SimdFInt32 simdoffset,
65 assert(std::size_t(base) % 16 == 0);
66 assert(align % 4 == 0);
68 // All instructions might be latency ~4 on MIC, so we use shifts where we
69 // only need a single instruction (since the shift parameter is an immediate),
70 // but multiplication otherwise.
73 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 2);
77 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 3);
81 simdoffset = simdoffset * SimdFInt32(align);
84 v0->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base, sizeof(float));
85 v1->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base + 1, sizeof(float));
86 v2->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base + 2, sizeof(float));
87 v3->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base + 3, sizeof(float));
91 static inline void gmx_simdcall
92 gatherLoadUBySimdIntTranspose(const float* base, SimdFInt32 simdoffset, SimdFloat* v0, SimdFloat* v1)
94 // All instructions might be latency ~4 on MIC, so we use shifts where we
95 // only need a single instruction (since the shift parameter is an immediate),
96 // but multiplication otherwise.
97 // For align == 2 we can merge the constant into the scale parameter,
98 // which can take constants up to 8 in total.
101 v0->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base, align * sizeof(float));
103 _mm512_i32gather_ps(simdoffset.simdInternal_, base + 1, align * sizeof(float));
109 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 2);
113 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 3);
117 simdoffset = simdoffset * SimdFInt32(align);
119 v0->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base, sizeof(float));
120 v1->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base + 1, sizeof(float));
125 static inline void gmx_simdcall
126 gatherLoadBySimdIntTranspose(const float* base, SimdFInt32 simdoffset, SimdFloat* v0, SimdFloat* v1)
128 assert(std::size_t(base) % 8 == 0);
129 assert(align % 2 == 0);
130 gatherLoadUBySimdIntTranspose<align>(base, simdoffset, v0, v1);
134 static inline void gmx_simdcall gatherLoadTranspose(const float* base,
135 const std::int32_t offset[],
141 gatherLoadBySimdIntTranspose<align>(base, simdLoad(offset, SimdFInt32Tag()), v0, v1, v2, v3);
145 static inline void gmx_simdcall
146 gatherLoadTranspose(const float* base, const std::int32_t offset[], SimdFloat* v0, SimdFloat* v1)
148 gatherLoadBySimdIntTranspose<align>(base, simdLoad(offset, SimdFInt32Tag()), v0, v1);
151 static const int c_simdBestPairAlignmentFloat = 2;
154 static inline void gmx_simdcall gatherLoadUTranspose(const float* base,
155 const std::int32_t offset[],
160 SimdFInt32 simdoffset;
162 assert(std::size_t(offset) % 64 == 0);
164 simdoffset = simdLoad(offset, SimdFInt32Tag());
166 // All instructions might be latency ~4 on MIC, so we use shifts where we
167 // only need a single instruction (since the shift parameter is an immediate),
168 // but multiplication otherwise.
171 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 2);
175 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 3);
179 simdoffset = simdoffset * SimdFInt32(align);
182 v0->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base, sizeof(float));
183 v1->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base + 1, sizeof(float));
184 v2->simdInternal_ = _mm512_i32gather_ps(simdoffset.simdInternal_, base + 2, sizeof(float));
189 static inline void gmx_simdcall
190 transposeScatterStoreU(float* base, const std::int32_t offset[], SimdFloat v0, SimdFloat v1, SimdFloat v2)
192 SimdFInt32 simdoffset;
194 assert(std::size_t(offset) % 64 == 0);
196 simdoffset = simdLoad(offset, SimdFInt32Tag());
198 // All instructions might be latency ~4 on MIC, so we use shifts where we
199 // only need a single instruction (since the shift parameter is an immediate),
200 // but multiplication otherwise.
203 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 2);
207 simdoffset.simdInternal_ = _mm512_slli_epi32(simdoffset.simdInternal_, 3);
211 simdoffset = simdoffset * SimdFInt32(align);
214 _mm512_i32scatter_ps(base, simdoffset.simdInternal_, v0.simdInternal_, sizeof(float));
215 _mm512_i32scatter_ps(base + 1, simdoffset.simdInternal_, v1.simdInternal_, sizeof(float));
216 _mm512_i32scatter_ps(base + 2, simdoffset.simdInternal_, v2.simdInternal_, sizeof(float));
221 static inline void gmx_simdcall
222 transposeScatterIncrU(float* base, const std::int32_t offset[], SimdFloat v0, SimdFloat v1, SimdFloat v2)
224 alignas(GMX_SIMD_ALIGNMENT) float rdata0[GMX_SIMD_FLOAT_WIDTH];
225 alignas(GMX_SIMD_ALIGNMENT) float rdata1[GMX_SIMD_FLOAT_WIDTH];
226 alignas(GMX_SIMD_ALIGNMENT) float rdata2[GMX_SIMD_FLOAT_WIDTH];
232 for (int i = 0; i < GMX_SIMD_FLOAT_WIDTH; i++)
234 base[align * offset[i] + 0] += rdata0[i];
235 base[align * offset[i] + 1] += rdata1[i];
236 base[align * offset[i] + 2] += rdata2[i];
241 static inline void gmx_simdcall
242 transposeScatterDecrU(float* base, const std::int32_t offset[], SimdFloat v0, SimdFloat v1, SimdFloat v2)
244 alignas(GMX_SIMD_ALIGNMENT) float rdata0[GMX_SIMD_FLOAT_WIDTH];
245 alignas(GMX_SIMD_ALIGNMENT) float rdata1[GMX_SIMD_FLOAT_WIDTH];
246 alignas(GMX_SIMD_ALIGNMENT) float rdata2[GMX_SIMD_FLOAT_WIDTH];
252 for (int i = 0; i < GMX_SIMD_FLOAT_WIDTH; i++)
254 base[align * offset[i] + 0] -= rdata0[i];
255 base[align * offset[i] + 1] -= rdata1[i];
256 base[align * offset[i] + 2] -= rdata2[i];
260 static inline void gmx_simdcall expandScalarsToTriplets(SimdFloat scalar,
261 SimdFloat* triplets0,
262 SimdFloat* triplets1,
263 SimdFloat* triplets2)
265 triplets0->simdInternal_ = _mm512_castsi512_ps(
266 _mm512_permutevar_epi32(_mm512_set_epi32(5, 4, 4, 4, 3, 3, 3, 2, 2, 2, 1, 1, 1, 0, 0, 0),
267 _mm512_castps_si512(scalar.simdInternal_)));
268 triplets1->simdInternal_ = _mm512_castsi512_ps(_mm512_permutevar_epi32(
269 _mm512_set_epi32(10, 10, 9, 9, 9, 8, 8, 8, 7, 7, 7, 6, 6, 6, 5, 5),
270 _mm512_castps_si512(scalar.simdInternal_)));
271 triplets2->simdInternal_ = _mm512_castsi512_ps(_mm512_permutevar_epi32(
272 _mm512_set_epi32(15, 15, 15, 14, 14, 14, 13, 13, 13, 12, 12, 12, 11, 11, 11, 10),
273 _mm512_castps_si512(scalar.simdInternal_)));
277 static inline float gmx_simdcall reduceIncr4ReturnSum(float* m, SimdFloat v0, SimdFloat v1, SimdFloat v2, SimdFloat v3)
280 __m512 t0, t1, t2, t3;
282 assert(std::size_t(m) % 16 == 0);
284 t0 = _mm512_add_ps(v0.simdInternal_, _mm512_swizzle_ps(v0.simdInternal_, _MM_SWIZ_REG_BADC));
285 t0 = _mm512_mask_add_ps(t0, _mm512_int2mask(0xCCCC), v2.simdInternal_,
286 _mm512_swizzle_ps(v2.simdInternal_, _MM_SWIZ_REG_BADC));
287 t1 = _mm512_add_ps(v1.simdInternal_, _mm512_swizzle_ps(v1.simdInternal_, _MM_SWIZ_REG_BADC));
288 t1 = _mm512_mask_add_ps(t1, _mm512_int2mask(0xCCCC), v3.simdInternal_,
289 _mm512_swizzle_ps(v3.simdInternal_, _MM_SWIZ_REG_BADC));
290 t2 = _mm512_add_ps(t0, _mm512_swizzle_ps(t0, _MM_SWIZ_REG_CDAB));
291 t2 = _mm512_mask_add_ps(t2, _mm512_int2mask(0xAAAA), t1, _mm512_swizzle_ps(t1, _MM_SWIZ_REG_CDAB));
293 t2 = _mm512_add_ps(t2, _mm512_permute4f128_ps(t2, _MM_PERM_BADC));
294 t2 = _mm512_add_ps(t2, _mm512_permute4f128_ps(t2, _MM_PERM_CDAB));
296 t0 = _mm512_mask_extload_ps(_mm512_undefined_ps(), _mm512_int2mask(0xF), m, _MM_UPCONV_PS_NONE,
297 _MM_BROADCAST_4X16, _MM_HINT_NONE);
298 t0 = _mm512_add_ps(t0, t2);
299 _mm512_mask_packstorelo_ps(m, _mm512_int2mask(0xF), t0);
301 t2 = _mm512_add_ps(t2, _mm512_swizzle_ps(t2, _MM_SWIZ_REG_BADC));
302 t2 = _mm512_add_ps(t2, _mm512_swizzle_ps(t2, _MM_SWIZ_REG_CDAB));
304 _mm512_mask_packstorelo_ps(&f, _mm512_mask2int(0x1), t2);
308 static inline SimdFloat gmx_simdcall loadDualHsimd(const float* m0, const float* m1)
310 assert(std::size_t(m0) % 32 == 0);
311 assert(std::size_t(m1) % 32 == 0);
313 return _mm512_castpd_ps(_mm512_mask_extload_pd(
314 _mm512_extload_pd(reinterpret_cast<const double*>(m0), _MM_UPCONV_PD_NONE,
315 _MM_BROADCAST_4X8, _MM_HINT_NONE),
316 _mm512_int2mask(0xF0), reinterpret_cast<const double*>(m1), _MM_UPCONV_PD_NONE,
317 _MM_BROADCAST_4X8, _MM_HINT_NONE));
320 static inline SimdFloat gmx_simdcall loadDuplicateHsimd(const float* m)
322 assert(std::size_t(m) % 32 == 0);
324 return _mm512_castpd_ps(_mm512_extload_pd(reinterpret_cast<const double*>(m),
325 _MM_UPCONV_PD_NONE, _MM_BROADCAST_4X8, _MM_HINT_NONE));
328 static inline SimdFloat gmx_simdcall loadU1DualHsimd(const float* m)
330 return _mm512_mask_extload_ps(
331 _mm512_extload_ps(m, _MM_UPCONV_PS_NONE, _MM_BROADCAST_1X16, _MM_HINT_NONE),
332 _mm512_int2mask(0xFF00), m + 1, _MM_UPCONV_PS_NONE, _MM_BROADCAST_1X16, _MM_HINT_NONE);
336 static inline void gmx_simdcall storeDualHsimd(float* m0, float* m1, SimdFloat a)
340 assert(std::size_t(m0) % 32 == 0);
341 assert(std::size_t(m1) % 32 == 0);
343 _mm512_mask_packstorelo_ps(m0, _mm512_int2mask(0x00FF), a.simdInternal_);
344 _mm512_mask_packstorelo_ps(m1, _mm512_int2mask(0xFF00), a.simdInternal_);
347 static inline void gmx_simdcall incrDualHsimd(float* m0, float* m1, SimdFloat a)
349 assert(std::size_t(m0) % 32 == 0);
350 assert(std::size_t(m1) % 32 == 0);
355 x = _mm512_castpd_ps(_mm512_extload_pd(reinterpret_cast<const double*>(m0), _MM_UPCONV_PD_NONE,
356 _MM_BROADCAST_4X8, _MM_HINT_NONE));
357 x = _mm512_add_ps(x, a.simdInternal_);
358 _mm512_mask_packstorelo_ps(m0, _mm512_int2mask(0x00FF), x);
361 x = _mm512_castpd_ps(_mm512_extload_pd(reinterpret_cast<const double*>(m1), _MM_UPCONV_PD_NONE,
362 _MM_BROADCAST_4X8, _MM_HINT_NONE));
363 x = _mm512_add_ps(x, a.simdInternal_);
364 _mm512_mask_packstorelo_ps(m1, _mm512_int2mask(0xFF00), x);
367 static inline void gmx_simdcall decrHsimd(float* m, SimdFloat a)
371 assert(std::size_t(m) % 32 == 0);
373 t = _mm512_castpd_ps(_mm512_extload_pd(reinterpret_cast<const double*>(m), _MM_UPCONV_PD_NONE,
374 _MM_BROADCAST_4X8, _MM_HINT_NONE));
375 a = _mm512_add_ps(a.simdInternal_, _mm512_permute4f128_ps(a.simdInternal_, _MM_PERM_BADC));
376 t = _mm512_sub_ps(t, a.simdInternal_);
377 _mm512_mask_packstorelo_ps(m, _mm512_int2mask(0x00FF), t);
382 static inline void gmx_simdcall gatherLoadTransposeHsimd(const float* base0,
384 const std::int32_t offset[],
388 __m512i idx0, idx1, idx;
391 assert(std::size_t(offset) % 32 == 0);
392 assert(std::size_t(base0) % 8 == 0);
393 assert(std::size_t(base1) % 8 == 0);
394 assert(std::size_t(align) % 2 == 0);
396 idx0 = _mm512_loadunpacklo_epi32(_mm512_undefined_epi32(), offset);
398 idx0 = _mm512_mullo_epi32(idx0, _mm512_set1_epi32(align));
399 idx1 = _mm512_add_epi32(idx0, _mm512_set1_epi32(1));
401 idx = _mm512_mask_permute4f128_epi32(idx0, _mm512_int2mask(0xFF00), idx1, _MM_PERM_BABA);
403 tmp1 = _mm512_i32gather_ps(idx, base0, sizeof(float));
404 tmp2 = _mm512_i32gather_ps(idx, base1, sizeof(float));
406 v0->simdInternal_ = _mm512_mask_permute4f128_ps(tmp1, _mm512_int2mask(0xFF00), tmp2, _MM_PERM_BABA);
407 v1->simdInternal_ = _mm512_mask_permute4f128_ps(tmp2, _mm512_int2mask(0x00FF), tmp1, _MM_PERM_DCDC);
410 static inline float gmx_simdcall reduceIncr4ReturnSumHsimd(float* m, SimdFloat v0, SimdFloat v1)
415 assert(std::size_t(m) % 32 == 0);
417 t0 = _mm512_add_ps(v0.simdInternal_, _mm512_swizzle_ps(v0.simdInternal_, _MM_SWIZ_REG_BADC));
418 t0 = _mm512_mask_add_ps(t0, _mm512_int2mask(0xCCCC), v1.simdInternal_,
419 _mm512_swizzle_ps(v1.simdInternal_, _MM_SWIZ_REG_BADC));
420 t0 = _mm512_add_ps(t0, _mm512_swizzle_ps(t0, _MM_SWIZ_REG_CDAB));
421 t0 = _mm512_add_ps(t0, _mm512_castpd_ps(_mm512_swizzle_pd(_mm512_castps_pd(t0), _MM_SWIZ_REG_BADC)));
422 t0 = _mm512_mask_permute4f128_ps(t0, _mm512_int2mask(0xAAAA), t0, _MM_PERM_BADC);
423 t1 = _mm512_mask_extload_ps(_mm512_undefined_ps(), _mm512_int2mask(0xF), m, _MM_UPCONV_PS_NONE,
424 _MM_BROADCAST_4X16, _MM_HINT_NONE);
425 t1 = _mm512_add_ps(t1, t0);
426 _mm512_mask_packstorelo_ps(m, _mm512_int2mask(0xF), t1);
428 t0 = _mm512_add_ps(t0, _mm512_swizzle_ps(t0, _MM_SWIZ_REG_BADC));
429 t0 = _mm512_add_ps(t0, _mm512_swizzle_ps(t0, _MM_SWIZ_REG_CDAB));
431 _mm512_mask_packstorelo_ps(&f, _mm512_mask2int(0x1), t0);
437 #endif // GMX_SIMD_IMPL_X86_MIC_UTIL_FLOAT_H