Add CPU reference and SIMD pruning kernels
[gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / simd_2xnn / nbnxn_kernel_simd_2xnn_prune.cpp
blob29e22c5615d5d6f0d1db1a251a4ab91ba17f549f
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017, 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.
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.
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.
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.
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.
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.
36 #include "gmxpre.h"
38 #include "nbnxn_kernel_simd_2xnn_prune.h"
40 #include "gromacs/mdlib/nbnxn_pairlist.h"
41 #include "gromacs/mdlib/nbnxn_simd.h"
42 #include "gromacs/utility/gmxassert.h"
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
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)
56 #ifdef GMX_NBNXN_SIMD_2XNN
57 const nbnxn_ci_t * gmx_restrict ciOuter = nbl->ciOuter;
58 nbnxn_ci_t * gmx_restrict ciInner = nbl->ci;
60 const nbnxn_cj_t * gmx_restrict cjOuter = nbl->cjOuter;
61 nbnxn_cj_t * gmx_restrict cjInner = nbl->cj;
63 const real * gmx_restrict shiftvec = shift_vec[0];
64 const real * gmx_restrict x = nbat->x;
66 const SimdReal rlist2_S(rlistInner*rlistInner);
68 /* Initialize the new list as empty and add pairs that are in range */
69 int nciInner = 0;
70 int ncjInner = 0;
71 for (int i = 0; i < nbl->nciOuter; i++)
73 const nbnxn_ci_t * gmx_restrict ciEntry = &ciOuter[i];
75 /* Copy the original list entry to the pruned entry */
76 ciInner[nciInner].ci = ciEntry->ci;
77 ciInner[nciInner].shift = ciEntry->shift;
78 ciInner[nciInner].cj_ind_start = ncjInner;
80 /* Extract shift data */
81 int ish = (ciEntry->shift & NBNXN_CI_SHIFT);
82 int ish3 = ish*3;
83 int ci = ciEntry->ci;
85 SimdReal shX_S = SimdReal(shiftvec[ish3 ]);
86 SimdReal shY_S = SimdReal(shiftvec[ish3 + 1]);
87 SimdReal shZ_S = SimdReal(shiftvec[ish3 + 2]);
89 #if UNROLLJ <= 4
90 int sci = ci*STRIDE;
91 int scix = sci*DIM;
92 #else
93 int sci = (ci >> 1)*STRIDE;
94 int scix = sci*DIM + (ci & 1)*(STRIDE >> 1);
95 sci += (ci & 1)*(STRIDE >> 1);
96 #endif
98 /* Load i atom data */
99 int sciy = scix + STRIDE;
100 int sciz = sciy + STRIDE;
101 SimdReal ix_S0 = load1DualHsimd(x + scix ) + shX_S;
102 SimdReal ix_S2 = load1DualHsimd(x + scix + 2) + shX_S;
103 SimdReal iy_S0 = load1DualHsimd(x + sciy ) + shY_S;
104 SimdReal iy_S2 = load1DualHsimd(x + sciy + 2) + shY_S;
105 SimdReal iz_S0 = load1DualHsimd(x + sciz ) + shZ_S;
106 SimdReal iz_S2 = load1DualHsimd(x + sciz + 2) + shZ_S;
108 for (int cjind = ciEntry->cj_ind_start; cjind < ciEntry->cj_ind_end; cjind++)
110 /* j-cluster index */
111 int cj = cjOuter[cjind].cj;
113 /* Atom indices (of the first atom in the cluster) */
114 int aj = cj*UNROLLJ;
115 #if UNROLLJ == STRIDE
116 int ajx = aj*DIM;
117 #else
118 int ajx = (cj >> 1)*DIM*STRIDE + (cj & 1)*UNROLLJ;
119 #endif
120 int ajy = ajx + STRIDE;
121 int ajz = ajy + STRIDE;
123 /* load j atom coordinates */
124 SimdReal jx_S = loadDuplicateHsimd(x + ajx);
125 SimdReal jy_S = loadDuplicateHsimd(x + ajy);
126 SimdReal jz_S = loadDuplicateHsimd(x + ajz);
128 /* Calculate distance */
129 SimdReal dx_S0 = ix_S0 - jx_S;
130 SimdReal dy_S0 = iy_S0 - jy_S;
131 SimdReal dz_S0 = iz_S0 - jz_S;
132 SimdReal dx_S2 = ix_S2 - jx_S;
133 SimdReal dy_S2 = iy_S2 - jy_S;
134 SimdReal dz_S2 = iz_S2 - jz_S;
136 /* rsq = dx*dx+dy*dy+dz*dz */
137 SimdReal rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
138 SimdReal rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
140 /* Do the cut-off check */
141 SimdBool wco_S0 = (rsq_S0 < rlist2_S);
142 SimdBool wco_S2 = (rsq_S2 < rlist2_S);
144 wco_S0 = wco_S0 || wco_S2;
146 /* Putting the assignment inside the conditional is slower */
147 cjInner[ncjInner] = cjOuter[cjind];
148 if (anyTrue(wco_S0))
150 ncjInner++;
154 if (ncjInner > ciInner[nciInner].cj_ind_start)
156 ciInner[nciInner].cj_ind_end = ncjInner;
157 nciInner++;
161 nbl->nci = nciInner;
163 #else /* GMX_NBNXN_SIMD_2XNN */
165 GMX_RELEASE_ASSERT(false, "2xNN kernel called without 2xNN support");
167 GMX_UNUSED_VALUE(nbl);
168 GMX_UNUSED_VALUE(nbat);
169 GMX_UNUSED_VALUE(shift_vec);
170 GMX_UNUSED_VALUE(rlistInner);
172 #endif /* GMX_NBNXN_SIMD_2XNN */