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_AVX_256_SIMD_FLOAT_H
37 #define GMX_SIMD_IMPL_X86_AVX_256_SIMD_FLOAT_H
45 #include <immintrin.h>
47 #include "gromacs/math/utilities.h"
57 SimdFloat(float f) : simdInternal_(_mm256_set1_ps(f)) {}
59 // Internal utility constructor to simplify return statements
60 SimdFloat(__m256 simd) : simdInternal_(simd) {}
70 SimdFInt32(std::int32_t i) : simdInternal_(_mm256_set1_epi32(i)) {}
72 // Internal utility constructor to simplify return statements
73 SimdFInt32(__m256i simd) : simdInternal_(simd) {}
75 __m256i simdInternal_;
83 SimdFBool(bool b) : simdInternal_(_mm256_castsi256_ps(_mm256_set1_epi32(b ? 0xFFFFFFFF : 0))) {}
85 // Internal utility constructor to simplify return statements
86 SimdFBool(__m256 simd) : simdInternal_(simd) {}
91 static inline SimdFloat gmx_simdcall simdLoad(const float* m, SimdFloatTag /*unused*/ = {})
93 assert(std::size_t(m) % 32 == 0);
94 return { _mm256_load_ps(m) };
97 static inline void gmx_simdcall store(float* m, SimdFloat a)
99 assert(std::size_t(m) % 32 == 0);
100 _mm256_store_ps(m, a.simdInternal_);
103 static inline SimdFloat gmx_simdcall simdLoadU(const float* m, SimdFloatTag /*unused*/ = {})
105 return { _mm256_loadu_ps(m) };
108 static inline void gmx_simdcall storeU(float* m, SimdFloat a)
110 _mm256_storeu_ps(m, a.simdInternal_);
113 static inline SimdFloat gmx_simdcall setZeroF()
115 return { _mm256_setzero_ps() };
118 static inline SimdFInt32 gmx_simdcall simdLoad(const std::int32_t* m, SimdFInt32Tag /*unused*/)
120 assert(std::size_t(m) % 32 == 0);
121 return { _mm256_load_si256(reinterpret_cast<const __m256i*>(m)) };
124 static inline void gmx_simdcall store(std::int32_t* m, SimdFInt32 a)
126 assert(std::size_t(m) % 32 == 0);
127 _mm256_store_si256(reinterpret_cast<__m256i*>(m), a.simdInternal_);
130 static inline SimdFInt32 gmx_simdcall simdLoadU(const std::int32_t* m, SimdFInt32Tag /*unused*/)
132 return { _mm256_loadu_si256(reinterpret_cast<const __m256i*>(m)) };
135 static inline void gmx_simdcall storeU(std::int32_t* m, SimdFInt32 a)
137 _mm256_storeu_si256(reinterpret_cast<__m256i*>(m), a.simdInternal_);
140 static inline SimdFInt32 gmx_simdcall setZeroFI()
142 return { _mm256_setzero_si256() };
146 static inline std::int32_t gmx_simdcall extract(SimdFInt32 a)
148 return _mm_extract_epi32(_mm256_extractf128_si256(a.simdInternal_, index >> 2), index & 0x3);
151 static inline SimdFloat gmx_simdcall operator&(SimdFloat a, SimdFloat b)
153 return { _mm256_and_ps(a.simdInternal_, b.simdInternal_) };
156 static inline SimdFloat gmx_simdcall andNot(SimdFloat a, SimdFloat b)
158 return { _mm256_andnot_ps(a.simdInternal_, b.simdInternal_) };
161 static inline SimdFloat gmx_simdcall operator|(SimdFloat a, SimdFloat b)
163 return { _mm256_or_ps(a.simdInternal_, b.simdInternal_) };
166 static inline SimdFloat gmx_simdcall operator^(SimdFloat a, SimdFloat b)
168 return { _mm256_xor_ps(a.simdInternal_, b.simdInternal_) };
171 static inline SimdFloat gmx_simdcall operator+(SimdFloat a, SimdFloat b)
173 return { _mm256_add_ps(a.simdInternal_, b.simdInternal_) };
176 static inline SimdFloat gmx_simdcall operator-(SimdFloat a, SimdFloat b)
178 return { _mm256_sub_ps(a.simdInternal_, b.simdInternal_) };
181 static inline SimdFloat gmx_simdcall operator-(SimdFloat x)
183 return { _mm256_xor_ps(x.simdInternal_, _mm256_set1_ps(GMX_FLOAT_NEGZERO)) };
186 static inline SimdFloat gmx_simdcall operator*(SimdFloat a, SimdFloat b)
188 return { _mm256_mul_ps(a.simdInternal_, b.simdInternal_) };
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)
195 return { _mm256_add_ps(_mm256_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_) };
198 static inline SimdFloat gmx_simdcall fms(SimdFloat a, SimdFloat b, SimdFloat c)
200 return { _mm256_sub_ps(_mm256_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_) };
203 static inline SimdFloat gmx_simdcall fnma(SimdFloat a, SimdFloat b, SimdFloat c)
205 return { _mm256_sub_ps(c.simdInternal_, _mm256_mul_ps(a.simdInternal_, b.simdInternal_)) };
208 static inline SimdFloat gmx_simdcall fnms(SimdFloat a, SimdFloat b, SimdFloat c)
210 return { _mm256_sub_ps(_mm256_setzero_ps(),
211 _mm256_add_ps(_mm256_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_)) };
215 static inline SimdFloat gmx_simdcall rsqrt(SimdFloat x)
217 return { _mm256_rsqrt_ps(x.simdInternal_) };
220 static inline SimdFloat gmx_simdcall rcp(SimdFloat x)
222 return { _mm256_rcp_ps(x.simdInternal_) };
225 static inline SimdFloat gmx_simdcall maskAdd(SimdFloat a, SimdFloat b, SimdFBool m)
227 return { _mm256_add_ps(a.simdInternal_, _mm256_and_ps(b.simdInternal_, m.simdInternal_)) };
230 static inline SimdFloat gmx_simdcall maskzMul(SimdFloat a, SimdFloat b, SimdFBool m)
232 return { _mm256_and_ps(_mm256_mul_ps(a.simdInternal_, b.simdInternal_), m.simdInternal_) };
235 static inline SimdFloat maskzFma(SimdFloat a, SimdFloat b, SimdFloat c, SimdFBool m)
237 return { _mm256_and_ps(_mm256_add_ps(_mm256_mul_ps(a.simdInternal_, b.simdInternal_), c.simdInternal_),
241 static inline SimdFloat maskzRsqrt(SimdFloat x, SimdFBool m)
244 x.simdInternal_ = _mm256_blendv_ps(_mm256_set1_ps(1.0F), x.simdInternal_, m.simdInternal_);
246 return { _mm256_and_ps(_mm256_rsqrt_ps(x.simdInternal_), m.simdInternal_) };
249 static inline SimdFloat maskzRcp(SimdFloat x, SimdFBool m)
252 x.simdInternal_ = _mm256_blendv_ps(_mm256_set1_ps(1.0F), x.simdInternal_, m.simdInternal_);
254 return { _mm256_and_ps(_mm256_rcp_ps(x.simdInternal_), m.simdInternal_) };
257 static inline SimdFloat gmx_simdcall abs(SimdFloat x)
259 return { _mm256_andnot_ps(_mm256_set1_ps(GMX_FLOAT_NEGZERO), x.simdInternal_) };
262 static inline SimdFloat gmx_simdcall max(SimdFloat a, SimdFloat b)
264 return { _mm256_max_ps(a.simdInternal_, b.simdInternal_) };
267 static inline SimdFloat gmx_simdcall min(SimdFloat a, SimdFloat b)
269 return { _mm256_min_ps(a.simdInternal_, b.simdInternal_) };
272 static inline SimdFloat gmx_simdcall round(SimdFloat x)
274 return { _mm256_round_ps(x.simdInternal_, _MM_FROUND_NINT) };
277 static inline SimdFloat gmx_simdcall trunc(SimdFloat x)
279 return { _mm256_round_ps(x.simdInternal_, _MM_FROUND_TRUNC) };
282 // Override for AVX2 and higher
283 #if GMX_SIMD_X86_AVX_256
284 static inline SimdFloat gmx_simdcall frexp(SimdFloat value, SimdFInt32* exponent)
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()
291 __m128i iExponentLow, iExponentHigh;
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);
303 return { _mm256_or_ps(_mm256_and_ps(value.simdInternal_, mantissaMask), half) };
306 template<MathOptimization opt = MathOptimization::Safe>
307 static inline SimdFloat gmx_simdcall ldexp(SimdFloat value, SimdFInt32 exponent)
309 const __m128i exponentBias = _mm_set1_epi32(127);
311 __m128i iExponentLow, iExponentHigh;
313 iExponentHigh = _mm256_extractf128_si256(exponent.simdInternal_, 0x1);
314 iExponentLow = _mm256_castsi256_si128(exponent.simdInternal_);
316 iExponentLow = _mm_add_epi32(iExponentLow, exponentBias);
317 iExponentHigh = _mm_add_epi32(iExponentHigh, exponentBias);
319 if (opt == MathOptimization::Safe)
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());
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)) };
334 static inline float gmx_simdcall reduce(SimdFloat a)
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);
343 static inline SimdFBool gmx_simdcall operator==(SimdFloat a, SimdFloat b)
345 return { _mm256_cmp_ps(a.simdInternal_, b.simdInternal_, _CMP_EQ_OQ) };
348 static inline SimdFBool gmx_simdcall operator!=(SimdFloat a, SimdFloat b)
350 return { _mm256_cmp_ps(a.simdInternal_, b.simdInternal_, _CMP_NEQ_OQ) };
353 static inline SimdFBool gmx_simdcall operator<(SimdFloat a, SimdFloat b)
355 return { _mm256_cmp_ps(a.simdInternal_, b.simdInternal_, _CMP_LT_OQ) };
358 static inline SimdFBool gmx_simdcall operator<=(SimdFloat a, SimdFloat b)
360 return { _mm256_cmp_ps(a.simdInternal_, b.simdInternal_, _CMP_LE_OQ) };
363 // Override for AVX2 and higher
364 #if GMX_SIMD_X86_AVX_256
365 static inline SimdFBool gmx_simdcall testBits(SimdFloat a)
367 __m256 tst = _mm256_cvtepi32_ps(_mm256_castps_si256(a.simdInternal_));
369 return { _mm256_cmp_ps(tst, _mm256_setzero_ps(), _CMP_NEQ_OQ) };
373 static inline SimdFBool gmx_simdcall operator&&(SimdFBool a, SimdFBool b)
375 return { _mm256_and_ps(a.simdInternal_, b.simdInternal_) };
378 static inline SimdFBool gmx_simdcall operator||(SimdFBool a, SimdFBool b)
380 return { _mm256_or_ps(a.simdInternal_, b.simdInternal_) };
383 static inline bool gmx_simdcall anyTrue(SimdFBool a)
385 return _mm256_movemask_ps(a.simdInternal_) != 0;
388 static inline SimdFloat gmx_simdcall selectByMask(SimdFloat a, SimdFBool mask)
390 return { _mm256_and_ps(a.simdInternal_, mask.simdInternal_) };
393 static inline SimdFloat gmx_simdcall selectByNotMask(SimdFloat a, SimdFBool mask)
395 return { _mm256_andnot_ps(mask.simdInternal_, a.simdInternal_) };
398 static inline SimdFloat gmx_simdcall blend(SimdFloat a, SimdFloat b, SimdFBool sel)
400 return { _mm256_blendv_ps(a.simdInternal_, b.simdInternal_, sel.simdInternal_) };
403 static inline SimdFInt32 gmx_simdcall cvtR2I(SimdFloat a)
405 return { _mm256_cvtps_epi32(a.simdInternal_) };
408 static inline SimdFInt32 gmx_simdcall cvttR2I(SimdFloat a)
410 return { _mm256_cvttps_epi32(a.simdInternal_) };
413 static inline SimdFloat gmx_simdcall cvtI2R(SimdFInt32 a)
415 return { _mm256_cvtepi32_ps(a.simdInternal_) };
420 #endif // GMX_SIMD_IMPL_X86_AVX_256_SIMD_FLOAT_H