Remove useless warnings from Doxygen logs.
[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  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38
39
40 #if GMX_SIMD_WIDTH_HERE >= NBNXN_CPU_CLUSTER_I_SIZE
41 #define STRIDE_S  (GMX_SIMD_WIDTH_HERE)
42 #else
43 #define STRIDE_S  NBNXN_CPU_CLUSTER_I_SIZE
44 #endif
45
46 /* Copies PBC shifted i-cell packed atom coordinates to working array */
47 static gmx_inline void
48 icell_set_x_simd_4xn(int ci,
49                      real shx, real shy, real shz,
50                      int gmx_unused na_c,
51                      int gmx_unused stride, const real *x,
52                      nbnxn_list_work_t *work)
53 {
54     int                    ia;
55     nbnxn_x_ci_simd_4xn_t *x_ci;
56
57     x_ci = work->x_ci_simd_4xn;
58
59     ia = X_IND_CI_SIMD_4XN(ci);
60
61     x_ci->ix_SSE0 = gmx_set1_pr(x[ia + 0*STRIDE_S    ] + shx);
62     x_ci->iy_SSE0 = gmx_set1_pr(x[ia + 1*STRIDE_S    ] + shy);
63     x_ci->iz_SSE0 = gmx_set1_pr(x[ia + 2*STRIDE_S    ] + shz);
64     x_ci->ix_SSE1 = gmx_set1_pr(x[ia + 0*STRIDE_S + 1] + shx);
65     x_ci->iy_SSE1 = gmx_set1_pr(x[ia + 1*STRIDE_S + 1] + shy);
66     x_ci->iz_SSE1 = gmx_set1_pr(x[ia + 2*STRIDE_S + 1] + shz);
67     x_ci->ix_SSE2 = gmx_set1_pr(x[ia + 0*STRIDE_S + 2] + shx);
68     x_ci->iy_SSE2 = gmx_set1_pr(x[ia + 1*STRIDE_S + 2] + shy);
69     x_ci->iz_SSE2 = gmx_set1_pr(x[ia + 2*STRIDE_S + 2] + shz);
70     x_ci->ix_SSE3 = gmx_set1_pr(x[ia + 0*STRIDE_S + 3] + shx);
71     x_ci->iy_SSE3 = gmx_set1_pr(x[ia + 1*STRIDE_S + 3] + shy);
72     x_ci->iz_SSE3 = gmx_set1_pr(x[ia + 2*STRIDE_S + 3] + shz);
73 }
74
75 #ifndef GMX_SIMD_HAVE_ANYTRUE
76 /* Fallback function in case gmx_anytrue_pr is not present */
77 static gmx_inline gmx_bool
78 gmx_anytrue_4xn_pb(gmx_mm_pb bool_S)
79 {
80     real     bools_array[2*GMX_SIMD_WIDTH_HERE], *bools;
81     gmx_bool any;
82     int      s;
83
84     bools = gmx_simd_align_real(bools_array);
85
86     gmx_store_pb(bools, bool_S);
87
88     any = FALSE;
89     for (s = 0; s < GMX_SIMD_WIDTH_HERE; s++)
90     {
91         if (GMX_SIMD_IS_TRUE(bools[s]))
92         {
93             any = TRUE;
94         }
95     }
96
97     return any;
98 }
99 #endif
100
101 /* SIMD code for making a pair list of cell ci vs cell cjf-cjl
102  * for coordinates in packed format.
103  * Checks bouding box distances and possibly atom pair distances.
104  * This is an accelerated version of make_cluster_list_simple.
105  */
106 static gmx_inline void
107 make_cluster_list_simd_4xn(const nbnxn_grid_t *gridj,
108                            nbnxn_pairlist_t *nbl,
109                            int ci, int cjf, int cjl,
110                            gmx_bool remove_sub_diag,
111                            const real *x_j,
112                            real rl2, float rbb2,
113                            int *ndistc)
114 {
115     const nbnxn_x_ci_simd_4xn_t *work;
116     const nbnxn_bb_t            *bb_ci;
117
118     gmx_mm_pr                    jx_SSE, jy_SSE, jz_SSE;
119
120     gmx_mm_pr                    dx_SSE0, dy_SSE0, dz_SSE0;
121     gmx_mm_pr                    dx_SSE1, dy_SSE1, dz_SSE1;
122     gmx_mm_pr                    dx_SSE2, dy_SSE2, dz_SSE2;
123     gmx_mm_pr                    dx_SSE3, dy_SSE3, dz_SSE3;
124
125     gmx_mm_pr                    rsq_SSE0;
126     gmx_mm_pr                    rsq_SSE1;
127     gmx_mm_pr                    rsq_SSE2;
128     gmx_mm_pr                    rsq_SSE3;
129
130     gmx_mm_pb                    wco_SSE0;
131     gmx_mm_pb                    wco_SSE1;
132     gmx_mm_pb                    wco_SSE2;
133     gmx_mm_pb                    wco_SSE3;
134     gmx_mm_pb                    wco_any_SSE01, wco_any_SSE23, wco_any_SSE;
135
136     gmx_mm_pr                    rc2_SSE;
137
138     gmx_bool                     InRange;
139     float                        d2;
140     int                          xind_f, xind_l, cj;
141
142     cjf = CI_TO_CJ_SIMD_4XN(cjf);
143     cjl = CI_TO_CJ_SIMD_4XN(cjl+1) - 1;
144
145     work = nbl->work->x_ci_simd_4xn;
146
147     bb_ci = nbl->work->bb_ci;
148
149     rc2_SSE   = gmx_set1_pr(rl2);
150
151     InRange = FALSE;
152     while (!InRange && cjf <= cjl)
153     {
154 #ifdef NBNXN_SEARCH_BB_SSE
155         d2 = subc_bb_dist2_sse(0, bb_ci, cjf, gridj->bbj);
156 #else
157         d2 = subc_bb_dist2(0, bb_ci, cjf, gridj->bbj);
158 #endif
159         *ndistc += 2;
160
161         /* Check if the distance is within the distance where
162          * we use only the bounding box distance rbb,
163          * or within the cut-off and there is at least one atom pair
164          * within the cut-off.
165          */
166         if (d2 < rbb2)
167         {
168             InRange = TRUE;
169         }
170         else if (d2 < rl2)
171         {
172             xind_f  = X_IND_CJ_SIMD_4XN(CI_TO_CJ_SIMD_4XN(gridj->cell0) + cjf);
173
174             jx_SSE  = gmx_load_pr(x_j+xind_f+0*STRIDE_S);
175             jy_SSE  = gmx_load_pr(x_j+xind_f+1*STRIDE_S);
176             jz_SSE  = gmx_load_pr(x_j+xind_f+2*STRIDE_S);
177
178
179             /* Calculate distance */
180             dx_SSE0            = gmx_sub_pr(work->ix_SSE0, jx_SSE);
181             dy_SSE0            = gmx_sub_pr(work->iy_SSE0, jy_SSE);
182             dz_SSE0            = gmx_sub_pr(work->iz_SSE0, jz_SSE);
183             dx_SSE1            = gmx_sub_pr(work->ix_SSE1, jx_SSE);
184             dy_SSE1            = gmx_sub_pr(work->iy_SSE1, jy_SSE);
185             dz_SSE1            = gmx_sub_pr(work->iz_SSE1, jz_SSE);
186             dx_SSE2            = gmx_sub_pr(work->ix_SSE2, jx_SSE);
187             dy_SSE2            = gmx_sub_pr(work->iy_SSE2, jy_SSE);
188             dz_SSE2            = gmx_sub_pr(work->iz_SSE2, jz_SSE);
189             dx_SSE3            = gmx_sub_pr(work->ix_SSE3, jx_SSE);
190             dy_SSE3            = gmx_sub_pr(work->iy_SSE3, jy_SSE);
191             dz_SSE3            = gmx_sub_pr(work->iz_SSE3, jz_SSE);
192
193             /* rsq = dx*dx+dy*dy+dz*dz */
194             rsq_SSE0           = gmx_calc_rsq_pr(dx_SSE0, dy_SSE0, dz_SSE0);
195             rsq_SSE1           = gmx_calc_rsq_pr(dx_SSE1, dy_SSE1, dz_SSE1);
196             rsq_SSE2           = gmx_calc_rsq_pr(dx_SSE2, dy_SSE2, dz_SSE2);
197             rsq_SSE3           = gmx_calc_rsq_pr(dx_SSE3, dy_SSE3, dz_SSE3);
198
199             wco_SSE0           = gmx_cmplt_pr(rsq_SSE0, rc2_SSE);
200             wco_SSE1           = gmx_cmplt_pr(rsq_SSE1, rc2_SSE);
201             wco_SSE2           = gmx_cmplt_pr(rsq_SSE2, rc2_SSE);
202             wco_SSE3           = gmx_cmplt_pr(rsq_SSE3, rc2_SSE);
203
204             wco_any_SSE01      = gmx_or_pb(wco_SSE0, wco_SSE1);
205             wco_any_SSE23      = gmx_or_pb(wco_SSE2, wco_SSE3);
206             wco_any_SSE        = gmx_or_pb(wco_any_SSE01, wco_any_SSE23);
207
208 #ifdef GMX_SIMD_HAVE_ANYTRUE
209             InRange            = gmx_anytrue_pb(wco_any_SSE);
210 #else
211             InRange            = gmx_anytrue_4xn_pb(wco_any_SSE);
212 #endif
213
214             *ndistc += 4*GMX_SIMD_WIDTH_HERE;
215         }
216         if (!InRange)
217         {
218             cjf++;
219         }
220     }
221     if (!InRange)
222     {
223         return;
224     }
225
226     InRange = FALSE;
227     while (!InRange && cjl > cjf)
228     {
229 #ifdef NBNXN_SEARCH_BB_SSE
230         d2 = subc_bb_dist2_sse(0, bb_ci, cjl, gridj->bbj);
231 #else
232         d2 = subc_bb_dist2(0, bb_ci, cjl, gridj->bbj);
233 #endif
234         *ndistc += 2;
235
236         /* Check if the distance is within the distance where
237          * we use only the bounding box distance rbb,
238          * or within the cut-off and there is at least one atom pair
239          * within the cut-off.
240          */
241         if (d2 < rbb2)
242         {
243             InRange = TRUE;
244         }
245         else if (d2 < rl2)
246         {
247             xind_l  = X_IND_CJ_SIMD_4XN(CI_TO_CJ_SIMD_4XN(gridj->cell0) + cjl);
248
249             jx_SSE  = gmx_load_pr(x_j+xind_l+0*STRIDE_S);
250             jy_SSE  = gmx_load_pr(x_j+xind_l+1*STRIDE_S);
251             jz_SSE  = gmx_load_pr(x_j+xind_l+2*STRIDE_S);
252
253             /* Calculate distance */
254             dx_SSE0            = gmx_sub_pr(work->ix_SSE0, jx_SSE);
255             dy_SSE0            = gmx_sub_pr(work->iy_SSE0, jy_SSE);
256             dz_SSE0            = gmx_sub_pr(work->iz_SSE0, jz_SSE);
257             dx_SSE1            = gmx_sub_pr(work->ix_SSE1, jx_SSE);
258             dy_SSE1            = gmx_sub_pr(work->iy_SSE1, jy_SSE);
259             dz_SSE1            = gmx_sub_pr(work->iz_SSE1, jz_SSE);
260             dx_SSE2            = gmx_sub_pr(work->ix_SSE2, jx_SSE);
261             dy_SSE2            = gmx_sub_pr(work->iy_SSE2, jy_SSE);
262             dz_SSE2            = gmx_sub_pr(work->iz_SSE2, jz_SSE);
263             dx_SSE3            = gmx_sub_pr(work->ix_SSE3, jx_SSE);
264             dy_SSE3            = gmx_sub_pr(work->iy_SSE3, jy_SSE);
265             dz_SSE3            = gmx_sub_pr(work->iz_SSE3, jz_SSE);
266
267             /* rsq = dx*dx+dy*dy+dz*dz */
268             rsq_SSE0           = gmx_calc_rsq_pr(dx_SSE0, dy_SSE0, dz_SSE0);
269             rsq_SSE1           = gmx_calc_rsq_pr(dx_SSE1, dy_SSE1, dz_SSE1);
270             rsq_SSE2           = gmx_calc_rsq_pr(dx_SSE2, dy_SSE2, dz_SSE2);
271             rsq_SSE3           = gmx_calc_rsq_pr(dx_SSE3, dy_SSE3, dz_SSE3);
272
273             wco_SSE0           = gmx_cmplt_pr(rsq_SSE0, rc2_SSE);
274             wco_SSE1           = gmx_cmplt_pr(rsq_SSE1, rc2_SSE);
275             wco_SSE2           = gmx_cmplt_pr(rsq_SSE2, rc2_SSE);
276             wco_SSE3           = gmx_cmplt_pr(rsq_SSE3, rc2_SSE);
277
278             wco_any_SSE01      = gmx_or_pb(wco_SSE0, wco_SSE1);
279             wco_any_SSE23      = gmx_or_pb(wco_SSE2, wco_SSE3);
280             wco_any_SSE        = gmx_or_pb(wco_any_SSE01, wco_any_SSE23);
281
282 #ifdef GMX_SIMD_HAVE_ANYTRUE
283             InRange            = gmx_anytrue_pb(wco_any_SSE);
284 #else
285             InRange            = gmx_anytrue_4xn_pb(wco_any_SSE);
286 #endif
287
288             *ndistc += 4*GMX_SIMD_WIDTH_HERE;
289         }
290         if (!InRange)
291         {
292             cjl--;
293         }
294     }
295
296     if (cjf <= cjl)
297     {
298         for (cj = cjf; cj <= cjl; cj++)
299         {
300             /* Store cj and the interaction mask */
301             nbl->cj[nbl->ncj].cj   = CI_TO_CJ_SIMD_4XN(gridj->cell0) + cj;
302             nbl->cj[nbl->ncj].excl = get_imask_simd_4xn(remove_sub_diag, ci, cj);
303             nbl->ncj++;
304         }
305         /* Increase the closing index in i super-cell list */
306         nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
307     }
308 }
309
310 #undef STRIDE_S
311