Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / simd / impl_x86_avx_256 / impl_x86_avx_256_simd_float.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_256_SIMD_FLOAT_H
37 #define GMX_SIMD_IMPL_X86_AVX_256_SIMD_FLOAT_H
38
39 #include "config.h"
40
41 #include <cassert>
42 #include <cstddef>
43 #include <cstdint>
44
45 #include <immintrin.h>
46
47 #include "gromacs/math/utilities.h"
48
49 namespace gmx
50 {
51
52 class SimdFloat
53 {
54 public:
55     SimdFloat() {}
56
57     SimdFloat(float f) : simdInternal_(_mm256_set1_ps(f)) {}
58
59     // Internal utility constructor to simplify return statements
60     SimdFloat(__m256 simd) : simdInternal_(simd) {}
61
62     __m256 simdInternal_;
63 };
64
65 class SimdFInt32
66 {
67 public:
68     SimdFInt32() {}
69
70     SimdFInt32(std::int32_t i) : simdInternal_(_mm256_set1_epi32(i)) {}
71
72     // Internal utility constructor to simplify return statements
73     SimdFInt32(__m256i simd) : simdInternal_(simd) {}
74
75     __m256i simdInternal_;
76 };
77
78 class SimdFBool
79 {
80 public:
81     SimdFBool() {}
82
83     SimdFBool(bool b) : simdInternal_(_mm256_castsi256_ps(_mm256_set1_epi32(b ? 0xFFFFFFFF : 0))) {}
84
85     // Internal utility constructor to simplify return statements
86     SimdFBool(__m256 simd) : simdInternal_(simd) {}
87
88     __m256 simdInternal_;
89 };
90
91 static inline SimdFloat gmx_simdcall simdLoad(const float* m, SimdFloatTag /*unused*/ = {})
92 {
93     assert(std::size_t(m) % 32 == 0);
94     return { _mm256_load_ps(m) };
95 }
96
97 static inline void gmx_simdcall store(float* m, SimdFloat a)
98 {
99     assert(std::size_t(m) % 32 == 0);
100     _mm256_store_ps(m, a.simdInternal_);
101 }
102
103 static inline SimdFloat gmx_simdcall simdLoadU(const float* m, SimdFloatTag /*unused*/ = {})
104 {
105     return { _mm256_loadu_ps(m) };
106 }
107
108 static inline void gmx_simdcall storeU(float* m, SimdFloat a)
109 {
110     _mm256_storeu_ps(m, a.simdInternal_);
111 }
112
113 static inline SimdFloat gmx_simdcall setZeroF()
114 {
115     return { _mm256_setzero_ps() };
116 }
117
118 static inline SimdFInt32 gmx_simdcall simdLoad(const std::int32_t* m, SimdFInt32Tag /*unused*/)
119 {
120     assert(std::size_t(m) % 32 == 0);
121     return { _mm256_load_si256(reinterpret_cast<const __m256i*>(m)) };
122 }
123
124 static inline void gmx_simdcall store(std::int32_t* m, SimdFInt32 a)
125 {
126     assert(std::size_t(m) % 32 == 0);
127     _mm256_store_si256(reinterpret_cast<__m256i*>(m), a.simdInternal_);
128 }
129
130 static inline SimdFInt32 gmx_simdcall simdLoadU(const std::int32_t* m, SimdFInt32Tag /*unused*/)
131 {
132     return { _mm256_loadu_si256(reinterpret_cast<const __m256i*>(m)) };
133 }
134
135 static inline void gmx_simdcall storeU(std::int32_t* m, SimdFInt32 a)
136 {
137     _mm256_storeu_si256(reinterpret_cast<__m256i*>(m), a.simdInternal_);
138 }
139
140 static inline SimdFInt32 gmx_simdcall setZeroFI()
141 {
142     return { _mm256_setzero_si256() };
143 }
144
145 template<int index>
146 static inline std::int32_t gmx_simdcall extract(SimdFInt32 a)
147 {
148     return _mm_extract_epi32(_mm256_extractf128_si256(a.simdInternal_, index >> 2), index & 0x3);
149 }
150
151 static inline SimdFloat gmx_simdcall operator&(SimdFloat a, SimdFloat b)
152 {
153     return { _mm256_and_ps(a.simdInternal_, b.simdInternal_) };
154 }
155
156 static inline SimdFloat gmx_simdcall andNot(SimdFloat a, SimdFloat b)
157 {
158     return { _mm256_andnot_ps(a.simdInternal_, b.simdInternal_) };
159 }
160
161 static inline SimdFloat gmx_simdcall operator|(SimdFloat a, SimdFloat b)
162 {
163     return { _mm256_or_ps(a.simdInternal_, b.simdInternal_) };
164 }
165
166 static inline SimdFloat gmx_simdcall operator^(SimdFloat a, SimdFloat b)
167 {
168     return { _mm256_xor_ps(a.simdInternal_, b.simdInternal_) };
169 }
170
171 static inline SimdFloat gmx_simdcall operator+(SimdFloat a, SimdFloat b)
172 {
173     return { _mm256_add_ps(a.simdInternal_, b.simdInternal_) };
174 }
175
176 static inline SimdFloat gmx_simdcall operator-(SimdFloat a, SimdFloat b)
177 {
178     return { _mm256_sub_ps(a.simdInternal_, b.simdInternal_) };
179 }
180
181 static inline SimdFloat gmx_simdcall operator-(SimdFloat x)
182 {
183     return { _mm256_xor_ps(x.simdInternal_, _mm256_set1_ps(GMX_FLOAT_NEGZERO)) };
184 }
185
186 static inline SimdFloat gmx_simdcall operator*(SimdFloat a, SimdFloat b)
187 {
188     return { _mm256_mul_ps(a.simdInternal_, b.simdInternal_) };
189 }
190
191 // Override for AVX2 and higher
192 #if GMX_SIMD_X86_AVX_256
193 static inline SimdFloat gmx_simdcall fma(SimdFloat a, SimdFloat b, SimdFloat c)
194 {
195     return { _mm256_add_ps(_mm256_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_) };
196 }
197
198 static inline SimdFloat gmx_simdcall fms(SimdFloat a, SimdFloat b, SimdFloat c)
199 {
200     return { _mm256_sub_ps(_mm256_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_) };
201 }
202
203 static inline SimdFloat gmx_simdcall fnma(SimdFloat a, SimdFloat b, SimdFloat c)
204 {
205     return { _mm256_sub_ps(c.simdInternal_, _mm256_mul_ps(a.simdInternal_, b.simdInternal_)) };
206 }
207
208 static inline SimdFloat gmx_simdcall fnms(SimdFloat a, SimdFloat b, SimdFloat c)
209 {
210     return { _mm256_sub_ps(_mm256_setzero_ps(),
211                            _mm256_add_ps(_mm256_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_)) };
212 }
213 #endif
214
215 static inline SimdFloat gmx_simdcall rsqrt(SimdFloat x)
216 {
217     return { _mm256_rsqrt_ps(x.simdInternal_) };
218 }
219
220 static inline SimdFloat gmx_simdcall rcp(SimdFloat x)
221 {
222     return { _mm256_rcp_ps(x.simdInternal_) };
223 }
224
225 static inline SimdFloat gmx_simdcall maskAdd(SimdFloat a, SimdFloat b, SimdFBool m)
226 {
227     return { _mm256_add_ps(a.simdInternal_, _mm256_and_ps(b.simdInternal_, m.simdInternal_)) };
228 }
229
230 static inline SimdFloat gmx_simdcall maskzMul(SimdFloat a, SimdFloat b, SimdFBool m)
231 {
232     return { _mm256_and_ps(_mm256_mul_ps(a.simdInternal_, b.simdInternal_), m.simdInternal_) };
233 }
234
235 static inline SimdFloat maskzFma(SimdFloat a, SimdFloat b, SimdFloat c, SimdFBool m)
236 {
237     return { _mm256_and_ps(_mm256_add_ps(_mm256_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_),
238                            m.simdInternal_) };
239 }
240
241 static inline SimdFloat maskzRsqrt(SimdFloat x, SimdFBool m)
242 {
243 #ifndef NDEBUG
244     x.simdInternal_ = _mm256_blendv_ps(_mm256_set1_ps(1.0F), x.simdInternal_, m.simdInternal_);
245 #endif
246     return { _mm256_and_ps(_mm256_rsqrt_ps(x.simdInternal_), m.simdInternal_) };
247 }
248
249 static inline SimdFloat maskzRcp(SimdFloat x, SimdFBool m)
250 {
251 #ifndef NDEBUG
252     x.simdInternal_ = _mm256_blendv_ps(_mm256_set1_ps(1.0F), x.simdInternal_, m.simdInternal_);
253 #endif
254     return { _mm256_and_ps(_mm256_rcp_ps(x.simdInternal_), m.simdInternal_) };
255 }
256
257 static inline SimdFloat gmx_simdcall abs(SimdFloat x)
258 {
259     return { _mm256_andnot_ps(_mm256_set1_ps(GMX_FLOAT_NEGZERO), x.simdInternal_) };
260 }
261
262 static inline SimdFloat gmx_simdcall max(SimdFloat a, SimdFloat b)
263 {
264     return { _mm256_max_ps(a.simdInternal_, b.simdInternal_) };
265 }
266
267 static inline SimdFloat gmx_simdcall min(SimdFloat a, SimdFloat b)
268 {
269     return { _mm256_min_ps(a.simdInternal_, b.simdInternal_) };
270 }
271
272 static inline SimdFloat gmx_simdcall round(SimdFloat x)
273 {
274     return { _mm256_round_ps(x.simdInternal_, _MM_FROUND_NINT) };
275 }
276
277 static inline SimdFloat gmx_simdcall trunc(SimdFloat x)
278 {
279     return { _mm256_round_ps(x.simdInternal_, _MM_FROUND_TRUNC) };
280 }
281
282 // Override for AVX2 and higher
283 #if GMX_SIMD_X86_AVX_256
284 static inline SimdFloat gmx_simdcall frexp(SimdFloat value, SimdFInt32* exponent)
285 {
286     const __m256 exponentMask = _mm256_castsi256_ps(_mm256_set1_epi32(0x7F800000));
287     const __m256 mantissaMask = _mm256_castsi256_ps(_mm256_set1_epi32(0x807FFFFF));
288     const __m256 half         = _mm256_set1_ps(0.5);
289     const __m128i exponentBias = _mm_set1_epi32(126); // add 1 to make our definition identical to frexp()
290     __m256i iExponent;
291     __m128i iExponentLow, iExponentHigh;
292
293     iExponent               = _mm256_castps_si256(_mm256_and_ps(value.simdInternal_, exponentMask));
294     iExponentHigh           = _mm256_extractf128_si256(iExponent, 0x1);
295     iExponentLow            = _mm256_castsi256_si128(iExponent);
296     iExponentLow            = _mm_srli_epi32(iExponentLow, 23);
297     iExponentHigh           = _mm_srli_epi32(iExponentHigh, 23);
298     iExponentLow            = _mm_sub_epi32(iExponentLow, exponentBias);
299     iExponentHigh           = _mm_sub_epi32(iExponentHigh, exponentBias);
300     iExponent               = _mm256_castsi128_si256(iExponentLow);
301     exponent->simdInternal_ = _mm256_insertf128_si256(iExponent, iExponentHigh, 0x1);
302
303     return { _mm256_or_ps(_mm256_and_ps(value.simdInternal_, mantissaMask), half) };
304 }
305
306 template<MathOptimization opt = MathOptimization::Safe>
307 static inline SimdFloat gmx_simdcall ldexp(SimdFloat value, SimdFInt32 exponent)
308 {
309     const __m128i exponentBias = _mm_set1_epi32(127);
310     __m256i       iExponent;
311     __m128i       iExponentLow, iExponentHigh;
312
313     iExponentHigh = _mm256_extractf128_si256(exponent.simdInternal_, 0x1);
314     iExponentLow  = _mm256_castsi256_si128(exponent.simdInternal_);
315
316     iExponentLow  = _mm_add_epi32(iExponentLow, exponentBias);
317     iExponentHigh = _mm_add_epi32(iExponentHigh, exponentBias);
318
319     if (opt == MathOptimization::Safe)
320     {
321         // Make sure biased argument is not negative
322         iExponentLow  = _mm_max_epi32(iExponentLow, _mm_setzero_si128());
323         iExponentHigh = _mm_max_epi32(iExponentHigh, _mm_setzero_si128());
324     }
325
326     iExponentLow  = _mm_slli_epi32(iExponentLow, 23);
327     iExponentHigh = _mm_slli_epi32(iExponentHigh, 23);
328     iExponent     = _mm256_castsi128_si256(iExponentLow);
329     iExponent     = _mm256_insertf128_si256(iExponent, iExponentHigh, 0x1);
330     return { _mm256_mul_ps(value.simdInternal_, _mm256_castsi256_ps(iExponent)) };
331 }
332 #endif
333
334 static inline float gmx_simdcall reduce(SimdFloat a)
335 {
336     __m128 t0;
337     t0 = _mm_add_ps(_mm256_castps256_ps128(a.simdInternal_), _mm256_extractf128_ps(a.simdInternal_, 0x1));
338     t0 = _mm_add_ps(t0, _mm_permute_ps(t0, _MM_SHUFFLE(1, 0, 3, 2)));
339     t0 = _mm_add_ss(t0, _mm_permute_ps(t0, _MM_SHUFFLE(0, 3, 2, 1)));
340     return *reinterpret_cast<float*>(&t0);
341 }
342
343 static inline SimdFBool gmx_simdcall operator==(SimdFloat a, SimdFloat b)
344 {
345     return { _mm256_cmp_ps(a.simdInternal_, b.simdInternal_, _CMP_EQ_OQ) };
346 }
347
348 static inline SimdFBool gmx_simdcall operator!=(SimdFloat a, SimdFloat b)
349 {
350     return { _mm256_cmp_ps(a.simdInternal_, b.simdInternal_, _CMP_NEQ_OQ) };
351 }
352
353 static inline SimdFBool gmx_simdcall operator<(SimdFloat a, SimdFloat b)
354 {
355     return { _mm256_cmp_ps(a.simdInternal_, b.simdInternal_, _CMP_LT_OQ) };
356 }
357
358 static inline SimdFBool gmx_simdcall operator<=(SimdFloat a, SimdFloat b)
359 {
360     return { _mm256_cmp_ps(a.simdInternal_, b.simdInternal_, _CMP_LE_OQ) };
361 }
362
363 // Override for AVX2 and higher
364 #if GMX_SIMD_X86_AVX_256
365 static inline SimdFBool gmx_simdcall testBits(SimdFloat a)
366 {
367     __m256 tst = _mm256_cvtepi32_ps(_mm256_castps_si256(a.simdInternal_));
368
369     return { _mm256_cmp_ps(tst, _mm256_setzero_ps(), _CMP_NEQ_OQ) };
370 }
371 #endif
372
373 static inline SimdFBool gmx_simdcall operator&&(SimdFBool a, SimdFBool b)
374 {
375     return { _mm256_and_ps(a.simdInternal_, b.simdInternal_) };
376 }
377
378 static inline SimdFBool gmx_simdcall operator||(SimdFBool a, SimdFBool b)
379 {
380     return { _mm256_or_ps(a.simdInternal_, b.simdInternal_) };
381 }
382
383 static inline bool gmx_simdcall anyTrue(SimdFBool a)
384 {
385     return _mm256_movemask_ps(a.simdInternal_) != 0;
386 }
387
388 static inline SimdFloat gmx_simdcall selectByMask(SimdFloat a, SimdFBool mask)
389 {
390     return { _mm256_and_ps(a.simdInternal_, mask.simdInternal_) };
391 }
392
393 static inline SimdFloat gmx_simdcall selectByNotMask(SimdFloat a, SimdFBool mask)
394 {
395     return { _mm256_andnot_ps(mask.simdInternal_, a.simdInternal_) };
396 }
397
398 static inline SimdFloat gmx_simdcall blend(SimdFloat a, SimdFloat b, SimdFBool sel)
399 {
400     return { _mm256_blendv_ps(a.simdInternal_, b.simdInternal_, sel.simdInternal_) };
401 }
402
403 static inline SimdFInt32 gmx_simdcall cvtR2I(SimdFloat a)
404 {
405     return { _mm256_cvtps_epi32(a.simdInternal_) };
406 }
407
408 static inline SimdFInt32 gmx_simdcall cvttR2I(SimdFloat a)
409 {
410     return { _mm256_cvttps_epi32(a.simdInternal_) };
411 }
412
413 static inline SimdFloat gmx_simdcall cvtI2R(SimdFInt32 a)
414 {
415     return { _mm256_cvtepi32_ps(a.simdInternal_) };
416 }
417
418 } // namespace gmx
419
420 #endif // GMX_SIMD_IMPL_X86_AVX_256_SIMD_FLOAT_H