Merge branch release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_search_simd_2xnn.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 /* Get the half-width SIMD stuff from the kernel utils files */
40 #include "nbnxn_kernels/nbnxn_kernel_simd_utils.h"
41
42
43 #if GMX_SIMD_WIDTH_HERE >= 2*NBNXN_CPU_CLUSTER_I_SIZE
44 #define STRIDE_S  (GMX_SIMD_WIDTH_HERE/2)
45 #else
46 #define STRIDE_S  NBNXN_CPU_CLUSTER_I_SIZE
47 #endif
48
49 static gmx_inline gmx_mm_pr gmx_load_hpr_hilo_pr(const real *a)
50 {
51     gmx_mm_hpr a_S;
52     gmx_mm_pr  a_a_S;
53
54     gmx_load_hpr(&a_S, a);
55
56     gmx_2hpr_to_pr(a_S, a_S, &a_a_S);
57
58     return a_a_S;
59 }
60
61 static gmx_inline gmx_mm_pr gmx_set_2real_shift_pr(const real *a, real shift)
62 {
63     gmx_mm_hpr a0_S, a1_S;
64     gmx_mm_pr  a0_a1_S;
65
66     gmx_set1_hpr(&a0_S, a[0] + shift);
67     gmx_set1_hpr(&a1_S, a[1] + shift);
68
69     gmx_2hpr_to_pr(a0_S, a1_S, &a0_a1_S);
70
71     return a0_a1_S;
72 }
73
74 /* Copies PBC shifted i-cell packed atom coordinates to working array */
75 static gmx_inline void
76 icell_set_x_simd_2xnn(int ci,
77                       real shx, real shy, real shz,
78                       int gmx_unused na_c,
79                       int gmx_unused stride, const real *x,
80                       nbnxn_list_work_t *work)
81 {
82     int                     ia;
83     nbnxn_x_ci_simd_2xnn_t *x_ci;
84
85     x_ci = work->x_ci_simd_2xnn;
86
87     ia = X_IND_CI_SIMD_2XNN(ci);
88
89     x_ci->ix_SSE0 = gmx_set_2real_shift_pr(x + ia + 0*STRIDE_S + 0, shx);
90     x_ci->iy_SSE0 = gmx_set_2real_shift_pr(x + ia + 1*STRIDE_S + 0, shy);
91     x_ci->iz_SSE0 = gmx_set_2real_shift_pr(x + ia + 2*STRIDE_S + 0, shz);
92     x_ci->ix_SSE2 = gmx_set_2real_shift_pr(x + ia + 0*STRIDE_S + 2, shx);
93     x_ci->iy_SSE2 = gmx_set_2real_shift_pr(x + ia + 1*STRIDE_S + 2, shy);
94     x_ci->iz_SSE2 = gmx_set_2real_shift_pr(x + ia + 2*STRIDE_S + 2, shz);
95 }
96
97 #ifndef GMX_SIMD_HAVE_ANYTRUE
98 /* Fallback function in case gmx_anytrue_pr is not present */
99 static gmx_inline gmx_bool
100 gmx_anytrue_2xn_pb(gmx_mm_pb bool_S)
101 {
102     real     bools_array[2*GMX_SIMD_WIDTH_HERE], *bools;
103     gmx_bool any;
104     int      s;
105
106     bools = gmx_simd_align_real(bools_array);
107
108     gmx_store_pb(bools, bool_S);
109
110     any = FALSE;
111     for (s = 0; s < GMX_SIMD_WIDTH_HERE; s++)
112     {
113         if (GMX_SIMD_IS_TRUE(s))
114         {
115             any = TRUE;
116         }
117     }
118
119     return any;
120 }
121 #endif
122
123 /* SIMD code for making a pair list of cell ci vs cell cjf-cjl
124  * for coordinates in packed format.
125  * Checks bouding box distances and possibly atom pair distances.
126  * This is an accelerated version of make_cluster_list_simple.
127  */
128 static gmx_inline void
129 make_cluster_list_simd_2xnn(const nbnxn_grid_t *gridj,
130                             nbnxn_pairlist_t *nbl,
131                             int ci, int cjf, int cjl,
132                             gmx_bool remove_sub_diag,
133                             const real *x_j,
134                             real rl2, float rbb2,
135                             int *ndistc)
136 {
137     const nbnxn_x_ci_simd_2xnn_t *work;
138     const float                  *bb_ci;
139
140     gmx_mm_pr                     jx_SSE, jy_SSE, jz_SSE;
141
142     gmx_mm_pr                     dx_SSE0, dy_SSE0, dz_SSE0;
143     gmx_mm_pr                     dx_SSE2, dy_SSE2, dz_SSE2;
144
145     gmx_mm_pr                     rsq_SSE0;
146     gmx_mm_pr                     rsq_SSE2;
147
148     gmx_mm_pb                     wco_SSE0;
149     gmx_mm_pb                     wco_SSE2;
150     gmx_mm_pb                     wco_any_SSE;
151
152     gmx_mm_pr                     rc2_SSE;
153
154     gmx_bool                      InRange;
155     float                         d2;
156     int                           xind_f, xind_l, cj;
157
158     cjf = CI_TO_CJ_SIMD_2XNN(cjf);
159     cjl = CI_TO_CJ_SIMD_2XNN(cjl+1) - 1;
160
161     work = nbl->work->x_ci_simd_2xnn;
162
163     bb_ci = nbl->work->bb_ci;
164
165     rc2_SSE   = gmx_set1_pr(rl2);
166
167     InRange = FALSE;
168     while (!InRange && cjf <= cjl)
169     {
170 #ifdef NBNXN_SEARCH_BB_SSE
171         d2 = subc_bb_dist2_sse(0, bb_ci, cjf, gridj->bbj);
172 #else
173         d2 = subc_bb_dist2(0, bb_ci, cjf, gridj->bbj);
174 #endif
175         *ndistc += 2;
176
177         /* Check if the distance is within the distance where
178          * we use only the bounding box distance rbb,
179          * or within the cut-off and there is at least one atom pair
180          * within the cut-off.
181          */
182         if (d2 < rbb2)
183         {
184             InRange = TRUE;
185         }
186         else if (d2 < rl2)
187         {
188             xind_f  = X_IND_CJ_SIMD_2XNN(CI_TO_CJ_SIMD_2XNN(gridj->cell0) + cjf);
189
190             jx_SSE  = gmx_load_hpr_hilo_pr(x_j+xind_f+0*STRIDE_S);
191             jy_SSE  = gmx_load_hpr_hilo_pr(x_j+xind_f+1*STRIDE_S);
192             jz_SSE  = gmx_load_hpr_hilo_pr(x_j+xind_f+2*STRIDE_S);
193
194             /* Calculate distance */
195             dx_SSE0            = gmx_sub_pr(work->ix_SSE0, jx_SSE);
196             dy_SSE0            = gmx_sub_pr(work->iy_SSE0, jy_SSE);
197             dz_SSE0            = gmx_sub_pr(work->iz_SSE0, jz_SSE);
198             dx_SSE2            = gmx_sub_pr(work->ix_SSE2, jx_SSE);
199             dy_SSE2            = gmx_sub_pr(work->iy_SSE2, jy_SSE);
200             dz_SSE2            = gmx_sub_pr(work->iz_SSE2, jz_SSE);
201
202             /* rsq = dx*dx+dy*dy+dz*dz */
203             rsq_SSE0           = gmx_calc_rsq_pr(dx_SSE0, dy_SSE0, dz_SSE0);
204             rsq_SSE2           = gmx_calc_rsq_pr(dx_SSE2, dy_SSE2, dz_SSE2);
205
206             wco_SSE0           = gmx_cmplt_pr(rsq_SSE0, rc2_SSE);
207             wco_SSE2           = gmx_cmplt_pr(rsq_SSE2, rc2_SSE);
208
209             wco_any_SSE        = gmx_or_pb(wco_SSE0, wco_SSE2);
210
211 #ifdef GMX_SIMD_HAVE_ANYTRUE
212             InRange            = gmx_anytrue_pb(wco_any_SSE);
213 #else
214             InRange            = gmx_anytrue_2xn_pb(wco_any_SSE);
215 #endif
216
217             *ndistc += 2*GMX_SIMD_WIDTH_HERE;
218         }
219         if (!InRange)
220         {
221             cjf++;
222         }
223     }
224     if (!InRange)
225     {
226         return;
227     }
228
229     InRange = FALSE;
230     while (!InRange && cjl > cjf)
231     {
232 #ifdef NBNXN_SEARCH_BB_SSE
233         d2 = subc_bb_dist2_sse(0, bb_ci, cjl, gridj->bbj);
234 #else
235         d2 = subc_bb_dist2(0, bb_ci, cjl, gridj->bbj);
236 #endif
237         *ndistc += 2;
238
239         /* Check if the distance is within the distance where
240          * we use only the bounding box distance rbb,
241          * or within the cut-off and there is at least one atom pair
242          * within the cut-off.
243          */
244         if (d2 < rbb2)
245         {
246             InRange = TRUE;
247         }
248         else if (d2 < rl2)
249         {
250             xind_l  = X_IND_CJ_SIMD_2XNN(CI_TO_CJ_SIMD_2XNN(gridj->cell0) + cjl);
251
252             jx_SSE  = gmx_load_hpr_hilo_pr(x_j+xind_l+0*STRIDE_S);
253             jy_SSE  = gmx_load_hpr_hilo_pr(x_j+xind_l+1*STRIDE_S);
254             jz_SSE  = gmx_load_hpr_hilo_pr(x_j+xind_l+2*STRIDE_S);
255
256             /* Calculate distance */
257             dx_SSE0            = gmx_sub_pr(work->ix_SSE0, jx_SSE);
258             dy_SSE0            = gmx_sub_pr(work->iy_SSE0, jy_SSE);
259             dz_SSE0            = gmx_sub_pr(work->iz_SSE0, 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
264             /* rsq = dx*dx+dy*dy+dz*dz */
265             rsq_SSE0           = gmx_calc_rsq_pr(dx_SSE0, dy_SSE0, dz_SSE0);
266             rsq_SSE2           = gmx_calc_rsq_pr(dx_SSE2, dy_SSE2, dz_SSE2);
267
268             wco_SSE0           = gmx_cmplt_pr(rsq_SSE0, rc2_SSE);
269             wco_SSE2           = gmx_cmplt_pr(rsq_SSE2, rc2_SSE);
270
271             wco_any_SSE        = gmx_or_pb(wco_SSE0, wco_SSE2);
272
273 #ifdef GMX_SIMD_HAVE_ANYTRUE
274             InRange            = gmx_anytrue_pb(wco_any_SSE);
275 #else
276             InRange            = gmx_anytrue_2xn_pb(wco_any_SSE);
277 #endif
278
279             *ndistc += 2*GMX_SIMD_WIDTH_HERE;
280         }
281         if (!InRange)
282         {
283             cjl--;
284         }
285     }
286
287     if (cjf <= cjl)
288     {
289         for (cj = cjf; cj <= cjl; cj++)
290         {
291             /* Store cj and the interaction mask */
292             nbl->cj[nbl->ncj].cj   = CI_TO_CJ_SIMD_2XNN(gridj->cell0) + cj;
293             nbl->cj[nbl->ncj].excl = get_imask_simd_2xnn(remove_sub_diag, ci, cj);
294             nbl->ncj++;
295         }
296         /* Increase the closing index in i super-cell list */
297         nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
298     }
299 }
300
301 #undef STRIDE_S