Reformat existing LGPL copyright notices.
[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) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2012, The GROMACS development team.
6  * Copyright (c) 2012, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38
39 #if GMX_SIMD_WIDTH_HERE >= NBNXN_CPU_CLUSTER_I_SIZE
40 #define STRIDE_S  (GMX_SIMD_WIDTH_HERE)
41 #else
42 #define STRIDE_S  NBNXN_CPU_CLUSTER_I_SIZE
43 #endif
44
45 /* Copies PBC shifted i-cell packed atom coordinates to working array */
46 static gmx_inline void
47 icell_set_x_simd_4xn(int ci,
48                      real shx, real shy, real shz,
49                      int gmx_unused na_c,
50                      int gmx_unused stride, const real *x,
51                      nbnxn_list_work_t *work)
52 {
53     int                    ia;
54     nbnxn_x_ci_simd_4xn_t *x_ci;
55
56     x_ci = work->x_ci_simd_4xn;
57
58     ia = X_IND_CI_SIMD_4XN(ci);
59
60     x_ci->ix_S0 = gmx_set1_pr(x[ia + 0*STRIDE_S    ] + shx);
61     x_ci->iy_S0 = gmx_set1_pr(x[ia + 1*STRIDE_S    ] + shy);
62     x_ci->iz_S0 = gmx_set1_pr(x[ia + 2*STRIDE_S    ] + shz);
63     x_ci->ix_S1 = gmx_set1_pr(x[ia + 0*STRIDE_S + 1] + shx);
64     x_ci->iy_S1 = gmx_set1_pr(x[ia + 1*STRIDE_S + 1] + shy);
65     x_ci->iz_S1 = gmx_set1_pr(x[ia + 2*STRIDE_S + 1] + shz);
66     x_ci->ix_S2 = gmx_set1_pr(x[ia + 0*STRIDE_S + 2] + shx);
67     x_ci->iy_S2 = gmx_set1_pr(x[ia + 1*STRIDE_S + 2] + shy);
68     x_ci->iz_S2 = gmx_set1_pr(x[ia + 2*STRIDE_S + 2] + shz);
69     x_ci->ix_S3 = gmx_set1_pr(x[ia + 0*STRIDE_S + 3] + shx);
70     x_ci->iy_S3 = gmx_set1_pr(x[ia + 1*STRIDE_S + 3] + shy);
71     x_ci->iz_S3 = gmx_set1_pr(x[ia + 2*STRIDE_S + 3] + shz);
72 }
73
74 /* SIMD code for making a pair list of cell ci vs cell cjf-cjl
75  * for coordinates in packed format.
76  * Checks bouding box distances and possibly atom pair distances.
77  * This is an accelerated version of make_cluster_list_simple.
78  */
79 static gmx_inline void
80 make_cluster_list_simd_4xn(const nbnxn_grid_t *gridj,
81                            nbnxn_pairlist_t *nbl,
82                            int ci, int cjf, int cjl,
83                            gmx_bool remove_sub_diag,
84                            const real *x_j,
85                            real rl2, float rbb2,
86                            int *ndistc)
87 {
88     const nbnxn_x_ci_simd_4xn_t *work;
89     const nbnxn_bb_t            *bb_ci;
90
91     gmx_mm_pr                    jx_S, jy_S, jz_S;
92
93     gmx_mm_pr                    dx_S0, dy_S0, dz_S0;
94     gmx_mm_pr                    dx_S1, dy_S1, dz_S1;
95     gmx_mm_pr                    dx_S2, dy_S2, dz_S2;
96     gmx_mm_pr                    dx_S3, dy_S3, dz_S3;
97
98     gmx_mm_pr                    rsq_S0;
99     gmx_mm_pr                    rsq_S1;
100     gmx_mm_pr                    rsq_S2;
101     gmx_mm_pr                    rsq_S3;
102
103     gmx_mm_pb                    wco_S0;
104     gmx_mm_pb                    wco_S1;
105     gmx_mm_pb                    wco_S2;
106     gmx_mm_pb                    wco_S3;
107     gmx_mm_pb                    wco_any_S01, wco_any_S23, wco_any_S;
108
109     gmx_mm_pr                    rc2_S;
110
111     gmx_bool                     InRange;
112     float                        d2;
113     int                          xind_f, xind_l, cj;
114
115     cjf = CI_TO_CJ_SIMD_4XN(cjf);
116     cjl = CI_TO_CJ_SIMD_4XN(cjl+1) - 1;
117
118     work = nbl->work->x_ci_simd_4xn;
119
120     bb_ci = nbl->work->bb_ci;
121
122     rc2_S   = gmx_set1_pr(rl2);
123
124     InRange = FALSE;
125     while (!InRange && cjf <= cjl)
126     {
127 #ifdef NBNXN_SEARCH_BB_SIMD4
128         d2 = subc_bb_dist2_simd4(0, bb_ci, cjf, gridj->bbj);
129 #else
130         d2 = subc_bb_dist2(0, bb_ci, cjf, gridj->bbj);
131 #endif
132         *ndistc += 2;
133
134         /* Check if the distance is within the distance where
135          * we use only the bounding box distance rbb,
136          * or within the cut-off and there is at least one atom pair
137          * within the cut-off.
138          */
139         if (d2 < rbb2)
140         {
141             InRange = TRUE;
142         }
143         else if (d2 < rl2)
144         {
145             xind_f  = X_IND_CJ_SIMD_4XN(CI_TO_CJ_SIMD_4XN(gridj->cell0) + cjf);
146
147             jx_S  = gmx_load_pr(x_j+xind_f+0*STRIDE_S);
148             jy_S  = gmx_load_pr(x_j+xind_f+1*STRIDE_S);
149             jz_S  = gmx_load_pr(x_j+xind_f+2*STRIDE_S);
150
151
152             /* Calculate distance */
153             dx_S0            = gmx_sub_pr(work->ix_S0, jx_S);
154             dy_S0            = gmx_sub_pr(work->iy_S0, jy_S);
155             dz_S0            = gmx_sub_pr(work->iz_S0, jz_S);
156             dx_S1            = gmx_sub_pr(work->ix_S1, jx_S);
157             dy_S1            = gmx_sub_pr(work->iy_S1, jy_S);
158             dz_S1            = gmx_sub_pr(work->iz_S1, jz_S);
159             dx_S2            = gmx_sub_pr(work->ix_S2, jx_S);
160             dy_S2            = gmx_sub_pr(work->iy_S2, jy_S);
161             dz_S2            = gmx_sub_pr(work->iz_S2, jz_S);
162             dx_S3            = gmx_sub_pr(work->ix_S3, jx_S);
163             dy_S3            = gmx_sub_pr(work->iy_S3, jy_S);
164             dz_S3            = gmx_sub_pr(work->iz_S3, jz_S);
165
166             /* rsq = dx*dx+dy*dy+dz*dz */
167             rsq_S0           = gmx_calc_rsq_pr(dx_S0, dy_S0, dz_S0);
168             rsq_S1           = gmx_calc_rsq_pr(dx_S1, dy_S1, dz_S1);
169             rsq_S2           = gmx_calc_rsq_pr(dx_S2, dy_S2, dz_S2);
170             rsq_S3           = gmx_calc_rsq_pr(dx_S3, dy_S3, dz_S3);
171
172             wco_S0           = gmx_cmplt_pr(rsq_S0, rc2_S);
173             wco_S1           = gmx_cmplt_pr(rsq_S1, rc2_S);
174             wco_S2           = gmx_cmplt_pr(rsq_S2, rc2_S);
175             wco_S3           = gmx_cmplt_pr(rsq_S3, rc2_S);
176
177             wco_any_S01      = gmx_or_pb(wco_S0, wco_S1);
178             wco_any_S23      = gmx_or_pb(wco_S2, wco_S3);
179             wco_any_S        = gmx_or_pb(wco_any_S01, wco_any_S23);
180
181             InRange          = gmx_anytrue_pb(wco_any_S);
182
183             *ndistc += 4*GMX_SIMD_WIDTH_HERE;
184         }
185         if (!InRange)
186         {
187             cjf++;
188         }
189     }
190     if (!InRange)
191     {
192         return;
193     }
194
195     InRange = FALSE;
196     while (!InRange && cjl > cjf)
197     {
198 #ifdef NBNXN_SEARCH_BB_SIMD4
199         d2 = subc_bb_dist2_simd4(0, bb_ci, cjl, gridj->bbj);
200 #else
201         d2 = subc_bb_dist2(0, bb_ci, cjl, gridj->bbj);
202 #endif
203         *ndistc += 2;
204
205         /* Check if the distance is within the distance where
206          * we use only the bounding box distance rbb,
207          * or within the cut-off and there is at least one atom pair
208          * within the cut-off.
209          */
210         if (d2 < rbb2)
211         {
212             InRange = TRUE;
213         }
214         else if (d2 < rl2)
215         {
216             xind_l  = X_IND_CJ_SIMD_4XN(CI_TO_CJ_SIMD_4XN(gridj->cell0) + cjl);
217
218             jx_S  = gmx_load_pr(x_j+xind_l+0*STRIDE_S);
219             jy_S  = gmx_load_pr(x_j+xind_l+1*STRIDE_S);
220             jz_S  = gmx_load_pr(x_j+xind_l+2*STRIDE_S);
221
222             /* Calculate distance */
223             dx_S0            = gmx_sub_pr(work->ix_S0, jx_S);
224             dy_S0            = gmx_sub_pr(work->iy_S0, jy_S);
225             dz_S0            = gmx_sub_pr(work->iz_S0, jz_S);
226             dx_S1            = gmx_sub_pr(work->ix_S1, jx_S);
227             dy_S1            = gmx_sub_pr(work->iy_S1, jy_S);
228             dz_S1            = gmx_sub_pr(work->iz_S1, jz_S);
229             dx_S2            = gmx_sub_pr(work->ix_S2, jx_S);
230             dy_S2            = gmx_sub_pr(work->iy_S2, jy_S);
231             dz_S2            = gmx_sub_pr(work->iz_S2, jz_S);
232             dx_S3            = gmx_sub_pr(work->ix_S3, jx_S);
233             dy_S3            = gmx_sub_pr(work->iy_S3, jy_S);
234             dz_S3            = gmx_sub_pr(work->iz_S3, jz_S);
235
236             /* rsq = dx*dx+dy*dy+dz*dz */
237             rsq_S0           = gmx_calc_rsq_pr(dx_S0, dy_S0, dz_S0);
238             rsq_S1           = gmx_calc_rsq_pr(dx_S1, dy_S1, dz_S1);
239             rsq_S2           = gmx_calc_rsq_pr(dx_S2, dy_S2, dz_S2);
240             rsq_S3           = gmx_calc_rsq_pr(dx_S3, dy_S3, dz_S3);
241
242             wco_S0           = gmx_cmplt_pr(rsq_S0, rc2_S);
243             wco_S1           = gmx_cmplt_pr(rsq_S1, rc2_S);
244             wco_S2           = gmx_cmplt_pr(rsq_S2, rc2_S);
245             wco_S3           = gmx_cmplt_pr(rsq_S3, rc2_S);
246
247             wco_any_S01      = gmx_or_pb(wco_S0, wco_S1);
248             wco_any_S23      = gmx_or_pb(wco_S2, wco_S3);
249             wco_any_S        = gmx_or_pb(wco_any_S01, wco_any_S23);
250
251             InRange          = gmx_anytrue_pb(wco_any_S);
252
253             *ndistc += 4*GMX_SIMD_WIDTH_HERE;
254         }
255         if (!InRange)
256         {
257             cjl--;
258         }
259     }
260
261     if (cjf <= cjl)
262     {
263         for (cj = cjf; cj <= cjl; cj++)
264         {
265             /* Store cj and the interaction mask */
266             nbl->cj[nbl->ncj].cj   = CI_TO_CJ_SIMD_4XN(gridj->cell0) + cj;
267             nbl->cj[nbl->ncj].excl = get_imask_simd_4xn(remove_sub_diag, ci, cj);
268 #ifdef GMX_CPU_ACCELERATION_IBM_QPX
269             nbl->cj[nbl->ncj].interaction_mask_indices[0] = (nbl->cj[nbl->ncj].excl & 0x000F) >> (0 * 4);
270             nbl->cj[nbl->ncj].interaction_mask_indices[1] = (nbl->cj[nbl->ncj].excl & 0x00F0) >> (1 * 4);
271             nbl->cj[nbl->ncj].interaction_mask_indices[2] = (nbl->cj[nbl->ncj].excl & 0x0F00) >> (2 * 4);
272             nbl->cj[nbl->ncj].interaction_mask_indices[3] = (nbl->cj[nbl->ncj].excl & 0xF000) >> (3 * 4);
273 #endif
274             nbl->ncj++;
275         }
276         /* Increase the closing index in i super-cell list */
277         nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
278     }
279 }
280
281 #undef STRIDE_S