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.
39 #if !(GMX_NBNXN_SIMD_BITWIDTH == 128 || GMX_NBNXN_SIMD_BITWIDTH == 256)
40 #error "unsupported GMX_NBNXN_SIMD_BITWIDTH"
43 #ifdef GMX_NBNXN_HALF_WIDTH_SIMD
44 #define GMX_USE_HALF_WIDTH_SIMD_HERE
46 #include "gmx_simd_macros.h"
48 #if GMX_SIMD_WIDTH_HERE >= NBNXN_CPU_CLUSTER_I_SIZE
49 #define STRIDE_S (GMX_SIMD_WIDTH_HERE)
51 #define STRIDE_S NBNXN_CPU_CLUSTER_I_SIZE
54 /* Copies PBC shifted i-cell packed atom coordinates to working array */
55 static gmx_inline void
56 icell_set_x_simd_4xn(int ci,
57 real shx, real shy, real shz,
59 int stride, const real *x,
60 nbnxn_list_work_t *work)
63 nbnxn_x_ci_simd_4xn_t *x_ci;
65 x_ci = work->x_ci_simd_4xn;
67 ia = X_IND_CI_SIMD_4XN(ci);
69 x_ci->ix_SSE0 = gmx_set1_pr(x[ia + 0*STRIDE_S ] + shx);
70 x_ci->iy_SSE0 = gmx_set1_pr(x[ia + 1*STRIDE_S ] + shy);
71 x_ci->iz_SSE0 = gmx_set1_pr(x[ia + 2*STRIDE_S ] + shz);
72 x_ci->ix_SSE1 = gmx_set1_pr(x[ia + 0*STRIDE_S + 1] + shx);
73 x_ci->iy_SSE1 = gmx_set1_pr(x[ia + 1*STRIDE_S + 1] + shy);
74 x_ci->iz_SSE1 = gmx_set1_pr(x[ia + 2*STRIDE_S + 1] + shz);
75 x_ci->ix_SSE2 = gmx_set1_pr(x[ia + 0*STRIDE_S + 2] + shx);
76 x_ci->iy_SSE2 = gmx_set1_pr(x[ia + 1*STRIDE_S + 2] + shy);
77 x_ci->iz_SSE2 = gmx_set1_pr(x[ia + 2*STRIDE_S + 2] + shz);
78 x_ci->ix_SSE3 = gmx_set1_pr(x[ia + 0*STRIDE_S + 3] + shx);
79 x_ci->iy_SSE3 = gmx_set1_pr(x[ia + 1*STRIDE_S + 3] + shy);
80 x_ci->iz_SSE3 = gmx_set1_pr(x[ia + 2*STRIDE_S + 3] + shz);
83 #ifndef GMX_HAVE_SIMD_ANYTRUE
84 /* Fallback function in case gmx_anytrue_pr is not present */
85 static gmx_inline gmx_bool
86 gmx_anytrue_4xn_pr(gmx_mm_pr bool_S)
88 real bools_array[2*GMX_SIMD_WIDTH_HERE], *bools;
92 bools = gmx_simd_align_real(bools_array);
94 gmx_store_pr(bools, bool_S);
97 for (s = 0; s < GMX_SIMD_WIDTH_HERE; s++)
99 if (GMX_SIMD_IS_TRUE(bools[s]))
109 /* SIMD code for making a pair list of cell ci vs cell cjf-cjl
110 * for coordinates in packed format.
111 * Checks bouding box distances and possibly atom pair distances.
112 * This is an accelerated version of make_cluster_list_simple.
114 static gmx_inline void
115 make_cluster_list_simd_4xn(const nbnxn_grid_t *gridj,
116 nbnxn_pairlist_t *nbl,
117 int ci, int cjf, int cjl,
118 gmx_bool remove_sub_diag,
120 real rl2, float rbb2,
123 const nbnxn_x_ci_simd_4xn_t *work;
126 gmx_mm_pr jx_SSE, jy_SSE, jz_SSE;
128 gmx_mm_pr dx_SSE0, dy_SSE0, dz_SSE0;
129 gmx_mm_pr dx_SSE1, dy_SSE1, dz_SSE1;
130 gmx_mm_pr dx_SSE2, dy_SSE2, dz_SSE2;
131 gmx_mm_pr dx_SSE3, dy_SSE3, dz_SSE3;
142 gmx_mm_pr wco_any_SSE01, wco_any_SSE23, wco_any_SSE;
148 int xind_f, xind_l, cj;
150 cjf = CI_TO_CJ_SIMD_4XN(cjf);
151 cjl = CI_TO_CJ_SIMD_4XN(cjl+1) - 1;
153 work = nbl->work->x_ci_simd_4xn;
155 bb_ci = nbl->work->bb_ci;
157 rc2_SSE = gmx_set1_pr(rl2);
160 while (!InRange && cjf <= cjl)
162 d2 = subc_bb_dist2_sse(4, 0, bb_ci, cjf, gridj->bbj);
165 /* Check if the distance is within the distance where
166 * we use only the bounding box distance rbb,
167 * or within the cut-off and there is at least one atom pair
168 * within the cut-off.
176 xind_f = X_IND_CJ_SIMD_4XN(CI_TO_CJ_SIMD_4XN(gridj->cell0) + cjf);
178 jx_SSE = gmx_load_pr(x_j+xind_f+0*STRIDE_S);
179 jy_SSE = gmx_load_pr(x_j+xind_f+1*STRIDE_S);
180 jz_SSE = gmx_load_pr(x_j+xind_f+2*STRIDE_S);
183 /* Calculate distance */
184 dx_SSE0 = gmx_sub_pr(work->ix_SSE0, jx_SSE);
185 dy_SSE0 = gmx_sub_pr(work->iy_SSE0, jy_SSE);
186 dz_SSE0 = gmx_sub_pr(work->iz_SSE0, jz_SSE);
187 dx_SSE1 = gmx_sub_pr(work->ix_SSE1, jx_SSE);
188 dy_SSE1 = gmx_sub_pr(work->iy_SSE1, jy_SSE);
189 dz_SSE1 = gmx_sub_pr(work->iz_SSE1, jz_SSE);
190 dx_SSE2 = gmx_sub_pr(work->ix_SSE2, jx_SSE);
191 dy_SSE2 = gmx_sub_pr(work->iy_SSE2, jy_SSE);
192 dz_SSE2 = gmx_sub_pr(work->iz_SSE2, jz_SSE);
193 dx_SSE3 = gmx_sub_pr(work->ix_SSE3, jx_SSE);
194 dy_SSE3 = gmx_sub_pr(work->iy_SSE3, jy_SSE);
195 dz_SSE3 = gmx_sub_pr(work->iz_SSE3, jz_SSE);
197 /* rsq = dx*dx+dy*dy+dz*dz */
198 rsq_SSE0 = gmx_calc_rsq_pr(dx_SSE0, dy_SSE0, dz_SSE0);
199 rsq_SSE1 = gmx_calc_rsq_pr(dx_SSE1, dy_SSE1, dz_SSE1);
200 rsq_SSE2 = gmx_calc_rsq_pr(dx_SSE2, dy_SSE2, dz_SSE2);
201 rsq_SSE3 = gmx_calc_rsq_pr(dx_SSE3, dy_SSE3, dz_SSE3);
203 wco_SSE0 = gmx_cmplt_pr(rsq_SSE0, rc2_SSE);
204 wco_SSE1 = gmx_cmplt_pr(rsq_SSE1, rc2_SSE);
205 wco_SSE2 = gmx_cmplt_pr(rsq_SSE2, rc2_SSE);
206 wco_SSE3 = gmx_cmplt_pr(rsq_SSE3, rc2_SSE);
208 wco_any_SSE01 = gmx_or_pr(wco_SSE0, wco_SSE1);
209 wco_any_SSE23 = gmx_or_pr(wco_SSE2, wco_SSE3);
210 wco_any_SSE = gmx_or_pr(wco_any_SSE01, wco_any_SSE23);
212 #ifdef GMX_HAVE_SIMD_ANYTRUE
213 InRange = gmx_anytrue_pr(wco_any_SSE);
215 InRange = gmx_anytrue_4xn_pr(wco_any_SSE);
218 *ndistc += 4*GMX_SIMD_WIDTH_HERE;
231 while (!InRange && cjl > cjf)
233 d2 = subc_bb_dist2_sse(4, 0, bb_ci, cjl, gridj->bbj);
236 /* Check if the distance is within the distance where
237 * we use only the bounding box distance rbb,
238 * or within the cut-off and there is at least one atom pair
239 * within the cut-off.
247 xind_l = X_IND_CJ_SIMD_4XN(CI_TO_CJ_SIMD_4XN(gridj->cell0) + cjl);
249 jx_SSE = gmx_load_pr(x_j+xind_l+0*STRIDE_S);
250 jy_SSE = gmx_load_pr(x_j+xind_l+1*STRIDE_S);
251 jz_SSE = gmx_load_pr(x_j+xind_l+2*STRIDE_S);
253 /* Calculate distance */
254 dx_SSE0 = gmx_sub_pr(work->ix_SSE0, jx_SSE);
255 dy_SSE0 = gmx_sub_pr(work->iy_SSE0, jy_SSE);
256 dz_SSE0 = gmx_sub_pr(work->iz_SSE0, jz_SSE);
257 dx_SSE1 = gmx_sub_pr(work->ix_SSE1, jx_SSE);
258 dy_SSE1 = gmx_sub_pr(work->iy_SSE1, jy_SSE);
259 dz_SSE1 = gmx_sub_pr(work->iz_SSE1, jz_SSE);
260 dx_SSE2 = gmx_sub_pr(work->ix_SSE2, jx_SSE);
261 dy_SSE2 = gmx_sub_pr(work->iy_SSE2, jy_SSE);
262 dz_SSE2 = gmx_sub_pr(work->iz_SSE2, jz_SSE);
263 dx_SSE3 = gmx_sub_pr(work->ix_SSE3, jx_SSE);
264 dy_SSE3 = gmx_sub_pr(work->iy_SSE3, jy_SSE);
265 dz_SSE3 = gmx_sub_pr(work->iz_SSE3, jz_SSE);
267 /* rsq = dx*dx+dy*dy+dz*dz */
268 rsq_SSE0 = gmx_calc_rsq_pr(dx_SSE0, dy_SSE0, dz_SSE0);
269 rsq_SSE1 = gmx_calc_rsq_pr(dx_SSE1, dy_SSE1, dz_SSE1);
270 rsq_SSE2 = gmx_calc_rsq_pr(dx_SSE2, dy_SSE2, dz_SSE2);
271 rsq_SSE3 = gmx_calc_rsq_pr(dx_SSE3, dy_SSE3, dz_SSE3);
273 wco_SSE0 = gmx_cmplt_pr(rsq_SSE0, rc2_SSE);
274 wco_SSE1 = gmx_cmplt_pr(rsq_SSE1, rc2_SSE);
275 wco_SSE2 = gmx_cmplt_pr(rsq_SSE2, rc2_SSE);
276 wco_SSE3 = gmx_cmplt_pr(rsq_SSE3, rc2_SSE);
278 wco_any_SSE01 = gmx_or_pr(wco_SSE0, wco_SSE1);
279 wco_any_SSE23 = gmx_or_pr(wco_SSE2, wco_SSE3);
280 wco_any_SSE = gmx_or_pr(wco_any_SSE01, wco_any_SSE23);
282 #ifdef GMX_HAVE_SIMD_ANYTRUE
283 InRange = gmx_anytrue_pr(wco_any_SSE);
285 InRange = gmx_anytrue_4xn_pr(wco_any_SSE);
288 *ndistc += 4*GMX_SIMD_WIDTH_HERE;
298 for (cj = cjf; cj <= cjl; cj++)
300 /* Store cj and the interaction mask */
301 nbl->cj[nbl->ncj].cj = CI_TO_CJ_SIMD_4XN(gridj->cell0) + cj;
302 nbl->cj[nbl->ncj].excl = get_imask_simd_4xn(remove_sub_diag, ci, cj);
305 /* Increase the closing index in i super-cell list */
306 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
311 #undef GMX_USE_HALF_WIDTH_SIMD_HERE