1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2012, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
36 /* GMX_MM128_HERE or GMX_MM256_HERE should be set before including this file.
37 * gmx_sse_or_avh.h should be included before including this file.
40 /* Copies PBC shifted i-cell packed atom coordinates to working array */
42 static void icell_set_x_x86_simd128
45 static void icell_set_x_x86_simd256
47 "error: GMX_MM128_HERE or GMX_MM256_HERE not defined"
51 real shx,real shy,real shz,
53 int stride,const real *x,
54 nbnxn_list_work_t *work)
58 nbnxn_x_ci_x86_simd128_t *x_ci;
60 x_ci = work->x_ci_x86_simd128;
62 ia = X_IND_CI_S128(ci);
64 nbnxn_x_ci_x86_simd256_t *x_ci;
66 x_ci = work->x_ci_x86_simd256;
68 ia = X_IND_CI_S256(ci);
71 x_ci->ix_SSE0 = gmx_set1_pr(x[ia + 0*STRIDE_S ] + shx);
72 x_ci->iy_SSE0 = gmx_set1_pr(x[ia + 1*STRIDE_S ] + shy);
73 x_ci->iz_SSE0 = gmx_set1_pr(x[ia + 2*STRIDE_S ] + shz);
74 x_ci->ix_SSE1 = gmx_set1_pr(x[ia + 0*STRIDE_S + 1] + shx);
75 x_ci->iy_SSE1 = gmx_set1_pr(x[ia + 1*STRIDE_S + 1] + shy);
76 x_ci->iz_SSE1 = gmx_set1_pr(x[ia + 2*STRIDE_S + 1] + shz);
77 x_ci->ix_SSE2 = gmx_set1_pr(x[ia + 0*STRIDE_S + 2] + shx);
78 x_ci->iy_SSE2 = gmx_set1_pr(x[ia + 1*STRIDE_S + 2] + shy);
79 x_ci->iz_SSE2 = gmx_set1_pr(x[ia + 2*STRIDE_S + 2] + shz);
80 x_ci->ix_SSE3 = gmx_set1_pr(x[ia + 0*STRIDE_S + 3] + shx);
81 x_ci->iy_SSE3 = gmx_set1_pr(x[ia + 1*STRIDE_S + 3] + shy);
82 x_ci->iz_SSE3 = gmx_set1_pr(x[ia + 2*STRIDE_S + 3] + shz);
85 /* SSE or AVX code for making a pair list of cell ci vs cell cjf-cjl
86 * for coordinates in packed format.
87 * Checks bouding box distances and possibly atom pair distances.
88 * This is an accelerated version of make_cluster_list_simple.
91 static void make_cluster_list_x86_simd128
94 static void make_cluster_list_x86_simd256
96 "error: GMX_MM128_HERE or GMX_MM256_HERE not defined"
99 (const nbnxn_grid_t *gridj,
100 nbnxn_pairlist_t *nbl,
101 int ci,int cjf,int cjl,
102 gmx_bool remove_sub_diag,
107 #ifdef GMX_MM128_HERE
108 const nbnxn_x_ci_x86_simd128_t *work;
110 const nbnxn_x_ci_x86_simd256_t *work;
115 gmx_mm_pr jx_SSE,jy_SSE,jz_SSE;
117 gmx_mm_pr dx_SSE0,dy_SSE0,dz_SSE0;
118 gmx_mm_pr dx_SSE1,dy_SSE1,dz_SSE1;
119 gmx_mm_pr dx_SSE2,dy_SSE2,dz_SSE2;
120 gmx_mm_pr dx_SSE3,dy_SSE3,dz_SSE3;
131 gmx_mm_pr wco_any_SSE01,wco_any_SSE23,wco_any_SSE;
137 int xind_f,xind_l,cj;
139 #ifdef GMX_MM128_HERE
140 cjf = CI_TO_CJ_S128(cjf);
141 cjl = CI_TO_CJ_S128(cjl+1) - 1;
143 work = nbl->work->x_ci_x86_simd128;
145 cjf = CI_TO_CJ_S256(cjf);
146 cjl = CI_TO_CJ_S256(cjl+1) - 1;
148 work = nbl->work->x_ci_x86_simd256;
151 bb_ci = nbl->work->bb_ci;
153 rc2_SSE = gmx_set1_pr(rl2);
156 while (!InRange && cjf <= cjl)
158 d2 = subc_bb_dist2_sse(4,0,bb_ci,cjf,gridj->bbj);
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.
172 #ifdef GMX_MM128_HERE
173 xind_f = X_IND_CJ_S128(CI_TO_CJ_S128(gridj->cell0) + cjf);
175 xind_f = X_IND_CJ_S256(CI_TO_CJ_S256(gridj->cell0) + cjf);
177 jx_SSE = gmx_load_pr(x_j+xind_f+0*STRIDE_S);
178 jy_SSE = gmx_load_pr(x_j+xind_f+1*STRIDE_S);
179 jz_SSE = gmx_load_pr(x_j+xind_f+2*STRIDE_S);
182 /* Calculate distance */
183 dx_SSE0 = gmx_sub_pr(work->ix_SSE0,jx_SSE);
184 dy_SSE0 = gmx_sub_pr(work->iy_SSE0,jy_SSE);
185 dz_SSE0 = gmx_sub_pr(work->iz_SSE0,jz_SSE);
186 dx_SSE1 = gmx_sub_pr(work->ix_SSE1,jx_SSE);
187 dy_SSE1 = gmx_sub_pr(work->iy_SSE1,jy_SSE);
188 dz_SSE1 = gmx_sub_pr(work->iz_SSE1,jz_SSE);
189 dx_SSE2 = gmx_sub_pr(work->ix_SSE2,jx_SSE);
190 dy_SSE2 = gmx_sub_pr(work->iy_SSE2,jy_SSE);
191 dz_SSE2 = gmx_sub_pr(work->iz_SSE2,jz_SSE);
192 dx_SSE3 = gmx_sub_pr(work->ix_SSE3,jx_SSE);
193 dy_SSE3 = gmx_sub_pr(work->iy_SSE3,jy_SSE);
194 dz_SSE3 = gmx_sub_pr(work->iz_SSE3,jz_SSE);
196 /* rsq = dx*dx+dy*dy+dz*dz */
197 rsq_SSE0 = gmx_calc_rsq_pr(dx_SSE0,dy_SSE0,dz_SSE0);
198 rsq_SSE1 = gmx_calc_rsq_pr(dx_SSE1,dy_SSE1,dz_SSE1);
199 rsq_SSE2 = gmx_calc_rsq_pr(dx_SSE2,dy_SSE2,dz_SSE2);
200 rsq_SSE3 = gmx_calc_rsq_pr(dx_SSE3,dy_SSE3,dz_SSE3);
202 wco_SSE0 = gmx_cmplt_pr(rsq_SSE0,rc2_SSE);
203 wco_SSE1 = gmx_cmplt_pr(rsq_SSE1,rc2_SSE);
204 wco_SSE2 = gmx_cmplt_pr(rsq_SSE2,rc2_SSE);
205 wco_SSE3 = gmx_cmplt_pr(rsq_SSE3,rc2_SSE);
207 wco_any_SSE01 = gmx_or_pr(wco_SSE0,wco_SSE1);
208 wco_any_SSE23 = gmx_or_pr(wco_SSE2,wco_SSE3);
209 wco_any_SSE = gmx_or_pr(wco_any_SSE01,wco_any_SSE23);
211 InRange = gmx_movemask_pr(wco_any_SSE);
213 *ndistc += 4*GMX_X86_SIMD_WIDTH_HERE;
226 while (!InRange && cjl > cjf)
228 d2 = subc_bb_dist2_sse(4,0,bb_ci,cjl,gridj->bbj);
231 /* Check if the distance is within the distance where
232 * we use only the bounding box distance rbb,
233 * or within the cut-off and there is at least one atom pair
234 * within the cut-off.
242 #ifdef GMX_MM128_HERE
243 xind_l = X_IND_CJ_S128(CI_TO_CJ_S128(gridj->cell0) + cjl);
245 xind_l = X_IND_CJ_S256(CI_TO_CJ_S256(gridj->cell0) + cjl);
247 jx_SSE = gmx_load_pr(x_j+xind_l+0*STRIDE_S);
248 jy_SSE = gmx_load_pr(x_j+xind_l+1*STRIDE_S);
249 jz_SSE = gmx_load_pr(x_j+xind_l+2*STRIDE_S);
251 /* Calculate distance */
252 dx_SSE0 = gmx_sub_pr(work->ix_SSE0,jx_SSE);
253 dy_SSE0 = gmx_sub_pr(work->iy_SSE0,jy_SSE);
254 dz_SSE0 = gmx_sub_pr(work->iz_SSE0,jz_SSE);
255 dx_SSE1 = gmx_sub_pr(work->ix_SSE1,jx_SSE);
256 dy_SSE1 = gmx_sub_pr(work->iy_SSE1,jy_SSE);
257 dz_SSE1 = gmx_sub_pr(work->iz_SSE1,jz_SSE);
258 dx_SSE2 = gmx_sub_pr(work->ix_SSE2,jx_SSE);
259 dy_SSE2 = gmx_sub_pr(work->iy_SSE2,jy_SSE);
260 dz_SSE2 = gmx_sub_pr(work->iz_SSE2,jz_SSE);
261 dx_SSE3 = gmx_sub_pr(work->ix_SSE3,jx_SSE);
262 dy_SSE3 = gmx_sub_pr(work->iy_SSE3,jy_SSE);
263 dz_SSE3 = gmx_sub_pr(work->iz_SSE3,jz_SSE);
265 /* rsq = dx*dx+dy*dy+dz*dz */
266 rsq_SSE0 = gmx_calc_rsq_pr(dx_SSE0,dy_SSE0,dz_SSE0);
267 rsq_SSE1 = gmx_calc_rsq_pr(dx_SSE1,dy_SSE1,dz_SSE1);
268 rsq_SSE2 = gmx_calc_rsq_pr(dx_SSE2,dy_SSE2,dz_SSE2);
269 rsq_SSE3 = gmx_calc_rsq_pr(dx_SSE3,dy_SSE3,dz_SSE3);
271 wco_SSE0 = gmx_cmplt_pr(rsq_SSE0,rc2_SSE);
272 wco_SSE1 = gmx_cmplt_pr(rsq_SSE1,rc2_SSE);
273 wco_SSE2 = gmx_cmplt_pr(rsq_SSE2,rc2_SSE);
274 wco_SSE3 = gmx_cmplt_pr(rsq_SSE3,rc2_SSE);
276 wco_any_SSE01 = gmx_or_pr(wco_SSE0,wco_SSE1);
277 wco_any_SSE23 = gmx_or_pr(wco_SSE2,wco_SSE3);
278 wco_any_SSE = gmx_or_pr(wco_any_SSE01,wco_any_SSE23);
280 InRange = gmx_movemask_pr(wco_any_SSE);
282 *ndistc += 4*GMX_X86_SIMD_WIDTH_HERE;
292 for(cj=cjf; cj<=cjl; cj++)
294 /* Store cj and the interaction mask */
295 #ifdef GMX_MM128_HERE
296 nbl->cj[nbl->ncj].cj = CI_TO_CJ_S128(gridj->cell0) + cj;
297 nbl->cj[nbl->ncj].excl = get_imask_x86_simd128(remove_sub_diag,ci,cj);
299 nbl->cj[nbl->ncj].cj = CI_TO_CJ_S256(gridj->cell0) + cj;
300 nbl->cj[nbl->ncj].excl = get_imask_x86_simd256(remove_sub_diag,ci,cj);
304 /* Increase the closing index in i super-cell list */
305 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;