2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
39 #include "kernel_prune.h"
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"
46 #ifdef GMX_NBNXN_SIMD_2XNN
47 # define GMX_SIMD_J_UNROLL_SIZE 2
48 # include "kernel_common.h"
51 /* Prune a single nbnxn_pairtlist_t entry with distance rlistInner */
52 void nbnxn_kernel_prune_2xnn(NbnxnPairlistCpu* nbl,
53 const nbnxn_atomdata_t* nbat,
54 gmx::ArrayRef<const gmx::RVec> shiftvec,
57 #ifdef GMX_NBNXN_SIMD_2XNN
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());
64 const nbnxn_ci_t* gmx_restrict ciOuter = nbl->ciOuter.data();
65 nbnxn_ci_t* gmx_restrict ciInner = nbl->ci.data();
67 const nbnxn_cj_t* gmx_restrict cjOuter = nbl->cjOuter.data();
68 nbnxn_cj_t* gmx_restrict cjInner = nbl->cj.data();
70 const real* gmx_restrict x = nbat->x().data();
72 const SimdReal rlist2_S(rlistInner * rlistInner);
74 /* Initialize the new list count as empty and add pairs that are in range */
77 const int nciOuter = nbl->ciOuter.size();
78 for (int i = 0; i < nciOuter; i++)
80 const nbnxn_ci_t* gmx_restrict ciEntry = &ciOuter[i];
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;
87 /* Extract shift data */
88 int ish = (ciEntry->shift & NBNXN_CI_SHIFT);
91 SimdReal shX_S = SimdReal(shiftvec[ish][XX]);
92 SimdReal shY_S = SimdReal(shiftvec[ish][YY]);
93 SimdReal shZ_S = SimdReal(shiftvec[ish][ZZ]);
96 int scix = ci * STRIDE * DIM;
98 int scix = (ci >> 1) * STRIDE * DIM + (ci & 1) * (STRIDE >> 1);
101 /* Load i atom data */
102 int sciy = scix + STRIDE;
103 int sciz = sciy + STRIDE;
104 SimdReal ix_S0 = loadU1DualHsimd(x + scix) + shX_S;
105 SimdReal ix_S2 = loadU1DualHsimd(x + scix + 2) + shX_S;
106 SimdReal iy_S0 = loadU1DualHsimd(x + sciy) + shY_S;
107 SimdReal iy_S2 = loadU1DualHsimd(x + sciy + 2) + shY_S;
108 SimdReal iz_S0 = loadU1DualHsimd(x + sciz) + shZ_S;
109 SimdReal iz_S2 = loadU1DualHsimd(x + sciz + 2) + shZ_S;
111 for (int cjind = ciEntry->cj_ind_start; cjind < ciEntry->cj_ind_end; cjind++)
113 /* j-cluster index */
114 int cj = cjOuter[cjind].cj;
116 /* Atom indices (of the first atom in the cluster) */
117 # if UNROLLJ == STRIDE
118 int aj = cj * UNROLLJ;
121 int ajx = (cj >> 1) * DIM * STRIDE + (cj & 1) * UNROLLJ;
123 int ajy = ajx + STRIDE;
124 int ajz = ajy + STRIDE;
126 /* load j atom coordinates */
127 SimdReal jx_S = loadDuplicateHsimd(x + ajx);
128 SimdReal jy_S = loadDuplicateHsimd(x + ajy);
129 SimdReal jz_S = loadDuplicateHsimd(x + ajz);
131 /* Calculate distance */
132 SimdReal dx_S0 = ix_S0 - jx_S;
133 SimdReal dy_S0 = iy_S0 - jy_S;
134 SimdReal dz_S0 = iz_S0 - jz_S;
135 SimdReal dx_S2 = ix_S2 - jx_S;
136 SimdReal dy_S2 = iy_S2 - jy_S;
137 SimdReal dz_S2 = iz_S2 - jz_S;
139 /* rsq = dx*dx+dy*dy+dz*dz */
140 SimdReal rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
141 SimdReal rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
143 /* Do the cut-off check */
144 SimdBool wco_S0 = (rsq_S0 < rlist2_S);
145 SimdBool wco_S2 = (rsq_S2 < rlist2_S);
147 wco_S0 = wco_S0 || wco_S2;
149 /* Putting the assignment inside the conditional is slower */
150 cjInner[ncjInner] = cjOuter[cjind];
157 if (ncjInner > ciInner[nciInner].cj_ind_start)
159 ciInner[nciInner].cj_ind_end = ncjInner;
164 nbl->ci.resize(nciInner);
165 nbl->cj.resize(ncjInner);
167 #else /* GMX_NBNXN_SIMD_2XNN */
169 GMX_RELEASE_ASSERT(false, "2xNN kernel called without 2xNN support");
171 GMX_UNUSED_VALUE(nbl);
172 GMX_UNUSED_VALUE(nbat);
173 GMX_UNUSED_VALUE(shiftvec);
174 GMX_UNUSED_VALUE(rlistInner);
176 #endif /* GMX_NBNXN_SIMD_2XNN */