2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2017,2019,2020, 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_AVX2_256_SIMD_FLOAT_H
37 #define GMX_SIMD_IMPL_X86_AVX2_256_SIMD_FLOAT_H
41 #include <immintrin.h>
43 #include "gromacs/math/utilities.h"
44 #include "gromacs/simd/impl_x86_avx_256/impl_x86_avx_256_simd_float.h"
52 MSVC_DIAGNOSTIC_IGNORE(26495) // simdInternal_ is not being initialized!
55 SimdFIBool(bool b) : simdInternal_(_mm256_set1_epi32(b ? 0xFFFFFFFF : 0)) {}
57 // Internal utility constructor to simplify return statements
58 SimdFIBool(__m256i simd) : simdInternal_(simd) {}
60 __m256i simdInternal_;
63 static inline SimdFloat gmx_simdcall fma(SimdFloat a, SimdFloat b, SimdFloat c)
65 return { _mm256_fmadd_ps(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
68 static inline SimdFloat gmx_simdcall fms(SimdFloat a, SimdFloat b, SimdFloat c)
70 return { _mm256_fmsub_ps(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
73 static inline SimdFloat gmx_simdcall fnma(SimdFloat a, SimdFloat b, SimdFloat c)
75 return { _mm256_fnmadd_ps(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
78 static inline SimdFloat gmx_simdcall fnms(SimdFloat a, SimdFloat b, SimdFloat c)
80 return { _mm256_fnmsub_ps(a.simdInternal_, b.simdInternal_, c.simdInternal_) };
83 static inline SimdFBool gmx_simdcall testBits(SimdFloat a)
85 __m256i ia = _mm256_castps_si256(a.simdInternal_);
86 __m256i res = _mm256_andnot_si256(_mm256_cmpeq_epi32(ia, _mm256_setzero_si256()),
87 _mm256_cmpeq_epi32(ia, ia));
89 return { _mm256_castsi256_ps(res) };
92 static inline SimdFloat gmx_simdcall frexp(SimdFloat value, SimdFInt32* exponent)
94 const __m256 exponentMask = _mm256_castsi256_ps(_mm256_set1_epi32(0x7F800000));
95 const __m256 mantissaMask = _mm256_castsi256_ps(_mm256_set1_epi32(0x807FFFFF));
96 const __m256i exponentBias = _mm256_set1_epi32(126); // add 1 to make our definition identical to frexp()
97 const __m256 half = _mm256_set1_ps(0.5);
100 iExponent = _mm256_castps_si256(_mm256_and_ps(value.simdInternal_, exponentMask));
101 exponent->simdInternal_ = _mm256_sub_epi32(_mm256_srli_epi32(iExponent, 23), exponentBias);
103 return { _mm256_or_ps(_mm256_and_ps(value.simdInternal_, mantissaMask), half) };
106 template<MathOptimization opt = MathOptimization::Safe>
107 static inline SimdFloat gmx_simdcall ldexp(SimdFloat value, SimdFInt32 exponent)
109 const __m256i exponentBias = _mm256_set1_epi32(127);
110 __m256i iExponent = _mm256_add_epi32(exponent.simdInternal_, exponentBias);
112 if (opt == MathOptimization::Safe)
114 // Make sure biased argument is not negative
115 iExponent = _mm256_max_epi32(iExponent, _mm256_setzero_si256());
118 iExponent = _mm256_slli_epi32(iExponent, 23);
119 return { _mm256_mul_ps(value.simdInternal_, _mm256_castsi256_ps(iExponent)) };
122 static inline SimdFInt32 gmx_simdcall operator&(SimdFInt32 a, SimdFInt32 b)
124 return { _mm256_and_si256(a.simdInternal_, b.simdInternal_) };
127 static inline SimdFInt32 gmx_simdcall andNot(SimdFInt32 a, SimdFInt32 b)
129 return { _mm256_andnot_si256(a.simdInternal_, b.simdInternal_) };
132 static inline SimdFInt32 gmx_simdcall operator|(SimdFInt32 a, SimdFInt32 b)
134 return { _mm256_or_si256(a.simdInternal_, b.simdInternal_) };
137 static inline SimdFInt32 gmx_simdcall operator^(SimdFInt32 a, SimdFInt32 b)
139 return { _mm256_xor_si256(a.simdInternal_, b.simdInternal_) };
142 static inline SimdFInt32 gmx_simdcall operator+(SimdFInt32 a, SimdFInt32 b)
144 return { _mm256_add_epi32(a.simdInternal_, b.simdInternal_) };
147 static inline SimdFInt32 gmx_simdcall operator-(SimdFInt32 a, SimdFInt32 b)
149 return { _mm256_sub_epi32(a.simdInternal_, b.simdInternal_) };
152 static inline SimdFInt32 gmx_simdcall operator*(SimdFInt32 a, SimdFInt32 b)
154 return { _mm256_mullo_epi32(a.simdInternal_, b.simdInternal_) };
157 static inline SimdFIBool gmx_simdcall operator==(SimdFInt32 a, SimdFInt32 b)
159 return { _mm256_cmpeq_epi32(a.simdInternal_, b.simdInternal_) };
162 static inline SimdFIBool gmx_simdcall testBits(SimdFInt32 a)
164 return { _mm256_andnot_si256(_mm256_cmpeq_epi32(a.simdInternal_, _mm256_setzero_si256()),
165 _mm256_cmpeq_epi32(a.simdInternal_, a.simdInternal_)) };
168 static inline SimdFIBool gmx_simdcall operator<(SimdFInt32 a, SimdFInt32 b)
170 return { _mm256_cmpgt_epi32(b.simdInternal_, a.simdInternal_) };
173 static inline SimdFIBool gmx_simdcall operator&&(SimdFIBool a, SimdFIBool b)
175 return { _mm256_and_si256(a.simdInternal_, b.simdInternal_) };
178 static inline SimdFIBool gmx_simdcall operator||(SimdFIBool a, SimdFIBool b)
180 return { _mm256_or_si256(a.simdInternal_, b.simdInternal_) };
183 static inline bool gmx_simdcall anyTrue(SimdFIBool a)
185 return _mm256_movemask_epi8(a.simdInternal_) != 0;
188 static inline SimdFInt32 gmx_simdcall selectByMask(SimdFInt32 a, SimdFIBool mask)
190 return { _mm256_and_si256(a.simdInternal_, mask.simdInternal_) };
193 static inline SimdFInt32 gmx_simdcall selectByNotMask(SimdFInt32 a, SimdFIBool mask)
195 return { _mm256_andnot_si256(mask.simdInternal_, a.simdInternal_) };
198 static inline SimdFInt32 gmx_simdcall blend(SimdFInt32 a, SimdFInt32 b, SimdFIBool sel)
200 return { _mm256_blendv_epi8(a.simdInternal_, b.simdInternal_, sel.simdInternal_) };
203 static inline SimdFIBool gmx_simdcall cvtB2IB(SimdFBool a)
205 return { _mm256_castps_si256(a.simdInternal_) };
208 static inline SimdFBool gmx_simdcall cvtIB2B(SimdFIBool a)
210 return { _mm256_castsi256_ps(a.simdInternal_) };
215 #endif // GMX_SIMD_IMPL_X86_AVX2_256_SIMD_FLOAT_H