Remove all unnecessary HAVE_CONFIG_H
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_search_simd_4xn.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,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
36
37 #if GMX_SIMD_REAL_WIDTH >= NBNXN_CPU_CLUSTER_I_SIZE
38 #define STRIDE_S  (GMX_SIMD_REAL_WIDTH)
39 #else
40 #define STRIDE_S  NBNXN_CPU_CLUSTER_I_SIZE
41 #endif
42
43 /* Copies PBC shifted i-cell packed atom coordinates to working array */
44 static gmx_inline void
45 icell_set_x_simd_4xn(int ci,
46                      real shx, real shy, real shz,
47                      int gmx_unused na_c,
48                      int gmx_unused stride, const real *x,
49                      nbnxn_list_work_t *work)
50 {
51     int                    ia;
52     nbnxn_x_ci_simd_4xn_t *x_ci;
53
54     x_ci = work->x_ci_simd_4xn;
55
56     ia = X_IND_CI_SIMD_4XN(ci);
57
58     x_ci->ix_S0 = gmx_simd_set1_r(x[ia + 0*STRIDE_S    ] + shx);
59     x_ci->iy_S0 = gmx_simd_set1_r(x[ia + 1*STRIDE_S    ] + shy);
60     x_ci->iz_S0 = gmx_simd_set1_r(x[ia + 2*STRIDE_S    ] + shz);
61     x_ci->ix_S1 = gmx_simd_set1_r(x[ia + 0*STRIDE_S + 1] + shx);
62     x_ci->iy_S1 = gmx_simd_set1_r(x[ia + 1*STRIDE_S + 1] + shy);
63     x_ci->iz_S1 = gmx_simd_set1_r(x[ia + 2*STRIDE_S + 1] + shz);
64     x_ci->ix_S2 = gmx_simd_set1_r(x[ia + 0*STRIDE_S + 2] + shx);
65     x_ci->iy_S2 = gmx_simd_set1_r(x[ia + 1*STRIDE_S + 2] + shy);
66     x_ci->iz_S2 = gmx_simd_set1_r(x[ia + 2*STRIDE_S + 2] + shz);
67     x_ci->ix_S3 = gmx_simd_set1_r(x[ia + 0*STRIDE_S + 3] + shx);
68     x_ci->iy_S3 = gmx_simd_set1_r(x[ia + 1*STRIDE_S + 3] + shy);
69     x_ci->iz_S3 = gmx_simd_set1_r(x[ia + 2*STRIDE_S + 3] + shz);
70 }
71
72 /* SIMD code for making a pair list of cell ci vs cell cjf-cjl
73  * for coordinates in packed format.
74  * Checks bouding box distances and possibly atom pair distances.
75  * This is an accelerated version of make_cluster_list_simple.
76  */
77 static gmx_inline void
78 make_cluster_list_simd_4xn(const nbnxn_grid_t *gridj,
79                            nbnxn_pairlist_t *nbl,
80                            int ci, int cjf, int cjl,
81                            gmx_bool remove_sub_diag,
82                            const real *x_j,
83                            real rl2, float rbb2,
84                            int *ndistc)
85 {
86     const nbnxn_x_ci_simd_4xn_t       *work;
87     const nbnxn_bb_t                  *bb_ci;
88
89     gmx_simd_real_t                    jx_S, jy_S, jz_S;
90
91     gmx_simd_real_t                    dx_S0, dy_S0, dz_S0;
92     gmx_simd_real_t                    dx_S1, dy_S1, dz_S1;
93     gmx_simd_real_t                    dx_S2, dy_S2, dz_S2;
94     gmx_simd_real_t                    dx_S3, dy_S3, dz_S3;
95
96     gmx_simd_real_t                    rsq_S0;
97     gmx_simd_real_t                    rsq_S1;
98     gmx_simd_real_t                    rsq_S2;
99     gmx_simd_real_t                    rsq_S3;
100
101     gmx_simd_bool_t                    wco_S0;
102     gmx_simd_bool_t                    wco_S1;
103     gmx_simd_bool_t                    wco_S2;
104     gmx_simd_bool_t                    wco_S3;
105     gmx_simd_bool_t                    wco_any_S01, wco_any_S23, wco_any_S;
106
107     gmx_simd_real_t                    rc2_S;
108
109     gmx_bool                           InRange;
110     float                              d2;
111     int                                xind_f, xind_l, cj;
112
113     /* cppcheck-suppress selfAssignment . selfAssignment for width 4.*/
114     cjf = CI_TO_CJ_SIMD_4XN(cjf);
115     cjl = CI_TO_CJ_SIMD_4XN(cjl+1) - 1;
116
117     work = nbl->work->x_ci_simd_4xn;
118
119     bb_ci = nbl->work->bb_ci;
120
121     rc2_S   = gmx_simd_set1_r(rl2);
122
123     InRange = FALSE;
124     while (!InRange && cjf <= cjl)
125     {
126 #ifdef NBNXN_SEARCH_BB_SIMD4
127         d2 = subc_bb_dist2_simd4(0, bb_ci, cjf, gridj->bbj);
128 #else
129         d2 = subc_bb_dist2(0, bb_ci, cjf, gridj->bbj);
130 #endif
131         *ndistc += 2;
132
133         /* Check if the distance is within the distance where
134          * we use only the bounding box distance rbb,
135          * or within the cut-off and there is at least one atom pair
136          * within the cut-off.
137          */
138         if (d2 < rbb2)
139         {
140             InRange = TRUE;
141         }
142         else if (d2 < rl2)
143         {
144             xind_f  = X_IND_CJ_SIMD_4XN(CI_TO_CJ_SIMD_4XN(gridj->cell0) + cjf);
145
146             jx_S  = gmx_simd_load_r(x_j+xind_f+0*STRIDE_S);
147             jy_S  = gmx_simd_load_r(x_j+xind_f+1*STRIDE_S);
148             jz_S  = gmx_simd_load_r(x_j+xind_f+2*STRIDE_S);
149
150
151             /* Calculate distance */
152             dx_S0            = gmx_simd_sub_r(work->ix_S0, jx_S);
153             dy_S0            = gmx_simd_sub_r(work->iy_S0, jy_S);
154             dz_S0            = gmx_simd_sub_r(work->iz_S0, jz_S);
155             dx_S1            = gmx_simd_sub_r(work->ix_S1, jx_S);
156             dy_S1            = gmx_simd_sub_r(work->iy_S1, jy_S);
157             dz_S1            = gmx_simd_sub_r(work->iz_S1, jz_S);
158             dx_S2            = gmx_simd_sub_r(work->ix_S2, jx_S);
159             dy_S2            = gmx_simd_sub_r(work->iy_S2, jy_S);
160             dz_S2            = gmx_simd_sub_r(work->iz_S2, jz_S);
161             dx_S3            = gmx_simd_sub_r(work->ix_S3, jx_S);
162             dy_S3            = gmx_simd_sub_r(work->iy_S3, jy_S);
163             dz_S3            = gmx_simd_sub_r(work->iz_S3, jz_S);
164
165             /* rsq = dx*dx+dy*dy+dz*dz */
166             rsq_S0           = gmx_simd_calc_rsq_r(dx_S0, dy_S0, dz_S0);
167             rsq_S1           = gmx_simd_calc_rsq_r(dx_S1, dy_S1, dz_S1);
168             rsq_S2           = gmx_simd_calc_rsq_r(dx_S2, dy_S2, dz_S2);
169             rsq_S3           = gmx_simd_calc_rsq_r(dx_S3, dy_S3, dz_S3);
170
171             wco_S0           = gmx_simd_cmplt_r(rsq_S0, rc2_S);
172             wco_S1           = gmx_simd_cmplt_r(rsq_S1, rc2_S);
173             wco_S2           = gmx_simd_cmplt_r(rsq_S2, rc2_S);
174             wco_S3           = gmx_simd_cmplt_r(rsq_S3, rc2_S);
175
176             wco_any_S01      = gmx_simd_or_b(wco_S0, wco_S1);
177             wco_any_S23      = gmx_simd_or_b(wco_S2, wco_S3);
178             wco_any_S        = gmx_simd_or_b(wco_any_S01, wco_any_S23);
179
180             InRange          = gmx_simd_anytrue_b(wco_any_S);
181
182             *ndistc += 4*GMX_SIMD_REAL_WIDTH;
183         }
184         if (!InRange)
185         {
186             cjf++;
187         }
188     }
189     if (!InRange)
190     {
191         return;
192     }
193
194     InRange = FALSE;
195     while (!InRange && cjl > cjf)
196     {
197 #ifdef NBNXN_SEARCH_BB_SIMD4
198         d2 = subc_bb_dist2_simd4(0, bb_ci, cjl, gridj->bbj);
199 #else
200         d2 = subc_bb_dist2(0, bb_ci, cjl, gridj->bbj);
201 #endif
202         *ndistc += 2;
203
204         /* Check if the distance is within the distance where
205          * we use only the bounding box distance rbb,
206          * or within the cut-off and there is at least one atom pair
207          * within the cut-off.
208          */
209         if (d2 < rbb2)
210         {
211             InRange = TRUE;
212         }
213         else if (d2 < rl2)
214         {
215             xind_l  = X_IND_CJ_SIMD_4XN(CI_TO_CJ_SIMD_4XN(gridj->cell0) + cjl);
216
217             jx_S  = gmx_simd_load_r(x_j+xind_l+0*STRIDE_S);
218             jy_S  = gmx_simd_load_r(x_j+xind_l+1*STRIDE_S);
219             jz_S  = gmx_simd_load_r(x_j+xind_l+2*STRIDE_S);
220
221             /* Calculate distance */
222             dx_S0            = gmx_simd_sub_r(work->ix_S0, jx_S);
223             dy_S0            = gmx_simd_sub_r(work->iy_S0, jy_S);
224             dz_S0            = gmx_simd_sub_r(work->iz_S0, jz_S);
225             dx_S1            = gmx_simd_sub_r(work->ix_S1, jx_S);
226             dy_S1            = gmx_simd_sub_r(work->iy_S1, jy_S);
227             dz_S1            = gmx_simd_sub_r(work->iz_S1, jz_S);
228             dx_S2            = gmx_simd_sub_r(work->ix_S2, jx_S);
229             dy_S2            = gmx_simd_sub_r(work->iy_S2, jy_S);
230             dz_S2            = gmx_simd_sub_r(work->iz_S2, jz_S);
231             dx_S3            = gmx_simd_sub_r(work->ix_S3, jx_S);
232             dy_S3            = gmx_simd_sub_r(work->iy_S3, jy_S);
233             dz_S3            = gmx_simd_sub_r(work->iz_S3, jz_S);
234
235             /* rsq = dx*dx+dy*dy+dz*dz */
236             rsq_S0           = gmx_simd_calc_rsq_r(dx_S0, dy_S0, dz_S0);
237             rsq_S1           = gmx_simd_calc_rsq_r(dx_S1, dy_S1, dz_S1);
238             rsq_S2           = gmx_simd_calc_rsq_r(dx_S2, dy_S2, dz_S2);
239             rsq_S3           = gmx_simd_calc_rsq_r(dx_S3, dy_S3, dz_S3);
240
241             wco_S0           = gmx_simd_cmplt_r(rsq_S0, rc2_S);
242             wco_S1           = gmx_simd_cmplt_r(rsq_S1, rc2_S);
243             wco_S2           = gmx_simd_cmplt_r(rsq_S2, rc2_S);
244             wco_S3           = gmx_simd_cmplt_r(rsq_S3, rc2_S);
245
246             wco_any_S01      = gmx_simd_or_b(wco_S0, wco_S1);
247             wco_any_S23      = gmx_simd_or_b(wco_S2, wco_S3);
248             wco_any_S        = gmx_simd_or_b(wco_any_S01, wco_any_S23);
249
250             InRange          = gmx_simd_anytrue_b(wco_any_S);
251
252             *ndistc += 4*GMX_SIMD_REAL_WIDTH;
253         }
254         if (!InRange)
255         {
256             cjl--;
257         }
258     }
259
260     if (cjf <= cjl)
261     {
262         for (cj = cjf; cj <= cjl; cj++)
263         {
264             /* Store cj and the interaction mask */
265             nbl->cj[nbl->ncj].cj   = CI_TO_CJ_SIMD_4XN(gridj->cell0) + cj;
266             nbl->cj[nbl->ncj].excl = get_imask_simd_4xn(remove_sub_diag, ci, cj);
267 #ifdef GMX_SIMD_IBM_QPX
268             nbl->cj[nbl->ncj].interaction_mask_indices[0] = (nbl->cj[nbl->ncj].excl & 0x000F) >> (0 * 4);
269             nbl->cj[nbl->ncj].interaction_mask_indices[1] = (nbl->cj[nbl->ncj].excl & 0x00F0) >> (1 * 4);
270             nbl->cj[nbl->ncj].interaction_mask_indices[2] = (nbl->cj[nbl->ncj].excl & 0x0F00) >> (2 * 4);
271             nbl->cj[nbl->ncj].interaction_mask_indices[3] = (nbl->cj[nbl->ncj].excl & 0xF000) >> (3 * 4);
272 #endif
273             nbl->ncj++;
274         }
275         /* Increase the closing index in i super-cell list */
276         nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
277     }
278 }
279
280 #undef STRIDE_S