2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 #if GMX_SIMD_REAL_WIDTH >= NBNXN_CPU_CLUSTER_I_SIZE
38 #define STRIDE_S (GMX_SIMD_REAL_WIDTH)
40 #define STRIDE_S NBNXN_CPU_CLUSTER_I_SIZE
43 /* Copies PBC shifted i-cell packed atom coordinates to working array */
44 static gmx_inline void
45 icell_set_x_simd_4xn(int ci,
46 real shx, real shy, real shz,
48 int gmx_unused stride, const real *x,
49 nbnxn_list_work_t *work)
52 nbnxn_x_ci_simd_4xn_t *x_ci;
54 x_ci = work->x_ci_simd_4xn;
56 ia = X_IND_CI_SIMD_4XN(ci);
58 x_ci->ix_S0 = gmx_simd_set1_r(x[ia + 0*STRIDE_S ] + shx);
59 x_ci->iy_S0 = gmx_simd_set1_r(x[ia + 1*STRIDE_S ] + shy);
60 x_ci->iz_S0 = gmx_simd_set1_r(x[ia + 2*STRIDE_S ] + shz);
61 x_ci->ix_S1 = gmx_simd_set1_r(x[ia + 0*STRIDE_S + 1] + shx);
62 x_ci->iy_S1 = gmx_simd_set1_r(x[ia + 1*STRIDE_S + 1] + shy);
63 x_ci->iz_S1 = gmx_simd_set1_r(x[ia + 2*STRIDE_S + 1] + shz);
64 x_ci->ix_S2 = gmx_simd_set1_r(x[ia + 0*STRIDE_S + 2] + shx);
65 x_ci->iy_S2 = gmx_simd_set1_r(x[ia + 1*STRIDE_S + 2] + shy);
66 x_ci->iz_S2 = gmx_simd_set1_r(x[ia + 2*STRIDE_S + 2] + shz);
67 x_ci->ix_S3 = gmx_simd_set1_r(x[ia + 0*STRIDE_S + 3] + shx);
68 x_ci->iy_S3 = gmx_simd_set1_r(x[ia + 1*STRIDE_S + 3] + shy);
69 x_ci->iz_S3 = gmx_simd_set1_r(x[ia + 2*STRIDE_S + 3] + shz);
72 /* SIMD code for making a pair list of cell ci vs cell cjf-cjl
73 * for coordinates in packed format.
74 * Checks bouding box distances and possibly atom pair distances.
75 * This is an accelerated version of make_cluster_list_simple.
77 static gmx_inline void
78 make_cluster_list_simd_4xn(const nbnxn_grid_t *gridj,
79 nbnxn_pairlist_t *nbl,
80 int ci, int cjf, int cjl,
81 gmx_bool remove_sub_diag,
86 const nbnxn_x_ci_simd_4xn_t *work;
87 const nbnxn_bb_t *bb_ci;
89 gmx_simd_real_t jx_S, jy_S, jz_S;
91 gmx_simd_real_t dx_S0, dy_S0, dz_S0;
92 gmx_simd_real_t dx_S1, dy_S1, dz_S1;
93 gmx_simd_real_t dx_S2, dy_S2, dz_S2;
94 gmx_simd_real_t dx_S3, dy_S3, dz_S3;
96 gmx_simd_real_t rsq_S0;
97 gmx_simd_real_t rsq_S1;
98 gmx_simd_real_t rsq_S2;
99 gmx_simd_real_t rsq_S3;
101 gmx_simd_bool_t wco_S0;
102 gmx_simd_bool_t wco_S1;
103 gmx_simd_bool_t wco_S2;
104 gmx_simd_bool_t wco_S3;
105 gmx_simd_bool_t wco_any_S01, wco_any_S23, wco_any_S;
107 gmx_simd_real_t rc2_S;
111 int xind_f, xind_l, cj;
113 cjf = CI_TO_CJ_SIMD_4XN(cjf);
114 cjl = CI_TO_CJ_SIMD_4XN(cjl+1) - 1;
116 work = nbl->work->x_ci_simd_4xn;
118 bb_ci = nbl->work->bb_ci;
120 rc2_S = gmx_simd_set1_r(rl2);
123 while (!InRange && cjf <= cjl)
125 #ifdef NBNXN_SEARCH_BB_SIMD4
126 d2 = subc_bb_dist2_simd4(0, bb_ci, cjf, gridj->bbj);
128 d2 = subc_bb_dist2(0, bb_ci, cjf, gridj->bbj);
132 /* Check if the distance is within the distance where
133 * we use only the bounding box distance rbb,
134 * or within the cut-off and there is at least one atom pair
135 * within the cut-off.
143 xind_f = X_IND_CJ_SIMD_4XN(CI_TO_CJ_SIMD_4XN(gridj->cell0) + cjf);
145 jx_S = gmx_simd_load_r(x_j+xind_f+0*STRIDE_S);
146 jy_S = gmx_simd_load_r(x_j+xind_f+1*STRIDE_S);
147 jz_S = gmx_simd_load_r(x_j+xind_f+2*STRIDE_S);
150 /* Calculate distance */
151 dx_S0 = gmx_simd_sub_r(work->ix_S0, jx_S);
152 dy_S0 = gmx_simd_sub_r(work->iy_S0, jy_S);
153 dz_S0 = gmx_simd_sub_r(work->iz_S0, jz_S);
154 dx_S1 = gmx_simd_sub_r(work->ix_S1, jx_S);
155 dy_S1 = gmx_simd_sub_r(work->iy_S1, jy_S);
156 dz_S1 = gmx_simd_sub_r(work->iz_S1, jz_S);
157 dx_S2 = gmx_simd_sub_r(work->ix_S2, jx_S);
158 dy_S2 = gmx_simd_sub_r(work->iy_S2, jy_S);
159 dz_S2 = gmx_simd_sub_r(work->iz_S2, jz_S);
160 dx_S3 = gmx_simd_sub_r(work->ix_S3, jx_S);
161 dy_S3 = gmx_simd_sub_r(work->iy_S3, jy_S);
162 dz_S3 = gmx_simd_sub_r(work->iz_S3, jz_S);
164 /* rsq = dx*dx+dy*dy+dz*dz */
165 rsq_S0 = gmx_simd_calc_rsq_r(dx_S0, dy_S0, dz_S0);
166 rsq_S1 = gmx_simd_calc_rsq_r(dx_S1, dy_S1, dz_S1);
167 rsq_S2 = gmx_simd_calc_rsq_r(dx_S2, dy_S2, dz_S2);
168 rsq_S3 = gmx_simd_calc_rsq_r(dx_S3, dy_S3, dz_S3);
170 wco_S0 = gmx_simd_cmplt_r(rsq_S0, rc2_S);
171 wco_S1 = gmx_simd_cmplt_r(rsq_S1, rc2_S);
172 wco_S2 = gmx_simd_cmplt_r(rsq_S2, rc2_S);
173 wco_S3 = gmx_simd_cmplt_r(rsq_S3, rc2_S);
175 wco_any_S01 = gmx_simd_or_b(wco_S0, wco_S1);
176 wco_any_S23 = gmx_simd_or_b(wco_S2, wco_S3);
177 wco_any_S = gmx_simd_or_b(wco_any_S01, wco_any_S23);
179 InRange = gmx_simd_anytrue_b(wco_any_S);
181 *ndistc += 4*GMX_SIMD_REAL_WIDTH;
194 while (!InRange && cjl > cjf)
196 #ifdef NBNXN_SEARCH_BB_SIMD4
197 d2 = subc_bb_dist2_simd4(0, bb_ci, cjl, gridj->bbj);
199 d2 = subc_bb_dist2(0, bb_ci, cjl, gridj->bbj);
203 /* Check if the distance is within the distance where
204 * we use only the bounding box distance rbb,
205 * or within the cut-off and there is at least one atom pair
206 * within the cut-off.
214 xind_l = X_IND_CJ_SIMD_4XN(CI_TO_CJ_SIMD_4XN(gridj->cell0) + cjl);
216 jx_S = gmx_simd_load_r(x_j+xind_l+0*STRIDE_S);
217 jy_S = gmx_simd_load_r(x_j+xind_l+1*STRIDE_S);
218 jz_S = gmx_simd_load_r(x_j+xind_l+2*STRIDE_S);
220 /* Calculate distance */
221 dx_S0 = gmx_simd_sub_r(work->ix_S0, jx_S);
222 dy_S0 = gmx_simd_sub_r(work->iy_S0, jy_S);
223 dz_S0 = gmx_simd_sub_r(work->iz_S0, jz_S);
224 dx_S1 = gmx_simd_sub_r(work->ix_S1, jx_S);
225 dy_S1 = gmx_simd_sub_r(work->iy_S1, jy_S);
226 dz_S1 = gmx_simd_sub_r(work->iz_S1, jz_S);
227 dx_S2 = gmx_simd_sub_r(work->ix_S2, jx_S);
228 dy_S2 = gmx_simd_sub_r(work->iy_S2, jy_S);
229 dz_S2 = gmx_simd_sub_r(work->iz_S2, jz_S);
230 dx_S3 = gmx_simd_sub_r(work->ix_S3, jx_S);
231 dy_S3 = gmx_simd_sub_r(work->iy_S3, jy_S);
232 dz_S3 = gmx_simd_sub_r(work->iz_S3, jz_S);
234 /* rsq = dx*dx+dy*dy+dz*dz */
235 rsq_S0 = gmx_simd_calc_rsq_r(dx_S0, dy_S0, dz_S0);
236 rsq_S1 = gmx_simd_calc_rsq_r(dx_S1, dy_S1, dz_S1);
237 rsq_S2 = gmx_simd_calc_rsq_r(dx_S2, dy_S2, dz_S2);
238 rsq_S3 = gmx_simd_calc_rsq_r(dx_S3, dy_S3, dz_S3);
240 wco_S0 = gmx_simd_cmplt_r(rsq_S0, rc2_S);
241 wco_S1 = gmx_simd_cmplt_r(rsq_S1, rc2_S);
242 wco_S2 = gmx_simd_cmplt_r(rsq_S2, rc2_S);
243 wco_S3 = gmx_simd_cmplt_r(rsq_S3, rc2_S);
245 wco_any_S01 = gmx_simd_or_b(wco_S0, wco_S1);
246 wco_any_S23 = gmx_simd_or_b(wco_S2, wco_S3);
247 wco_any_S = gmx_simd_or_b(wco_any_S01, wco_any_S23);
249 InRange = gmx_simd_anytrue_b(wco_any_S);
251 *ndistc += 4*GMX_SIMD_REAL_WIDTH;
261 for (cj = cjf; cj <= cjl; cj++)
263 /* Store cj and the interaction mask */
264 nbl->cj[nbl->ncj].cj = CI_TO_CJ_SIMD_4XN(gridj->cell0) + cj;
265 nbl->cj[nbl->ncj].excl = get_imask_simd_4xn(remove_sub_diag, ci, cj);
266 #ifdef GMX_SIMD_IBM_QPX
267 nbl->cj[nbl->ncj].interaction_mask_indices[0] = (nbl->cj[nbl->ncj].excl & 0x000F) >> (0 * 4);
268 nbl->cj[nbl->ncj].interaction_mask_indices[1] = (nbl->cj[nbl->ncj].excl & 0x00F0) >> (1 * 4);
269 nbl->cj[nbl->ncj].interaction_mask_indices[2] = (nbl->cj[nbl->ncj].excl & 0x0F00) >> (2 * 4);
270 nbl->cj[nbl->ncj].interaction_mask_indices[3] = (nbl->cj[nbl->ncj].excl & 0xF000) >> (3 * 4);
274 /* Increase the closing index in i super-cell list */
275 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;