Improve accuracy of SIMD exp for small args
[gromacs.git] / src / gromacs / simd / impl_x86_avx2_256 / impl_x86_avx2_256_simd_float.h
blobe549657a0c1aa4242a4a6d692cf131d8ebc85caa
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2017, 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
39 #include "config.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"
46 namespace gmx
49 class SimdFIBool
51 public:
52 SimdFIBool() {}
54 SimdFIBool(bool b) : simdInternal_(_mm256_set1_epi32( b ? 0xFFFFFFFF : 0)) {}
56 // Internal utility constructor to simplify return statements
57 SimdFIBool(__m256i simd) : simdInternal_(simd) {}
59 __m256i simdInternal_;
62 static inline SimdFloat gmx_simdcall
63 fma(SimdFloat a, SimdFloat b, SimdFloat c)
65 return {
66 _mm256_fmadd_ps(a.simdInternal_, b.simdInternal_, c.simdInternal_)
70 static inline SimdFloat gmx_simdcall
71 fms(SimdFloat a, SimdFloat b, SimdFloat c)
73 return {
74 _mm256_fmsub_ps(a.simdInternal_, b.simdInternal_, c.simdInternal_)
78 static inline SimdFloat gmx_simdcall
79 fnma(SimdFloat a, SimdFloat b, SimdFloat c)
81 return {
82 _mm256_fnmadd_ps(a.simdInternal_, b.simdInternal_, c.simdInternal_)
86 static inline SimdFloat gmx_simdcall
87 fnms(SimdFloat a, SimdFloat b, SimdFloat c)
89 return {
90 _mm256_fnmsub_ps(a.simdInternal_, b.simdInternal_, c.simdInternal_)
94 static inline SimdFBool gmx_simdcall
95 testBits(SimdFloat a)
97 __m256i ia = _mm256_castps_si256(a.simdInternal_);
98 __m256i res = _mm256_andnot_si256( _mm256_cmpeq_epi32(ia, _mm256_setzero_si256()), _mm256_cmpeq_epi32(ia, ia));
100 return {
101 _mm256_castsi256_ps(res)
105 static inline SimdFloat gmx_simdcall
106 frexp(SimdFloat value, SimdFInt32 * exponent)
108 const __m256 exponentMask = _mm256_castsi256_ps(_mm256_set1_epi32(0x7F800000));
109 const __m256 mantissaMask = _mm256_castsi256_ps(_mm256_set1_epi32(0x807FFFFF));
110 const __m256i exponentBias = _mm256_set1_epi32(126); // add 1 to make our definition identical to frexp()
111 const __m256 half = _mm256_set1_ps(0.5);
112 __m256i iExponent;
114 iExponent = _mm256_castps_si256(_mm256_and_ps(value.simdInternal_, exponentMask));
115 exponent->simdInternal_ = _mm256_sub_epi32(_mm256_srli_epi32(iExponent, 23), exponentBias);
117 return {
118 _mm256_or_ps(_mm256_and_ps(value.simdInternal_, mantissaMask), half)
122 template <MathOptimization opt = MathOptimization::Safe>
123 static inline SimdFloat gmx_simdcall
124 ldexp(SimdFloat value, SimdFInt32 exponent)
126 const __m256i exponentBias = _mm256_set1_epi32(127);
127 __m256i iExponent = _mm256_add_epi32(exponent.simdInternal_, exponentBias);
129 if (opt == MathOptimization::Safe)
131 // Make sure biased argument is not negative
132 iExponent = _mm256_max_epi32(iExponent, _mm256_setzero_si256());
135 iExponent = _mm256_slli_epi32(iExponent, 23);
136 return {
137 _mm256_mul_ps(value.simdInternal_, _mm256_castsi256_ps(iExponent))
141 static inline SimdFInt32 gmx_simdcall
142 operator<<(SimdFInt32 a, int n)
144 return {
145 _mm256_slli_epi32(a.simdInternal_, n)
149 static inline SimdFInt32 gmx_simdcall
150 operator>>(SimdFInt32 a, int n)
152 return {
153 _mm256_srli_epi32(a.simdInternal_, n)
157 static inline SimdFInt32 gmx_simdcall
158 operator&(SimdFInt32 a, SimdFInt32 b)
160 return {
161 _mm256_and_si256(a.simdInternal_, b.simdInternal_)
165 static inline SimdFInt32 gmx_simdcall
166 andNot(SimdFInt32 a, SimdFInt32 b)
168 return {
169 _mm256_andnot_si256(a.simdInternal_, b.simdInternal_)
173 static inline SimdFInt32 gmx_simdcall
174 operator|(SimdFInt32 a, SimdFInt32 b)
176 return {
177 _mm256_or_si256(a.simdInternal_, b.simdInternal_)
181 static inline SimdFInt32 gmx_simdcall
182 operator^(SimdFInt32 a, SimdFInt32 b)
184 return {
185 _mm256_xor_si256(a.simdInternal_, b.simdInternal_)
189 static inline SimdFInt32 gmx_simdcall
190 operator+(SimdFInt32 a, SimdFInt32 b)
192 return {
193 _mm256_add_epi32(a.simdInternal_, b.simdInternal_)
197 static inline SimdFInt32 gmx_simdcall
198 operator-(SimdFInt32 a, SimdFInt32 b)
200 return {
201 _mm256_sub_epi32(a.simdInternal_, b.simdInternal_)
205 static inline SimdFInt32 gmx_simdcall
206 operator*(SimdFInt32 a, SimdFInt32 b)
208 return {
209 _mm256_mullo_epi32(a.simdInternal_, b.simdInternal_)
213 static inline SimdFIBool gmx_simdcall
214 operator==(SimdFInt32 a, SimdFInt32 b)
216 return {
217 _mm256_cmpeq_epi32(a.simdInternal_, b.simdInternal_)
221 static inline SimdFIBool gmx_simdcall
222 testBits(SimdFInt32 a)
224 return {
225 _mm256_andnot_si256(_mm256_cmpeq_epi32(a.simdInternal_, _mm256_setzero_si256()),
226 _mm256_cmpeq_epi32(a.simdInternal_, a.simdInternal_))
230 static inline SimdFIBool gmx_simdcall
231 operator<(SimdFInt32 a, SimdFInt32 b)
233 return {
234 _mm256_cmpgt_epi32(b.simdInternal_, a.simdInternal_)
238 static inline SimdFIBool gmx_simdcall
239 operator&&(SimdFIBool a, SimdFIBool b)
241 return {
242 _mm256_and_si256(a.simdInternal_, b.simdInternal_)
246 static inline SimdFIBool gmx_simdcall
247 operator||(SimdFIBool a, SimdFIBool b)
249 return {
250 _mm256_or_si256(a.simdInternal_, b.simdInternal_)
254 static inline bool gmx_simdcall
255 anyTrue(SimdFIBool a) { return _mm256_movemask_epi8(a.simdInternal_) != 0; }
257 static inline SimdFInt32 gmx_simdcall
258 selectByMask(SimdFInt32 a, SimdFIBool mask)
260 return {
261 _mm256_and_si256(a.simdInternal_, mask.simdInternal_)
265 static inline SimdFInt32 gmx_simdcall
266 selectByNotMask(SimdFInt32 a, SimdFIBool mask)
268 return {
269 _mm256_andnot_si256(mask.simdInternal_, a.simdInternal_)
273 static inline SimdFInt32 gmx_simdcall
274 blend(SimdFInt32 a, SimdFInt32 b, SimdFIBool sel)
276 return {
277 _mm256_blendv_epi8(a.simdInternal_, b.simdInternal_, sel.simdInternal_)
281 static inline SimdFIBool gmx_simdcall
282 cvtB2IB(SimdFBool a)
284 return {
285 _mm256_castps_si256(a.simdInternal_)
289 static inline SimdFBool gmx_simdcall
290 cvtIB2B(SimdFIBool a)
292 return {
293 _mm256_castsi256_ps(a.simdInternal_)
297 } // namespace gmx
299 #endif // GMX_SIMD_IMPL_X86_AVX2_256_SIMD_FLOAT_H