nbnxn utils performance improvement for Phi
[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,2015, 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 /* Half-width SIMD real type */
48 typedef __m512 gmx_mm_hpr; /* high half is ignored */
49
50 /* Half-width SIMD operations */
51
52 /* Load reals at half-width aligned pointer b into half-width SIMD register a */
53 static gmx_inline void
54 gmx_load_hpr(gmx_mm_hpr *a, const real *b)
55 {
56     *a = _mm512_castpd_ps(_mm512_extload_pd((const double*)b, _MM_UPCONV_PD_NONE, _MM_BROADCAST_4X8, _MM_HINT_NONE));
57 }
58
59 /* Set all entries in half-width SIMD register *a to b */
60 static gmx_inline void
61 gmx_set1_hpr(gmx_mm_hpr *a, real b)
62 {
63     *a = _mm512_set1_ps(b);
64 }
65
66 /* Load one real at b and one real at b+1 into halves of a, respectively */
67 static gmx_inline void
68 gmx_load1p1_pr(gmx_simd_float_t *a, const real *b)
69 {
70
71     *a = _mm512_mask_extload_ps(_mm512_extload_ps(b, _MM_UPCONV_PS_NONE, _MM_BROADCAST_1X16, _MM_HINT_NONE), mask_hih,
72                                 b+1, _MM_UPCONV_PS_NONE, _MM_BROADCAST_1X16, _MM_HINT_NONE);
73 }
74
75 /* Load reals at half-width aligned pointer b into two halves of a */
76 static gmx_inline void
77 gmx_loaddh_pr(gmx_simd_float_t *a, const real *b)
78 {
79     *a = _mm512_castpd_ps(_mm512_extload_pd((const double*)b, _MM_UPCONV_PD_NONE, _MM_BROADCAST_4X8, _MM_HINT_NONE));
80 }
81
82 /* Store half-width SIMD register b into half width aligned memory a */
83 static gmx_inline void
84 gmx_store_hpr(real *a, gmx_mm_hpr b)
85 {
86     assert((size_t)a%32 == 0);
87     _mm512_mask_packstorelo_ps(a, mask_loh, b);
88 }
89
90 #define gmx_add_hpr _mm512_add_ps
91 #define gmx_sub_hpr _mm512_sub_ps
92
93 /* Sum over 4 half SIMD registers */
94 static gmx_inline gmx_mm_hpr
95 gmx_sum4_hpr(gmx_simd_float_t a, gmx_simd_float_t b)
96 {
97     a = _mm512_add_ps(a, b);
98     b = _mm512_permute4f128_ps(a, PERM_HIGH2LOW);
99     return _mm512_add_ps(a, b);
100 }
101
102 /* Sum the elements of halfs of each input register and store sums in out */
103 static gmx_inline __m512
104 gmx_mm_transpose_sum4h_pr(gmx_simd_float_t a, gmx_simd_float_t b)
105 {
106     a = _mm512_add_ps(a, _mm512_swizzle_ps(a, _MM_SWIZ_REG_CDAB));
107     a = _mm512_add_ps(a, _mm512_swizzle_ps(a, _MM_SWIZ_REG_BADC));
108     a = _mm512_add_ps(a, _mm512_permute4f128_ps(a, _MM_PERM_CDAB));
109
110     b = _mm512_add_ps(b, _mm512_swizzle_ps(b, _MM_SWIZ_REG_CDAB));
111     b = _mm512_add_ps(b, _mm512_swizzle_ps(b, _MM_SWIZ_REG_BADC));
112     a = _mm512_mask_add_ps(a, _mm512_int2mask(0xF0F0), b, _mm512_permute4f128_ps(b, _MM_PERM_CDAB));
113
114     return _mm512_castsi512_ps(_mm512_permutevar_epi32(_mm512_setr4_epi32(0, 8, 4, 12), _mm512_castps_si512(a)));
115 }
116
117 static gmx_inline void
118 gmx_pr_to_2hpr(gmx_simd_float_t a, gmx_mm_hpr *b, gmx_mm_hpr *c)
119 {
120     *b = a;
121     *c = _mm512_permute4f128_ps(a, PERM_HIGH2LOW);
122 }
123
124 static gmx_inline void
125 gmx_2hpr_to_pr(gmx_mm_hpr a, gmx_mm_hpr b, gmx_simd_float_t *c)
126 {
127     *c = _mm512_mask_permute4f128_ps(a, mask_hih, b, PERM_LOW2HIGH);
128 }
129
130 /* recombine the 2 high half into c */
131 static gmx_inline void
132 gmx_2hpr_high_to_pr(gmx_mm_hpr a, gmx_mm_hpr b, gmx_simd_float_t *c)
133 {
134     *c = _mm512_mask_permute4f128_ps(b, mask_loh, a, PERM_HIGH2LOW);
135 }
136
137 static gmx_inline void
138 gmx_2hepi_to_epi(gmx_simd_int32_t a, gmx_simd_int32_t b, gmx_simd_int32_t *c)
139 {
140     *c = _mm512_mask_permute4f128_epi32(a, mask_hih, b, PERM_LOW2HIGH);
141 }
142
143 /* recombine the 2 high half into c */
144 static gmx_inline void
145 gmx_2hepi_high_to_epi(gmx_simd_int32_t a, gmx_simd_int32_t b, gmx_simd_int32_t *c)
146 {
147     *c = _mm512_mask_permute4f128_epi32(b, mask_loh, a, PERM_HIGH2LOW);
148 }
149
150 /* Align a stack-based thread-local working array. work-array (currently) not used by load_table_f*/
151 static gmx_inline int *
152 prepare_table_load_buffer(const int *array)
153 {
154     return NULL;
155 }
156
157 /* Using TAB_FDV0 is slower (for non-analytical PME). For the _mm512_i32gather_ps it helps
158    to have the 16 elements in fewer cache lines. This is why it is faster to do F&D toghether
159    and low/high half after each other, then simply doing a gather for tab_coul_F and tab_coul_F+1.
160    The ording of the 16 elements doesn't matter, so it doesn't help to get FD sorted as odd/even
161    instead of low/high.
162  */
163 static gmx_inline void
164 load_table_f(const real *tab_coul_F, gmx_simd_int32_t ti_S, int *ti,
165              gmx_simd_float_t *ctab0_S, gmx_simd_float_t *ctab1_S)
166 {
167     __m512i idx;
168     __m512i ti1 = _mm512_add_epi32(ti_S, _mm512_set1_epi32(1)); /* incr by 1 for tab1 */
169     gmx_2hepi_to_epi(ti_S, ti1, &idx);
170     __m512  tmp1 = _mm512_i32gather_ps(idx, tab_coul_F, sizeof(float));
171     gmx_2hepi_high_to_epi(ti_S, ti1, &idx);
172     __m512  tmp2 = _mm512_i32gather_ps(idx, tab_coul_F, sizeof(float));
173
174     gmx_2hpr_to_pr(tmp1, tmp2, ctab0_S);
175     gmx_2hpr_high_to_pr(tmp1, tmp2, ctab1_S);
176
177     *ctab1_S  = gmx_simd_sub_r(*ctab1_S, *ctab0_S);
178 }
179
180 static gmx_inline void
181 load_table_f_v(const real *tab_coul_F, const real *tab_coul_V,
182                gmx_simd_int32_t ti_S, int *ti,
183                gmx_simd_float_t *ctab0_S, gmx_simd_float_t *ctab1_S,
184                gmx_simd_float_t *ctabv_S)
185 {
186     load_table_f(tab_coul_F, ti_S, ti, ctab0_S, ctab1_S);
187     *ctabv_S = _mm512_i32gather_ps(ti_S, tab_coul_V, sizeof(float));
188 }
189
190 static gmx_inline __m512
191 gmx_mm_transpose_sum4_pr(gmx_simd_float_t in0, gmx_simd_float_t in1,
192                          gmx_simd_float_t in2, gmx_simd_float_t in3)
193 {
194     return _mm512_setr4_ps(_mm512_reduce_add_ps(in0),
195                            _mm512_reduce_add_ps(in1),
196                            _mm512_reduce_add_ps(in2),
197                            _mm512_reduce_add_ps(in3));
198 }
199
200 static gmx_inline void
201 load_lj_pair_params2(const real *nbfp0, const real *nbfp1,
202                      const int *type, int aj,
203                      gmx_simd_float_t *c6_S, gmx_simd_float_t *c12_S)
204 {
205     __m512i idx0, idx1, idx;
206
207     /* load all 8 unaligned requires 2 load. */
208     idx0 = _mm512_loadunpacklo_epi32(_mm512_undefined_epi32(), type+aj);
209     idx0 = _mm512_loadunpackhi_epi32(idx0, type+aj+16);
210
211     idx0 = _mm512_mullo_epi32(idx0, _mm512_set1_epi32(nbfp_stride));
212     idx1 = _mm512_add_epi32(idx0, _mm512_set1_epi32(1)); /* incr by 1 for c12 */
213
214     gmx_2hepi_to_epi(idx0, idx1, &idx);
215     __m512 tmp1 = _mm512_i32gather_ps(idx, nbfp0, sizeof(float));
216     __m512 tmp2 = _mm512_i32gather_ps(idx, nbfp1, sizeof(float));
217
218     gmx_2hpr_to_pr(tmp1, tmp2, c6_S);
219     gmx_2hpr_high_to_pr(tmp1, tmp2, c12_S);
220 }
221
222 /* Code for handling loading exclusions and converting them into
223    interactions. */
224 #define gmx_load1_exclfilter _mm512_set1_epi32
225 #define gmx_load_exclusion_filter _mm512_load_epi32
226 #define gmx_checkbitmask_pb _mm512_test_epi32_mask
227
228 #endif /* _nbnxn_kernel_simd_utils_ref_h_ */