Tests of restrained listed potentials.
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / nbnxn_kernel_ref_prune.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 #include "gmxpre.h"
37
38 #include "nbnxn_kernel_ref_prune.h"
39
40 #include "gromacs/mdlib/nbnxn_consts.h"
41 #include "gromacs/mdlib/nbnxn_pairlist.h"
42 #include "gromacs/utility/gmxassert.h"
43
44
45 /* Prune a single NbnxnPairlistCpu entry with distance rlistInner */
46 void
47 nbnxn_kernel_prune_ref(NbnxnPairlistCpu *         nbl,
48                        const nbnxn_atomdata_t *   nbat,
49                        const rvec * gmx_restrict  shift_vec,
50                        real                       rlistInner)
51 {
52     /* We avoid push_back() for efficiency reasons and resize after filling */
53     nbl->ci.resize(nbl->ciOuter.size());
54     nbl->cj.resize(nbl->cjOuter.size());
55
56     const nbnxn_ci_t * gmx_restrict ciOuter  = nbl->ciOuter.data();
57     nbnxn_ci_t       * gmx_restrict ciInner  = nbl->ci.data();
58
59     const nbnxn_cj_t * gmx_restrict cjOuter   = nbl->cjOuter.data();
60     nbnxn_cj_t       * gmx_restrict cjInner   = nbl->cj.data();
61
62     const real       * gmx_restrict shiftvec = shift_vec[0];
63     const real       * gmx_restrict x        = nbat->x().data();
64
65     const real                      rlist2   = rlistInner*rlistInner;
66
67     /* Use compile time constants to speed up the code */
68     constexpr int c_xStride  = 3;
69     GMX_ASSERT(c_xStride == nbat->xstride, "xStride should match nbat->xstride");
70     constexpr int c_xiStride = 3;
71
72     constexpr int c_iUnroll  = c_nbnxnCpuIClusterSize;
73     constexpr int c_jUnroll  = c_nbnxnCpuIClusterSize;
74
75     /* Initialize the new list as empty and add pairs that are in range */
76     int       nciInner = 0;
77     int       ncjInner = 0;
78     const int nciOuter = nbl->ciOuter.size();
79     for (int ciIndex = 0; ciIndex < nciOuter; ciIndex++)
80     {
81         const nbnxn_ci_t * gmx_restrict ciEntry = &ciOuter[ciIndex];
82
83         /* Copy the original list entry to the pruned entry */
84         ciInner[nciInner].ci           = ciEntry->ci;
85         ciInner[nciInner].shift        = ciEntry->shift;
86         ciInner[nciInner].cj_ind_start = ncjInner;
87
88         /* Extract shift data */
89         int ish = (ciEntry->shift & NBNXN_CI_SHIFT);
90         int ci  = ciEntry->ci;
91
92         /* Load i atom coordinates */
93         real xi[c_iUnroll*c_xiStride];
94         for (int i = 0; i < c_iUnroll; i++)
95         {
96             for (int d = 0; d < DIM; d++)
97             {
98                 xi[i*c_xiStride + d] = x[(ci*c_iUnroll + i)*c_xStride + d] + shiftvec[ish*DIM + d];
99             }
100         }
101
102         for (int cjind = ciEntry->cj_ind_start; cjind < ciEntry->cj_ind_end; cjind++)
103         {
104             /* j-cluster index */
105             int  cj        = cjOuter[cjind].cj;
106
107             bool isInRange = false;
108             for (int i = 0; i < c_iUnroll && !isInRange; i++)
109             {
110                 for (int j = 0; j < c_jUnroll; j++)
111                 {
112                     int  aj  = cj*c_jUnroll + j;
113
114                     real dx  = xi[i*c_xiStride + XX] - x[aj*c_xStride + XX];
115                     real dy  = xi[i*c_xiStride + YY] - x[aj*c_xStride + YY];
116                     real dz  = xi[i*c_xiStride + ZZ] - x[aj*c_xStride + ZZ];
117
118                     real rsq = dx*dx + dy*dy + dz*dz;
119
120                     if (rsq < rlist2)
121                     {
122                         isInRange = true;
123                     }
124                 }
125             }
126
127             if (isInRange)
128             {
129                 /* This cluster is in range, put it in the pruned list */
130                 cjInner[ncjInner++] = cjOuter[cjind];
131             }
132         }
133
134         /* Check if there are any j's in the list, if so, add the i-entry */
135         if (ncjInner > ciInner[nciInner].cj_ind_start)
136         {
137             ciInner[nciInner].cj_ind_end = ncjInner;
138             nciInner++;
139         }
140     }
141
142     nbl->ci.resize(nciInner);
143     nbl->cj.resize(ncjInner);
144 }