5b3096bee322f5fe367c39da2165a030ea6b1db7
[alexxy/gromacs.git] / src / gromacs / simd / impl_x86_avx_128_fma / impl_x86_avx_128_fma.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_128_FMA_H
37 #define GMX_SIMD_IMPL_X86_AVX_128_FMA_H
38
39 #include <math.h>
40 #include <immintrin.h>
41 #include <x86intrin.h>
42
43 /* x86 128-bit AVX with FMA SIMD instruction wrappers
44  *
45  * Please see documentation in gromacs/simd/simd.h for details
46  */
47
48 /* Inherit parts of AVX_128_FMA from SSE4.1 */
49 #include "gromacs/simd/impl_x86_sse4_1/impl_x86_sse4_1.h"
50 /* Increment over SSE4.1 capabilities */
51 #define GMX_SIMD_X86_AVX_128_FMA_OR_HIGHER
52
53 /* Override some capability definitions for things added in AVX over SSE4.1 */
54 #define GMX_SIMD_HAVE_FMA
55 #define GMX_SIMD_HAVE_FRACTION
56 #define GMX_SIMD4_HAVE_DOUBLE  /* We can use 256-bit operations for this */
57
58 /* SINGLE */
59 #undef  gmx_simd_fmadd_ps
60 #define gmx_simd_fmadd_ps           _mm_macc_ps
61 #undef  gmx_simd_fmsub_ps
62 #define gmx_simd_fmsub_ps(a, b, c)    _mm_msub_ps
63 #undef  gmx_simd_fnmadd_ps
64 #define gmx_simd_fnmadd_ps(a, b, c)   _mm_nmacc_ps
65 #undef  gmx_simd_fnmsub_ps
66 #define gmx_simd_fnmsub_ps(a, b, c)   _mm_nmsub_ps
67 #undef  gmx_simd_fraction_f
68 #define gmx_simd_fraction_f         _mm_frcz_ps
69
70 /* DOUBLE */
71 #undef  gmx_simd_fmadd_pd
72 #define gmx_simd_fmadd_pd            _mm_macc_pd
73 #undef  gmx_simd_fmsub_pd
74 #define gmx_simd_fmsub_pd(a, b, c)     _mm_msub_pd
75 #undef  gmx_simd_fnmadd_pd
76 #define gmx_simd_fnmadd_pd(a, b, c)    _mm_nmacc_pd
77 #undef  gmx_simd_fnmsub_pd
78 #define gmx_simd_fnmsub_pd(a, b, c)    _mm_nmsub_pd
79 #undef  gmx_simd_fraction_d
80 #define gmx_simd_fraction_d          _mm_frcz_pd
81
82 /* Even if the _main_ SIMD implementation for this architecture file corresponds
83  * to 128-bit AVX (since it will be faster), the 256-bit operations will always
84  * be available in AVX, so we can use them for double precision SIMD4!
85  */
86 /* SIMD4 Double precision floating point */
87 #define gmx_simd4_double_t               __m256d
88 #define gmx_simd4_load_d                 _mm256_load_pd
89 #define gmx_simd4_load1_d                _mm256_broadcast_sd
90 #define gmx_simd4_set1_d                 _mm256_set1_pd
91 #define gmx_simd4_store_d                _mm256_store_pd
92 #define gmx_simd4_loadu_d                _mm256_loadu_pd
93 #define gmx_simd4_storeu_d               _mm256_storeu_pd
94 #define gmx_simd4_setzero_d              _mm256_setzero_pd
95 #define gmx_simd4_add_d                  _mm256_add_pd
96 #define gmx_simd4_sub_d                  _mm256_sub_pd
97 #define gmx_simd4_mul_d                  _mm256_mul_pd
98 #define gmx_simd4_fmadd_d                _mm256_macc_pd
99 #define gmx_simd4_fmsub_d                _mm256_msub_pd
100 #define gmx_simd4_fnmadd_d               _mm256_nmacc_pd
101 #define gmx_simd4_fnmsub_d               _mm256_nmsub_pd
102 #define gmx_simd4_and_d                  _mm256_and_pd
103 #define gmx_simd4_andnot_d               _mm256_andnot_pd
104 #define gmx_simd4_or_d                   _mm256_or_pd
105 #define gmx_simd4_xor_d                  _mm256_xor_pd
106 #define gmx_simd4_rsqrt_d(x)             _mm256_cvtps_pd(_mm_rsqrt_ps(_mm256_cvtpd_ps(x)))
107 #define gmx_simd4_fabs_d(x)              _mm256_andnot_pd(_mm256_set1_pd(GMX_DOUBLE_NEGZERO), x)
108 #define gmx_simd4_fneg_d(x)              _mm256_xor_pd(x, _mm256_set1_pd(GMX_DOUBLE_NEGZERO))
109 #define gmx_simd4_max_d                  _mm256_max_pd
110 #define gmx_simd4_min_d                  _mm256_min_pd
111 #define gmx_simd4_round_d(x)             _mm256_round_pd(x, _MM_FROUND_NINT)
112 #define gmx_simd4_trunc_d(x)             _mm256_round_pd(x, _MM_FROUND_TRUNC)
113 #define gmx_simd4_dotproduct3_d          gmx_simd4_dotproduct3_d_avx_128_fma
114 /* SIMD4 booleans corresponding to double */
115 #define gmx_simd4_dbool_t                __m256d
116 #define gmx_simd4_cmpeq_d(a, b)           _mm256_cmp_pd(a, b, _CMP_EQ_OQ)
117 #define gmx_simd4_cmplt_d(a, b)           _mm256_cmp_pd(a, b, _CMP_LT_OQ)
118 #define gmx_simd4_cmple_d(a, b)           _mm256_cmp_pd(a, b, _CMP_LE_OQ)
119 #define gmx_simd4_and_db                 _mm256_and_pd
120 #define gmx_simd4_or_db                  _mm256_or_pd
121 #define gmx_simd4_anytrue_db             _mm256_movemask_pd
122 #define gmx_simd4_blendzero_d            _mm256_and_pd
123 #define gmx_simd4_blendnotzero_d(a, sel)  _mm256_andnot_pd(sel, a)
124 #define gmx_simd4_blendv_d               _mm256_blendv_pd
125 #define gmx_simd4_reduce_d               gmx_simd4_reduce_d_avx_128_fma
126 /* SIMD4 float/double conversion */
127 #define gmx_simd4_cvt_f2d                _mm256_cvtps_pd
128 #define gmx_simd4_cvt_d2f                _mm256_cvtpd_ps
129
130 static gmx_inline double gmx_simdcall
131 gmx_simd4_reduce_d_avx_128_fma(__m256d a)
132 {
133     double  f;
134     __m128d a0, a1;
135     a  = _mm256_hadd_pd(a, a);
136     a0 = _mm256_castpd256_pd128(a);
137     a1 = _mm256_extractf128_pd(a, 0x1);
138     a0 = _mm_add_sd(a0, a1);
139     _mm_store_sd(&f, a0);
140     return f;
141 }
142
143 static gmx_inline double gmx_simdcall
144 gmx_simd4_dotproduct3_d_avx_128_fma(__m256d a, __m256d b)
145 {
146     double  d;
147     __m128d tmp1, tmp2;
148     a    = _mm256_mul_pd(a, b);
149     tmp1 = _mm256_castpd256_pd128(a);
150     tmp2 = _mm256_extractf128_pd(a, 0x1);
151
152     tmp1 = _mm_add_pd(tmp1, _mm_permute_pd(tmp1, _MM_SHUFFLE2(0, 1)));
153     tmp1 = _mm_add_pd(tmp1, tmp2);
154     _mm_store_sd(&d, tmp1);
155     return d;
156 }
157
158 #endif /* GMX_SIMD_IMPL_X86_AVX_128_FMA_H */