Split lines with many copyright years
[alexxy/gromacs.git] / src / gromacs / nbnxm / pairlist_simd_2xmm.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
9  *
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.
14  *
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.
19  *
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.
24  *
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.
32  *
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.
35  */
36
37 /* Stride of the packed x coordinate array */
38 static constexpr int c_xStride2xNN = (GMX_SIMD_REAL_WIDTH >= 2 * c_nbnxnCpuIClusterSize)
39                                              ? GMX_SIMD_REAL_WIDTH / 2
40                                              : c_nbnxnCpuIClusterSize;
41
42 /* Copies PBC shifted i-cell packed atom coordinates to working array */
43 static inline void icell_set_x_simd_2xnn(int  ci,
44                                          real shx,
45                                          real shy,
46                                          real shz,
47                                          int gmx_unused        stride,
48                                          const real*           x,
49                                          NbnxnPairlistCpuWork* work)
50 {
51     int   ia;
52     real* x_ci_simd = work->iClusterData.xSimd.data();
53
54     ia = xIndexFromCi<NbnxnLayout::Simd2xNN>(ci);
55
56     store(x_ci_simd + 0 * GMX_SIMD_REAL_WIDTH,
57           loadU1DualHsimd(x + ia + 0 * c_xStride2xNN + 0) + SimdReal(shx));
58     store(x_ci_simd + 1 * GMX_SIMD_REAL_WIDTH,
59           loadU1DualHsimd(x + ia + 1 * c_xStride2xNN + 0) + SimdReal(shy));
60     store(x_ci_simd + 2 * GMX_SIMD_REAL_WIDTH,
61           loadU1DualHsimd(x + ia + 2 * c_xStride2xNN + 0) + SimdReal(shz));
62     store(x_ci_simd + 3 * GMX_SIMD_REAL_WIDTH,
63           loadU1DualHsimd(x + ia + 0 * c_xStride2xNN + 2) + SimdReal(shx));
64     store(x_ci_simd + 4 * GMX_SIMD_REAL_WIDTH,
65           loadU1DualHsimd(x + ia + 1 * c_xStride2xNN + 2) + SimdReal(shy));
66     store(x_ci_simd + 5 * GMX_SIMD_REAL_WIDTH,
67           loadU1DualHsimd(x + ia + 2 * c_xStride2xNN + 2) + SimdReal(shz));
68 }
69
70 /* SIMD code for checking and adding cluster-pairs to the list using coordinates in packed format.
71  *
72  * Checks bouding box distances and possibly atom pair distances.
73  * This is an accelerated version of make_cluster_list_simple.
74  *
75  * \param[in]     jGrid               The j-grid
76  * \param[in,out] nbl                 The pair-list to store the cluster pairs in
77  * \param[in]     icluster            The index of the i-cluster
78  * \param[in]     firstCell           The first cluster in the j-range, using i-cluster size indexing
79  * \param[in]     lastCell            The last cluster in the j-range, using i-cluster size indexing
80  * \param[in]     excludeSubDiagonal  Exclude atom pairs with i-index > j-index
81  * \param[in]     x_j                 Coordinates for the j-atom, in SIMD packed format
82  * \param[in]     rlist2              The squared list cut-off
83  * \param[in]     rbb2                The squared cut-off for putting cluster-pairs in the list based on bounding box distance only
84  * \param[in,out] numDistanceChecks   The number of distance checks performed
85  */
86 static inline void makeClusterListSimd2xnn(const Grid&       jGrid,
87                                            NbnxnPairlistCpu* nbl,
88                                            int               icluster,
89                                            int               firstCell,
90                                            int               lastCell,
91                                            bool              excludeSubDiagonal,
92                                            const real* gmx_restrict x_j,
93                                            real                     rlist2,
94                                            float                    rbb2,
95                                            int* gmx_restrict numDistanceChecks)
96 {
97     using namespace gmx;
98     const real* gmx_restrict x_ci_simd    = nbl->work->iClusterData.xSimd.data();
99     const BoundingBox* gmx_restrict bb_ci = nbl->work->iClusterData.bb.data();
100
101     SimdReal jx_S, jy_S, jz_S;
102
103     SimdReal dx_S0, dy_S0, dz_S0;
104     SimdReal dx_S2, dy_S2, dz_S2;
105
106     SimdReal rsq_S0;
107     SimdReal rsq_S2;
108
109     SimdBool wco_S0;
110     SimdBool wco_S2;
111     SimdBool wco_any_S;
112
113     SimdReal rc2_S;
114
115     gmx_bool InRange;
116     float    d2;
117     int      xind_f, xind_l;
118
119     int jclusterFirst = cjFromCi<NbnxnLayout::Simd2xNN, 0>(firstCell);
120     int jclusterLast  = cjFromCi<NbnxnLayout::Simd2xNN, 1>(lastCell);
121     GMX_ASSERT(jclusterLast >= jclusterFirst,
122                "We should have a non-empty j-cluster range, since the calling code should have "
123                "ensured a non-empty cell range");
124
125     rc2_S = SimdReal(rlist2);
126
127     InRange = FALSE;
128     while (!InRange && jclusterFirst <= jclusterLast)
129     {
130         d2 = clusterBoundingBoxDistance2(bb_ci[0], jGrid.jBoundingBoxes()[jclusterFirst]);
131         *numDistanceChecks += 2;
132
133         /* Check if the distance is within the distance where
134          * we use only the bounding box distance rbb,
135          * or within the cut-off and there is at least one atom pair
136          * within the cut-off.
137          */
138         if (d2 < rbb2)
139         {
140             InRange = TRUE;
141         }
142         else if (d2 < rlist2)
143         {
144             xind_f = xIndexFromCj<NbnxnLayout::Simd2xNN>(
145                     cjFromCi<NbnxnLayout::Simd2xNN, 0>(jGrid.cellOffset()) + jclusterFirst);
146
147             jx_S = loadDuplicateHsimd(x_j + xind_f + 0 * c_xStride2xNN);
148             jy_S = loadDuplicateHsimd(x_j + xind_f + 1 * c_xStride2xNN);
149             jz_S = loadDuplicateHsimd(x_j + xind_f + 2 * c_xStride2xNN);
150
151             /* Calculate distance */
152             dx_S0 = load<SimdReal>(x_ci_simd + 0 * GMX_SIMD_REAL_WIDTH) - jx_S;
153             dy_S0 = load<SimdReal>(x_ci_simd + 1 * GMX_SIMD_REAL_WIDTH) - jy_S;
154             dz_S0 = load<SimdReal>(x_ci_simd + 2 * GMX_SIMD_REAL_WIDTH) - jz_S;
155             dx_S2 = load<SimdReal>(x_ci_simd + 3 * GMX_SIMD_REAL_WIDTH) - jx_S;
156             dy_S2 = load<SimdReal>(x_ci_simd + 4 * GMX_SIMD_REAL_WIDTH) - jy_S;
157             dz_S2 = load<SimdReal>(x_ci_simd + 5 * GMX_SIMD_REAL_WIDTH) - jz_S;
158
159             /* rsq = dx*dx+dy*dy+dz*dz */
160             rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
161             rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
162
163             wco_S0 = (rsq_S0 < rc2_S);
164             wco_S2 = (rsq_S2 < rc2_S);
165
166             wco_any_S = wco_S0 || wco_S2;
167
168             InRange = anyTrue(wco_any_S);
169
170             *numDistanceChecks += 2 * GMX_SIMD_REAL_WIDTH;
171         }
172         if (!InRange)
173         {
174             jclusterFirst++;
175         }
176     }
177     if (!InRange)
178     {
179         return;
180     }
181
182     InRange = FALSE;
183     while (!InRange && jclusterLast > jclusterFirst)
184     {
185         d2 = clusterBoundingBoxDistance2(bb_ci[0], jGrid.jBoundingBoxes()[jclusterLast]);
186         *numDistanceChecks += 2;
187
188         /* Check if the distance is within the distance where
189          * we use only the bounding box distance rbb,
190          * or within the cut-off and there is at least one atom pair
191          * within the cut-off.
192          */
193         if (d2 < rbb2)
194         {
195             InRange = TRUE;
196         }
197         else if (d2 < rlist2)
198         {
199             xind_l = xIndexFromCj<NbnxnLayout::Simd2xNN>(
200                     cjFromCi<NbnxnLayout::Simd2xNN, 0>(jGrid.cellOffset()) + jclusterLast);
201
202             jx_S = loadDuplicateHsimd(x_j + xind_l + 0 * c_xStride2xNN);
203             jy_S = loadDuplicateHsimd(x_j + xind_l + 1 * c_xStride2xNN);
204             jz_S = loadDuplicateHsimd(x_j + xind_l + 2 * c_xStride2xNN);
205
206             /* Calculate distance */
207             dx_S0 = load<SimdReal>(x_ci_simd + 0 * GMX_SIMD_REAL_WIDTH) - jx_S;
208             dy_S0 = load<SimdReal>(x_ci_simd + 1 * GMX_SIMD_REAL_WIDTH) - jy_S;
209             dz_S0 = load<SimdReal>(x_ci_simd + 2 * GMX_SIMD_REAL_WIDTH) - jz_S;
210             dx_S2 = load<SimdReal>(x_ci_simd + 3 * GMX_SIMD_REAL_WIDTH) - jx_S;
211             dy_S2 = load<SimdReal>(x_ci_simd + 4 * GMX_SIMD_REAL_WIDTH) - jy_S;
212             dz_S2 = load<SimdReal>(x_ci_simd + 5 * GMX_SIMD_REAL_WIDTH) - jz_S;
213
214             /* rsq = dx*dx+dy*dy+dz*dz */
215             rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
216             rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
217
218             wco_S0 = (rsq_S0 < rc2_S);
219             wco_S2 = (rsq_S2 < rc2_S);
220
221             wco_any_S = wco_S0 || wco_S2;
222
223             InRange = anyTrue(wco_any_S);
224
225             *numDistanceChecks += 2 * GMX_SIMD_REAL_WIDTH;
226         }
227         if (!InRange)
228         {
229             jclusterLast--;
230         }
231     }
232
233     if (jclusterFirst <= jclusterLast)
234     {
235         for (int jcluster = jclusterFirst; jcluster <= jclusterLast; jcluster++)
236         {
237             /* Store cj and the interaction mask */
238             nbnxn_cj_t cjEntry;
239             cjEntry.cj   = cjFromCi<NbnxnLayout::Simd2xNN, 0>(jGrid.cellOffset()) + jcluster;
240             cjEntry.excl = get_imask_simd_2xnn(excludeSubDiagonal, icluster, jcluster);
241             nbl->cj.push_back(cjEntry);
242         }
243         /* Increase the closing index in the i list */
244         nbl->ci.back().cj_ind_end = nbl->cj.size();
245     }
246 }