2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, 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.
36 /* Stride of the packed x coordinate array */
37 static constexpr int c_xStride2xNN
= c_nbnxnCpuIClusterSize
;
39 /* Copies PBC shifted i-cell packed atom coordinates to working array */
41 icell_set_x_simd_2xnn(int ci
,
42 real shx
, real shy
, real shz
,
43 int gmx_unused stride
, const real
*x
,
44 NbnxnPairlistCpuWork
*work
)
47 real
*x_ci_simd
= work
->iClusterData
.xSimd
.data();
49 ia
= xIndexFromCi
<NbnxnLayout::Simd2xNN
>(ci
);
51 store(x_ci_simd
+ 0*GMX_SIMD_REAL_WIDTH
, loadU1DualHsimd(x
+ ia
+ 0*c_xStride2xNN
+ 0) + SimdReal(shx
) );
52 store(x_ci_simd
+ 1*GMX_SIMD_REAL_WIDTH
, loadU1DualHsimd(x
+ ia
+ 1*c_xStride2xNN
+ 0) + SimdReal(shy
) );
53 store(x_ci_simd
+ 2*GMX_SIMD_REAL_WIDTH
, loadU1DualHsimd(x
+ ia
+ 2*c_xStride2xNN
+ 0) + SimdReal(shz
) );
54 store(x_ci_simd
+ 3*GMX_SIMD_REAL_WIDTH
, loadU1DualHsimd(x
+ ia
+ 0*c_xStride2xNN
+ 2) + SimdReal(shx
) );
55 store(x_ci_simd
+ 4*GMX_SIMD_REAL_WIDTH
, loadU1DualHsimd(x
+ ia
+ 1*c_xStride2xNN
+ 2) + SimdReal(shy
) );
56 store(x_ci_simd
+ 5*GMX_SIMD_REAL_WIDTH
, loadU1DualHsimd(x
+ ia
+ 2*c_xStride2xNN
+ 2) + SimdReal(shz
) );
59 /* SIMD code for checking and adding cluster-pairs to the list using coordinates in packed format.
61 * Checks bouding box distances and possibly atom pair distances.
62 * This is an accelerated version of make_cluster_list_simple.
64 * \param[in] jGrid The j-grid
65 * \param[in,out] nbl The pair-list to store the cluster pairs in
66 * \param[in] icluster The index of the i-cluster
67 * \param[in] firstCell The first cluster in the j-range, using i-cluster size indexing
68 * \param[in] lastCell The last cluster in the j-range, using i-cluster size indexing
69 * \param[in] excludeSubDiagonal Exclude atom pairs with i-index > j-index
70 * \param[in] x_j Coordinates for the j-atom, in SIMD packed format
71 * \param[in] rlist2 The squared list cut-off
72 * \param[in] rbb2 The squared cut-off for putting cluster-pairs in the list based on bounding box distance only
73 * \param[in,out] numDistanceChecks The number of distance checks performed
76 makeClusterListSimd2xnn(const Grid
&jGrid
,
77 NbnxnPairlistCpu
* nbl
,
81 bool excludeSubDiagonal
,
82 const real
* gmx_restrict x_j
,
85 int * gmx_restrict numDistanceChecks
)
88 const real
* gmx_restrict x_ci_simd
= nbl
->work
->iClusterData
.xSimd
.data();
89 const BoundingBox
* gmx_restrict bb_ci
= nbl
->work
->iClusterData
.bb
.data();
91 SimdReal jx_S
, jy_S
, jz_S
;
93 SimdReal dx_S0
, dy_S0
, dz_S0
;
94 SimdReal dx_S2
, dy_S2
, dz_S2
;
109 int jclusterFirst
= cjFromCi
<NbnxnLayout::Simd2xNN
, 0>(firstCell
);
110 int jclusterLast
= cjFromCi
<NbnxnLayout::Simd2xNN
, 1>(lastCell
);
111 GMX_ASSERT(jclusterLast
>= jclusterFirst
, "We should have a non-empty j-cluster range, since the calling code should have ensured a non-empty cell range");
113 rc2_S
= SimdReal(rlist2
);
116 while (!InRange
&& jclusterFirst
<= jclusterLast
)
118 d2
= clusterBoundingBoxDistance2(bb_ci
[0], jGrid
.jBoundingBoxes()[jclusterFirst
]);
119 *numDistanceChecks
+= 2;
121 /* Check if the distance is within the distance where
122 * we use only the bounding box distance rbb,
123 * or within the cut-off and there is at least one atom pair
124 * within the cut-off.
130 else if (d2
< rlist2
)
132 xind_f
= xIndexFromCj
<NbnxnLayout::Simd2xNN
>(cjFromCi
<NbnxnLayout::Simd2xNN
, 0>(jGrid
.cellOffset()) + jclusterFirst
);
134 jx_S
= loadDuplicateHsimd(x_j
+ xind_f
+ 0*c_xStride2xNN
);
135 jy_S
= loadDuplicateHsimd(x_j
+ xind_f
+ 1*c_xStride2xNN
);
136 jz_S
= loadDuplicateHsimd(x_j
+ xind_f
+ 2*c_xStride2xNN
);
138 /* Calculate distance */
139 dx_S0
= load
<SimdReal
>(x_ci_simd
+ 0*GMX_SIMD_REAL_WIDTH
) - jx_S
;
140 dy_S0
= load
<SimdReal
>(x_ci_simd
+ 1*GMX_SIMD_REAL_WIDTH
) - jy_S
;
141 dz_S0
= load
<SimdReal
>(x_ci_simd
+ 2*GMX_SIMD_REAL_WIDTH
) - jz_S
;
142 dx_S2
= load
<SimdReal
>(x_ci_simd
+ 3*GMX_SIMD_REAL_WIDTH
) - jx_S
;
143 dy_S2
= load
<SimdReal
>(x_ci_simd
+ 4*GMX_SIMD_REAL_WIDTH
) - jy_S
;
144 dz_S2
= load
<SimdReal
>(x_ci_simd
+ 5*GMX_SIMD_REAL_WIDTH
) - jz_S
;
146 /* rsq = dx*dx+dy*dy+dz*dz */
147 rsq_S0
= norm2(dx_S0
, dy_S0
, dz_S0
);
148 rsq_S2
= norm2(dx_S2
, dy_S2
, dz_S2
);
150 wco_S0
= (rsq_S0
< rc2_S
);
151 wco_S2
= (rsq_S2
< rc2_S
);
153 wco_any_S
= wco_S0
|| wco_S2
;
155 InRange
= anyTrue(wco_any_S
);
157 *numDistanceChecks
+= 2*GMX_SIMD_REAL_WIDTH
;
170 while (!InRange
&& jclusterLast
> jclusterFirst
)
172 d2
= clusterBoundingBoxDistance2(bb_ci
[0], jGrid
.jBoundingBoxes()[jclusterLast
]);
173 *numDistanceChecks
+= 2;
175 /* Check if the distance is within the distance where
176 * we use only the bounding box distance rbb,
177 * or within the cut-off and there is at least one atom pair
178 * within the cut-off.
184 else if (d2
< rlist2
)
186 xind_l
= xIndexFromCj
<NbnxnLayout::Simd2xNN
>(cjFromCi
<NbnxnLayout::Simd2xNN
, 0>(jGrid
.cellOffset()) + jclusterLast
);
188 jx_S
= loadDuplicateHsimd(x_j
+ xind_l
+ 0*c_xStride2xNN
);
189 jy_S
= loadDuplicateHsimd(x_j
+ xind_l
+ 1*c_xStride2xNN
);
190 jz_S
= loadDuplicateHsimd(x_j
+ xind_l
+ 2*c_xStride2xNN
);
192 /* Calculate distance */
193 dx_S0
= load
<SimdReal
>(x_ci_simd
+ 0*GMX_SIMD_REAL_WIDTH
) - jx_S
;
194 dy_S0
= load
<SimdReal
>(x_ci_simd
+ 1*GMX_SIMD_REAL_WIDTH
) - jy_S
;
195 dz_S0
= load
<SimdReal
>(x_ci_simd
+ 2*GMX_SIMD_REAL_WIDTH
) - jz_S
;
196 dx_S2
= load
<SimdReal
>(x_ci_simd
+ 3*GMX_SIMD_REAL_WIDTH
) - jx_S
;
197 dy_S2
= load
<SimdReal
>(x_ci_simd
+ 4*GMX_SIMD_REAL_WIDTH
) - jy_S
;
198 dz_S2
= load
<SimdReal
>(x_ci_simd
+ 5*GMX_SIMD_REAL_WIDTH
) - jz_S
;
200 /* rsq = dx*dx+dy*dy+dz*dz */
201 rsq_S0
= norm2(dx_S0
, dy_S0
, dz_S0
);
202 rsq_S2
= norm2(dx_S2
, dy_S2
, dz_S2
);
204 wco_S0
= (rsq_S0
< rc2_S
);
205 wco_S2
= (rsq_S2
< rc2_S
);
207 wco_any_S
= wco_S0
|| wco_S2
;
209 InRange
= anyTrue(wco_any_S
);
211 *numDistanceChecks
+= 2*GMX_SIMD_REAL_WIDTH
;
219 if (jclusterFirst
<= jclusterLast
)
221 for (int jcluster
= jclusterFirst
; jcluster
<= jclusterLast
; jcluster
++)
223 /* Store cj and the interaction mask */
225 cjEntry
.cj
= cjFromCi
<NbnxnLayout::Simd2xNN
, 0>(jGrid
.cellOffset()) + jcluster
;
226 cjEntry
.excl
= get_imask_simd_2xnn(excludeSubDiagonal
, icluster
, jcluster
);
227 nbl
->cj
.push_back(cjEntry
);
229 /* Increase the closing index in the i list */
230 nbl
->ci
.back().cj_ind_end
= nbl
->cj
.size();