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