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