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,2021, 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.
39 #include "kernel_ref_prune.h"
41 #include "gromacs/nbnxm/atomdata.h"
42 #include "gromacs/nbnxm/pairlist.h"
43 #include "gromacs/utility/gmxassert.h"
45 /* Prune a single NbnxnPairlistCpu entry with distance rlistInner */
46 void nbnxn_kernel_prune_ref(NbnxnPairlistCpu* nbl,
47 const nbnxn_atomdata_t* nbat,
48 gmx::ArrayRef<const gmx::RVec> shiftvec,
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 x = nbat->x().data();
63 const real rlist2 = rlistInner * rlistInner;
65 /* Use compile time constants to speed up the code */
66 constexpr int c_xStride = 3;
67 GMX_ASSERT(c_xStride == nbat->xstride, "xStride should match nbat->xstride");
68 constexpr int c_xiStride = 3;
70 constexpr int c_iUnroll = c_nbnxnCpuIClusterSize;
71 constexpr int c_jUnroll = c_nbnxnCpuIClusterSize;
73 /* Initialize the new list as empty and add pairs that are in range */
76 const int nciOuter = nbl->ciOuter.size();
77 for (int ciIndex = 0; ciIndex < nciOuter; ciIndex++)
79 const nbnxn_ci_t* gmx_restrict ciEntry = &ciOuter[ciIndex];
81 /* Copy the original list entry to the pruned entry */
82 ciInner[nciInner].ci = ciEntry->ci;
83 ciInner[nciInner].shift = ciEntry->shift;
84 ciInner[nciInner].cj_ind_start = ncjInner;
86 /* Extract shift data */
87 int ish = (ciEntry->shift & NBNXN_CI_SHIFT);
90 /* Load i atom coordinates */
91 real xi[c_iUnroll * c_xiStride];
92 for (int i = 0; i < c_iUnroll; i++)
94 for (int d = 0; d < DIM; d++)
96 xi[i * c_xiStride + d] = x[(ci * c_iUnroll + i) * c_xStride + d] + shiftvec[ish][d];
100 for (int cjind = ciEntry->cj_ind_start; cjind < ciEntry->cj_ind_end; cjind++)
102 /* j-cluster index */
103 int cj = cjOuter[cjind].cj;
105 bool isInRange = false;
106 for (int i = 0; i < c_iUnroll && !isInRange; i++)
108 for (int j = 0; j < c_jUnroll; j++)
110 int aj = cj * c_jUnroll + j;
112 real dx = xi[i * c_xiStride + XX] - x[aj * c_xStride + XX];
113 real dy = xi[i * c_xiStride + YY] - x[aj * c_xStride + YY];
114 real dz = xi[i * c_xiStride + ZZ] - x[aj * c_xStride + ZZ];
116 real rsq = dx * dx + dy * dy + dz * dz;
127 /* This cluster is in range, put it in the pruned list */
128 cjInner[ncjInner++] = cjOuter[cjind];
132 /* Check if there are any j's in the list, if so, add the i-entry */
133 if (ncjInner > ciInner[nciInner].cj_ind_start)
135 ciInner[nciInner].cj_ind_end = ncjInner;
140 nbl->ci.resize(nciInner);
141 nbl->cj.resize(ncjInner);