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