2 * This file is part of the GROMACS molecular simulation package.
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,2013, 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.
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.
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.
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.
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.
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.
40 #if GMX_SIMD_WIDTH_HERE >= NBNXN_CPU_CLUSTER_I_SIZE
41 #define STRIDE_S (GMX_SIMD_WIDTH_HERE)
43 #define STRIDE_S NBNXN_CPU_CLUSTER_I_SIZE
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,
51 int stride, const real *x,
52 nbnxn_list_work_t *work)
55 nbnxn_x_ci_simd_4xn_t *x_ci;
57 x_ci = work->x_ci_simd_4xn;
59 ia = X_IND_CI_SIMD_4XN(ci);
61 x_ci->ix_S0 = gmx_set1_pr(x[ia + 0*STRIDE_S ] + shx);
62 x_ci->iy_S0 = gmx_set1_pr(x[ia + 1*STRIDE_S ] + shy);
63 x_ci->iz_S0 = gmx_set1_pr(x[ia + 2*STRIDE_S ] + shz);
64 x_ci->ix_S1 = gmx_set1_pr(x[ia + 0*STRIDE_S + 1] + shx);
65 x_ci->iy_S1 = gmx_set1_pr(x[ia + 1*STRIDE_S + 1] + shy);
66 x_ci->iz_S1 = gmx_set1_pr(x[ia + 2*STRIDE_S + 1] + shz);
67 x_ci->ix_S2 = gmx_set1_pr(x[ia + 0*STRIDE_S + 2] + shx);
68 x_ci->iy_S2 = gmx_set1_pr(x[ia + 1*STRIDE_S + 2] + shy);
69 x_ci->iz_S2 = gmx_set1_pr(x[ia + 2*STRIDE_S + 2] + shz);
70 x_ci->ix_S3 = gmx_set1_pr(x[ia + 0*STRIDE_S + 3] + shx);
71 x_ci->iy_S3 = gmx_set1_pr(x[ia + 1*STRIDE_S + 3] + shy);
72 x_ci->iz_S3 = gmx_set1_pr(x[ia + 2*STRIDE_S + 3] + shz);
75 /* SIMD code for making a pair list of cell ci vs cell cjf-cjl
76 * for coordinates in packed format.
77 * Checks bouding box distances and possibly atom pair distances.
78 * This is an accelerated version of make_cluster_list_simple.
80 static gmx_inline void
81 make_cluster_list_simd_4xn(const nbnxn_grid_t *gridj,
82 nbnxn_pairlist_t *nbl,
83 int ci, int cjf, int cjl,
84 gmx_bool remove_sub_diag,
89 const nbnxn_x_ci_simd_4xn_t *work;
90 const nbnxn_bb_t *bb_ci;
92 gmx_mm_pr jx_S, jy_S, jz_S;
94 gmx_mm_pr dx_S0, dy_S0, dz_S0;
95 gmx_mm_pr dx_S1, dy_S1, dz_S1;
96 gmx_mm_pr dx_S2, dy_S2, dz_S2;
97 gmx_mm_pr dx_S3, dy_S3, dz_S3;
108 gmx_mm_pb wco_any_S01, wco_any_S23, wco_any_S;
114 int xind_f, xind_l, cj;
116 cjf = CI_TO_CJ_SIMD_4XN(cjf);
117 cjl = CI_TO_CJ_SIMD_4XN(cjl+1) - 1;
119 work = nbl->work->x_ci_simd_4xn;
121 bb_ci = nbl->work->bb_ci;
123 rc2_S = gmx_set1_pr(rl2);
126 while (!InRange && cjf <= cjl)
128 #ifdef NBNXN_SEARCH_BB_SIMD4
129 d2 = subc_bb_dist2_simd4(0, bb_ci, cjf, gridj->bbj);
131 d2 = subc_bb_dist2(0, bb_ci, cjf, gridj->bbj);
135 /* Check if the distance is within the distance where
136 * we use only the bounding box distance rbb,
137 * or within the cut-off and there is at least one atom pair
138 * within the cut-off.
146 xind_f = X_IND_CJ_SIMD_4XN(CI_TO_CJ_SIMD_4XN(gridj->cell0) + cjf);
148 jx_S = gmx_load_pr(x_j+xind_f+0*STRIDE_S);
149 jy_S = gmx_load_pr(x_j+xind_f+1*STRIDE_S);
150 jz_S = gmx_load_pr(x_j+xind_f+2*STRIDE_S);
153 /* Calculate distance */
154 dx_S0 = gmx_sub_pr(work->ix_S0, jx_S);
155 dy_S0 = gmx_sub_pr(work->iy_S0, jy_S);
156 dz_S0 = gmx_sub_pr(work->iz_S0, jz_S);
157 dx_S1 = gmx_sub_pr(work->ix_S1, jx_S);
158 dy_S1 = gmx_sub_pr(work->iy_S1, jy_S);
159 dz_S1 = gmx_sub_pr(work->iz_S1, jz_S);
160 dx_S2 = gmx_sub_pr(work->ix_S2, jx_S);
161 dy_S2 = gmx_sub_pr(work->iy_S2, jy_S);
162 dz_S2 = gmx_sub_pr(work->iz_S2, jz_S);
163 dx_S3 = gmx_sub_pr(work->ix_S3, jx_S);
164 dy_S3 = gmx_sub_pr(work->iy_S3, jy_S);
165 dz_S3 = gmx_sub_pr(work->iz_S3, jz_S);
167 /* rsq = dx*dx+dy*dy+dz*dz */
168 rsq_S0 = gmx_calc_rsq_pr(dx_S0, dy_S0, dz_S0);
169 rsq_S1 = gmx_calc_rsq_pr(dx_S1, dy_S1, dz_S1);
170 rsq_S2 = gmx_calc_rsq_pr(dx_S2, dy_S2, dz_S2);
171 rsq_S3 = gmx_calc_rsq_pr(dx_S3, dy_S3, dz_S3);
173 wco_S0 = gmx_cmplt_pr(rsq_S0, rc2_S);
174 wco_S1 = gmx_cmplt_pr(rsq_S1, rc2_S);
175 wco_S2 = gmx_cmplt_pr(rsq_S2, rc2_S);
176 wco_S3 = gmx_cmplt_pr(rsq_S3, rc2_S);
178 wco_any_S01 = gmx_or_pb(wco_S0, wco_S1);
179 wco_any_S23 = gmx_or_pb(wco_S2, wco_S3);
180 wco_any_S = gmx_or_pb(wco_any_S01, wco_any_S23);
182 InRange = gmx_anytrue_pb(wco_any_S);
184 *ndistc += 4*GMX_SIMD_WIDTH_HERE;
197 while (!InRange && cjl > cjf)
199 #ifdef NBNXN_SEARCH_BB_SIMD4
200 d2 = subc_bb_dist2_simd4(0, bb_ci, cjl, gridj->bbj);
202 d2 = subc_bb_dist2(0, bb_ci, cjl, gridj->bbj);
206 /* Check if the distance is within the distance where
207 * we use only the bounding box distance rbb,
208 * or within the cut-off and there is at least one atom pair
209 * within the cut-off.
217 xind_l = X_IND_CJ_SIMD_4XN(CI_TO_CJ_SIMD_4XN(gridj->cell0) + cjl);
219 jx_S = gmx_load_pr(x_j+xind_l+0*STRIDE_S);
220 jy_S = gmx_load_pr(x_j+xind_l+1*STRIDE_S);
221 jz_S = gmx_load_pr(x_j+xind_l+2*STRIDE_S);
223 /* Calculate distance */
224 dx_S0 = gmx_sub_pr(work->ix_S0, jx_S);
225 dy_S0 = gmx_sub_pr(work->iy_S0, jy_S);
226 dz_S0 = gmx_sub_pr(work->iz_S0, jz_S);
227 dx_S1 = gmx_sub_pr(work->ix_S1, jx_S);
228 dy_S1 = gmx_sub_pr(work->iy_S1, jy_S);
229 dz_S1 = gmx_sub_pr(work->iz_S1, jz_S);
230 dx_S2 = gmx_sub_pr(work->ix_S2, jx_S);
231 dy_S2 = gmx_sub_pr(work->iy_S2, jy_S);
232 dz_S2 = gmx_sub_pr(work->iz_S2, jz_S);
233 dx_S3 = gmx_sub_pr(work->ix_S3, jx_S);
234 dy_S3 = gmx_sub_pr(work->iy_S3, jy_S);
235 dz_S3 = gmx_sub_pr(work->iz_S3, jz_S);
237 /* rsq = dx*dx+dy*dy+dz*dz */
238 rsq_S0 = gmx_calc_rsq_pr(dx_S0, dy_S0, dz_S0);
239 rsq_S1 = gmx_calc_rsq_pr(dx_S1, dy_S1, dz_S1);
240 rsq_S2 = gmx_calc_rsq_pr(dx_S2, dy_S2, dz_S2);
241 rsq_S3 = gmx_calc_rsq_pr(dx_S3, dy_S3, dz_S3);
243 wco_S0 = gmx_cmplt_pr(rsq_S0, rc2_S);
244 wco_S1 = gmx_cmplt_pr(rsq_S1, rc2_S);
245 wco_S2 = gmx_cmplt_pr(rsq_S2, rc2_S);
246 wco_S3 = gmx_cmplt_pr(rsq_S3, rc2_S);
248 wco_any_S01 = gmx_or_pb(wco_S0, wco_S1);
249 wco_any_S23 = gmx_or_pb(wco_S2, wco_S3);
250 wco_any_S = gmx_or_pb(wco_any_S01, wco_any_S23);
252 InRange = gmx_anytrue_pb(wco_any_S);
254 *ndistc += 4*GMX_SIMD_WIDTH_HERE;
264 for (cj = cjf; cj <= cjl; cj++)
266 /* Store cj and the interaction mask */
267 nbl->cj[nbl->ncj].cj = CI_TO_CJ_SIMD_4XN(gridj->cell0) + cj;
268 nbl->cj[nbl->ncj].excl = get_imask_simd_4xn(remove_sub_diag, ci, cj);
269 #ifdef GMX_CPU_ACCELERATION_IBM_QPX
270 nbl->cj[nbl->ncj].interaction_mask_indices[0] = (nbl->cj[nbl->ncj].excl & 0x000F) >> (0 * 4);
271 nbl->cj[nbl->ncj].interaction_mask_indices[1] = (nbl->cj[nbl->ncj].excl & 0x00F0) >> (1 * 4);
272 nbl->cj[nbl->ncj].interaction_mask_indices[2] = (nbl->cj[nbl->ncj].excl & 0x0F00) >> (2 * 4);
273 nbl->cj[nbl->ncj].interaction_mask_indices[3] = (nbl->cj[nbl->ncj].excl & 0xF000) >> (3 * 4);
277 /* Increase the closing index in i super-cell list */
278 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;