Redesigned SIMD module and unit tests.
[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 /* 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_loadunpacklo_ps(_mm512_undefined_ps(), b);
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_permute4f128_ps(_mm512_loadunpacklo_ps(_mm512_undefined_ps(), b), PERM_LOW2HIGH);
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     _mm512_mask_packstorelo_ps(a, mask_loh, b);
87 }
88
89 #define gmx_add_hpr _mm512_add_ps
90 #define gmx_sub_hpr _mm512_sub_ps
91
92 /* Sum over 4 half SIMD registers */
93 static gmx_inline gmx_mm_hpr
94 gmx_sum4_hpr(gmx_simd_float_t a, gmx_simd_float_t b)
95 {
96     a = _mm512_add_ps(a, b);
97     b = _mm512_permute4f128_ps(a, PERM_HIGH2LOW);
98     return _mm512_add_ps(a, b);
99 }
100
101 /* Sum the elements of halfs of each input register and store sums in out */
102 static gmx_inline __m512
103 gmx_mm_transpose_sum4h_pr(gmx_simd_float_t a, gmx_simd_float_t b)
104 {
105     return _mm512_setr4_ps(_mm512_mask_reduce_add_ps(mask_loh, a),
106                            _mm512_mask_reduce_add_ps(mask_hih, a),
107                            _mm512_mask_reduce_add_ps(mask_loh, b),
108                            _mm512_mask_reduce_add_ps(mask_hih, b));
109 }
110
111 static gmx_inline void
112 gmx_pr_to_2hpr(gmx_simd_float_t a, gmx_mm_hpr *b, gmx_mm_hpr *c)
113 {
114     *b = a;
115     *c = _mm512_permute4f128_ps(a, PERM_HIGH2LOW);
116 }
117
118 static gmx_inline void
119 gmx_2hpr_to_pr(gmx_mm_hpr a, gmx_mm_hpr b, gmx_simd_float_t *c)
120 {
121     *c = _mm512_mask_permute4f128_ps(a, mask_hih, b, PERM_LOW2HIGH);
122 }
123
124 /* recombine the 2 high half into c */
125 static gmx_inline void
126 gmx_2hpr_high_to_pr(gmx_mm_hpr a, gmx_mm_hpr b, gmx_simd_float_t *c)
127 {
128     *c = _mm512_mask_permute4f128_ps(b, mask_loh, a, PERM_HIGH2LOW);
129 }
130
131 static gmx_inline void
132 gmx_2hepi_to_epi(gmx_simd_int32_t a, gmx_simd_int32_t b, gmx_simd_int32_t *c)
133 {
134     *c = _mm512_mask_permute4f128_epi32(a, mask_hih, b, PERM_LOW2HIGH);
135 }
136
137 /* recombine the 2 high half into c */
138 static gmx_inline void
139 gmx_2hepi_high_to_epi(gmx_simd_int32_t a, gmx_simd_int32_t b, gmx_simd_int32_t *c)
140 {
141     *c = _mm512_mask_permute4f128_epi32(b, mask_loh, a, PERM_HIGH2LOW);
142 }
143
144 /* Align a stack-based thread-local working array. work-array (currently) not used by load_table_f*/
145 static gmx_inline int *
146 prepare_table_load_buffer(const int *array)
147 {
148     return NULL;
149 }
150
151 /* Using TAB_FDV0 is slower (for non-analytical PME). For the _mm512_i32gather_ps it helps
152    to have the 16 elements in fewer cache lines. This is why it is faster to do F&D toghether
153    and low/high half after each other, then simply doing a gather for tab_coul_F and tab_coul_F+1.
154    The ording of the 16 elements doesn't matter, so it doesn't help to get FD sorted as odd/even
155    instead of low/high.
156  */
157 static gmx_inline void
158 load_table_f(const real *tab_coul_F, gmx_simd_int32_t ti_S, int *ti,
159              gmx_simd_float_t *ctab0_S, gmx_simd_float_t *ctab1_S)
160 {
161     __m512i idx;
162     __m512i ti1 = _mm512_add_epi32(ti_S, _mm512_set1_epi32(1)); /* incr by 1 for tab1 */
163     gmx_2hepi_to_epi(ti_S, ti1, &idx);
164     __m512  tmp1 = _mm512_i32gather_ps(idx, tab_coul_F, sizeof(float));
165     gmx_2hepi_high_to_epi(ti_S, ti1, &idx);
166     __m512  tmp2 = _mm512_i32gather_ps(idx, tab_coul_F, sizeof(float));
167
168     gmx_2hpr_to_pr(tmp1, tmp2, ctab0_S);
169     gmx_2hpr_high_to_pr(tmp1, tmp2, ctab1_S);
170
171     *ctab1_S  = gmx_simd_sub_r(*ctab1_S, *ctab0_S);
172 }
173
174 static gmx_inline void
175 load_table_f_v(const real *tab_coul_F, const real *tab_coul_V,
176                gmx_simd_int32_t ti_S, int *ti,
177                gmx_simd_float_t *ctab0_S, gmx_simd_float_t *ctab1_S,
178                gmx_simd_float_t *ctabv_S)
179 {
180     load_table_f(tab_coul_F, ti_S, ti, ctab0_S, ctab1_S);
181     *ctabv_S = _mm512_i32gather_ps(ti_S, tab_coul_V, sizeof(float));
182 }
183
184 static gmx_inline __m512
185 gmx_mm_transpose_sum4_pr(gmx_simd_float_t in0, gmx_simd_float_t in1,
186                          gmx_simd_float_t in2, gmx_simd_float_t in3)
187 {
188     return _mm512_setr4_ps(_mm512_reduce_add_ps(in0),
189                            _mm512_reduce_add_ps(in1),
190                            _mm512_reduce_add_ps(in2),
191                            _mm512_reduce_add_ps(in3));
192 }
193
194 static gmx_inline void
195 load_lj_pair_params2(const real *nbfp0, const real *nbfp1,
196                      const int *type, int aj,
197                      gmx_simd_float_t *c6_S, gmx_simd_float_t *c12_S)
198 {
199     __m512i idx0, idx1, idx;
200
201     /* load all 8 unaligned requires 2 load. */
202     idx0 = _mm512_loadunpacklo_epi32(_mm512_undefined_epi32(), type+aj);
203     idx0 = _mm512_loadunpackhi_epi32(idx0, type+aj+16);
204
205     idx0 = _mm512_mullo_epi32(idx0, _mm512_set1_epi32(nbfp_stride));
206     idx1 = _mm512_add_epi32(idx0, _mm512_set1_epi32(1)); /* incr by 1 for c12 */
207
208     gmx_2hepi_to_epi(idx0, idx1, &idx);
209     __m512 tmp1 = _mm512_i32gather_ps(idx, nbfp0, sizeof(float));
210     __m512 tmp2 = _mm512_i32gather_ps(idx, nbfp1, sizeof(float));
211
212     gmx_2hpr_to_pr(tmp1, tmp2, c6_S);
213     gmx_2hpr_high_to_pr(tmp1, tmp2, c12_S);
214 }
215
216 /* Code for handling loading exclusions and converting them into
217    interactions. */
218 #define gmx_load1_exclfilter _mm512_set1_epi32
219 #define gmx_load_exclusion_filter _mm512_load_epi32
220 #define gmx_checkbitmask_pb _mm512_test_epi32_mask
221
222 #endif /* _nbnxn_kernel_simd_utils_ref_h_ */