Update copyright statements and change license to LGPL
[alexxy/gromacs.git] / src / mdlib / nbnxn_search_x86_simd.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 /* GMX_MM128_HERE or GMX_MM256_HERE should be set before including this file.
40  * gmx_sse_or_avh.h should be included before including this file.
41  */
42
43 /* Copies PBC shifted i-cell packed atom coordinates to working array */
44 #ifdef GMX_MM128_HERE
45 static void icell_set_x_x86_simd128
46 #else
47 #ifdef GMX_MM256_HERE
48 static void icell_set_x_x86_simd256
49 #else
50 "error: GMX_MM128_HERE or GMX_MM256_HERE not defined"
51 #endif
52 #endif
53                                    (int ci,
54                                     real shx,real shy,real shz,
55                                     int na_c,
56                                     int stride,const real *x,
57                                     nbnxn_list_work_t *work)
58 {
59     int  ia;
60 #ifdef GMX_MM128_HERE
61     nbnxn_x_ci_x86_simd128_t *x_ci;
62
63     x_ci = work->x_ci_x86_simd128;
64
65     ia = X_IND_CI_S128(ci);
66 #else
67     nbnxn_x_ci_x86_simd256_t *x_ci;
68
69     x_ci = work->x_ci_x86_simd256;
70
71     ia = X_IND_CI_S256(ci);
72 #endif
73
74     x_ci->ix_SSE0 = gmx_set1_pr(x[ia + 0*STRIDE_S    ] + shx);
75     x_ci->iy_SSE0 = gmx_set1_pr(x[ia + 1*STRIDE_S    ] + shy);
76     x_ci->iz_SSE0 = gmx_set1_pr(x[ia + 2*STRIDE_S    ] + shz);
77     x_ci->ix_SSE1 = gmx_set1_pr(x[ia + 0*STRIDE_S + 1] + shx);
78     x_ci->iy_SSE1 = gmx_set1_pr(x[ia + 1*STRIDE_S + 1] + shy);
79     x_ci->iz_SSE1 = gmx_set1_pr(x[ia + 2*STRIDE_S + 1] + shz);
80     x_ci->ix_SSE2 = gmx_set1_pr(x[ia + 0*STRIDE_S + 2] + shx);
81     x_ci->iy_SSE2 = gmx_set1_pr(x[ia + 1*STRIDE_S + 2] + shy);
82     x_ci->iz_SSE2 = gmx_set1_pr(x[ia + 2*STRIDE_S + 2] + shz);
83     x_ci->ix_SSE3 = gmx_set1_pr(x[ia + 0*STRIDE_S + 3] + shx);
84     x_ci->iy_SSE3 = gmx_set1_pr(x[ia + 1*STRIDE_S + 3] + shy);
85     x_ci->iz_SSE3 = gmx_set1_pr(x[ia + 2*STRIDE_S + 3] + shz);
86 }
87
88 /* SSE or AVX code for making a pair list of cell ci vs cell cjf-cjl
89  * for coordinates in packed format.
90  * Checks bouding box distances and possibly atom pair distances.
91  * This is an accelerated version of make_cluster_list_simple.
92  */
93 #ifdef GMX_MM128_HERE
94 static void make_cluster_list_x86_simd128
95 #else
96 #ifdef GMX_MM256_HERE
97 static void make_cluster_list_x86_simd256
98 #else
99 "error: GMX_MM128_HERE or GMX_MM256_HERE not defined"
100 #endif
101 #endif
102                                          (const nbnxn_grid_t *gridj,
103                                           nbnxn_pairlist_t *nbl,
104                                           int ci,int cjf,int cjl,
105                                           gmx_bool remove_sub_diag,
106                                           const real *x_j,
107                                           real rl2,float rbb2,
108                                           int *ndistc)
109 {
110 #ifdef GMX_MM128_HERE
111     const nbnxn_x_ci_x86_simd128_t *work;
112 #else
113     const nbnxn_x_ci_x86_simd256_t *work;
114 #endif
115
116     const float *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_pr  wco_SSE0;
131     gmx_mm_pr  wco_SSE1;
132     gmx_mm_pr  wco_SSE2;
133     gmx_mm_pr  wco_SSE3;
134     gmx_mm_pr  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 #ifdef GMX_MM128_HERE
143     cjf = CI_TO_CJ_S128(cjf);
144     cjl = CI_TO_CJ_S128(cjl+1) - 1;
145
146     work = nbl->work->x_ci_x86_simd128;
147 #else
148     cjf = CI_TO_CJ_S256(cjf);
149     cjl = CI_TO_CJ_S256(cjl+1) - 1;
150
151     work = nbl->work->x_ci_x86_simd256;
152 #endif
153
154     bb_ci = nbl->work->bb_ci;
155
156     rc2_SSE   = gmx_set1_pr(rl2);
157
158     InRange = FALSE;
159     while (!InRange && cjf <= cjl)
160     {
161         d2 = subc_bb_dist2_sse(4,0,bb_ci,cjf,gridj->bbj);
162         *ndistc += 2;
163         
164         /* Check if the distance is within the distance where
165          * we use only the bounding box distance rbb,
166          * or within the cut-off and there is at least one atom pair
167          * within the cut-off.
168          */
169         if (d2 < rbb2)
170         {
171             InRange = TRUE;
172         }
173         else if (d2 < rl2)
174         {
175 #ifdef GMX_MM128_HERE
176             xind_f  = X_IND_CJ_S128(CI_TO_CJ_S128(gridj->cell0) + cjf);
177 #else
178             xind_f  = X_IND_CJ_S256(CI_TO_CJ_S256(gridj->cell0) + cjf);
179 #endif
180             jx_SSE  = gmx_load_pr(x_j+xind_f+0*STRIDE_S);
181             jy_SSE  = gmx_load_pr(x_j+xind_f+1*STRIDE_S);
182             jz_SSE  = gmx_load_pr(x_j+xind_f+2*STRIDE_S);
183
184             
185             /* Calculate distance */
186             dx_SSE0            = gmx_sub_pr(work->ix_SSE0,jx_SSE);
187             dy_SSE0            = gmx_sub_pr(work->iy_SSE0,jy_SSE);
188             dz_SSE0            = gmx_sub_pr(work->iz_SSE0,jz_SSE);
189             dx_SSE1            = gmx_sub_pr(work->ix_SSE1,jx_SSE);
190             dy_SSE1            = gmx_sub_pr(work->iy_SSE1,jy_SSE);
191             dz_SSE1            = gmx_sub_pr(work->iz_SSE1,jz_SSE);
192             dx_SSE2            = gmx_sub_pr(work->ix_SSE2,jx_SSE);
193             dy_SSE2            = gmx_sub_pr(work->iy_SSE2,jy_SSE);
194             dz_SSE2            = gmx_sub_pr(work->iz_SSE2,jz_SSE);
195             dx_SSE3            = gmx_sub_pr(work->ix_SSE3,jx_SSE);
196             dy_SSE3            = gmx_sub_pr(work->iy_SSE3,jy_SSE);
197             dz_SSE3            = gmx_sub_pr(work->iz_SSE3,jz_SSE);
198             
199             /* rsq = dx*dx+dy*dy+dz*dz */
200             rsq_SSE0           = gmx_calc_rsq_pr(dx_SSE0,dy_SSE0,dz_SSE0);
201             rsq_SSE1           = gmx_calc_rsq_pr(dx_SSE1,dy_SSE1,dz_SSE1);
202             rsq_SSE2           = gmx_calc_rsq_pr(dx_SSE2,dy_SSE2,dz_SSE2);
203             rsq_SSE3           = gmx_calc_rsq_pr(dx_SSE3,dy_SSE3,dz_SSE3);
204             
205             wco_SSE0           = gmx_cmplt_pr(rsq_SSE0,rc2_SSE);
206             wco_SSE1           = gmx_cmplt_pr(rsq_SSE1,rc2_SSE);
207             wco_SSE2           = gmx_cmplt_pr(rsq_SSE2,rc2_SSE);
208             wco_SSE3           = gmx_cmplt_pr(rsq_SSE3,rc2_SSE);
209             
210             wco_any_SSE01      = gmx_or_pr(wco_SSE0,wco_SSE1);
211             wco_any_SSE23      = gmx_or_pr(wco_SSE2,wco_SSE3);
212             wco_any_SSE        = gmx_or_pr(wco_any_SSE01,wco_any_SSE23);
213             
214             InRange            = gmx_movemask_pr(wco_any_SSE);
215
216             *ndistc += 4*GMX_X86_SIMD_WIDTH_HERE;
217         }
218         if (!InRange)
219         {
220             cjf++;
221         }
222     }
223     if (!InRange)
224     {
225         return;
226     }
227
228     InRange = FALSE;
229     while (!InRange && cjl > cjf)
230     {
231         d2 = subc_bb_dist2_sse(4,0,bb_ci,cjl,gridj->bbj);
232         *ndistc += 2;
233         
234         /* Check if the distance is within the distance where
235          * we use only the bounding box distance rbb,
236          * or within the cut-off and there is at least one atom pair
237          * within the cut-off.
238          */
239         if (d2 < rbb2)
240         {
241             InRange = TRUE;
242         }
243         else if (d2 < rl2)
244         {
245 #ifdef GMX_MM128_HERE
246             xind_l  = X_IND_CJ_S128(CI_TO_CJ_S128(gridj->cell0) + cjl);
247 #else
248             xind_l  = X_IND_CJ_S256(CI_TO_CJ_S256(gridj->cell0) + cjl);
249 #endif
250             jx_SSE  = gmx_load_pr(x_j+xind_l+0*STRIDE_S);
251             jy_SSE  = gmx_load_pr(x_j+xind_l+1*STRIDE_S);
252             jz_SSE  = gmx_load_pr(x_j+xind_l+2*STRIDE_S);
253             
254             /* Calculate distance */
255             dx_SSE0            = gmx_sub_pr(work->ix_SSE0,jx_SSE);
256             dy_SSE0            = gmx_sub_pr(work->iy_SSE0,jy_SSE);
257             dz_SSE0            = gmx_sub_pr(work->iz_SSE0,jz_SSE);
258             dx_SSE1            = gmx_sub_pr(work->ix_SSE1,jx_SSE);
259             dy_SSE1            = gmx_sub_pr(work->iy_SSE1,jy_SSE);
260             dz_SSE1            = gmx_sub_pr(work->iz_SSE1,jz_SSE);
261             dx_SSE2            = gmx_sub_pr(work->ix_SSE2,jx_SSE);
262             dy_SSE2            = gmx_sub_pr(work->iy_SSE2,jy_SSE);
263             dz_SSE2            = gmx_sub_pr(work->iz_SSE2,jz_SSE);
264             dx_SSE3            = gmx_sub_pr(work->ix_SSE3,jx_SSE);
265             dy_SSE3            = gmx_sub_pr(work->iy_SSE3,jy_SSE);
266             dz_SSE3            = gmx_sub_pr(work->iz_SSE3,jz_SSE);
267             
268             /* rsq = dx*dx+dy*dy+dz*dz */
269             rsq_SSE0           = gmx_calc_rsq_pr(dx_SSE0,dy_SSE0,dz_SSE0);
270             rsq_SSE1           = gmx_calc_rsq_pr(dx_SSE1,dy_SSE1,dz_SSE1);
271             rsq_SSE2           = gmx_calc_rsq_pr(dx_SSE2,dy_SSE2,dz_SSE2);
272             rsq_SSE3           = gmx_calc_rsq_pr(dx_SSE3,dy_SSE3,dz_SSE3);
273             
274             wco_SSE0           = gmx_cmplt_pr(rsq_SSE0,rc2_SSE);
275             wco_SSE1           = gmx_cmplt_pr(rsq_SSE1,rc2_SSE);
276             wco_SSE2           = gmx_cmplt_pr(rsq_SSE2,rc2_SSE);
277             wco_SSE3           = gmx_cmplt_pr(rsq_SSE3,rc2_SSE);
278             
279             wco_any_SSE01      = gmx_or_pr(wco_SSE0,wco_SSE1);
280             wco_any_SSE23      = gmx_or_pr(wco_SSE2,wco_SSE3);
281             wco_any_SSE        = gmx_or_pr(wco_any_SSE01,wco_any_SSE23);
282             
283             InRange            = gmx_movemask_pr(wco_any_SSE);
284
285             *ndistc += 4*GMX_X86_SIMD_WIDTH_HERE;
286         }
287         if (!InRange)
288         {
289             cjl--;
290         }
291     }
292
293     if (cjf <= cjl)
294     {
295         for(cj=cjf; cj<=cjl; cj++)
296         {
297             /* Store cj and the interaction mask */
298 #ifdef GMX_MM128_HERE
299             nbl->cj[nbl->ncj].cj   = CI_TO_CJ_S128(gridj->cell0) + cj;
300             nbl->cj[nbl->ncj].excl = get_imask_x86_simd128(remove_sub_diag,ci,cj);
301 #else
302             nbl->cj[nbl->ncj].cj   = CI_TO_CJ_S256(gridj->cell0) + cj;
303             nbl->cj[nbl->ncj].excl = get_imask_x86_simd256(remove_sub_diag,ci,cj);
304 #endif
305             nbl->ncj++;
306         }
307         /* Increase the closing index in i super-cell list */
308         nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
309     }
310 }