cb2059a713230ad8ca8b28200fd6a5201ae7a77b
[alexxy/gromacs.git] / src / gromacs / nbnxm / kernels_simd_4xm / kernel_prune.cpp
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,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.
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 #include "gmxpre.h"
38
39 #include "kernel_prune.h"
40
41 #include "gromacs/nbnxm/atomdata.h"
42 #include "gromacs/nbnxm/nbnxm_simd.h"
43 #include "gromacs/nbnxm/pairlist.h"
44 #include "gromacs/utility/gmxassert.h"
45
46 #ifdef GMX_NBNXN_SIMD_4XN
47 #    define GMX_SIMD_J_UNROLL_SIZE 1
48 #    include "kernel_common.h"
49 #endif
50
51 /* Prune a single nbnxn_pairtlist_t entry with distance rlistInner */
52 void nbnxn_kernel_prune_4xn(NbnxnPairlistCpu*              nbl,
53                             const nbnxn_atomdata_t*        nbat,
54                             gmx::ArrayRef<const gmx::RVec> shiftvec,
55                             real                           rlistInner)
56 {
57 #ifdef GMX_NBNXN_SIMD_4XN
58     using namespace gmx;
59
60     /* We avoid push_back() for efficiency reasons and resize after filling */
61     nbl->ci.resize(nbl->ciOuter.size());
62     nbl->cj.resize(nbl->cjOuter.size());
63
64     const nbnxn_ci_t* gmx_restrict ciOuter = nbl->ciOuter.data();
65     nbnxn_ci_t* gmx_restrict ciInner       = nbl->ci.data();
66
67     const nbnxn_cj_t* gmx_restrict cjOuter = nbl->cjOuter.data();
68     nbnxn_cj_t* gmx_restrict cjInner       = nbl->cj.data();
69
70     const real* gmx_restrict x = nbat->x().data();
71
72     const SimdReal rlist2_S(rlistInner * rlistInner);
73
74     /* Initialize the new list count as empty and add pairs that are in range */
75     int       nciInner = 0;
76     int       ncjInner = 0;
77     const int nciOuter = nbl->ciOuter.size();
78     for (int i = 0; i < nciOuter; i++)
79     {
80         const nbnxn_ci_t* gmx_restrict ciEntry = &ciOuter[i];
81
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;
86
87         /* Extract shift data */
88         int ish = (ciEntry->shift & NBNXN_CI_SHIFT);
89         int ci  = ciEntry->ci;
90
91         SimdReal shX_S = SimdReal(shiftvec[ish][XX]);
92         SimdReal shY_S = SimdReal(shiftvec[ish][YY]);
93         SimdReal shZ_S = SimdReal(shiftvec[ish][ZZ]);
94
95 #    if UNROLLJ <= 4
96         int scix = ci * STRIDE * DIM;
97 #    else
98         int scix = (ci >> 1) * STRIDE * DIM + (ci & 1) * (STRIDE >> 1);
99 #    endif
100
101         /* Load i atom data */
102         int      sciy  = scix + STRIDE;
103         int      sciz  = sciy + STRIDE;
104         SimdReal ix_S0 = SimdReal(x[scix]) + shX_S;
105         SimdReal ix_S1 = SimdReal(x[scix + 1]) + shX_S;
106         SimdReal ix_S2 = SimdReal(x[scix + 2]) + shX_S;
107         SimdReal ix_S3 = SimdReal(x[scix + 3]) + shX_S;
108         SimdReal iy_S0 = SimdReal(x[sciy]) + shY_S;
109         SimdReal iy_S1 = SimdReal(x[sciy + 1]) + shY_S;
110         SimdReal iy_S2 = SimdReal(x[sciy + 2]) + shY_S;
111         SimdReal iy_S3 = SimdReal(x[sciy + 3]) + shY_S;
112         SimdReal iz_S0 = SimdReal(x[sciz]) + shZ_S;
113         SimdReal iz_S1 = SimdReal(x[sciz + 1]) + shZ_S;
114         SimdReal iz_S2 = SimdReal(x[sciz + 2]) + shZ_S;
115         SimdReal iz_S3 = SimdReal(x[sciz + 3]) + shZ_S;
116
117         for (int cjind = ciEntry->cj_ind_start; cjind < ciEntry->cj_ind_end; cjind++)
118         {
119             /* j-cluster index */
120             int cj = cjOuter[cjind].cj;
121
122             /* Atom indices (of the first atom in the cluster) */
123 #    if UNROLLJ == STRIDE
124             int aj  = cj * UNROLLJ;
125             int ajx = aj * DIM;
126 #    else
127             int ajx = (cj >> 1) * DIM * STRIDE + (cj & 1) * UNROLLJ;
128 #    endif
129             int ajy = ajx + STRIDE;
130             int ajz = ajy + STRIDE;
131
132             /* load j atom coordinates */
133             SimdReal jx_S = load<SimdReal>(x + ajx);
134             SimdReal jy_S = load<SimdReal>(x + ajy);
135             SimdReal jz_S = load<SimdReal>(x + ajz);
136
137             /* Calculate distance */
138             SimdReal dx_S0 = ix_S0 - jx_S;
139             SimdReal dy_S0 = iy_S0 - jy_S;
140             SimdReal dz_S0 = iz_S0 - jz_S;
141             SimdReal dx_S1 = ix_S1 - jx_S;
142             SimdReal dy_S1 = iy_S1 - jy_S;
143             SimdReal dz_S1 = iz_S1 - jz_S;
144             SimdReal dx_S2 = ix_S2 - jx_S;
145             SimdReal dy_S2 = iy_S2 - jy_S;
146             SimdReal dz_S2 = iz_S2 - jz_S;
147             SimdReal dx_S3 = ix_S3 - jx_S;
148             SimdReal dy_S3 = iy_S3 - jy_S;
149             SimdReal dz_S3 = iz_S3 - jz_S;
150
151             /* rsq = dx*dx+dy*dy+dz*dz */
152             SimdReal rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
153             SimdReal rsq_S1 = norm2(dx_S1, dy_S1, dz_S1);
154             SimdReal rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
155             SimdReal rsq_S3 = norm2(dx_S3, dy_S3, dz_S3);
156
157             /* Do the cut-off check */
158             SimdBool wco_S0 = (rsq_S0 < rlist2_S);
159             SimdBool wco_S1 = (rsq_S1 < rlist2_S);
160             SimdBool wco_S2 = (rsq_S2 < rlist2_S);
161             SimdBool wco_S3 = (rsq_S3 < rlist2_S);
162
163             wco_S0 = wco_S0 || wco_S1;
164             wco_S2 = wco_S2 || wco_S3;
165             wco_S0 = wco_S0 || wco_S2;
166
167             /* Putting the assignment inside the conditional is slower */
168             cjInner[ncjInner] = cjOuter[cjind];
169             if (anyTrue(wco_S0))
170             {
171                 ncjInner++;
172             }
173         }
174
175         if (ncjInner > ciInner[nciInner].cj_ind_start)
176         {
177             ciInner[nciInner].cj_ind_end = ncjInner;
178             nciInner++;
179         }
180     }
181
182     nbl->ci.resize(nciInner);
183     nbl->cj.resize(ncjInner);
184
185 #else /* GMX_NBNXN_SIMD_4XN */
186
187     GMX_RELEASE_ASSERT(false, "4xN kernel called without 4xN support");
188
189     GMX_UNUSED_VALUE(nbl);
190     GMX_UNUSED_VALUE(nbat);
191     GMX_UNUSED_VALUE(shiftvec);
192     GMX_UNUSED_VALUE(rlistInner);
193
194 #endif /* GMX_NBNXN_SIMD_4XN */
195 }