e688d14c94d7beb4678777131968fea812f595de
[alexxy/gromacs.git] / src / gromacs / simd / impl_x86_avx_256 / impl_x86_avx_256.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014, 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_H
37 #define GMX_SIMD_IMPL_X86_AVX_256_H
38
39 #include <math.h>
40 #include <immintrin.h>
41
42 #include "config.h"
43
44 /* It is cleaner to start the AVX implementation from scratch rather than
45  * first inheriting from SSE4.1, which in turn inherits from SSE2. However,
46  * the capabilities still form a superset.
47  */
48 #define GMX_SIMD_X86_SSE2_OR_HIGHER
49 #define GMX_SIMD_X86_SSE4_1_OR_HIGHER
50 #define GMX_SIMD_X86_AVX_256_OR_HIGHER
51
52
53 /* x86 256-bit AVX SIMD instruction wrappers
54  *
55  * Please see documentation in gromacs/simd/simd.h for defines.
56  */
57
58 /* Capability definitions for 256-bit AVX - no inheritance from SSE */
59 #define GMX_SIMD_HAVE_FLOAT
60 #define GMX_SIMD_HAVE_DOUBLE
61 #define GMX_SIMD_HAVE_SIMD_HARDWARE
62 #define GMX_SIMD_HAVE_LOADU
63 #define GMX_SIMD_HAVE_STOREU
64 #define GMX_SIMD_HAVE_LOGICAL
65 #undef  GMX_SIMD_HAVE_FMA
66 #undef  GMX_SIMD_HAVE_FRACTION
67 #define GMX_SIMD_HAVE_FINT32
68 #define GMX_SIMD_HAVE_FINT32_EXTRACT     /* Emulated */
69 #undef  GMX_SIMD_HAVE_FINT32_LOGICAL     /* AVX1 cannot do 256-bit int shifts */
70 #undef  GMX_SIMD_HAVE_FINT32_ARITHMETICS /* AVX1 cannot do 256-bit int +,-,*  */
71 #define GMX_SIMD_HAVE_DINT32
72 #define GMX_SIMD_HAVE_DINT32_EXTRACT     /* Native, dint uses 128-bit SIMD    */
73 #define GMX_SIMD_HAVE_DINT32_LOGICAL
74 #define GMX_SIMD_HAVE_DINT32_ARITHMETICS
75 #define GMX_SIMD4_HAVE_FLOAT
76 #define GMX_SIMD4_HAVE_DOUBLE
77
78 /* Implementation details */
79 #define GMX_SIMD_FLOAT_WIDTH         8
80 #define GMX_SIMD_DOUBLE_WIDTH        4
81 #define GMX_SIMD_FINT32_WIDTH        8
82 #define GMX_SIMD_DINT32_WIDTH        4
83 #define GMX_SIMD_RSQRT_BITS         11
84 #define GMX_SIMD_RCP_BITS           11
85
86 /****************************************************
87  *      SINGLE PRECISION SIMD IMPLEMENTATION        *
88  ****************************************************/
89 #define gmx_simd_float_t           __m256
90 #define gmx_simd_load_f            _mm256_load_ps
91 #define gmx_simd_load1_f           _mm256_broadcast_ss
92 #define gmx_simd_set1_f            _mm256_set1_ps
93 #define gmx_simd_store_f           _mm256_store_ps
94 #define gmx_simd_loadu_f           _mm256_loadu_ps
95 #define gmx_simd_storeu_f          _mm256_storeu_ps
96 #define gmx_simd_setzero_f         _mm256_setzero_ps
97 #define gmx_simd_add_f             _mm256_add_ps
98 #define gmx_simd_sub_f             _mm256_sub_ps
99 #define gmx_simd_mul_f             _mm256_mul_ps
100 #define gmx_simd_fmadd_f(a, b, c)    _mm256_add_ps(_mm256_mul_ps(a, b), c)
101 #define gmx_simd_fmsub_f(a, b, c)    _mm256_sub_ps(_mm256_mul_ps(a, b), c)
102 #define gmx_simd_fnmadd_f(a, b, c)   _mm256_sub_ps(c, _mm256_mul_ps(a, b))
103 #define gmx_simd_fnmsub_f(a, b, c)   _mm256_sub_ps(_mm256_setzero_ps(), gmx_simd_fmadd_f(a, b, c))
104 #define gmx_simd_and_f             _mm256_and_ps
105 #define gmx_simd_andnot_f          _mm256_andnot_ps
106 #define gmx_simd_or_f              _mm256_or_ps
107 #define gmx_simd_xor_f             _mm256_xor_ps
108 #define gmx_simd_rsqrt_f           _mm256_rsqrt_ps
109 #define gmx_simd_rcp_f             _mm256_rcp_ps
110 #define gmx_simd_fabs_f(x)         _mm256_andnot_ps(_mm256_set1_ps(GMX_FLOAT_NEGZERO), x)
111 #define gmx_simd_fneg_f(x)         _mm256_xor_ps(x, _mm256_set1_ps(GMX_FLOAT_NEGZERO))
112 #define gmx_simd_max_f             _mm256_max_ps
113 #define gmx_simd_min_f             _mm256_min_ps
114 #define gmx_simd_round_f(x)        _mm256_round_ps(x, _MM_FROUND_NINT)
115 #define gmx_simd_trunc_f(x)        _mm256_round_ps(x, _MM_FROUND_TRUNC)
116 #define gmx_simd_fraction_f(x)     _mm256_sub_ps(x, gmx_simd_trunc_f(x))
117 #define gmx_simd_get_exponent_f    gmx_simd_get_exponent_f_avx_256
118 #define gmx_simd_get_mantissa_f    gmx_simd_get_mantissa_f_avx_256
119 #define gmx_simd_set_exponent_f    gmx_simd_set_exponent_f_avx_256
120 /* integer datatype corresponding to float: gmx_simd_fint32_t */
121 #define gmx_simd_fint32_t          __m256i
122 #define gmx_simd_load_fi(m)        _mm256_load_si256((__m256i const*)m)
123 #define gmx_simd_set1_fi           _mm256_set1_epi32
124 #define gmx_simd_store_fi(m, x)    _mm256_store_si256((__m256i *)m, x)
125 #define gmx_simd_loadu_fi(m)       _mm256_loadu_si256((__m256i const*)m)
126 #define gmx_simd_storeu_fi(m, x)   _mm256_storeu_si256((__m256i *)m, x)
127 #define gmx_simd_setzero_fi        _mm256_setzero_si256
128 #define gmx_simd_cvt_f2i           _mm256_cvtps_epi32
129 #define gmx_simd_cvtt_f2i          _mm256_cvttps_epi32
130 #define gmx_simd_cvt_i2f           _mm256_cvtepi32_ps
131 #define gmx_simd_extract_fi(x, i)   _mm_extract_epi32(_mm256_extractf128_si256(x, (i)>>2), (i)&0x3)
132 /* Integer logical ops on gmx_simd_fint32_t */
133 /* gmx_simd_add_fi not supported     */
134 /* gmx_simd_sub_fi not supported     */
135 /* gmx_simd_mul_fi not supported     */
136 /* gmx_simd_slli_fi not supported    */
137 /* gmx_simd_srli_fi not supported    */
138 /* gmx_simd_and_fi not supported     */
139 /* gmx_simd_andnot_fi not supported  */
140 /* gmx_simd_or_fi not supported      */
141 /* gmx_simd_xor_fi not supported     */
142 /* Integer arithmetic ops on gmx_simd_fint32_t */
143 /* gmx_simd_add_fi not supported     */
144 /* gmx_simd_sub_fi not supported     */
145 /* gmx_simd_mul_fi not supported     */
146 /* Boolean & comparison operations on gmx_simd_float_t */
147 #define gmx_simd_fbool_t           __m256
148 #define gmx_simd_cmpeq_f(a, b)      _mm256_cmp_ps(a, b, _CMP_EQ_OQ)
149 #define gmx_simd_cmplt_f(a, b)      _mm256_cmp_ps(a, b, _CMP_LT_OQ)
150 #define gmx_simd_cmple_f(a, b)      _mm256_cmp_ps(a, b, _CMP_LE_OQ)
151 #define gmx_simd_and_fb            _mm256_and_ps
152 #define gmx_simd_or_fb             _mm256_or_ps
153 #define gmx_simd_anytrue_fb        _mm256_movemask_ps
154 #define gmx_simd_blendzero_f       _mm256_and_ps
155 #define gmx_simd_blendnotzero_f(a, sel)  _mm256_andnot_ps(sel, a)
156 #define gmx_simd_blendv_f          _mm256_blendv_ps
157 #define gmx_simd_reduce_f          gmx_simd_reduce_f_avx_256
158 /* Boolean & comparison operations on gmx_simd_fint32_t */
159 #define gmx_simd_fibool_t          __m256i
160 /* gmx_simd_cmpeq_fi not supported        */
161 /* gmx_simd_cmplt_fi not supported        */
162 /* gmx_simd_and_fib not supported         */
163 /* gmx_simd_or_fib not supported          */
164 /* gmx_simd_anytrue_fib not supported     */
165 /* gmx_simd_blendzero_fi not supported    */
166 /* gmx_simd_blendnotzero_fi not supported    */
167 /* gmx_simd_blendv_fi not supported       */
168 /* Conversions between different booleans */
169 #define gmx_simd_cvt_fb2fib        _mm256_castps_si256
170 #define gmx_simd_cvt_fib2fb        _mm256_castsi256_ps
171
172 /****************************************************
173  *      DOUBLE PRECISION SIMD IMPLEMENTATION        *
174  ****************************************************/
175 #define gmx_simd_double_t          __m256d
176 #define gmx_simd_load_d            _mm256_load_pd
177 #define gmx_simd_load1_d           _mm256_broadcast_sd
178 #define gmx_simd_set1_d            _mm256_set1_pd
179 #define gmx_simd_store_d           _mm256_store_pd
180 #define gmx_simd_loadu_d           _mm256_loadu_pd
181 #define gmx_simd_storeu_d          _mm256_storeu_pd
182 #define gmx_simd_setzero_d         _mm256_setzero_pd
183 #define gmx_simd_add_d             _mm256_add_pd
184 #define gmx_simd_sub_d             _mm256_sub_pd
185 #define gmx_simd_mul_d             _mm256_mul_pd
186 #define gmx_simd_fmadd_d(a, b, c)    _mm256_add_pd(_mm256_mul_pd(a, b), c)
187 #define gmx_simd_fmsub_d(a, b, c)    _mm256_sub_pd(_mm256_mul_pd(a, b), c)
188 #define gmx_simd_fnmadd_d(a, b, c)   _mm256_sub_pd(c, _mm256_mul_pd(a, b))
189 #define gmx_simd_fnmsub_d(a, b, c)   _mm256_sub_pd(_mm256_setzero_pd(), gmx_simd_fmadd_d(a, b, c))
190 #define gmx_simd_and_d             _mm256_and_pd
191 #define gmx_simd_andnot_d          _mm256_andnot_pd
192 #define gmx_simd_or_d              _mm256_or_pd
193 #define gmx_simd_xor_d             _mm256_xor_pd
194 #define gmx_simd_rsqrt_d(x)        _mm256_cvtps_pd(_mm_rsqrt_ps(_mm256_cvtpd_ps(x)))
195 #define gmx_simd_rcp_d(x)          _mm256_cvtps_pd(_mm_rcp_ps(_mm256_cvtpd_ps(x)))
196 #define gmx_simd_fabs_d(x)         _mm256_andnot_pd(_mm256_set1_pd(-0.0), x)
197 #define gmx_simd_fneg_d(x)         _mm256_xor_pd(x, _mm256_set1_pd(-0.0))
198 #define gmx_simd_max_d             _mm256_max_pd
199 #define gmx_simd_min_d             _mm256_min_pd
200 #define gmx_simd_round_d(x)        _mm256_round_pd(x, _MM_FROUND_NINT)
201 #define gmx_simd_trunc_d(x)        _mm256_round_pd(x, _MM_FROUND_TRUNC)
202 #define gmx_simd_fraction_d(x)     _mm256_sub_pd(x, gmx_simd_trunc_d(x))
203 #define gmx_simd_get_exponent_d    gmx_simd_get_exponent_d_avx_256
204 #define gmx_simd_get_mantissa_d    gmx_simd_get_mantissa_d_avx_256
205 #define gmx_simd_set_exponent_d    gmx_simd_set_exponent_d_avx_256
206 /* integer datatype corresponding to double: gmx_simd_dint32_t */
207 #define gmx_simd_dint32_t          __m128i
208 #define gmx_simd_load_di(m)        _mm_load_si128((const __m128i *)m)
209 #define gmx_simd_set1_di           _mm_set1_epi32
210 #define gmx_simd_store_di(m, x)     _mm_store_si128((__m128i *)m, x)
211 #define gmx_simd_loadu_di(m)       _mm_loadu_si128((const __m128i *)m)
212 #define gmx_simd_storeu_di(m, x)    _mm_storeu_si128((__m128i *)m, x)
213 #define gmx_simd_setzero_di        _mm_setzero_si128
214 #define gmx_simd_cvt_d2i           _mm256_cvtpd_epi32
215 #define gmx_simd_cvtt_d2i          _mm256_cvttpd_epi32
216 #define gmx_simd_cvt_i2d           _mm256_cvtepi32_pd
217 #define gmx_simd_extract_di        _mm_extract_epi32
218 /* Integer logical ops on gmx_simd_dint32_t */
219 #define gmx_simd_slli_di           _mm_slli_epi32
220 #define gmx_simd_srli_di           _mm_srli_epi32
221 #define gmx_simd_and_di            _mm_and_si128
222 #define gmx_simd_andnot_di         _mm_andnot_si128
223 #define gmx_simd_or_di             _mm_or_si128
224 #define gmx_simd_xor_di            _mm_xor_si128
225 /* Integer arithmetic ops on integer datatype corresponding to double */
226 #define gmx_simd_add_di            _mm_add_epi32
227 #define gmx_simd_sub_di            _mm_sub_epi32
228 #define gmx_simd_mul_di            _mm_mullo_epi32
229 /* Boolean & comparison operations on gmx_simd_double_t */
230 #define gmx_simd_dbool_t           __m256d
231 #define gmx_simd_cmpeq_d(a, b)      _mm256_cmp_pd(a, b, _CMP_EQ_OQ)
232 #define gmx_simd_cmplt_d(a, b)      _mm256_cmp_pd(a, b, _CMP_LT_OQ)
233 #define gmx_simd_cmple_d(a, b)      _mm256_cmp_pd(a, b, _CMP_LE_OQ)
234 #define gmx_simd_and_db            _mm256_and_pd
235 #define gmx_simd_or_db             _mm256_or_pd
236 #define gmx_simd_anytrue_db        _mm256_movemask_pd
237 #define gmx_simd_blendzero_d       _mm256_and_pd
238 #define gmx_simd_blendnotzero_d(a, sel)  _mm256_andnot_pd(sel, a)
239 #define gmx_simd_blendv_d          _mm256_blendv_pd
240 #define gmx_simd_reduce_d          gmx_simd_reduce_d_avx_256
241 /* Boolean & comparison operations on gmx_simd_dint32_t */
242 #define gmx_simd_dibool_t          __m128i
243 #define gmx_simd_cmpeq_di          _mm_cmpeq_epi32
244 #define gmx_simd_cmplt_di          _mm_cmplt_epi32
245 #define gmx_simd_and_dib           _mm_and_si128
246 #define gmx_simd_or_dib            _mm_or_si128
247 #define gmx_simd_anytrue_dib       _mm_movemask_epi8
248 #define gmx_simd_blendzero_di      _mm_and_si128
249 #define gmx_simd_blendnotzero_di(a, sel)  _mm_andnot_si128(sel, a)
250 #define gmx_simd_blendv_di         _mm_blendv_epi8
251 /* Conversions between different booleans */
252 #define gmx_simd_cvt_db2dib        gmx_simd_cvt_db2dib_avx_256
253 #define gmx_simd_cvt_dib2db        gmx_simd_cvt_dib2db_avx_256
254 /* Float/double conversion */
255 #define gmx_simd_cvt_f2dd          gmx_simd_cvt_f2dd_avx_256
256 #define gmx_simd_cvt_dd2f          gmx_simd_cvt_dd2f_avx_256
257
258 /****************************************************
259  *      SINGLE PRECISION SIMD4 IMPLEMENTATION       *
260  ****************************************************/
261 #define gmx_simd4_float_t          __m128
262 #define gmx_simd4_load_f           _mm_load_ps
263 #define gmx_simd4_load1_f          _mm_broadcast_ss
264 #define gmx_simd4_set1_f           _mm_set1_ps
265 #define gmx_simd4_store_f          _mm_store_ps
266 #define gmx_simd4_loadu_f          _mm_loadu_ps
267 #define gmx_simd4_storeu_f         _mm_storeu_ps
268 #define gmx_simd4_setzero_f        _mm_setzero_ps
269 #define gmx_simd4_add_f            _mm_add_ps
270 #define gmx_simd4_sub_f            _mm_sub_ps
271 #define gmx_simd4_mul_f            _mm_mul_ps
272 #define gmx_simd4_fmadd_f(a, b, c)   _mm_add_ps(_mm_mul_ps(a, b), c)
273 #define gmx_simd4_fmsub_f(a, b, c)   _mm_sub_ps(_mm_mul_ps(a, b), c)
274 #define gmx_simd4_fnmadd_f(a, b, c)  _mm_sub_ps(c, _mm_mul_ps(a, b))
275 #define gmx_simd4_fnmsub_f(a, b, c)  _mm_sub_ps(_mm_setzero_ps(), gmx_simd4_fmadd_f(a, b, c))
276 #define gmx_simd4_and_f            _mm_and_ps
277 #define gmx_simd4_andnot_f         _mm_andnot_ps
278 #define gmx_simd4_or_f             _mm_or_ps
279 #define gmx_simd4_xor_f            _mm_xor_ps
280 #define gmx_simd4_rsqrt_f          _mm_rsqrt_ps
281 #define gmx_simd4_fabs_f(x)        _mm_andnot_ps(_mm_set1_ps(-0.0), x)
282 #define gmx_simd4_fneg_f(x)        _mm_xor_ps(x, _mm_set1_ps(-0.0))
283 #define gmx_simd4_max_f            _mm_max_ps
284 #define gmx_simd4_min_f            _mm_min_ps
285 #define gmx_simd4_round_f(x)       _mm_round_ps(x, _MM_FROUND_NINT)
286 #define gmx_simd4_trunc_f(x)       _mm_round_ps(x, _MM_FROUND_TRUNC)
287 #define gmx_simd4_dotproduct3_f    gmx_simd4_dotproduct3_f_avx_256
288 #define gmx_simd4_fbool_t          __m128
289 #define gmx_simd4_cmpeq_f          _mm_cmpeq_ps
290 #define gmx_simd4_cmplt_f          _mm_cmplt_ps
291 #define gmx_simd4_cmple_f          _mm_cmple_ps
292 #define gmx_simd4_and_fb           _mm_and_ps
293 #define gmx_simd4_or_fb            _mm_or_ps
294 #define gmx_simd4_anytrue_fb       _mm_movemask_ps
295 #define gmx_simd4_blendzero_f      _mm_and_ps
296 #define gmx_simd4_blendnotzero_f(a, sel)  _mm_andnot_ps(sel, a)
297 #define gmx_simd4_blendv_f         _mm_blendv_ps
298 #define gmx_simd4_reduce_f         gmx_simd4_reduce_f_avx_256
299
300 /****************************************************
301  *      DOUBLE PRECISION SIMD4 IMPLEMENTATION       *
302  ****************************************************/
303 #define gmx_simd4_double_t          gmx_simd_double_t
304 #define gmx_simd4_load_d            gmx_simd_load_d
305 #define gmx_simd4_load1_d           gmx_simd_load1_d
306 #define gmx_simd4_set1_d            gmx_simd_set1_d
307 #define gmx_simd4_store_d           gmx_simd_store_d
308 #define gmx_simd4_loadu_d           gmx_simd_loadu_d
309 #define gmx_simd4_storeu_d          gmx_simd_storeu_d
310 #define gmx_simd4_setzero_d         gmx_simd_setzero_d
311 #define gmx_simd4_add_d             gmx_simd_add_d
312 #define gmx_simd4_sub_d             gmx_simd_sub_d
313 #define gmx_simd4_mul_d             gmx_simd_mul_d
314 #define gmx_simd4_fmadd_d           gmx_simd_fmadd_d
315 #define gmx_simd4_fmsub_d           gmx_simd_fmsub_d
316 #define gmx_simd4_fnmadd_d          gmx_simd_fnmadd_d
317 #define gmx_simd4_fnmsub_d          gmx_simd_fnmsub_d
318 #define gmx_simd4_and_d             gmx_simd_and_d
319 #define gmx_simd4_andnot_d          gmx_simd_andnot_d
320 #define gmx_simd4_or_d              gmx_simd_or_d
321 #define gmx_simd4_xor_d             gmx_simd_xor_d
322 #define gmx_simd4_rsqrt_d           gmx_simd_rsqrt_d
323 #define gmx_simd4_fabs_d            gmx_simd_fabs_d
324 #define gmx_simd4_fneg_d            gmx_simd_fneg_d
325 #define gmx_simd4_max_d             gmx_simd_max_d
326 #define gmx_simd4_min_d             gmx_simd_min_d
327 #define gmx_simd4_round_d           gmx_simd_round_d
328 #define gmx_simd4_trunc_d           gmx_simd_trunc_d
329 #define gmx_simd4_dotproduct3_d     gmx_simd4_dotproduct3_d_avx_256
330 #define gmx_simd4_dbool_t           gmx_simd_dbool_t
331 #define gmx_simd4_cmpeq_d           gmx_simd_cmpeq_d
332 #define gmx_simd4_cmplt_d           gmx_simd_cmplt_d
333 #define gmx_simd4_cmple_d           gmx_simd_cmple_d
334 #define gmx_simd4_and_db            gmx_simd_and_db
335 #define gmx_simd4_or_db             gmx_simd_or_db
336 #define gmx_simd4_anytrue_db        gmx_simd_anytrue_db
337 #define gmx_simd4_blendzero_d       gmx_simd_blendzero_d
338 #define gmx_simd4_blendnotzero_d    gmx_simd_blendnotzero_d
339 #define gmx_simd4_blendv_d          gmx_simd_blendv_d
340 #define gmx_simd4_reduce_d          gmx_simd_reduce_d
341 /* SIMD4 float/double conversion */
342 #define gmx_simd4_cvt_f2d           _mm256_cvtps_pd
343 #define gmx_simd4_cvt_d2f           _mm256_cvtpd_ps
344
345 /*********************************************************
346  * SIMD SINGLE PRECISION IMPLEMENTATION HELPER FUNCTIONS *
347  *********************************************************/
348 static gmx_inline __m256 gmx_simdcall
349 gmx_simd_get_exponent_f_avx_256(__m256 x)
350 {
351     const __m256  expmask      = _mm256_castsi256_ps(_mm256_set1_epi32(0x7F800000));
352     const __m128i expbias      = _mm_set1_epi32(127);
353     __m256i       iexp256;
354     __m128i       iexp128a, iexp128b;
355
356     iexp256   = _mm256_castps_si256(_mm256_and_ps(x, expmask));
357     iexp128b  = _mm256_extractf128_si256(iexp256, 0x1);
358     iexp128a  = _mm256_castsi256_si128(iexp256);
359     iexp128a  = _mm_srli_epi32(iexp128a, 23);
360     iexp128b  = _mm_srli_epi32(iexp128b, 23);
361     iexp128a  = _mm_sub_epi32(iexp128a, expbias);
362     iexp128b  = _mm_sub_epi32(iexp128b, expbias);
363     iexp256   = _mm256_castsi128_si256(iexp128a);
364     iexp256   = _mm256_insertf128_si256(iexp256, iexp128b, 0x1);
365     return _mm256_cvtepi32_ps(iexp256);
366 }
367
368 static gmx_inline __m256 gmx_simdcall
369 gmx_simd_get_mantissa_f_avx_256(__m256 x)
370 {
371     const __m256 mantmask   = _mm256_castsi256_ps(_mm256_set1_epi32(0x007FFFFF));
372     const __m256 one        = _mm256_set1_ps(1.0);
373
374     x = _mm256_and_ps(x, mantmask);
375     return _mm256_or_ps(x, one);
376 }
377
378 static gmx_inline __m256 gmx_simdcall
379 gmx_simd_set_exponent_f_avx_256(__m256 x)
380 {
381     const __m128i expbias      = _mm_set1_epi32(127);
382     __m256i       iexp256;
383     __m128i       iexp128a, iexp128b;
384
385     iexp256   = _mm256_cvtps_epi32(x);
386     iexp128b  = _mm256_extractf128_si256(iexp256, 0x1);
387     iexp128a  = _mm256_castsi256_si128(iexp256);
388     iexp128a  = _mm_slli_epi32(_mm_add_epi32(iexp128a, expbias), 23);
389     iexp128b  = _mm_slli_epi32(_mm_add_epi32(iexp128b, expbias), 23);
390     iexp256   = _mm256_castsi128_si256(iexp128a);
391     iexp256   = _mm256_insertf128_si256(iexp256, iexp128b, 0x1);
392     return _mm256_castsi256_ps(iexp256);
393 }
394
395 static gmx_inline float gmx_simdcall
396 gmx_simd_reduce_f_avx_256(__m256 a)
397 {
398     float  f;
399
400     __m128 a0, a1;
401     a  = _mm256_hadd_ps(a, a);
402     a  = _mm256_hadd_ps(a, a);
403     a0 = _mm256_castps256_ps128(a);
404     a1 = _mm256_extractf128_ps(a, 0x1);
405     a0 = _mm_add_ss(a0, a1);
406     _mm_store_ss(&f, a0);
407     return f;
408 }
409
410 /*********************************************************
411  * SIMD DOUBLE PRECISION IMPLEMENTATION HELPER FUNCTIONS *
412  *********************************************************/
413 static gmx_inline __m256d gmx_simdcall
414 gmx_simd_get_exponent_d_avx_256(__m256d x)
415 {
416     const __m256d expmask      = _mm256_castsi256_pd( _mm256_set1_epi64x(0x7FF0000000000000LL));
417     const __m128i expbias      = _mm_set1_epi32(1023);
418     __m256i       iexp256;
419     __m128i       iexp128a, iexp128b;
420
421     iexp256   = _mm256_castpd_si256(_mm256_and_pd(x, expmask));
422     iexp128b  = _mm256_extractf128_si256(iexp256, 0x1);
423     iexp128a  = _mm256_castsi256_si128(iexp256);
424     iexp128a  = _mm_srli_epi64(iexp128a, 52);
425     iexp128b  = _mm_srli_epi64(iexp128b, 52);
426     iexp128a  = _mm_shuffle_epi32(iexp128a, _MM_SHUFFLE(1, 1, 2, 0));
427     iexp128b  = _mm_shuffle_epi32(iexp128b, _MM_SHUFFLE(2, 0, 1, 1));
428     iexp128a  = _mm_or_si128(iexp128a, iexp128b);
429     iexp128a  = _mm_sub_epi32(iexp128a, expbias);
430     return _mm256_cvtepi32_pd(iexp128a);
431 }
432
433 static gmx_inline __m256d gmx_simdcall
434 gmx_simd_get_mantissa_d_avx_256(__m256d x)
435 {
436     const __m256d mantmask   = _mm256_castsi256_pd(_mm256_set1_epi64x(0x000FFFFFFFFFFFFFLL));
437     const __m256d one        = _mm256_set1_pd(1.0);
438
439     x = _mm256_and_pd(x, mantmask);
440     return _mm256_or_pd(x, one);
441 }
442
443 static gmx_inline __m256d gmx_simdcall
444 gmx_simd_set_exponent_d_avx_256(__m256d x)
445 {
446     const __m128i expbias      = _mm_set1_epi32(1023);
447     __m128i       iexp128a, iexp128b;
448
449     iexp128a = _mm256_cvtpd_epi32(x);
450     iexp128a = _mm_add_epi32(iexp128a, expbias);
451     iexp128b = _mm_shuffle_epi32(iexp128a, _MM_SHUFFLE(3, 3, 2, 2));
452     iexp128a = _mm_shuffle_epi32(iexp128a, _MM_SHUFFLE(1, 1, 0, 0));
453     iexp128b = _mm_slli_epi64(iexp128b, 52);
454     iexp128a = _mm_slli_epi64(iexp128a, 52);
455     return _mm256_castsi256_pd(_mm256_insertf128_si256(_mm256_castsi128_si256(iexp128a), iexp128b, 0x1));
456 }
457
458 static gmx_inline double gmx_simdcall
459 gmx_simd_reduce_d_avx_256(__m256d a)
460 {
461     double  f;
462     __m128d a0, a1;
463     a  = _mm256_hadd_pd(a, a);
464     a0 = _mm256_castpd256_pd128(a);
465     a1 = _mm256_extractf128_pd(a, 0x1);
466     a0 = _mm_add_sd(a0, a1);
467     _mm_store_sd(&f, a0);
468     return f;
469 }
470
471 static gmx_inline gmx_simd_dibool_t gmx_simdcall
472 gmx_simd_cvt_db2dib_avx_256(gmx_simd_dbool_t a)
473 {
474     __m128i a1 = _mm256_extractf128_si256(_mm256_castpd_si256(a), 0x1);
475     __m128i a0 = _mm256_castsi256_si128(_mm256_castpd_si256(a));
476     a0 = _mm_shuffle_epi32(a0, _MM_SHUFFLE(2, 0, 2, 0));
477     a1 = _mm_shuffle_epi32(a1, _MM_SHUFFLE(2, 0, 2, 0));
478     return _mm_blend_epi16(a0, a1, 0xF0);
479 }
480
481 static gmx_inline gmx_simd_dbool_t gmx_simdcall
482 gmx_simd_cvt_dib2db_avx_256(gmx_simd_dibool_t a)
483 {
484     __m128i a1 = _mm_shuffle_epi32(a, _MM_SHUFFLE(3, 3, 2, 2));
485     __m128i a0 = _mm_shuffle_epi32(a, _MM_SHUFFLE(1, 1, 0, 0));
486     return _mm256_castsi256_pd(_mm256_insertf128_si256(_mm256_castsi128_si256(a0), a1, 0x1));
487 }
488
489 static gmx_inline void gmx_simdcall
490 gmx_simd_cvt_f2dd_avx_256(__m256 f, __m256d *d0, __m256d *d1)
491 {
492     *d0 = _mm256_cvtps_pd(_mm256_castps256_ps128(f));
493     *d1 = _mm256_cvtps_pd(_mm256_extractf128_ps(f, 0x1));
494 }
495
496 static gmx_inline __m256 gmx_simdcall
497 gmx_simd_cvt_dd2f_avx_256(__m256d d0, __m256d d1)
498 {
499     __m128 f0 = _mm256_cvtpd_ps(d0);
500     __m128 f1 = _mm256_cvtpd_ps(d1);
501     return _mm256_insertf128_ps(_mm256_castps128_ps256(f0), f1, 0x1);
502 }
503
504 /* SIMD4 reduce helper */
505 static gmx_inline float gmx_simdcall
506 gmx_simd4_reduce_f_avx_256(__m128 a)
507 {
508     float f;
509     a = _mm_hadd_ps(a, a);
510     a = _mm_hadd_ps(a, a);
511     _mm_store_ss(&f, a);
512     return f;
513 }
514
515 /* SIMD4 Dotproduct helper function */
516 static gmx_inline float gmx_simdcall
517 gmx_simd4_dotproduct3_f_avx_256(__m128 a, __m128 b)
518 {
519     float  f;
520     __m128 c;
521     a = _mm_mul_ps(a, b);
522     c = _mm_add_ps(a, _mm_permute_ps(a, _MM_SHUFFLE(0, 3, 2, 1)));
523     c = _mm_add_ps(c, _mm_permute_ps(a, _MM_SHUFFLE(1, 0, 3, 2)));
524     _mm_store_ss(&f, c);
525     return f;
526 }
527
528 static gmx_inline double gmx_simdcall
529 gmx_simd4_dotproduct3_d_avx_256(__m256d a, __m256d b)
530 {
531     double  d;
532     __m128d tmp1, tmp2;
533     a    = _mm256_mul_pd(a, b);
534     tmp1 = _mm256_castpd256_pd128(a);
535     tmp2 = _mm256_extractf128_pd(a, 0x1);
536
537     tmp1 = _mm_add_pd(tmp1, _mm_permute_pd(tmp1, _MM_SHUFFLE2(0, 1)));
538     tmp1 = _mm_add_pd(tmp1, tmp2);
539     _mm_store_sd(&d, tmp1);
540     return d;
541 }
542
543 /* Function to check whether SIMD operations have resulted in overflow */
544 static int
545 gmx_simd_check_and_reset_overflow(void)
546 {
547     int MXCSR;
548     int sse_overflow;
549
550     MXCSR = _mm_getcsr();
551     /* The overflow flag is bit 3 in the register */
552     if (MXCSR & 0x0008)
553     {
554         sse_overflow = 1;
555         /* Set the overflow flag to zero */
556         MXCSR = MXCSR & 0xFFF7;
557         _mm_setcsr(MXCSR);
558     }
559     else
560     {
561         sse_overflow = 0;
562     }
563     return sse_overflow;
564 }
565
566
567 #endif /* GMX_SIMD_IMPL_X86_AVX_256_H */