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