Update instructions in containers.rst
[gromacs.git] / src / gromacs / nbnxm / pairlist_simd_4xm.h
blobce8cf3f28f446204653cefa837bbaf5e7f8e8a96
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5 * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
37 /*! \internal \file
39 * \brief
40 * Declares inline-friendly code for making 4xN pairlists
42 * \author Berk Hess <hess@kth.se>
43 * \ingroup module_nbnxm
46 //! Stride of the packed x coordinate array
47 static constexpr int c_xStride4xN =
48 (GMX_SIMD_REAL_WIDTH > c_nbnxnCpuIClusterSize ? GMX_SIMD_REAL_WIDTH : c_nbnxnCpuIClusterSize);
50 //! Copies PBC shifted i-cell packed atom coordinates to working array
51 static inline void icell_set_x_simd_4xn(int ci,
52 real shx,
53 real shy,
54 real shz,
55 int gmx_unused stride,
56 const real* x,
57 NbnxnPairlistCpuWork* work)
59 int ia;
60 real* x_ci_simd = work->iClusterData.xSimd.data();
62 ia = xIndexFromCi<NbnxnLayout::Simd4xN>(ci);
64 store(x_ci_simd + 0 * GMX_SIMD_REAL_WIDTH, SimdReal(x[ia + 0 * c_xStride4xN] + shx));
65 store(x_ci_simd + 1 * GMX_SIMD_REAL_WIDTH, SimdReal(x[ia + 1 * c_xStride4xN] + shy));
66 store(x_ci_simd + 2 * GMX_SIMD_REAL_WIDTH, SimdReal(x[ia + 2 * c_xStride4xN] + shz));
67 store(x_ci_simd + 3 * GMX_SIMD_REAL_WIDTH, SimdReal(x[ia + 0 * c_xStride4xN + 1] + shx));
68 store(x_ci_simd + 4 * GMX_SIMD_REAL_WIDTH, SimdReal(x[ia + 1 * c_xStride4xN + 1] + shy));
69 store(x_ci_simd + 5 * GMX_SIMD_REAL_WIDTH, SimdReal(x[ia + 2 * c_xStride4xN + 1] + shz));
70 store(x_ci_simd + 6 * GMX_SIMD_REAL_WIDTH, SimdReal(x[ia + 0 * c_xStride4xN + 2] + shx));
71 store(x_ci_simd + 7 * GMX_SIMD_REAL_WIDTH, SimdReal(x[ia + 1 * c_xStride4xN + 2] + shy));
72 store(x_ci_simd + 8 * GMX_SIMD_REAL_WIDTH, SimdReal(x[ia + 2 * c_xStride4xN + 2] + shz));
73 store(x_ci_simd + 9 * GMX_SIMD_REAL_WIDTH, SimdReal(x[ia + 0 * c_xStride4xN + 3] + shx));
74 store(x_ci_simd + 10 * GMX_SIMD_REAL_WIDTH, SimdReal(x[ia + 1 * c_xStride4xN + 3] + shy));
75 store(x_ci_simd + 11 * GMX_SIMD_REAL_WIDTH, SimdReal(x[ia + 2 * c_xStride4xN + 3] + shz));
78 /*! \brief SIMD code for checking and adding cluster-pairs to the list using coordinates in packed format.
80 * Checks bouding box distances and possibly atom pair distances.
81 * This is an accelerated version of make_cluster_list_simple.
83 * \param[in] jGrid The j-grid
84 * \param[in,out] nbl The pair-list to store the cluster pairs in
85 * \param[in] icluster The index of the i-cluster
86 * \param[in] firstCell The first cluster in the j-range, using i-cluster size indexing
87 * \param[in] lastCell The last cluster in the j-range, using i-cluster size indexing
88 * \param[in] excludeSubDiagonal Exclude atom pairs with i-index > j-index
89 * \param[in] x_j Coordinates for the j-atom, in SIMD packed format
90 * \param[in] rlist2 The squared list cut-off
91 * \param[in] rbb2 The squared cut-off for putting cluster-pairs in the list based on bounding box distance only
92 * \param[in,out] numDistanceChecks The number of distance checks performed
94 static inline void makeClusterListSimd4xn(const Grid& jGrid,
95 NbnxnPairlistCpu* nbl,
96 int icluster,
97 int firstCell,
98 int lastCell,
99 bool excludeSubDiagonal,
100 const real* gmx_restrict x_j,
101 real rlist2,
102 float rbb2,
103 int* gmx_restrict numDistanceChecks)
105 using namespace gmx;
106 const real* gmx_restrict x_ci_simd = nbl->work->iClusterData.xSimd.data();
107 const BoundingBox* gmx_restrict bb_ci = nbl->work->iClusterData.bb.data();
109 SimdReal jx_S, jy_S, jz_S;
111 SimdReal dx_S0, dy_S0, dz_S0;
112 SimdReal dx_S1, dy_S1, dz_S1;
113 SimdReal dx_S2, dy_S2, dz_S2;
114 SimdReal dx_S3, dy_S3, dz_S3;
116 SimdReal rsq_S0;
117 SimdReal rsq_S1;
118 SimdReal rsq_S2;
119 SimdReal rsq_S3;
121 SimdBool wco_S0;
122 SimdBool wco_S1;
123 SimdBool wco_S2;
124 SimdBool wco_S3;
125 SimdBool wco_any_S01, wco_any_S23, wco_any_S;
127 SimdReal rc2_S;
129 gmx_bool InRange;
130 float d2;
131 int xind_f, xind_l;
133 /* Convert the j-range from i-cluster size indexing to j-cluster indexing */
134 int jclusterFirst = cjFromCi<NbnxnLayout::Simd4xN, 0>(firstCell);
135 int jclusterLast = cjFromCi<NbnxnLayout::Simd4xN, 1>(lastCell);
136 GMX_ASSERT(jclusterLast >= jclusterFirst,
137 "We should have a non-empty j-cluster range, since the calling code should have "
138 "ensured a non-empty cell range");
140 rc2_S = SimdReal(rlist2);
142 InRange = FALSE;
143 while (!InRange && jclusterFirst <= jclusterLast)
145 d2 = clusterBoundingBoxDistance2(bb_ci[0], jGrid.jBoundingBoxes()[jclusterFirst]);
146 *numDistanceChecks += 2;
148 /* Check if the distance is within the distance where
149 * we use only the bounding box distance rbb,
150 * or within the cut-off and there is at least one atom pair
151 * within the cut-off.
153 if (d2 < rbb2)
155 InRange = TRUE;
157 else if (d2 < rlist2)
159 xind_f = xIndexFromCj<NbnxnLayout::Simd4xN>(
160 cjFromCi<NbnxnLayout::Simd4xN, 0>(jGrid.cellOffset()) + jclusterFirst);
162 jx_S = load<SimdReal>(x_j + xind_f + 0 * c_xStride4xN);
163 jy_S = load<SimdReal>(x_j + xind_f + 1 * c_xStride4xN);
164 jz_S = load<SimdReal>(x_j + xind_f + 2 * c_xStride4xN);
167 /* Calculate distance */
168 dx_S0 = load<SimdReal>(x_ci_simd + 0 * GMX_SIMD_REAL_WIDTH) - jx_S;
169 dy_S0 = load<SimdReal>(x_ci_simd + 1 * GMX_SIMD_REAL_WIDTH) - jy_S;
170 dz_S0 = load<SimdReal>(x_ci_simd + 2 * GMX_SIMD_REAL_WIDTH) - jz_S;
171 dx_S1 = load<SimdReal>(x_ci_simd + 3 * GMX_SIMD_REAL_WIDTH) - jx_S;
172 dy_S1 = load<SimdReal>(x_ci_simd + 4 * GMX_SIMD_REAL_WIDTH) - jy_S;
173 dz_S1 = load<SimdReal>(x_ci_simd + 5 * GMX_SIMD_REAL_WIDTH) - jz_S;
174 dx_S2 = load<SimdReal>(x_ci_simd + 6 * GMX_SIMD_REAL_WIDTH) - jx_S;
175 dy_S2 = load<SimdReal>(x_ci_simd + 7 * GMX_SIMD_REAL_WIDTH) - jy_S;
176 dz_S2 = load<SimdReal>(x_ci_simd + 8 * GMX_SIMD_REAL_WIDTH) - jz_S;
177 dx_S3 = load<SimdReal>(x_ci_simd + 9 * GMX_SIMD_REAL_WIDTH) - jx_S;
178 dy_S3 = load<SimdReal>(x_ci_simd + 10 * GMX_SIMD_REAL_WIDTH) - jy_S;
179 dz_S3 = load<SimdReal>(x_ci_simd + 11 * GMX_SIMD_REAL_WIDTH) - jz_S;
181 /* rsq = dx*dx+dy*dy+dz*dz */
182 rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
183 rsq_S1 = norm2(dx_S1, dy_S1, dz_S1);
184 rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
185 rsq_S3 = norm2(dx_S3, dy_S3, dz_S3);
187 wco_S0 = (rsq_S0 < rc2_S);
188 wco_S1 = (rsq_S1 < rc2_S);
189 wco_S2 = (rsq_S2 < rc2_S);
190 wco_S3 = (rsq_S3 < rc2_S);
192 wco_any_S01 = wco_S0 || wco_S1;
193 wco_any_S23 = wco_S2 || wco_S3;
194 wco_any_S = wco_any_S01 || wco_any_S23;
196 InRange = anyTrue(wco_any_S);
198 *numDistanceChecks += 4 * GMX_SIMD_REAL_WIDTH;
200 if (!InRange)
202 jclusterFirst++;
205 if (!InRange)
207 return;
210 InRange = FALSE;
211 while (!InRange && jclusterLast > jclusterFirst)
213 d2 = clusterBoundingBoxDistance2(bb_ci[0], jGrid.jBoundingBoxes()[jclusterLast]);
214 *numDistanceChecks += 2;
216 /* Check if the distance is within the distance where
217 * we use only the bounding box distance rbb,
218 * or within the cut-off and there is at least one atom pair
219 * within the cut-off.
221 if (d2 < rbb2)
223 InRange = TRUE;
225 else if (d2 < rlist2)
227 xind_l = xIndexFromCj<NbnxnLayout::Simd4xN>(
228 cjFromCi<NbnxnLayout::Simd4xN, 0>(jGrid.cellOffset()) + jclusterLast);
230 jx_S = load<SimdReal>(x_j + xind_l + 0 * c_xStride4xN);
231 jy_S = load<SimdReal>(x_j + xind_l + 1 * c_xStride4xN);
232 jz_S = load<SimdReal>(x_j + xind_l + 2 * c_xStride4xN);
234 /* Calculate distance */
235 dx_S0 = load<SimdReal>(x_ci_simd + 0 * GMX_SIMD_REAL_WIDTH) - jx_S;
236 dy_S0 = load<SimdReal>(x_ci_simd + 1 * GMX_SIMD_REAL_WIDTH) - jy_S;
237 dz_S0 = load<SimdReal>(x_ci_simd + 2 * GMX_SIMD_REAL_WIDTH) - jz_S;
238 dx_S1 = load<SimdReal>(x_ci_simd + 3 * GMX_SIMD_REAL_WIDTH) - jx_S;
239 dy_S1 = load<SimdReal>(x_ci_simd + 4 * GMX_SIMD_REAL_WIDTH) - jy_S;
240 dz_S1 = load<SimdReal>(x_ci_simd + 5 * GMX_SIMD_REAL_WIDTH) - jz_S;
241 dx_S2 = load<SimdReal>(x_ci_simd + 6 * GMX_SIMD_REAL_WIDTH) - jx_S;
242 dy_S2 = load<SimdReal>(x_ci_simd + 7 * GMX_SIMD_REAL_WIDTH) - jy_S;
243 dz_S2 = load<SimdReal>(x_ci_simd + 8 * GMX_SIMD_REAL_WIDTH) - jz_S;
244 dx_S3 = load<SimdReal>(x_ci_simd + 9 * GMX_SIMD_REAL_WIDTH) - jx_S;
245 dy_S3 = load<SimdReal>(x_ci_simd + 10 * GMX_SIMD_REAL_WIDTH) - jy_S;
246 dz_S3 = load<SimdReal>(x_ci_simd + 11 * GMX_SIMD_REAL_WIDTH) - jz_S;
248 /* rsq = dx*dx+dy*dy+dz*dz */
249 rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
250 rsq_S1 = norm2(dx_S1, dy_S1, dz_S1);
251 rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
252 rsq_S3 = norm2(dx_S3, dy_S3, dz_S3);
254 wco_S0 = (rsq_S0 < rc2_S);
255 wco_S1 = (rsq_S1 < rc2_S);
256 wco_S2 = (rsq_S2 < rc2_S);
257 wco_S3 = (rsq_S3 < rc2_S);
259 wco_any_S01 = wco_S0 || wco_S1;
260 wco_any_S23 = wco_S2 || wco_S3;
261 wco_any_S = wco_any_S01 || wco_any_S23;
263 InRange = anyTrue(wco_any_S);
265 *numDistanceChecks += 4 * GMX_SIMD_REAL_WIDTH;
267 if (!InRange)
269 jclusterLast--;
273 if (jclusterFirst <= jclusterLast)
275 for (int jcluster = jclusterFirst; jcluster <= jclusterLast; jcluster++)
277 /* Store cj and the interaction mask */
278 nbnxn_cj_t cjEntry;
279 cjEntry.cj = cjFromCi<NbnxnLayout::Simd4xN, 0>(jGrid.cellOffset()) + jcluster;
280 cjEntry.excl = get_imask_simd_4xn(excludeSubDiagonal, icluster, jcluster);
281 nbl->cj.push_back(cjEntry);
283 /* Increase the closing index in the i list */
284 nbl->ci.back().cj_ind_end = nbl->cj.size();