Enable more warnings for Clang 6
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / simd_2xnn / nbnxn_kernel_simd_2xnn_prune.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018, 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_2xnn_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_2XNN
45 #define GMX_SIMD_J_UNROLL_SIZE 2
46 #include "gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn_common.h"
47 #endif
48
49 /* Prune a single nbnxn_pairtlist_t entry with distance rlistInner */
50 void
51 nbnxn_kernel_prune_2xnn(nbnxn_pairlist_t *         nbl,
52                         const nbnxn_atomdata_t *   nbat,
53                         const rvec * gmx_restrict  shift_vec,
54                         real                       rlistInner)
55 {
56 #ifdef GMX_NBNXN_SIMD_2XNN
57     using namespace gmx;
58     const nbnxn_ci_t * gmx_restrict ciOuter  = nbl->ciOuter;
59     nbnxn_ci_t       * gmx_restrict ciInner  = nbl->ci;
60
61     const nbnxn_cj_t * gmx_restrict cjOuter  = nbl->cjOuter;
62     nbnxn_cj_t       * gmx_restrict cjInner  = nbl->cj;
63
64     const real       * gmx_restrict shiftvec = shift_vec[0];
65     const real       * gmx_restrict x        = nbat->x;
66
67     const SimdReal                  rlist2_S(rlistInner*rlistInner);
68
69     /* Initialize the new list as empty and add pairs that are in range */
70     int nciInner = 0;
71     int ncjInner = 0;
72     for (int i = 0; i < nbl->nciOuter; i++)
73     {
74         const nbnxn_ci_t * gmx_restrict ciEntry = &ciOuter[i];
75
76         /* Copy the original list entry to the pruned entry */
77         ciInner[nciInner].ci           = ciEntry->ci;
78         ciInner[nciInner].shift        = ciEntry->shift;
79         ciInner[nciInner].cj_ind_start = ncjInner;
80
81         /* Extract shift data */
82         int      ish     = (ciEntry->shift & NBNXN_CI_SHIFT);
83         int      ish3    = ish*3;
84         int      ci      = ciEntry->ci;
85
86         SimdReal shX_S   = SimdReal(shiftvec[ish3    ]);
87         SimdReal shY_S   = SimdReal(shiftvec[ish3 + 1]);
88         SimdReal shZ_S   = SimdReal(shiftvec[ish3 + 2]);
89
90 #if UNROLLJ <= 4
91         int      scix    = ci*STRIDE*DIM;
92 #else
93         int      scix    = (ci >> 1)*STRIDE*DIM + (ci & 1)*(STRIDE >> 1);
94 #endif
95
96         /* Load i atom data */
97         int      sciy    = scix + STRIDE;
98         int      sciz    = sciy + STRIDE;
99         SimdReal ix_S0   = loadU1DualHsimd(x + scix    ) + shX_S;
100         SimdReal ix_S2   = loadU1DualHsimd(x + scix + 2) + shX_S;
101         SimdReal iy_S0   = loadU1DualHsimd(x + sciy    ) + shY_S;
102         SimdReal iy_S2   = loadU1DualHsimd(x + sciy + 2) + shY_S;
103         SimdReal iz_S0   = loadU1DualHsimd(x + sciz    ) + shZ_S;
104         SimdReal iz_S2   = loadU1DualHsimd(x + sciz + 2) + shZ_S;
105
106         for (int cjind = ciEntry->cj_ind_start; cjind < ciEntry->cj_ind_end; cjind++)
107         {
108             /* j-cluster index */
109             int cj      = cjOuter[cjind].cj;
110
111             /* Atom indices (of the first atom in the cluster) */
112 #if UNROLLJ == STRIDE
113             int aj      = cj*UNROLLJ;
114             int ajx     = aj*DIM;
115 #else
116             int ajx     = (cj >> 1)*DIM*STRIDE + (cj & 1)*UNROLLJ;
117 #endif
118             int ajy     = ajx + STRIDE;
119             int ajz     = ajy + STRIDE;
120
121             /* load j atom coordinates */
122             SimdReal jx_S   = loadDuplicateHsimd(x + ajx);
123             SimdReal jy_S   = loadDuplicateHsimd(x + ajy);
124             SimdReal jz_S   = loadDuplicateHsimd(x + ajz);
125
126             /* Calculate distance */
127             SimdReal dx_S0  = ix_S0 - jx_S;
128             SimdReal dy_S0  = iy_S0 - jy_S;
129             SimdReal dz_S0  = iz_S0 - jz_S;
130             SimdReal dx_S2  = ix_S2 - jx_S;
131             SimdReal dy_S2  = iy_S2 - jy_S;
132             SimdReal dz_S2  = iz_S2 - jz_S;
133
134             /* rsq = dx*dx+dy*dy+dz*dz */
135             SimdReal rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
136             SimdReal rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
137
138             /* Do the cut-off check */
139             SimdBool wco_S0 = (rsq_S0 < rlist2_S);
140             SimdBool wco_S2 = (rsq_S2 < rlist2_S);
141
142             wco_S0          = wco_S0 || wco_S2;
143
144             /* Putting the assignment inside the conditional is slower */
145             cjInner[ncjInner] = cjOuter[cjind];
146             if (anyTrue(wco_S0))
147             {
148                 ncjInner++;
149             }
150         }
151
152         if (ncjInner > ciInner[nciInner].cj_ind_start)
153         {
154             ciInner[nciInner].cj_ind_end = ncjInner;
155             nciInner++;
156         }
157     }
158
159     nbl->nci = nciInner;
160
161 #else  /* GMX_NBNXN_SIMD_2XNN */
162
163     GMX_RELEASE_ASSERT(false, "2xNN kernel called without 2xNN support");
164
165     GMX_UNUSED_VALUE(nbl);
166     GMX_UNUSED_VALUE(nbat);
167     GMX_UNUSED_VALUE(shift_vec);
168     GMX_UNUSED_VALUE(rlistInner);
169
170 #endif /* GMX_NBNXN_SIMD_2XNN */
171 }