2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
38 #include "kernel_ref_prune.h"
40 #include "gromacs/nbnxm/atomdata.h"
41 #include "gromacs/nbnxm/pairlist.h"
42 #include "gromacs/utility/gmxassert.h"
44 /* Prune a single NbnxnPairlistCpu entry with distance rlistInner */
46 nbnxn_kernel_prune_ref(NbnxnPairlistCpu * nbl,
47 const nbnxn_atomdata_t * nbat,
48 const rvec * gmx_restrict shift_vec,
51 /* We avoid push_back() for efficiency reasons and resize after filling */
52 nbl->ci.resize(nbl->ciOuter.size());
53 nbl->cj.resize(nbl->cjOuter.size());
55 const nbnxn_ci_t * gmx_restrict ciOuter = nbl->ciOuter.data();
56 nbnxn_ci_t * gmx_restrict ciInner = nbl->ci.data();
58 const nbnxn_cj_t * gmx_restrict cjOuter = nbl->cjOuter.data();
59 nbnxn_cj_t * gmx_restrict cjInner = nbl->cj.data();
61 const real * gmx_restrict shiftvec = shift_vec[0];
62 const real * gmx_restrict x = nbat->x().data();
64 const real rlist2 = rlistInner*rlistInner;
66 /* Use compile time constants to speed up the code */
67 constexpr int c_xStride = 3;
68 GMX_ASSERT(c_xStride == nbat->xstride, "xStride should match nbat->xstride");
69 constexpr int c_xiStride = 3;
71 constexpr int c_iUnroll = c_nbnxnCpuIClusterSize;
72 constexpr int c_jUnroll = c_nbnxnCpuIClusterSize;
74 /* Initialize the new list as empty and add pairs that are in range */
77 const int nciOuter = nbl->ciOuter.size();
78 for (int ciIndex = 0; ciIndex < nciOuter; ciIndex++)
80 const nbnxn_ci_t * gmx_restrict ciEntry = &ciOuter[ciIndex];
82 /* Copy the original list entry to the pruned entry */
83 ciInner[nciInner].ci = ciEntry->ci;
84 ciInner[nciInner].shift = ciEntry->shift;
85 ciInner[nciInner].cj_ind_start = ncjInner;
87 /* Extract shift data */
88 int ish = (ciEntry->shift & NBNXN_CI_SHIFT);
91 /* Load i atom coordinates */
92 real xi[c_iUnroll*c_xiStride];
93 for (int i = 0; i < c_iUnroll; i++)
95 for (int d = 0; d < DIM; d++)
97 xi[i*c_xiStride + d] = x[(ci*c_iUnroll + i)*c_xStride + d] + shiftvec[ish*DIM + d];
101 for (int cjind = ciEntry->cj_ind_start; cjind < ciEntry->cj_ind_end; cjind++)
103 /* j-cluster index */
104 int cj = cjOuter[cjind].cj;
106 bool isInRange = false;
107 for (int i = 0; i < c_iUnroll && !isInRange; i++)
109 for (int j = 0; j < c_jUnroll; j++)
111 int aj = cj*c_jUnroll + j;
113 real dx = xi[i*c_xiStride + XX] - x[aj*c_xStride + XX];
114 real dy = xi[i*c_xiStride + YY] - x[aj*c_xStride + YY];
115 real dz = xi[i*c_xiStride + ZZ] - x[aj*c_xStride + ZZ];
117 real rsq = dx*dx + dy*dy + dz*dz;
128 /* This cluster is in range, put it in the pruned list */
129 cjInner[ncjInner++] = cjOuter[cjind];
133 /* Check if there are any j's in the list, if so, add the i-entry */
134 if (ncjInner > ciInner[nciInner].cj_ind_start)
136 ciInner[nciInner].cj_ind_end = ncjInner;
141 nbl->ci.resize(nciInner);
142 nbl->cj.resize(ncjInner);