Redesigned SIMD module and unit tests.
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / nbnxn_kernel_simd_utils_ibm_qpx.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_ibm_qpx_h_
36 #define _nbnxn_kernel_simd_utils_ibm_qpx_h_
37
38 typedef gmx_simd_real_t gmx_exclfilter;
39 static const int filter_stride = 1;
40
41 /* This files contains all functions/macros for the SIMD kernels
42  * which have explicit dependencies on the j-cluster size and/or SIMD-width.
43  * The functionality which depends on the j-cluster size is:
44  *   LJ-parameter lookup
45  *   force table lookup
46  *   energy group pair energy storage
47  */
48
49 /* Collect all [0123] elements of the 4 inputs to out[0123], respectively */
50 static gmx_inline void
51 gmx_transpose_4_ps(gmx_simd_real_t a, gmx_simd_real_t b,
52                    gmx_simd_real_t c, gmx_simd_real_t d,
53                    gmx_simd_real_t *out0, gmx_simd_real_t *out1,
54                    gmx_simd_real_t *out2, gmx_simd_real_t *out3)
55 {
56     /* Prepare control vectors for swizzling. In its third input,
57        vec_perm accepts indices into the effective 8-wide SIMD vector
58        created by concatenating its first two inputs. Those indices
59        map data from the input vectors to the output vector.
60
61        vec_gpci() converts an octal literal of the indices into the
62        correct form for vec_perm() to use. That form is an octal digit
63        in bits 0-2 of the mantissa of each double. */
64     gmx_simd_real_t p6420 = vec_gpci(06420);
65     gmx_simd_real_t p7531 = vec_gpci(07531);
66
67     /* Four-way swizzle (i.e. transpose) of vectors a = a0a1a2a3, etc. */
68     gmx_simd_real_t b2b0a2a0 = vec_perm(a, b, p6420);
69     gmx_simd_real_t b3b1a3a1 = vec_perm(a, b, p7531);
70     gmx_simd_real_t d2d0c2c0 = vec_perm(c, d, p6420);
71     gmx_simd_real_t d3d1c3c1 = vec_perm(c, d, p7531);
72     *out0 = vec_perm(d2d0c2c0, b2b0a2a0, p7531);
73     *out1 = vec_perm(d3d1c3c1, b3b1a3a1, p7531);
74     *out2 = vec_perm(d2d0c2c0, b2b0a2a0, p6420);
75     *out3 = vec_perm(d3d1c3c1, b3b1a3a1, p6420);
76 }
77
78 /* Collect element 0 and 1 of the 4 inputs to out0 and out1, respectively */
79 static gmx_inline void
80 gmx_shuffle_4_ps_fil01_to_2_ps(gmx_simd_real_t a, gmx_simd_real_t b,
81                                gmx_simd_real_t c, gmx_simd_real_t d,
82                                gmx_simd_real_t *out0, gmx_simd_real_t *out1)
83 {
84     gmx_simd_real_t p6420 = vec_gpci(06420);
85     gmx_simd_real_t p7531 = vec_gpci(07531);
86
87     /* Partial four-way swizzle of vectors a = a0a1a2a3, etc. */
88     gmx_simd_real_t b2b0a2a0 = vec_perm(a, b, p6420);
89     gmx_simd_real_t b3b1a3a1 = vec_perm(a, b, p7531);
90     gmx_simd_real_t d2d0c2c0 = vec_perm(c, d, p6420);
91     gmx_simd_real_t d3d1c3c1 = vec_perm(c, d, p7531);
92     *out0 = vec_perm(d2d0c2c0, b2b0a2a0, p7531);
93     *out1 = vec_perm(d3d1c3c1, b3b1a3a1, p7531);
94 }
95
96 /* Collect element 2 of the 4 inputs to out */
97 static gmx_inline gmx_simd_real_t
98 gmx_shuffle_4_ps_fil2_to_1_ps(gmx_simd_real_t a, gmx_simd_real_t b,
99                               gmx_simd_real_t c, gmx_simd_real_t d)
100 {
101     gmx_simd_real_t p6420 = vec_gpci(06420);
102
103     /* Partial four-way swizzle of vectors a = a0a1a2a3, etc. */
104     gmx_simd_real_t b2b0a2a0 = vec_perm(a, b, p6420);
105     gmx_simd_real_t d2d0c2c0 = vec_perm(c, d, p6420);
106     return vec_perm(d2d0c2c0, b2b0a2a0, p6420);
107 }
108
109 #ifdef TAB_FDV0
110 /* Align a stack-based thread-local working array. Table loads on QPX
111  * use the array, but most other implementations do not. */
112 static gmx_inline int *
113 prepare_table_load_buffer(int *array)
114 {
115     return gmx_simd_align_i(array);
116 }
117
118 static gmx_inline void
119 load_table_f(const real *tab_coul_FDV0, gmx_simd_int32_t ti_S, int *ti,
120              gmx_simd_real_t *ctab0_S, gmx_simd_real_t *ctab1_S)
121 {
122 #ifdef NDEBUG
123     /* Just like 256-bit AVX, we need to use memory to get indices
124        into integer registers efficiently. */
125     vec_st(ti_S, 0, ti);
126 #else
127     vec_sta(ti_S, 0, ti);
128 #endif
129
130     /* Here we load 4 aligned reals, but we need just 2 elements of each */
131     gmx_simd_real_t a = gmx_simd_load_r(tab_coul_FDV0 + ti[0] * nbfp_stride);
132     gmx_simd_real_t b = gmx_simd_load_r(tab_coul_FDV0 + ti[1] * nbfp_stride);
133     gmx_simd_real_t c = gmx_simd_load_r(tab_coul_FDV0 + ti[2] * nbfp_stride);
134     gmx_simd_real_t d = gmx_simd_load_r(tab_coul_FDV0 + ti[3] * nbfp_stride);
135
136     gmx_shuffle_4_ps_fil01_to_2_ps(a, b, c, d, ctab0_S, ctab1_S);
137 }
138
139 static gmx_inline void
140 load_table_f_v(const real *tab_coul_FDV0,
141                gmx_simd_int32_t ti_S, int *ti,
142                gmx_simd_real_t *ctab0_S, gmx_simd_real_t *ctab1_S,
143                gmx_simd_real_t *ctabv_S)
144 {
145 #ifdef NDEBUG
146     /* Just like 256-bit AVX, we need to use memory to get indices
147        into integer registers efficiently. */
148     vec_st(ti_S, 0, ti);
149 #else
150     vec_sta(ti_S, 0, ti);
151 #endif
152
153     /* Here we load 4 aligned reals, but we need just 3 elements of each. */
154     gmx_simd_real_t a = gmx_simd_load_r(tab_coul_FDV0 + ti[0] * nbfp_stride);
155     gmx_simd_real_t b = gmx_simd_load_r(tab_coul_FDV0 + ti[1] * nbfp_stride);
156     gmx_simd_real_t c = gmx_simd_load_r(tab_coul_FDV0 + ti[2] * nbfp_stride);
157     gmx_simd_real_t d = gmx_simd_load_r(tab_coul_FDV0 + ti[3] * nbfp_stride);
158
159     gmx_shuffle_4_ps_fil01_to_2_ps(a, b, c, d, ctab0_S, ctab1_S);
160     *ctabv_S = gmx_shuffle_4_ps_fil2_to_1_ps(a, b, c, d);
161 }
162 #else
163
164 /* Not required for BlueGene/Q */
165
166 #endif
167
168 /* Sum the elements within each input register and store the sums in out.
169  */
170 static gmx_inline gmx_simd_real_t
171 gmx_mm_transpose_sum4_pr(gmx_simd_real_t a, gmx_simd_real_t b,
172                          gmx_simd_real_t c, gmx_simd_real_t d)
173 {
174     gmx_simd_real_t a0b0c0d0, a1b1c1d1, a2b2c2d2, a3b3c3d3;
175     gmx_transpose_4_ps(a, b, c, d,
176                        &a0b0c0d0,
177                        &a1b1c1d1,
178                        &a2b2c2d2,
179                        &a3b3c3d3);
180     /* Now reduce the transposed vectors */
181     gmx_simd_real_t sum01 = gmx_simd_add_r(a0b0c0d0, a1b1c1d1);
182     gmx_simd_real_t sim23 = gmx_simd_add_r(a2b2c2d2, a3b3c3d3);
183     return gmx_simd_add_r(sum01, sim23);
184 }
185
186 static gmx_inline void
187 load_lj_pair_params(const real *nbfp, const int *type, int aj,
188                     gmx_simd_real_t *c6_S, gmx_simd_real_t *c12_S)
189 {
190     /* Here we load 4 aligned reals, but we need just 2 elemnts of each. */
191     gmx_simd_real_t a = gmx_simd_load_r(nbfp + type[aj+0] * nbfp_stride);
192     gmx_simd_real_t b = gmx_simd_load_r(nbfp + type[aj+1] * nbfp_stride);
193     gmx_simd_real_t c = gmx_simd_load_r(nbfp + type[aj+2] * nbfp_stride);
194     gmx_simd_real_t d = gmx_simd_load_r(nbfp + type[aj+3] * nbfp_stride);
195
196     gmx_shuffle_4_ps_fil01_to_2_ps(a, b, c, d, c6_S, c12_S);
197 }
198
199 /* Define USE_FUNCTIONS_FOR_QPX to get the static inline functions
200  * that seem to exhaust xlC 12.1 during kernel compilation */
201
202 #ifndef USE_FUNCTIONS_FOR_QPX
203
204 #define gmx_load_exclusion_filter(a) vec_ldia(0, (int *) a)
205 #define gmx_load_interaction_mask_pb(a, b) vec_ld(a, (real *) b)
206
207 #else /* USE_FUNCTIONS_FOR_QPX */
208
209 static gmx_inline gmx_exclfilter gmx_load_exclusion_filter(const unsigned *a)
210 {
211 #ifdef NDEBUG
212     return vec_ldia(0, (int *) a);
213 #else
214     return vec_ldiaa(0, (int *) a);
215 #endif
216 }
217
218 /* Code for handling loading and applying exclusion masks. Note that
219    parameter a is not treated like an array index; it is naively added
220    to b, so should be in bytes. */
221 static gmx_inline gmx_simd_bool_t gmx_load_interaction_mask_pb(long a, const real *b)
222 {
223 #ifdef NDEBUG
224     return vec_ld(a, (real *) b);
225 #else
226     return vec_lda(a, (real *) b);
227 #endif
228 }
229
230 #endif /* USE_FUNCTIONS_FOR_QPX */
231
232 #endif /* _nbnxn_kernel_simd_utils_ibm_qpx_h_ */