Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / nbnxm / kernels_simd_4xm / kernel_prune.cpp
blob6d46acff1e60518c3411b348f9720020a2cfcc54
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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 "kernel_prune.h"
40 #include "gromacs/nbnxm/atomdata.h"
41 #include "gromacs/nbnxm/nbnxm_simd.h"
42 #include "gromacs/nbnxm/pairlist.h"
43 #include "gromacs/utility/gmxassert.h"
45 #ifdef GMX_NBNXN_SIMD_4XN
46 #define GMX_SIMD_J_UNROLL_SIZE 1
47 #include "gromacs/nbnxm/kernels_simd_4xm/kernel_common.h"
48 #endif
50 /* Prune a single nbnxn_pairtlist_t entry with distance rlistInner */
51 void
52 nbnxn_kernel_prune_4xn(NbnxnPairlistCpu * nbl,
53 const nbnxn_atomdata_t * nbat,
54 const rvec * gmx_restrict shift_vec,
55 real rlistInner)
57 #ifdef GMX_NBNXN_SIMD_4XN
58 using namespace gmx;
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 shiftvec = shift_vec[0];
71 const real * gmx_restrict x = nbat->x().data();
73 const SimdReal rlist2_S(rlistInner*rlistInner);
75 /* Initialize the new list count as empty and add pairs that are in range */
76 int nciInner = 0;
77 int ncjInner = 0;
78 const int nciOuter = nbl->ciOuter.size();
79 for (int i = 0; i < nciOuter; i++)
81 const nbnxn_ci_t * gmx_restrict ciEntry = &ciOuter[i];
83 /* Copy the original list entry to the pruned entry */
84 ciInner[nciInner].ci = ciEntry->ci;
85 ciInner[nciInner].shift = ciEntry->shift;
86 ciInner[nciInner].cj_ind_start = ncjInner;
88 /* Extract shift data */
89 int ish = (ciEntry->shift & NBNXN_CI_SHIFT);
90 int ish3 = ish*3;
91 int ci = ciEntry->ci;
93 SimdReal shX_S = SimdReal(shiftvec[ish3 ]);
94 SimdReal shY_S = SimdReal(shiftvec[ish3 + 1]);
95 SimdReal shZ_S = SimdReal(shiftvec[ish3 + 2]);
97 #if UNROLLJ <= 4
98 int scix = ci*STRIDE*DIM;
99 #else
100 int scix = (ci >> 1)*STRIDE*DIM + (ci & 1)*(STRIDE >> 1);
101 #endif
103 /* Load i atom data */
104 int sciy = scix + STRIDE;
105 int sciz = sciy + STRIDE;
106 SimdReal ix_S0 = SimdReal(x[scix ]) + shX_S;
107 SimdReal ix_S1 = SimdReal(x[scix + 1]) + shX_S;
108 SimdReal ix_S2 = SimdReal(x[scix + 2]) + shX_S;
109 SimdReal ix_S3 = SimdReal(x[scix + 3]) + shX_S;
110 SimdReal iy_S0 = SimdReal(x[sciy ]) + shY_S;
111 SimdReal iy_S1 = SimdReal(x[sciy + 1]) + shY_S;
112 SimdReal iy_S2 = SimdReal(x[sciy + 2]) + shY_S;
113 SimdReal iy_S3 = SimdReal(x[sciy + 3]) + shY_S;
114 SimdReal iz_S0 = SimdReal(x[sciz ]) + shZ_S;
115 SimdReal iz_S1 = SimdReal(x[sciz + 1]) + shZ_S;
116 SimdReal iz_S2 = SimdReal(x[sciz + 2]) + shZ_S;
117 SimdReal iz_S3 = SimdReal(x[sciz + 3]) + shZ_S;
119 for (int cjind = ciEntry->cj_ind_start; cjind < ciEntry->cj_ind_end; cjind++)
121 /* j-cluster index */
122 int cj = cjOuter[cjind].cj;
124 /* Atom indices (of the first atom in the cluster) */
125 #if UNROLLJ == STRIDE
126 int aj = cj*UNROLLJ;
127 int ajx = aj*DIM;
128 #else
129 int ajx = (cj >> 1)*DIM*STRIDE + (cj & 1)*UNROLLJ;
130 #endif
131 int ajy = ajx + STRIDE;
132 int ajz = ajy + STRIDE;
134 /* load j atom coordinates */
135 SimdReal jx_S = load<SimdReal>(x + ajx);
136 SimdReal jy_S = load<SimdReal>(x + ajy);
137 SimdReal jz_S = load<SimdReal>(x + ajz);
139 /* Calculate distance */
140 SimdReal dx_S0 = ix_S0 - jx_S;
141 SimdReal dy_S0 = iy_S0 - jy_S;
142 SimdReal dz_S0 = iz_S0 - jz_S;
143 SimdReal dx_S1 = ix_S1 - jx_S;
144 SimdReal dy_S1 = iy_S1 - jy_S;
145 SimdReal dz_S1 = iz_S1 - jz_S;
146 SimdReal dx_S2 = ix_S2 - jx_S;
147 SimdReal dy_S2 = iy_S2 - jy_S;
148 SimdReal dz_S2 = iz_S2 - jz_S;
149 SimdReal dx_S3 = ix_S3 - jx_S;
150 SimdReal dy_S3 = iy_S3 - jy_S;
151 SimdReal dz_S3 = iz_S3 - jz_S;
153 /* rsq = dx*dx+dy*dy+dz*dz */
154 SimdReal rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
155 SimdReal rsq_S1 = norm2(dx_S1, dy_S1, dz_S1);
156 SimdReal rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
157 SimdReal rsq_S3 = norm2(dx_S3, dy_S3, dz_S3);
159 /* Do the cut-off check */
160 SimdBool wco_S0 = (rsq_S0 < rlist2_S);
161 SimdBool wco_S1 = (rsq_S1 < rlist2_S);
162 SimdBool wco_S2 = (rsq_S2 < rlist2_S);
163 SimdBool wco_S3 = (rsq_S3 < rlist2_S);
165 wco_S0 = wco_S0 || wco_S1;
166 wco_S2 = wco_S2 || wco_S3;
167 wco_S0 = wco_S0 || wco_S2;
169 /* Putting the assignment inside the conditional is slower */
170 cjInner[ncjInner] = cjOuter[cjind];
171 if (anyTrue(wco_S0))
173 ncjInner++;
177 if (ncjInner > ciInner[nciInner].cj_ind_start)
179 ciInner[nciInner].cj_ind_end = ncjInner;
180 nciInner++;
184 nbl->ci.resize(nciInner);
185 nbl->cj.resize(ncjInner);
187 #else /* GMX_NBNXN_SIMD_4XN */
189 GMX_RELEASE_ASSERT(false, "4xN kernel called without 4xN support");
191 GMX_UNUSED_VALUE(nbl);
192 GMX_UNUSED_VALUE(nbat);
193 GMX_UNUSED_VALUE(shift_vec);
194 GMX_UNUSED_VALUE(rlistInner);
196 #endif /* GMX_NBNXN_SIMD_4XN */