96ec6800c42ff9da11a7d4b2ebd9dbb8437f26d4
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / nbnxn_kernel_simd_utils_x86_mic.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,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 #ifndef _nbnxn_kernel_simd_utils_x86_mic_h_
36 #define _nbnxn_kernel_simd_utils_x86_mic_h_
37
38 typedef gmx_simd_int32_t      gmx_exclfilter;
39 static const int filter_stride = GMX_SIMD_INT32_WIDTH/GMX_SIMD_REAL_WIDTH;
40
41 #define PERM_LOW2HIGH _MM_PERM_BABA
42 #define PERM_HIGH2LOW _MM_PERM_DCDC
43
44 #define mask_loh _mm512_int2mask(0x00FF) /* would be better a constant - but can't initialize with a function call. */
45 #define mask_hih _mm512_int2mask(0xFF00)
46
47 /* float/double SIMD register type */
48 typedef __m512 gmx_mm_pr4;
49
50 static gmx_inline gmx_mm_pr4
51 gmx_load_pr4(const real *r)
52 {
53     return _mm512_loadunpacklo_ps(_mm512_undefined_ps(), r);
54 }
55
56 static gmx_inline void
57 gmx_store_pr4(real *dest, gmx_mm_pr4 src)
58 {
59     _mm512_mask_packstorelo_ps(dest, _mm512_int2mask(0xF), src);
60 }
61
62 static gmx_inline gmx_mm_pr4
63 gmx_add_pr4(gmx_mm_pr4 a, gmx_mm_pr4 b)
64 {
65     return _mm512_add_ps(a, b);
66 }
67
68 /* Half-width SIMD real type */
69 typedef __m512 gmx_mm_hpr; /* high half is ignored */
70
71 /* Half-width SIMD operations */
72
73 /* Load reals at half-width aligned pointer b into half-width SIMD register a */
74 static gmx_inline void
75 gmx_load_hpr(gmx_mm_hpr *a, const real *b)
76 {
77     *a = _mm512_loadunpacklo_ps(_mm512_undefined_ps(), b);
78 }
79
80 /* Set all entries in half-width SIMD register *a to b */
81 static gmx_inline void
82 gmx_set1_hpr(gmx_mm_hpr *a, real b)
83 {
84     *a = _mm512_set1_ps(b);
85 }
86
87 /* Load one real at b and one real at b+1 into halves of a, respectively */
88 static gmx_inline void
89 gmx_load1p1_pr(gmx_mm_ps *a, const real *b)
90 {
91
92     *a = _mm512_mask_extload_ps(_mm512_extload_ps(b, _MM_UPCONV_PS_NONE, _MM_BROADCAST_1X16, _MM_HINT_NONE), mask_hih,
93                                 b+1, _MM_UPCONV_PS_NONE, _MM_BROADCAST_1X16, _MM_HINT_NONE);
94 }
95
96 /* Load reals at half-width aligned pointer b into two halves of a */
97 static gmx_inline void
98 gmx_loaddh_pr(gmx_mm_ps *a, const real *b)
99 {
100     *a = _mm512_permute4f128_ps(_mm512_loadunpacklo_ps(_mm512_undefined_ps(), b), PERM_LOW2HIGH);
101 }
102
103 /* Store half-width SIMD register b into half width aligned memory a */
104 static gmx_inline void
105 gmx_store_hpr(real *a, gmx_mm_hpr b)
106 {
107     _mm512_mask_packstorelo_ps(a, mask_loh, b);
108 }
109
110 #define gmx_add_hpr _mm512_add_ps
111 #define gmx_sub_hpr _mm512_sub_ps
112
113 /* Sum over 4 half SIMD registers */
114 static gmx_inline gmx_mm_hpr
115 gmx_sum4_hpr(gmx_mm_ps a, gmx_mm_ps b)
116 {
117     a = _mm512_add_ps(a, b);
118     b = _mm512_permute4f128_ps(a, PERM_HIGH2LOW);
119     return _mm512_add_ps(a, b);
120 }
121
122 /* Sum the elements of halfs of each input register and store sums in out */
123 static gmx_inline gmx_mm_pr4
124 gmx_mm_transpose_sum4h_pr(gmx_mm_ps a, gmx_mm_ps b)
125 {
126     return _mm512_setr4_ps(_mm512_mask_reduce_add_ps(mask_loh, a),
127                            _mm512_mask_reduce_add_ps(mask_hih, a),
128                            _mm512_mask_reduce_add_ps(mask_loh, b),
129                            _mm512_mask_reduce_add_ps(mask_hih, b));
130 }
131
132 static gmx_inline void
133 gmx_pr_to_2hpr(gmx_mm_ps a, gmx_mm_hpr *b, gmx_mm_hpr *c)
134 {
135     *b = a;
136     *c = _mm512_permute4f128_ps(a, PERM_HIGH2LOW);
137 }
138
139 static gmx_inline void
140 gmx_2hpr_to_pr(gmx_mm_hpr a, gmx_mm_hpr b, gmx_mm_ps *c)
141 {
142     *c = _mm512_mask_permute4f128_ps(a, mask_hih, b, PERM_LOW2HIGH);
143 }
144
145 /* recombine the 2 high half into c */
146 static gmx_inline void
147 gmx_2hpr_high_to_pr(gmx_mm_hpr a, gmx_mm_hpr b, gmx_mm_ps *c)
148 {
149     *c = _mm512_mask_permute4f128_ps(b, mask_loh, a, PERM_HIGH2LOW);
150 }
151
152 static gmx_inline void
153 gmx_2hepi_to_epi(gmx_simd_int32_t a, gmx_simd_int32_t b, gmx_simd_int32_t *c)
154 {
155     *c = _mm512_mask_permute4f128_epi32(a, mask_hih, b, PERM_LOW2HIGH);
156 }
157
158 /* recombine the 2 high half into c */
159 static gmx_inline void
160 gmx_2hepi_high_to_epi(gmx_simd_int32_t a, gmx_simd_int32_t b, gmx_simd_int32_t *c)
161 {
162     *c = _mm512_mask_permute4f128_epi32(b, mask_loh, a, PERM_HIGH2LOW);
163 }
164
165 /* Align a stack-based thread-local working array. work-array (currently) not used by load_table_f*/
166 static gmx_inline int *
167 prepare_table_load_buffer(const int *array)
168 {
169     return NULL;
170 }
171
172 /* Using TAB_FDV0 is slower (for non-analytical PME). For the _mm512_i32gather_ps it helps
173    to have the 16 elements in fewer cache lines. This is why it is faster to do F&D toghether
174    and low/high half after each other, then simply doing a gather for tab_coul_F and tab_coul_F+1.
175    The ording of the 16 elements doesn't matter, so it doesn't help to get FD sorted as odd/even
176    instead of low/high.
177  */
178 static gmx_inline void
179 load_table_f(const real *tab_coul_F, gmx_simd_int32_t ti_S, int *ti,
180              gmx_mm_ps *ctab0_S, gmx_mm_ps *ctab1_S)
181 {
182     __m512i idx;
183     __m512i ti1 = _mm512_add_epi32(ti_S, _mm512_set1_epi32(1)); /* incr by 1 for tab1 */
184     gmx_2hepi_to_epi(ti_S, ti1, &idx);
185     __m512  tmp1 = _mm512_i32gather_ps(idx, tab_coul_F, sizeof(float));
186     gmx_2hepi_high_to_epi(ti_S, ti1, &idx);
187     __m512  tmp2 = _mm512_i32gather_ps(idx, tab_coul_F, sizeof(float));
188
189     gmx_2hpr_to_pr(tmp1, tmp2, ctab0_S);
190     gmx_2hpr_high_to_pr(tmp1, tmp2, ctab1_S);
191     *ctab1_S  = gmx_simd_sub_r(*ctab1_S, *ctab0_S);
192 }
193
194 static gmx_inline void
195 load_table_f_v(const real *tab_coul_F, const real *tab_coul_V,
196                gmx_simd_int32_t ti_S, int *ti,
197                gmx_mm_ps *ctab0_S, gmx_mm_ps *ctab1_S,
198                gmx_mm_ps *ctabv_S)
199 {
200     load_table_f(tab_coul_F, ti_S, ti, ctab0_S, ctab1_S);
201     *ctabv_S = _mm512_i32gather_ps(ti_S, tab_coul_V, sizeof(float));
202 }
203
204 static gmx_inline gmx_mm_pr4
205 gmx_mm_transpose_sum4_pr(gmx_mm_ps in0, gmx_mm_ps in1,
206                          gmx_mm_ps in2, gmx_mm_ps in3)
207 {
208     return _mm512_setr4_ps(_mm512_reduce_add_ps(in0),
209                            _mm512_reduce_add_ps(in1),
210                            _mm512_reduce_add_ps(in2),
211                            _mm512_reduce_add_ps(in3));
212 }
213
214 static gmx_inline void
215 load_lj_pair_params2(const real *nbfp0, const real *nbfp1,
216                      const int *type, int aj,
217                      gmx_mm_ps *c6_S, gmx_mm_ps *c12_S)
218 {
219     __m512i idx0, idx1, idx;
220
221     /* load all 8 unaligned requires 2 load. */
222     idx0 = _mm512_loadunpacklo_epi32(_mm512_undefined_epi32(), type+aj);
223     idx0 = _mm512_loadunpackhi_epi32(idx0, type+aj+16);
224
225     idx0 = _mm512_mullo_epi32(idx0, _mm512_set1_epi32(nbfp_stride));
226     idx1 = _mm512_add_epi32(idx0, _mm512_set1_epi32(1)); /* incr by 1 for c12 */
227
228     gmx_2hepi_to_epi(idx0, idx1, &idx);
229     __m512 tmp1 = _mm512_i32gather_ps(idx, nbfp0, sizeof(float));
230     __m512 tmp2 = _mm512_i32gather_ps(idx, nbfp1, sizeof(float));
231
232     gmx_2hpr_to_pr(tmp1, tmp2, c6_S);
233     gmx_2hpr_high_to_pr(tmp1, tmp2, c12_S);
234 }
235
236 #define HAVE_GMX_SUM_SIMD
237 static gmx_inline real
238 gmx_sum_simd(gmx_simd_real_t x, real* b)
239 {
240     return _mm512_reduce_add_ps(x);
241 }
242 static gmx_inline real
243 gmx_sum_simd4(gmx_simd_real_t x, real* b)
244 {
245     return _mm512_mask_reduce_add_ps(_mm512_int2mask(0xF), x);
246 }
247
248 /* Code for handling loading exclusions and converting them into
249    interactions. */
250 #define gmx_load1_exclfilter _mm512_set1_epi32
251 #define gmx_load_exclusion_filter _mm512_load_epi32
252 #define gmx_checkbitmask_pb _mm512_test_epi32_mask
253
254 #endif /* _nbnxn_kernel_simd_utils_ref_h_ */