Merge 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 #if GMX_NBNXN_SIMD_BITWIDTH == 128
40 #define GMX_MM128_HERE
41 #else
42 #if GMX_NBNXN_SIMD_BITWIDTH == 256
43 #define GMX_MM256_HERE
44 #else
45 #error "unsupported GMX_NBNXN_SIMD_BITWIDTH"
46 #endif
47 #endif
48 #include "gmx_simd_macros.h"
49
50 #if GMX_SIMD_WIDTH_HERE >= 2*NBNXN_CPU_CLUSTER_I_SIZE
51 #define STRIDE_S  (GMX_SIMD_WIDTH_HERE/2)
52 #else
53 #define STRIDE_S  NBNXN_CPU_CLUSTER_I_SIZE
54 #endif
55
56 static gmx_inline gmx_mm_pr gmx_load_hpr_hilo_pr(const real *a)
57 {
58     gmx_mm_hpr a_SSE;
59
60     a_SSE = _mm_load_ps(a);
61
62     return gmx_2hpr_to_pr(a_SSE, a_SSE);
63 }
64
65 static gmx_inline gmx_mm_pr gmx_set_2real_shift_pr(const real *a, real shift)
66 {
67     gmx_mm_hpr a0, a1;
68
69     a0 = _mm_set1_ps(a[0] + shift);
70     a1 = _mm_set1_ps(a[1] + shift);
71
72     return gmx_2hpr_to_pr(a1, a0);
73 }
74
75 /* Copies PBC shifted i-cell packed atom coordinates to working array */
76 static gmx_inline void
77 icell_set_x_simd_2xnn(int ci,
78                       real shx, real shy, real shz,
79                       int na_c,
80                       int stride, const real *x,
81                       nbnxn_list_work_t *work)
82 {
83     int                     ia;
84     nbnxn_x_ci_simd_2xnn_t *x_ci;
85
86     x_ci = work->x_ci_simd_2xnn;
87
88     ia = X_IND_CI_SIMD_2XNN(ci);
89
90     x_ci->ix_SSE0 = gmx_set_2real_shift_pr(x + ia + 0*STRIDE_S + 0, shx);
91     x_ci->iy_SSE0 = gmx_set_2real_shift_pr(x + ia + 1*STRIDE_S + 0, shy);
92     x_ci->iz_SSE0 = gmx_set_2real_shift_pr(x + ia + 2*STRIDE_S + 0, shz);
93     x_ci->ix_SSE2 = gmx_set_2real_shift_pr(x + ia + 0*STRIDE_S + 2, shx);
94     x_ci->iy_SSE2 = gmx_set_2real_shift_pr(x + ia + 1*STRIDE_S + 2, shy);
95     x_ci->iz_SSE2 = gmx_set_2real_shift_pr(x + ia + 2*STRIDE_S + 2, shz);
96 }
97
98 /* SIMD code for making a pair list of cell ci vs cell cjf-cjl
99  * for coordinates in packed format.
100  * Checks bouding box distances and possibly atom pair distances.
101  * This is an accelerated version of make_cluster_list_simple.
102  */
103 static gmx_inline void
104 make_cluster_list_simd_2xnn(const nbnxn_grid_t *gridj,
105                             nbnxn_pairlist_t *nbl,
106                             int ci, int cjf, int cjl,
107                             gmx_bool remove_sub_diag,
108                             const real *x_j,
109                             real rl2, float rbb2,
110                             int *ndistc)
111 {
112     const nbnxn_x_ci_simd_2xnn_t *work;
113     const float                  *bb_ci;
114
115     gmx_mm_pr                     jx_SSE, jy_SSE, jz_SSE;
116
117     gmx_mm_pr                     dx_SSE0, dy_SSE0, dz_SSE0;
118     gmx_mm_pr                     dx_SSE2, dy_SSE2, dz_SSE2;
119
120     gmx_mm_pr                     rsq_SSE0;
121     gmx_mm_pr                     rsq_SSE2;
122
123     gmx_mm_pr                     wco_SSE0;
124     gmx_mm_pr                     wco_SSE2;
125     gmx_mm_pr                     wco_any_SSE;
126
127     gmx_mm_pr                     rc2_SSE;
128
129     gmx_bool                      InRange;
130     float                         d2;
131     int                           xind_f, xind_l, cj;
132
133     cjf = CI_TO_CJ_SIMD_2XNN(cjf);
134     cjl = CI_TO_CJ_SIMD_2XNN(cjl+1) - 1;
135
136     work = nbl->work->x_ci_simd_2xnn;
137
138     bb_ci = nbl->work->bb_ci;
139
140     rc2_SSE   = gmx_set1_pr(rl2);
141
142     InRange = FALSE;
143     while (!InRange && cjf <= cjl)
144     {
145         d2       = subc_bb_dist2_sse(4, 0, bb_ci, cjf, gridj->bbj);
146         *ndistc += 2;
147
148         /* Check if the distance is within the distance where
149          * we use only the bounding box distance rbb,
150          * or within the cut-off and there is at least one atom pair
151          * within the cut-off.
152          */
153         if (d2 < rbb2)
154         {
155             InRange = TRUE;
156         }
157         else if (d2 < rl2)
158         {
159             xind_f  = X_IND_CJ_SIMD_2XNN(CI_TO_CJ_SIMD_2XNN(gridj->cell0) + cjf);
160
161             jx_SSE  = gmx_load_hpr_hilo_pr(x_j+xind_f+0*STRIDE_S);
162             jy_SSE  = gmx_load_hpr_hilo_pr(x_j+xind_f+1*STRIDE_S);
163             jz_SSE  = gmx_load_hpr_hilo_pr(x_j+xind_f+2*STRIDE_S);
164
165             /* Calculate distance */
166             dx_SSE0            = gmx_sub_pr(work->ix_SSE0, jx_SSE);
167             dy_SSE0            = gmx_sub_pr(work->iy_SSE0, jy_SSE);
168             dz_SSE0            = gmx_sub_pr(work->iz_SSE0, jz_SSE);
169             dx_SSE2            = gmx_sub_pr(work->ix_SSE2, jx_SSE);
170             dy_SSE2            = gmx_sub_pr(work->iy_SSE2, jy_SSE);
171             dz_SSE2            = gmx_sub_pr(work->iz_SSE2, jz_SSE);
172
173             /* rsq = dx*dx+dy*dy+dz*dz */
174             rsq_SSE0           = gmx_calc_rsq_pr(dx_SSE0, dy_SSE0, dz_SSE0);
175             rsq_SSE2           = gmx_calc_rsq_pr(dx_SSE2, dy_SSE2, dz_SSE2);
176
177             wco_SSE0           = gmx_cmplt_pr(rsq_SSE0, rc2_SSE);
178             wco_SSE2           = gmx_cmplt_pr(rsq_SSE2, rc2_SSE);
179
180             wco_any_SSE        = gmx_or_pr(wco_SSE0, wco_SSE2);
181
182             InRange            = gmx_movemask_pr(wco_any_SSE);
183
184             *ndistc += 2*GMX_SIMD_WIDTH_HERE;
185         }
186         if (!InRange)
187         {
188             cjf++;
189         }
190     }
191     if (!InRange)
192     {
193         return;
194     }
195
196     InRange = FALSE;
197     while (!InRange && cjl > cjf)
198     {
199         d2       = subc_bb_dist2_sse(4, 0, bb_ci, cjl, gridj->bbj);
200         *ndistc += 2;
201
202         /* Check if the distance is within the distance where
203          * we use only the bounding box distance rbb,
204          * or within the cut-off and there is at least one atom pair
205          * within the cut-off.
206          */
207         if (d2 < rbb2)
208         {
209             InRange = TRUE;
210         }
211         else if (d2 < rl2)
212         {
213             xind_l  = X_IND_CJ_SIMD_2XNN(CI_TO_CJ_SIMD_2XNN(gridj->cell0) + cjl);
214
215             jx_SSE  = gmx_load_hpr_hilo_pr(x_j+xind_l+0*STRIDE_S);
216             jy_SSE  = gmx_load_hpr_hilo_pr(x_j+xind_l+1*STRIDE_S);
217             jz_SSE  = gmx_load_hpr_hilo_pr(x_j+xind_l+2*STRIDE_S);
218
219             /* Calculate distance */
220             dx_SSE0            = gmx_sub_pr(work->ix_SSE0, jx_SSE);
221             dy_SSE0            = gmx_sub_pr(work->iy_SSE0, jy_SSE);
222             dz_SSE0            = gmx_sub_pr(work->iz_SSE0, jz_SSE);
223             dx_SSE2            = gmx_sub_pr(work->ix_SSE2, jx_SSE);
224             dy_SSE2            = gmx_sub_pr(work->iy_SSE2, jy_SSE);
225             dz_SSE2            = gmx_sub_pr(work->iz_SSE2, jz_SSE);
226
227             /* rsq = dx*dx+dy*dy+dz*dz */
228             rsq_SSE0           = gmx_calc_rsq_pr(dx_SSE0, dy_SSE0, dz_SSE0);
229             rsq_SSE2           = gmx_calc_rsq_pr(dx_SSE2, dy_SSE2, dz_SSE2);
230
231             wco_SSE0           = gmx_cmplt_pr(rsq_SSE0, rc2_SSE);
232             wco_SSE2           = gmx_cmplt_pr(rsq_SSE2, rc2_SSE);
233
234             wco_any_SSE        = gmx_or_pr(wco_SSE0, wco_SSE2);
235
236             InRange            = gmx_movemask_pr(wco_any_SSE);
237
238             *ndistc += 2*GMX_SIMD_WIDTH_HERE;
239         }
240         if (!InRange)
241         {
242             cjl--;
243         }
244     }
245
246     if (cjf <= cjl)
247     {
248         for (cj = cjf; cj <= cjl; cj++)
249         {
250             /* Store cj and the interaction mask */
251             nbl->cj[nbl->ncj].cj   = CI_TO_CJ_SIMD_2XNN(gridj->cell0) + cj;
252             nbl->cj[nbl->ncj].excl = get_imask_x86_simd_2xnn(remove_sub_diag, ci, cj);
253             nbl->ncj++;
254         }
255         /* Increase the closing index in i super-cell list */
256         nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
257     }
258 }
259
260 #undef STRIDE_S
261 #undef GMX_MM128_HERE
262 #undef GMX_MM256_HERE