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.
49 #include "gromacs/domdec/domdec_struct.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/mdlib/gmx_omp_nthreads.h"
56 #include "gromacs/mdtypes/group.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/mdtypes/nblist.h"
59 #include "gromacs/nbnxm/atomdata.h"
60 #include "gromacs/nbnxm/gpu_data_mgmt.h"
61 #include "gromacs/pbcutil/ishift.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/simd/simd.h"
64 #include "gromacs/simd/vector_operations.h"
65 #include "gromacs/utility/exceptions.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/gmxomp.h"
68 #include "gromacs/utility/listoflists.h"
69 #include "gromacs/utility/smalloc.h"
71 #include "boundingboxes.h"
72 #include "clusterdistancekerneltype.h"
74 #include "nbnxm_geometry.h"
75 #include "nbnxm_simd.h"
76 #include "pairlistset.h"
77 #include "pairlistsets.h"
78 #include "pairlistwork.h"
79 #include "pairsearch.h"
81 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
83 using BoundingBox = Nbnxm::BoundingBox; // TODO: Remove when refactoring this file
84 using BoundingBox1D = Nbnxm::BoundingBox1D; // TODO: Remove when refactoring this file
86 using Grid = Nbnxm::Grid; // TODO: Remove when refactoring this file
88 // Convenience alias for partial Nbnxn namespace usage
89 using InteractionLocality = gmx::InteractionLocality;
91 /* We shift the i-particles backward for PBC.
92 * This leads to more conditionals than shifting forward.
93 * We do this to get more balanced pair lists.
95 constexpr bool c_pbcShiftBackward = true;
97 /* Layout for the nonbonded NxN pair lists */
98 enum class NbnxnLayout
100 NoSimd4x4, // i-cluster size 4, j-cluster size 4
101 Simd4xN, // i-cluster size 4, j-cluster size SIMD width
102 Simd2xNN, // i-cluster size 4, j-cluster size half SIMD width
103 Gpu8x8x8 // i-cluster size 8, j-cluster size 8 + super-clustering
106 #if defined(GMX_NBNXN_SIMD_4XN) || defined(GMX_NBNXN_SIMD_2XNN)
107 /* Returns the j-cluster size */
108 template<NbnxnLayout layout>
109 static constexpr int jClusterSize()
111 static_assert(layout == NbnxnLayout::NoSimd4x4 || layout == NbnxnLayout::Simd4xN
112 || layout == NbnxnLayout::Simd2xNN,
113 "Currently jClusterSize only supports CPU layouts");
115 return layout == NbnxnLayout::Simd4xN
116 ? GMX_SIMD_REAL_WIDTH
117 : (layout == NbnxnLayout::Simd2xNN ? GMX_SIMD_REAL_WIDTH / 2 : c_nbnxnCpuIClusterSize);
120 /*! \brief Returns the j-cluster index given the i-cluster index.
122 * \tparam jClusterSize The number of atoms in a j-cluster
123 * \tparam jSubClusterIndex The j-sub-cluster index (0/1), used when size(j-cluster) <
125 * \param[in] ci The i-cluster index
127 template<int jClusterSize, int jSubClusterIndex>
128 static inline int cjFromCi(int ci)
130 static_assert(jClusterSize == c_nbnxnCpuIClusterSize / 2 || jClusterSize == c_nbnxnCpuIClusterSize
131 || jClusterSize == c_nbnxnCpuIClusterSize * 2,
132 "Only j-cluster sizes 2, 4 and 8 are currently implemented");
134 static_assert(jSubClusterIndex == 0 || jSubClusterIndex == 1,
135 "Only sub-cluster indices 0 and 1 are supported");
137 if (jClusterSize == c_nbnxnCpuIClusterSize / 2)
139 if (jSubClusterIndex == 0)
145 return ((ci + 1) << 1) - 1;
148 else if (jClusterSize == c_nbnxnCpuIClusterSize)
158 /*! \brief Returns the j-cluster index given the i-cluster index.
160 * \tparam layout The pair-list layout
161 * \tparam jSubClusterIndex The j-sub-cluster index (0/1), used when size(j-cluster) <
163 * \param[in] ci The i-cluster index
165 template<NbnxnLayout layout, int jSubClusterIndex>
166 static inline int cjFromCi(int ci)
168 constexpr int clusterSize = jClusterSize<layout>();
170 return cjFromCi<clusterSize, jSubClusterIndex>(ci);
173 /* Returns the nbnxn coordinate data index given the i-cluster index */
174 template<NbnxnLayout layout>
175 static inline int xIndexFromCi(int ci)
177 constexpr int clusterSize = jClusterSize<layout>();
179 static_assert(clusterSize == c_nbnxnCpuIClusterSize / 2 || clusterSize == c_nbnxnCpuIClusterSize
180 || clusterSize == c_nbnxnCpuIClusterSize * 2,
181 "Only j-cluster sizes 2, 4 and 8 are currently implemented");
183 if (clusterSize <= c_nbnxnCpuIClusterSize)
185 /* Coordinates are stored packed in groups of 4 */
186 return ci * STRIDE_P4;
190 /* Coordinates packed in 8, i-cluster size is half the packing width */
191 return (ci >> 1) * STRIDE_P8 + (ci & 1) * (c_packX8 >> 1);
195 /* Returns the nbnxn coordinate data index given the j-cluster index */
196 template<NbnxnLayout layout>
197 static inline int xIndexFromCj(int cj)
199 constexpr int clusterSize = jClusterSize<layout>();
201 static_assert(clusterSize == c_nbnxnCpuIClusterSize / 2 || clusterSize == c_nbnxnCpuIClusterSize
202 || clusterSize == c_nbnxnCpuIClusterSize * 2,
203 "Only j-cluster sizes 2, 4 and 8 are currently implemented");
205 if (clusterSize == c_nbnxnCpuIClusterSize / 2)
207 /* Coordinates are stored packed in groups of 4 */
208 return (cj >> 1) * STRIDE_P4 + (cj & 1) * (c_packX4 >> 1);
210 else if (clusterSize == c_nbnxnCpuIClusterSize)
212 /* Coordinates are stored packed in groups of 4 */
213 return cj * STRIDE_P4;
217 /* Coordinates are stored packed in groups of 8 */
218 return cj * STRIDE_P8;
221 #endif // defined(GMX_NBNXN_SIMD_4XN) || defined(GMX_NBNXN_SIMD_2XNN)
223 static constexpr int sizeNeededForBufferFlags(const int numAtoms)
225 return (numAtoms + NBNXN_BUFFERFLAG_SIZE - 1) / NBNXN_BUFFERFLAG_SIZE;
228 // Resets current flags to 0 and adds more flags if needed.
229 static void resizeAndZeroBufferFlags(std::vector<gmx_bitmask_t>* flags, const int numAtoms)
232 flags->resize(sizeNeededForBufferFlags(numAtoms), gmx_bitmask_t{ 0 });
236 /* Returns the pair-list cutoff between a bounding box and a grid cell given an atom-to-atom pair-list cutoff
238 * Given a cutoff distance between atoms, this functions returns the cutoff
239 * distance2 between a bounding box of a group of atoms and a grid cell.
240 * Since atoms can be geometrically outside of the cell they have been
241 * assigned to (when atom groups instead of individual atoms are assigned
242 * to cells), this distance returned can be larger than the input.
244 static real listRangeForBoundingBoxToGridCell(real rlist, const Grid::Dimensions& gridDims)
246 return rlist + gridDims.maxAtomGroupRadius;
248 /* Returns the pair-list cutoff between a grid cells given an atom-to-atom pair-list cutoff
250 * Given a cutoff distance between atoms, this functions returns the cutoff
251 * distance2 between two grid cells.
252 * Since atoms can be geometrically outside of the cell they have been
253 * assigned to (when atom groups instead of individual atoms are assigned
254 * to cells), this distance returned can be larger than the input.
256 static real listRangeForGridCellToGridCell(real rlist,
257 const Grid::Dimensions& iGridDims,
258 const Grid::Dimensions& jGridDims)
260 return rlist + iGridDims.maxAtomGroupRadius + jGridDims.maxAtomGroupRadius;
263 /* Determines the cell range along one dimension that
264 * the bounding box b0 - b1 sees.
268 get_cell_range(real b0, real b1, const Grid::Dimensions& jGridDims, real d2, real rlist, int* cf, int* cl)
270 real listRangeBBToCell2 = gmx::square(listRangeForBoundingBoxToGridCell(rlist, jGridDims));
271 real distanceInCells = (b0 - jGridDims.lowerCorner[dim]) * jGridDims.invCellSize[dim];
272 *cf = std::max(static_cast<int>(distanceInCells), 0);
275 && d2 + gmx::square((b0 - jGridDims.lowerCorner[dim]) - (*cf - 1 + 1) * jGridDims.cellSize[dim])
276 < listRangeBBToCell2)
281 *cl = std::min(static_cast<int>((b1 - jGridDims.lowerCorner[dim]) * jGridDims.invCellSize[dim]),
282 jGridDims.numCells[dim] - 1);
283 while (*cl < jGridDims.numCells[dim] - 1
284 && d2 + gmx::square((*cl + 1) * jGridDims.cellSize[dim] - (b1 - jGridDims.lowerCorner[dim]))
285 < listRangeBBToCell2)
291 /* Reference code calculating the distance^2 between two bounding boxes */
293 static float box_dist2(float bx0, float bx1, float by0,
294 float by1, float bz0, float bz1,
295 const BoundingBox *bb)
298 float dl, dh, dm, dm0;
302 dl = bx0 - bb->upper.x;
303 dh = bb->lower.x - bx1;
304 dm = std::max(dl, dh);
305 dm0 = std::max(dm, 0.0f);
308 dl = by0 - bb->upper.y;
309 dh = bb->lower.y - by1;
310 dm = std::max(dl, dh);
311 dm0 = std::max(dm, 0.0f);
314 dl = bz0 - bb->upper.z;
315 dh = bb->lower.z - bz1;
316 dm = std::max(dl, dh);
317 dm0 = std::max(dm, 0.0f);
324 #if !NBNXN_SEARCH_BB_SIMD4
326 /*! \brief Plain C code calculating the distance^2 between two bounding boxes in xyz0 format
328 * \param[in] bb_i First bounding box
329 * \param[in] bb_j Second bounding box
331 static float clusterBoundingBoxDistance2(const BoundingBox& bb_i, const BoundingBox& bb_j)
333 float dl = bb_i.lower.x - bb_j.upper.x;
334 float dh = bb_j.lower.x - bb_i.upper.x;
335 float dm = std::max(dl, dh);
336 float dm0 = std::max(dm, 0.0F);
337 float d2 = dm0 * dm0;
339 dl = bb_i.lower.y - bb_j.upper.y;
340 dh = bb_j.lower.y - bb_i.upper.y;
341 dm = std::max(dl, dh);
342 dm0 = std::max(dm, 0.0F);
345 dl = bb_i.lower.z - bb_j.upper.z;
346 dh = bb_j.lower.z - bb_i.upper.z;
347 dm = std::max(dl, dh);
348 dm0 = std::max(dm, 0.0F);
354 #else /* NBNXN_SEARCH_BB_SIMD4 */
356 /*! \brief 4-wide SIMD code calculating the distance^2 between two bounding boxes in xyz0 format
358 * \param[in] bb_i First bounding box, should be aligned for 4-wide SIMD
359 * \param[in] bb_j Second bounding box, should be aligned for 4-wide SIMD
361 static float clusterBoundingBoxDistance2(const BoundingBox& bb_i, const BoundingBox& bb_j)
363 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
366 const Simd4Float bb_i_S0 = load4(bb_i.lower.ptr());
367 const Simd4Float bb_i_S1 = load4(bb_i.upper.ptr());
368 const Simd4Float bb_j_S0 = load4(bb_j.lower.ptr());
369 const Simd4Float bb_j_S1 = load4(bb_j.upper.ptr());
371 const Simd4Float dl_S = bb_i_S0 - bb_j_S1;
372 const Simd4Float dh_S = bb_j_S0 - bb_i_S1;
374 const Simd4Float dm_S = max(dl_S, dh_S);
375 const Simd4Float dm0_S = max(dm_S, simd4SetZeroF());
377 return dotProduct(dm0_S, dm0_S);
380 /* Calculate bb bounding distances of bb_i[si,...,si+3] and store them in d2 */
381 template<int boundingBoxStart>
382 static inline void gmx_simdcall clusterBoundingBoxDistance2_xxxx_simd4_inner(const float* bb_i,
384 const Simd4Float xj_l,
385 const Simd4Float yj_l,
386 const Simd4Float zj_l,
387 const Simd4Float xj_h,
388 const Simd4Float yj_h,
389 const Simd4Float zj_h)
391 constexpr int stride = c_packedBoundingBoxesDimSize;
393 const int shi = boundingBoxStart * Nbnxm::c_numBoundingBoxBounds1D * DIM;
395 const Simd4Float zero = setZero();
397 const Simd4Float xi_l = load4(bb_i + shi + 0 * stride);
398 const Simd4Float yi_l = load4(bb_i + shi + 1 * stride);
399 const Simd4Float zi_l = load4(bb_i + shi + 2 * stride);
400 const Simd4Float xi_h = load4(bb_i + shi + 3 * stride);
401 const Simd4Float yi_h = load4(bb_i + shi + 4 * stride);
402 const Simd4Float zi_h = load4(bb_i + shi + 5 * stride);
404 const Simd4Float dx_0 = xi_l - xj_h;
405 const Simd4Float dy_0 = yi_l - yj_h;
406 const Simd4Float dz_0 = zi_l - zj_h;
408 const Simd4Float dx_1 = xj_l - xi_h;
409 const Simd4Float dy_1 = yj_l - yi_h;
410 const Simd4Float dz_1 = zj_l - zi_h;
412 const Simd4Float mx = max(dx_0, dx_1);
413 const Simd4Float my = max(dy_0, dy_1);
414 const Simd4Float mz = max(dz_0, dz_1);
416 const Simd4Float m0x = max(mx, zero);
417 const Simd4Float m0y = max(my, zero);
418 const Simd4Float m0z = max(mz, zero);
420 const Simd4Float d2x = m0x * m0x;
421 const Simd4Float d2y = m0y * m0y;
422 const Simd4Float d2z = m0z * m0z;
424 const Simd4Float d2s = d2x + d2y;
425 const Simd4Float d2t = d2s + d2z;
427 store4(d2 + boundingBoxStart, d2t);
430 /* 4-wide SIMD code for nsi bb distances for bb format xxxxyyyyzzzz */
431 static void clusterBoundingBoxDistance2_xxxx_simd4(const float* bb_j, const int nsi, const float* bb_i, float* d2)
433 constexpr int stride = c_packedBoundingBoxesDimSize;
435 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
438 const Simd4Float xj_l = Simd4Float(bb_j[0 * stride]);
439 const Simd4Float yj_l = Simd4Float(bb_j[1 * stride]);
440 const Simd4Float zj_l = Simd4Float(bb_j[2 * stride]);
441 const Simd4Float xj_h = Simd4Float(bb_j[3 * stride]);
442 const Simd4Float yj_h = Simd4Float(bb_j[4 * stride]);
443 const Simd4Float zj_h = Simd4Float(bb_j[5 * stride]);
445 /* Here we "loop" over si (0,stride) from 0 to nsi with step stride.
446 * But as we know the number of iterations is 1 or 2, we unroll manually.
448 clusterBoundingBoxDistance2_xxxx_simd4_inner<0>(bb_i, d2, xj_l, yj_l, zj_l, xj_h, yj_h, zj_h);
451 clusterBoundingBoxDistance2_xxxx_simd4_inner<stride>(bb_i, d2, xj_l, yj_l, zj_l, xj_h, yj_h, zj_h);
455 #endif /* NBNXN_SEARCH_BB_SIMD4 */
458 /* Returns if any atom pair from two clusters is within distance sqrt(rlist2) */
459 static inline gmx_bool
460 clusterpair_in_range(const NbnxnPairlistGpuWork& work, int si, int csj, int stride, const real* x_j, real rlist2)
462 #if !GMX_SIMD4_HAVE_REAL
465 * All coordinates are stored as xyzxyz...
468 const real* x_i = work.iSuperClusterData.x.data();
470 for (int i = 0; i < c_nbnxnGpuClusterSize; i++)
472 int i0 = (si * c_nbnxnGpuClusterSize + i) * DIM;
473 for (int j = 0; j < c_nbnxnGpuClusterSize; j++)
475 int j0 = (csj * c_nbnxnGpuClusterSize + j) * stride;
477 real d2 = gmx::square(x_i[i0] - x_j[j0]) + gmx::square(x_i[i0 + 1] - x_j[j0 + 1])
478 + gmx::square(x_i[i0 + 2] - x_j[j0 + 2]);
489 #else /* !GMX_SIMD4_HAVE_REAL */
491 /* 4-wide SIMD version.
492 * The coordinates x_i are stored as xxxxyyyy..., x_j is stored xyzxyz...
493 * Using 8-wide AVX(2) is not faster on Intel Sandy Bridge and Haswell.
495 static_assert(c_nbnxnGpuClusterSize == 8 || c_nbnxnGpuClusterSize == 4,
496 "A cluster is hard-coded to 4/8 atoms.");
498 Simd4Real rc2_S{ rlist2 };
500 const real* x_i = work.iSuperClusterData.xSimd.data();
502 int dim_stride = c_nbnxnGpuClusterSize * DIM;
503 Simd4Real ix_S0 = load4(x_i + si * dim_stride + 0 * GMX_SIMD4_WIDTH);
504 Simd4Real iy_S0 = load4(x_i + si * dim_stride + 1 * GMX_SIMD4_WIDTH);
505 Simd4Real iz_S0 = load4(x_i + si * dim_stride + 2 * GMX_SIMD4_WIDTH);
507 Simd4Real ix_S1, iy_S1, iz_S1;
508 if (c_nbnxnGpuClusterSize == 8)
510 ix_S1 = load4(x_i + si * dim_stride + 3 * GMX_SIMD4_WIDTH);
511 iy_S1 = load4(x_i + si * dim_stride + 4 * GMX_SIMD4_WIDTH);
512 iz_S1 = load4(x_i + si * dim_stride + 5 * GMX_SIMD4_WIDTH);
514 /* We loop from the outer to the inner particles to maximize
515 * the chance that we find a pair in range quickly and return.
517 int j0 = csj * c_nbnxnGpuClusterSize;
518 int j1 = j0 + c_nbnxnGpuClusterSize - 1;
521 Simd4Real jx0_S, jy0_S, jz0_S;
522 Simd4Real jx1_S, jy1_S, jz1_S;
524 Simd4Real dx_S0, dy_S0, dz_S0;
525 Simd4Real dx_S1, dy_S1, dz_S1;
526 Simd4Real dx_S2, dy_S2, dz_S2;
527 Simd4Real dx_S3, dy_S3, dz_S3;
538 Simd4Bool wco_any_S01, wco_any_S23, wco_any_S;
540 jx0_S = Simd4Real(x_j[j0 * stride + 0]);
541 jy0_S = Simd4Real(x_j[j0 * stride + 1]);
542 jz0_S = Simd4Real(x_j[j0 * stride + 2]);
544 jx1_S = Simd4Real(x_j[j1 * stride + 0]);
545 jy1_S = Simd4Real(x_j[j1 * stride + 1]);
546 jz1_S = Simd4Real(x_j[j1 * stride + 2]);
548 /* Calculate distance */
549 dx_S0 = ix_S0 - jx0_S;
550 dy_S0 = iy_S0 - jy0_S;
551 dz_S0 = iz_S0 - jz0_S;
552 dx_S2 = ix_S0 - jx1_S;
553 dy_S2 = iy_S0 - jy1_S;
554 dz_S2 = iz_S0 - jz1_S;
555 if (c_nbnxnGpuClusterSize == 8)
557 dx_S1 = ix_S1 - jx0_S;
558 dy_S1 = iy_S1 - jy0_S;
559 dz_S1 = iz_S1 - jz0_S;
560 dx_S3 = ix_S1 - jx1_S;
561 dy_S3 = iy_S1 - jy1_S;
562 dz_S3 = iz_S1 - jz1_S;
565 /* rsq = dx*dx+dy*dy+dz*dz */
566 rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
567 rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
568 if (c_nbnxnGpuClusterSize == 8)
570 rsq_S1 = norm2(dx_S1, dy_S1, dz_S1);
571 rsq_S3 = norm2(dx_S3, dy_S3, dz_S3);
574 wco_S0 = (rsq_S0 < rc2_S);
575 wco_S2 = (rsq_S2 < rc2_S);
576 if (c_nbnxnGpuClusterSize == 8)
578 wco_S1 = (rsq_S1 < rc2_S);
579 wco_S3 = (rsq_S3 < rc2_S);
581 if (c_nbnxnGpuClusterSize == 8)
583 wco_any_S01 = wco_S0 || wco_S1;
584 wco_any_S23 = wco_S2 || wco_S3;
585 wco_any_S = wco_any_S01 || wco_any_S23;
589 wco_any_S = wco_S0 || wco_S2;
592 if (anyTrue(wco_any_S))
603 #endif /* !GMX_SIMD4_HAVE_REAL */
606 /* Returns the j-cluster index for index cjIndex in a cj list */
607 static inline int nblCj(gmx::ArrayRef<const nbnxn_cj_t> cjList, int cjIndex)
609 return cjList[cjIndex].cj;
612 /* Returns the j-cluster index for index cjIndex in a cj4 list */
613 static inline int nblCj(gmx::ArrayRef<const nbnxn_cj4_t> cj4List, int cjIndex)
615 return cj4List[cjIndex / c_nbnxnGpuJgroupSize].cj[cjIndex & (c_nbnxnGpuJgroupSize - 1)];
618 /* Returns the i-interaction mask of the j sub-cell for index cj_ind */
619 static unsigned int nbl_imask0(const NbnxnPairlistGpu* nbl, int cj_ind)
621 return nbl->cj4[cj_ind / c_nbnxnGpuJgroupSize].imei[0].imask;
624 NbnxnPairlistCpu::NbnxnPairlistCpu() :
625 na_ci(c_nbnxnCpuIClusterSize),
630 work(std::make_unique<NbnxnPairlistCpuWork>())
634 NbnxnPairlistGpu::NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy) :
635 na_ci(c_nbnxnGpuClusterSize),
636 na_cj(c_nbnxnGpuClusterSize),
637 na_sc(c_gpuNumClusterPerCell * c_nbnxnGpuClusterSize),
639 sci({}, { pinningPolicy }),
640 cj4({}, { pinningPolicy }),
641 excl({}, { pinningPolicy }),
643 work(std::make_unique<NbnxnPairlistGpuWork>())
645 static_assert(c_nbnxnGpuNumClusterPerSupercluster == c_gpuNumClusterPerCell,
646 "The search code assumes that the a super-cluster matches a search grid cell");
648 static_assert(sizeof(cj4[0].imei[0].imask) * 8 >= c_nbnxnGpuJgroupSize * c_gpuNumClusterPerCell,
649 "The i super-cluster cluster interaction mask does not contain a sufficient "
652 static_assert(sizeof(excl[0]) * 8 >= c_nbnxnGpuJgroupSize * c_gpuNumClusterPerCell,
653 "The GPU exclusion mask does not contain a sufficient number of bits");
655 // We always want a first entry without any exclusions
659 // TODO: Move to pairlistset.cpp
660 PairlistSet::PairlistSet(const PairlistParams& pairlistParams) :
661 params_(pairlistParams),
662 combineLists_(sc_isGpuPairListType[pairlistParams.pairlistType]), // Currently GPU lists are always combined
663 isCpuType_(!sc_isGpuPairListType[pairlistParams.pairlistType])
666 const int numLists = gmx_omp_nthreads_get(ModuleMultiThread::Nonbonded);
668 if (!combineLists_ && numLists > NBNXN_BUFFERFLAG_MAX_THREADS)
671 "%d OpenMP threads were requested. Since the non-bonded force buffer reduction "
672 "is prohibitively slow with more than %d threads, we do not allow this. Use %d "
673 "or less OpenMP threads.",
675 NBNXN_BUFFERFLAG_MAX_THREADS,
676 NBNXN_BUFFERFLAG_MAX_THREADS);
681 cpuLists_.resize(numLists);
684 cpuListsWork_.resize(numLists);
689 /* Only list 0 is used on the GPU, use normal allocation for i>0 */
690 gpuLists_.emplace_back(gmx::PinningPolicy::PinnedIfSupported);
691 /* Lists 0 to numLists are use for constructing lists in parallel
692 * on the CPU using numLists threads (and then merged into list 0).
694 for (int i = 1; i < numLists; i++)
696 gpuLists_.emplace_back(gmx::PinningPolicy::CannotBePinned);
701 fepLists_.resize(numLists);
703 /* Execute in order to avoid memory interleaving between threads */
704 #pragma omp parallel for num_threads(numLists) schedule(static)
705 for (int i = 0; i < numLists; i++)
709 /* We used to allocate all normal lists locally on each thread
710 * as well. The question is if allocating the object on the
711 * master thread (but all contained list memory thread local)
712 * impacts performance.
714 fepLists_[i] = std::make_unique<t_nblist>();
716 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
721 /* Print statistics of a pair list, used for debug output */
722 static void print_nblist_statistics(FILE* fp,
723 const NbnxnPairlistCpu& nbl,
724 const Nbnxm::GridSet& gridSet,
727 const Grid& grid = gridSet.grids()[0];
728 const Grid::Dimensions& dims = grid.dimensions();
730 fprintf(fp, "nbl nci %zu ncj %d\n", nbl.ci.size(), nbl.ncjInUse);
731 const int numAtomsJCluster = grid.geometry().numAtomsJCluster;
733 if (grid.numCells() == 0)
738 const double numAtomsPerCell = nbl.ncjInUse / static_cast<double>(grid.numCells()) * numAtomsJCluster;
740 "nbl na_cj %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
744 nbl.ncjInUse / static_cast<double>(grid.numCells()),
747 / (0.5 * 4.0 / 3.0 * M_PI * rl * rl * rl * grid.numCells() * numAtomsJCluster
748 / (dims.gridSize[XX] * dims.gridSize[YY] * dims.gridSize[ZZ])));
751 "nbl average j cell list length %.1f\n",
752 0.25 * nbl.ncjInUse / std::max(static_cast<double>(nbl.ci.size()), 1.0));
754 int cs[gmx::c_numShiftVectors] = { 0 };
756 for (const nbnxn_ci_t& ciEntry : nbl.ci)
758 cs[ciEntry.shift & NBNXN_CI_SHIFT] += ciEntry.cj_ind_end - ciEntry.cj_ind_start;
760 int j = ciEntry.cj_ind_start;
761 while (j < ciEntry.cj_ind_end && nbl.cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
768 "nbl cell pairs, total: %zu excl: %d %.1f%%\n",
771 100 * npexcl / std::max(static_cast<double>(nbl.cj.size()), 1.0));
772 for (int s = 0; s < gmx::c_numShiftVectors; s++)
776 fprintf(fp, "nbl shift %2d ncj %3d\n", s, cs[s]);
781 /* Print statistics of a pair lists, used for debug output */
782 static void print_nblist_statistics(FILE* fp,
783 const NbnxnPairlistGpu& nbl,
784 const Nbnxm::GridSet& gridSet,
787 const Grid& grid = gridSet.grids()[0];
788 const Grid::Dimensions& dims = grid.dimensions();
791 "nbl nsci %zu ncj4 %zu nsi %d excl4 %zu\n",
796 const int numAtomsCluster = grid.geometry().numAtomsICluster;
797 const double numAtomsPerCell = nbl.nci_tot / static_cast<double>(grid.numClusters()) * numAtomsCluster;
799 "nbl na_c %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
803 nbl.nci_tot / static_cast<double>(grid.numClusters()),
806 / (0.5 * 4.0 / 3.0 * M_PI * rl * rl * rl * grid.numClusters() * numAtomsCluster
807 / (dims.gridSize[XX] * dims.gridSize[YY] * dims.gridSize[ZZ])));
812 int c[c_gpuNumClusterPerCell + 1] = { 0 };
813 for (const nbnxn_sci_t& sci : nbl.sci)
816 for (int j4 = sci.cj4_ind_start; j4 < sci.cj4_ind_end; j4++)
818 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
821 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
823 if (nbl.cj4[j4].imei[0].imask & (1U << (j * c_gpuNumClusterPerCell + si)))
833 sum_nsp2 += nsp * nsp;
834 nsp_max = std::max(nsp_max, nsp);
836 if (!nbl.sci.empty())
838 sum_nsp /= nbl.sci.size();
839 sum_nsp2 /= nbl.sci.size();
842 "nbl #cluster-pairs: av %.1f stddev %.1f max %d\n",
844 std::sqrt(sum_nsp2 - sum_nsp * sum_nsp),
847 if (!nbl.cj4.empty())
849 for (int b = 0; b <= c_gpuNumClusterPerCell; b++)
852 "nbl j-list #i-subcell %d %7d %4.1f\n",
855 100.0 * c[b] / size_t{ nbl.cj4.size() * c_nbnxnGpuJgroupSize });
860 /* Returns a reference to the exclusion mask for j-cluster-group \p cj4 and warp \p warp
861 * Generates a new exclusion entry when the j-cluster-group uses
862 * the default all-interaction mask at call time, so the returned mask
863 * can be modified when needed.
865 static nbnxn_excl_t& get_exclusion_mask(NbnxnPairlistGpu* nbl, int cj4, int warp)
867 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
869 /* No exclusions set, make a new list entry */
870 const size_t oldSize = nbl->excl.size();
871 GMX_ASSERT(oldSize >= 1, "We should always have entry [0]");
872 /* Add entry with default values: no exclusions */
873 nbl->excl.resize(oldSize + 1);
874 nbl->cj4[cj4].imei[warp].excl_ind = oldSize;
877 return nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
880 /* Sets self exclusions and excludes half of the double pairs in the self cluster-pair \p nbl->cj4[cj4Index].cj[jOffsetInGroup]
882 * \param[in,out] nbl The cluster pair list
883 * \param[in] cj4Index The j-cluster group index into \p nbl->cj4
884 * \param[in] jOffsetInGroup The j-entry offset in \p nbl->cj4[cj4Index]
885 * \param[in] iClusterInCell The i-cluster index in the cell
887 static void setSelfAndNewtonExclusionsGpu(NbnxnPairlistGpu* nbl,
889 const int jOffsetInGroup,
890 const int iClusterInCell)
892 constexpr int numJatomsPerPart = c_nbnxnGpuClusterSize / c_nbnxnGpuClusterpairSplit;
894 /* The exclusions are stored separately for each part of the split */
895 for (int part = 0; part < c_nbnxnGpuClusterpairSplit; part++)
897 const int jOffset = part * numJatomsPerPart;
898 /* Make a new exclusion mask entry for each part, if we don't already have one yet */
899 nbnxn_excl_t& excl = get_exclusion_mask(nbl, cj4Index, part);
901 /* Set all bits with j-index <= i-index */
902 for (int jIndexInPart = 0; jIndexInPart < numJatomsPerPart; jIndexInPart++)
904 for (int i = jOffset + jIndexInPart; i < c_nbnxnGpuClusterSize; i++)
906 excl.pair[jIndexInPart * c_nbnxnGpuClusterSize + i] &=
907 ~(1U << (jOffsetInGroup * c_gpuNumClusterPerCell + iClusterInCell));
913 /* Returns a diagonal or off-diagonal interaction mask for plain C lists */
914 static unsigned int get_imask(gmx_bool rdiag, int ci, int cj)
916 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
919 /* Returns a diagonal or off-diagonal interaction mask for cj-size=2 */
920 gmx_unused static unsigned int get_imask_simd_j2(gmx_bool rdiag, int ci, int cj)
922 return (rdiag && ci * 2 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_0
923 : (rdiag && ci * 2 + 1 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_1
924 : NBNXN_INTERACTION_MASK_ALL));
927 /* Returns a diagonal or off-diagonal interaction mask for cj-size=4 */
928 gmx_unused static unsigned int get_imask_simd_j4(gmx_bool rdiag, int ci, int cj)
930 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
933 /* Returns a diagonal or off-diagonal interaction mask for cj-size=8 */
934 gmx_unused static unsigned int get_imask_simd_j8(gmx_bool rdiag, int ci, int cj)
936 return (rdiag && ci == cj * 2 ? NBNXN_INTERACTION_MASK_DIAG_J8_0
937 : (rdiag && ci == cj * 2 + 1 ? NBNXN_INTERACTION_MASK_DIAG_J8_1
938 : NBNXN_INTERACTION_MASK_ALL));
942 # if GMX_SIMD_REAL_WIDTH == 2
943 # define get_imask_simd_4xn get_imask_simd_j2
945 # if GMX_SIMD_REAL_WIDTH == 4
946 # define get_imask_simd_4xn get_imask_simd_j4
948 # if GMX_SIMD_REAL_WIDTH == 8
949 # define get_imask_simd_4xn get_imask_simd_j8
950 # define get_imask_simd_2xnn get_imask_simd_j4
952 # if GMX_SIMD_REAL_WIDTH == 16
953 # define get_imask_simd_2xnn get_imask_simd_j8
957 /* Plain C code for checking and adding cluster-pairs to the list.
959 * \param[in] gridj The j-grid
960 * \param[in,out] nbl The pair-list to store the cluster pairs in
961 * \param[in] icluster The index of the i-cluster
962 * \param[in] jclusterFirst The first cluster in the j-range
963 * \param[in] jclusterLast The last cluster in the j-range
964 * \param[in] excludeSubDiagonal Exclude atom pairs with i-index > j-index
965 * \param[in] x_j Coordinates for the j-atom, in xyz format
966 * \param[in] rlist2 The squared list cut-off
967 * \param[in] rbb2 The squared cut-off for putting cluster-pairs in the list based on bounding box distance only
968 * \param[in,out] numDistanceChecks The number of distance checks performed
970 static void makeClusterListSimple(const Grid& jGrid,
971 NbnxnPairlistCpu* nbl,
975 bool excludeSubDiagonal,
976 const real* gmx_restrict x_j,
979 int* gmx_restrict numDistanceChecks)
981 const BoundingBox* gmx_restrict bb_ci = nbl->work->iClusterData.bb.data();
982 const real* gmx_restrict x_ci = nbl->work->iClusterData.x.data();
984 bool InRange = false;
985 while (!InRange && jclusterFirst <= jclusterLast)
987 real d2 = clusterBoundingBoxDistance2(bb_ci[0], jGrid.jBoundingBoxes()[jclusterFirst]);
988 *numDistanceChecks += 2;
990 /* Check if the distance is within the distance where
991 * we use only the bounding box distance rbb,
992 * or within the cut-off and there is at least one atom pair
993 * within the cut-off.
999 else if (d2 < rlist2)
1001 int cjf_gl = jGrid.cellOffset() + jclusterFirst;
1002 for (int i = 0; i < c_nbnxnCpuIClusterSize && !InRange; i++)
1004 for (int j = 0; j < c_nbnxnCpuIClusterSize; j++)
1008 || (gmx::square(x_ci[i * STRIDE_XYZ + XX]
1009 - x_j[(cjf_gl * c_nbnxnCpuIClusterSize + j) * STRIDE_XYZ + XX])
1010 + gmx::square(x_ci[i * STRIDE_XYZ + YY]
1011 - x_j[(cjf_gl * c_nbnxnCpuIClusterSize + j) * STRIDE_XYZ + YY])
1012 + gmx::square(x_ci[i * STRIDE_XYZ + ZZ]
1013 - x_j[(cjf_gl * c_nbnxnCpuIClusterSize + j) * STRIDE_XYZ + ZZ])
1017 *numDistanceChecks += c_nbnxnCpuIClusterSize * c_nbnxnCpuIClusterSize;
1030 while (!InRange && jclusterLast > jclusterFirst)
1032 real d2 = clusterBoundingBoxDistance2(bb_ci[0], jGrid.jBoundingBoxes()[jclusterLast]);
1033 *numDistanceChecks += 2;
1035 /* Check if the distance is within the distance where
1036 * we use only the bounding box distance rbb,
1037 * or within the cut-off and there is at least one atom pair
1038 * within the cut-off.
1044 else if (d2 < rlist2)
1046 int cjl_gl = jGrid.cellOffset() + jclusterLast;
1047 for (int i = 0; i < c_nbnxnCpuIClusterSize && !InRange; i++)
1049 for (int j = 0; j < c_nbnxnCpuIClusterSize; j++)
1053 || (gmx::square(x_ci[i * STRIDE_XYZ + XX]
1054 - x_j[(cjl_gl * c_nbnxnCpuIClusterSize + j) * STRIDE_XYZ + XX])
1055 + gmx::square(x_ci[i * STRIDE_XYZ + YY]
1056 - x_j[(cjl_gl * c_nbnxnCpuIClusterSize + j) * STRIDE_XYZ + YY])
1057 + gmx::square(x_ci[i * STRIDE_XYZ + ZZ]
1058 - x_j[(cjl_gl * c_nbnxnCpuIClusterSize + j) * STRIDE_XYZ + ZZ])
1062 *numDistanceChecks += c_nbnxnCpuIClusterSize * c_nbnxnCpuIClusterSize;
1070 if (jclusterFirst <= jclusterLast)
1072 for (int jcluster = jclusterFirst; jcluster <= jclusterLast; jcluster++)
1074 /* Store cj and the interaction mask */
1076 cjEntry.cj = jGrid.cellOffset() + jcluster;
1077 cjEntry.excl = get_imask(excludeSubDiagonal, icluster, jcluster);
1078 nbl->cj.push_back(cjEntry);
1080 /* Increase the closing index in the i list */
1081 nbl->ci.back().cj_ind_end = nbl->cj.size();
1085 #ifdef GMX_NBNXN_SIMD_4XN
1086 # include "pairlist_simd_4xm.h"
1088 #ifdef GMX_NBNXN_SIMD_2XNN
1089 # include "pairlist_simd_2xmm.h"
1092 /* Plain C or SIMD4 code for making a pair list of super-cell sci vs scj.
1093 * Checks bounding box distances and possibly atom pair distances.
1095 static void make_cluster_list_supersub(const Grid& iGrid,
1097 NbnxnPairlistGpu* nbl,
1100 const bool excludeSubDiagonal,
1105 int* numDistanceChecks)
1107 NbnxnPairlistGpuWork& work = *nbl->work;
1110 const float* pbb_ci = work.iSuperClusterData.bbPacked.data();
1112 const BoundingBox* bb_ci = work.iSuperClusterData.bb.data();
1115 assert(c_nbnxnGpuClusterSize == iGrid.geometry().numAtomsICluster);
1116 assert(c_nbnxnGpuClusterSize == jGrid.geometry().numAtomsICluster);
1118 /* We generate the pairlist mainly based on bounding-box distances
1119 * and do atom pair distance based pruning on the GPU.
1120 * Only if a j-group contains a single cluster-pair, we try to prune
1121 * that pair based on atom distances on the CPU to avoid empty j-groups.
1123 #define PRUNE_LIST_CPU_ONE 1
1124 #define PRUNE_LIST_CPU_ALL 0
1126 #if PRUNE_LIST_CPU_ONE
1130 float* d2l = work.distanceBuffer.data();
1132 for (int subc = 0; subc < jGrid.numClustersPerCell()[scj]; subc++)
1134 const int cj4_ind = work.cj_ind / c_nbnxnGpuJgroupSize;
1135 const int cj_offset = work.cj_ind - cj4_ind * c_nbnxnGpuJgroupSize;
1136 const int cj = scj * c_gpuNumClusterPerCell + subc;
1138 const int cj_gl = jGrid.cellOffset() * c_gpuNumClusterPerCell + cj;
1140 int ci1 = (excludeSubDiagonal && sci == scj) ? subc + 1 : iGrid.numClustersPerCell()[sci];
1144 /* Determine all ci1 bb distances in one call with SIMD4 */
1145 const int offset = packedBoundingBoxesIndex(cj) + (cj & (c_packedBoundingBoxesDimSize - 1));
1146 clusterBoundingBoxDistance2_xxxx_simd4(
1147 jGrid.packedBoundingBoxes().data() + offset, ci1, pbb_ci, d2l);
1148 *numDistanceChecks += c_nbnxnGpuClusterSize * 2;
1152 unsigned int imask = 0;
1153 /* We use a fixed upper-bound instead of ci1 to help optimization */
1154 for (int ci = 0; ci < c_gpuNumClusterPerCell; ci++)
1162 /* Determine the bb distance between ci and cj */
1163 d2l[ci] = clusterBoundingBoxDistance2(bb_ci[ci], jGrid.jBoundingBoxes()[cj]);
1164 *numDistanceChecks += 2;
1168 #if PRUNE_LIST_CPU_ALL
1169 /* Check if the distance is within the distance where
1170 * we use only the bounding box distance rbb,
1171 * or within the cut-off and there is at least one atom pair
1172 * within the cut-off. This check is very costly.
1174 *numDistanceChecks += c_nbnxnGpuClusterSize * c_nbnxnGpuClusterSize;
1175 if (d2 < rbb2 || (d2 < rlist2 && clusterpair_in_range(work, ci, cj_gl, stride, x, rlist2)))
1177 /* Check if the distance between the two bounding boxes
1178 * in within the pair-list cut-off.
1183 /* Flag this i-subcell to be taken into account */
1184 imask |= (1U << (cj_offset * c_gpuNumClusterPerCell + ci));
1186 #if PRUNE_LIST_CPU_ONE
1194 #if PRUNE_LIST_CPU_ONE
1195 /* If we only found 1 pair, check if any atoms are actually
1196 * within the cut-off, so we could get rid of it.
1198 if (npair == 1 && d2l[ci_last] >= rbb2
1199 && !clusterpair_in_range(work, ci_last, cj_gl, stride, x, rlist2))
1201 imask &= ~(1U << (cj_offset * c_gpuNumClusterPerCell + ci_last));
1208 /* We have at least one cluster pair: add a j-entry */
1209 if (static_cast<size_t>(cj4_ind) == nbl->cj4.size())
1211 nbl->cj4.resize(nbl->cj4.size() + 1);
1213 nbnxn_cj4_t* cj4 = &nbl->cj4[cj4_ind];
1215 cj4->cj[cj_offset] = cj_gl;
1217 /* Set the exclusions for the ci==sj entry.
1218 * Here we don't bother to check if this entry is actually flagged,
1219 * as it will nearly always be in the list.
1221 if (excludeSubDiagonal && sci == scj)
1223 setSelfAndNewtonExclusionsGpu(nbl, cj4_ind, cj_offset, subc);
1226 /* Copy the cluster interaction mask to the list */
1227 for (int w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
1229 cj4->imei[w].imask |= imask;
1232 nbl->work->cj_ind++;
1234 /* Keep the count */
1235 nbl->nci_tot += npair;
1237 /* Increase the closing index in i super-cell list */
1238 nbl->sci.back().cj4_ind_end =
1239 (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1) / c_nbnxnGpuJgroupSize;
1244 /* Returns how many contiguous j-clusters we have starting in the i-list */
1245 template<typename CjListType>
1246 static int numContiguousJClusters(const int cjIndexStart,
1247 const int cjIndexEnd,
1248 gmx::ArrayRef<const CjListType> cjList)
1250 const int firstJCluster = nblCj(cjList, cjIndexStart);
1252 int numContiguous = 0;
1254 while (cjIndexStart + numContiguous < cjIndexEnd
1255 && nblCj(cjList, cjIndexStart + numContiguous) == firstJCluster + numContiguous)
1260 return numContiguous;
1264 * \brief Helper struct for efficient searching for excluded atoms in a j-list
1268 /*! \brief Constructs a j-list range from \p cjList with the given index range */
1269 template<typename CjListType>
1270 JListRanges(int cjIndexStart, int cjIndexEnd, gmx::ArrayRef<const CjListType> cjList);
1272 int cjIndexStart; //!< The start index in the j-list
1273 int cjIndexEnd; //!< The end index in the j-list
1274 int cjFirst; //!< The j-cluster with index cjIndexStart
1275 int cjLast; //!< The j-cluster with index cjIndexEnd-1
1276 int numDirect; //!< Up to cjIndexStart+numDirect the j-clusters are cjFirst + the index offset
1280 template<typename CjListType>
1281 JListRanges::JListRanges(int cjIndexStart, int cjIndexEnd, gmx::ArrayRef<const CjListType> cjList) :
1282 cjIndexStart(cjIndexStart), cjIndexEnd(cjIndexEnd)
1284 GMX_ASSERT(cjIndexEnd > cjIndexStart, "JListRanges should only be called with non-empty lists");
1286 cjFirst = nblCj(cjList, cjIndexStart);
1287 cjLast = nblCj(cjList, cjIndexEnd - 1);
1289 /* Determine how many contiguous j-cells we have starting
1290 * from the first i-cell. This number can be used to directly
1291 * calculate j-cell indices for excluded atoms.
1293 numDirect = numContiguousJClusters(cjIndexStart, cjIndexEnd, cjList);
1297 /* Return the index of \p jCluster in the given range or -1 when not present
1299 * Note: This code is executed very often and therefore performance is
1300 * important. It should be inlined and fully optimized.
1302 template<typename CjListType>
1303 static inline int findJClusterInJList(int jCluster,
1304 const JListRanges& ranges,
1305 gmx::ArrayRef<const CjListType> cjList)
1307 if (jCluster < ranges.cjFirst + ranges.numDirect)
1309 /* We can calculate the index directly using the offset */
1310 return ranges.cjIndexStart + jCluster - ranges.cjFirst;
1314 /* Search for jCluster using bisection */
1316 int rangeStart = ranges.cjIndexStart + ranges.numDirect;
1317 int rangeEnd = ranges.cjIndexEnd;
1318 while (index == -1 && rangeStart < rangeEnd)
1320 int rangeMiddle = (rangeStart + rangeEnd) >> 1;
1322 const int clusterMiddle = nblCj(cjList, rangeMiddle);
1324 if (jCluster == clusterMiddle)
1326 index = rangeMiddle;
1328 else if (jCluster < clusterMiddle)
1330 rangeEnd = rangeMiddle;
1334 rangeStart = rangeMiddle + 1;
1341 // TODO: Get rid of the two functions below by renaming sci to ci (or something better)
1343 /* Return the i-entry in the list we are currently operating on */
1344 static nbnxn_ci_t* getOpenIEntry(NbnxnPairlistCpu* nbl)
1346 return &nbl->ci.back();
1349 /* Return the i-entry in the list we are currently operating on */
1350 static nbnxn_sci_t* getOpenIEntry(NbnxnPairlistGpu* nbl)
1352 return &nbl->sci.back();
1355 /* Set all atom-pair exclusions for a simple type list i-entry
1357 * Set all atom-pair exclusions from the topology stored in exclusions
1358 * as masks in the pair-list for simple list entry iEntry.
1360 static void setExclusionsForIEntry(const Nbnxm::GridSet& gridSet,
1361 NbnxnPairlistCpu* nbl,
1362 gmx_bool diagRemoved,
1364 const nbnxn_ci_t& iEntry,
1365 const ListOfLists<int>& exclusions)
1367 if (iEntry.cj_ind_end == iEntry.cj_ind_start)
1369 /* Empty list: no exclusions */
1373 const JListRanges ranges(iEntry.cj_ind_start, iEntry.cj_ind_end, gmx::makeConstArrayRef(nbl->cj));
1375 const int iCluster = iEntry.ci;
1377 gmx::ArrayRef<const int> cell = gridSet.cells();
1378 gmx::ArrayRef<const int> atomIndices = gridSet.atomIndices();
1380 /* Loop over the atoms in the i-cluster */
1381 for (int i = 0; i < nbl->na_ci; i++)
1383 const int iIndex = iCluster * nbl->na_ci + i;
1384 const int iAtom = atomIndices[iIndex];
1387 /* Loop over the topology-based exclusions for this i-atom */
1388 for (const int jAtom : exclusions[iAtom])
1392 /* The self exclusion are already set, save some time */
1396 /* Get the index of the j-atom in the nbnxn atom data */
1397 const int jIndex = cell[jAtom];
1399 /* Without shifts we only calculate interactions j>i
1400 * for one-way pair-lists.
1402 if (diagRemoved && jIndex <= iIndex)
1407 const int jCluster = (jIndex >> na_cj_2log);
1409 /* Could the cluster se be in our list? */
1410 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
1413 findJClusterInJList(jCluster, ranges, gmx::makeConstArrayRef(nbl->cj));
1417 /* We found an exclusion, clear the corresponding
1420 const int innerJ = jIndex - (jCluster << na_cj_2log);
1422 nbl->cj[index].excl &= ~(1U << ((i << na_cj_2log) + innerJ));
1430 /* Add a new i-entry to the FEP list and copy the i-properties */
1431 static inline void fep_list_new_nri_copy(t_nblist* nlist)
1433 /* Add a new i-entry */
1436 assert(nlist->nri < nlist->maxnri);
1438 /* Duplicate the last i-entry, except for jindex, which continues */
1439 nlist->iinr[nlist->nri] = nlist->iinr[nlist->nri - 1];
1440 nlist->shift[nlist->nri] = nlist->shift[nlist->nri - 1];
1441 nlist->gid[nlist->nri] = nlist->gid[nlist->nri - 1];
1442 nlist->jindex[nlist->nri] = nlist->nrj;
1445 /* Rellocate FEP list for size nl->maxnri, TODO: replace by C++ */
1446 static void reallocate_nblist(t_nblist* nl)
1448 nl->iinr.resize(nl->maxnri);
1449 nl->gid.resize(nl->maxnri);
1450 nl->shift.resize(nl->maxnri);
1451 nl->jindex.resize(nl->maxnri + 1);
1454 /* For load balancing of the free-energy lists over threads, we set
1455 * the maximum nrj size of an i-entry to 40. This leads to good
1456 * load balancing in the worst case scenario of a single perturbed
1457 * particle on 16 threads, while not introducing significant overhead.
1458 * Note that half of the perturbed pairs will anyhow end up in very small lists,
1459 * since non perturbed i-particles will see few perturbed j-particles).
1461 const int max_nrj_fep = 40;
1463 /* Exclude the perturbed pairs from the Verlet list. This is only done to avoid
1464 * singularities for overlapping particles (0/0), since the charges and
1465 * LJ parameters have been zeroed in the nbnxn data structure.
1466 * Simultaneously make a group pair list for the perturbed pairs.
1468 static void make_fep_list(gmx::ArrayRef<const int> atomIndices,
1469 const nbnxn_atomdata_t* nbat,
1470 NbnxnPairlistCpu* nbl,
1471 gmx_bool bDiagRemoved,
1473 real gmx_unused shx,
1474 real gmx_unused shy,
1475 real gmx_unused shz,
1476 real gmx_unused rlist_fep2,
1484 if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
1490 const int ci = nbl_ci->ci;
1492 const int cj_ind_start = nbl_ci->cj_ind_start;
1493 const int cj_ind_end = nbl_ci->cj_ind_end;
1495 /* In worst case we have alternating energy groups
1496 * and create #atom-pair lists, which means we need the size
1497 * of a cluster pair (na_ci*na_cj) times the number of cj's.
1499 const int nri_max = nbl->na_ci * nbl->na_cj * (cj_ind_end - cj_ind_start);
1500 if (nlist->nri + nri_max > nlist->maxnri)
1502 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1503 reallocate_nblist(nlist);
1506 const int numAtomsJCluster = jGrid.geometry().numAtomsJCluster;
1508 const nbnxn_atomdata_t::Params& nbatParams = nbat->params();
1510 const int ngid = nbatParams.nenergrp;
1512 /* TODO: Consider adding a check in grompp and changing this to an assert */
1513 const int numBitsInEnergyGroupIdsForAtomsInJCluster = sizeof(gid_cj) * 8;
1514 if (ngid * numAtomsJCluster > numBitsInEnergyGroupIdsForAtomsInJCluster)
1517 "The Verlet scheme with %dx%d kernels and free-energy only supports up to %zu "
1519 iGrid.geometry().numAtomsICluster,
1521 (sizeof(gid_cj) * 8) / numAtomsJCluster);
1524 const int egp_shift = nbatParams.neg_2log;
1525 const int egp_mask = (1 << egp_shift) - 1;
1527 /* Loop over the atoms in the i sub-cell */
1528 bool bFEP_i_all = true;
1529 for (int i = 0; i < nbl->na_ci; i++)
1531 const int ind_i = ci * nbl->na_ci + i;
1532 const int ai = atomIndices[ind_i];
1535 int nri = nlist->nri;
1536 nlist->jindex[nri + 1] = nlist->jindex[nri];
1537 nlist->iinr[nri] = ai;
1538 /* The actual energy group pair index is set later */
1539 nlist->gid[nri] = 0;
1540 nlist->shift[nri] = nbl_ci->shift & NBNXN_CI_SHIFT;
1542 bool bFEP_i = iGrid.atomIsPerturbed(ci - iGrid.cellOffset(), i);
1544 bFEP_i_all = bFEP_i_all && bFEP_i;
1546 if (nlist->nrj + (cj_ind_end - cj_ind_start) * nbl->na_cj > nlist->maxnrj)
1548 nlist->maxnrj = over_alloc_small(nlist->nrj + (cj_ind_end - cj_ind_start) * nbl->na_cj);
1549 nlist->jjnr.resize(nlist->maxnrj);
1550 nlist->excl_fep.resize(nlist->maxnrj);
1555 gid_i = (nbatParams.energrp[ci] >> (egp_shift * i)) & egp_mask;
1558 for (int cj_ind = cj_ind_start; cj_ind < cj_ind_end; cj_ind++)
1560 unsigned int fep_cj = 0U;
1563 const int cja = nbl->cj[cj_ind].cj;
1565 if (numAtomsJCluster == jGrid.geometry().numAtomsICluster)
1567 const int cjr = cja - jGrid.cellOffset();
1568 fep_cj = jGrid.fepBits(cjr);
1571 gid_cj = nbatParams.energrp[cja];
1574 else if (2 * numAtomsJCluster == jGrid.geometry().numAtomsICluster)
1576 const int cjr = cja - jGrid.cellOffset() * 2;
1577 /* Extract half of the ci fep/energrp mask */
1578 fep_cj = (jGrid.fepBits(cjr >> 1) >> ((cjr & 1) * numAtomsJCluster))
1579 & ((1 << numAtomsJCluster) - 1);
1582 gid_cj = nbatParams.energrp[cja >> 1] >> ((cja & 1) * numAtomsJCluster * egp_shift)
1583 & ((1 << (numAtomsJCluster * egp_shift)) - 1);
1588 const int cjr = cja - (jGrid.cellOffset() >> 1);
1589 /* Combine two ci fep masks/energrp */
1590 fep_cj = jGrid.fepBits(cjr * 2)
1591 + (jGrid.fepBits(cjr * 2 + 1) << jGrid.geometry().numAtomsICluster);
1594 gid_cj = nbatParams.energrp[cja * 2]
1595 + (nbatParams.energrp[cja * 2 + 1]
1596 << (jGrid.geometry().numAtomsICluster * egp_shift));
1600 if (bFEP_i || fep_cj != 0)
1602 for (int j = 0; j < nbl->na_cj; j++)
1604 /* Is this interaction perturbed and not excluded? */
1605 const int ind_j = cja * nbl->na_cj + j;
1606 const int aj = atomIndices[ind_j];
1607 if (aj >= 0 && (bFEP_i || (fep_cj & (1 << j))) && (!bDiagRemoved || ind_j >= ind_i))
1611 const int gid_j = (gid_cj >> (j * egp_shift)) & egp_mask;
1612 const int gid = GID(gid_i, gid_j, ngid);
1614 if (nlist->nrj > nlist->jindex[nri] && nlist->gid[nri] != gid)
1616 /* Energy group pair changed: new list */
1617 fep_list_new_nri_copy(nlist);
1620 nlist->gid[nri] = gid;
1623 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1625 fep_list_new_nri_copy(nlist);
1629 /* Add it to the FEP list */
1630 nlist->jjnr[nlist->nrj] = aj;
1631 nlist->excl_fep[nlist->nrj] =
1632 (nbl->cj[cj_ind].excl >> (i * nbl->na_cj + j)) & 1;
1635 /* Exclude it from the normal list.
1636 * Note that the charge has been set to zero,
1637 * but we need to avoid 0/0, as perturbed atoms
1638 * can be on top of each other.
1640 nbl->cj[cj_ind].excl &= ~(1U << (i * nbl->na_cj + j));
1646 if (nlist->nrj > nlist->jindex[nri])
1648 /* Actually add this new, non-empty, list */
1650 nlist->jindex[nlist->nri] = nlist->nrj;
1657 /* All interactions are perturbed, we can skip this entry */
1658 nbl_ci->cj_ind_end = cj_ind_start;
1659 nbl->ncjInUse -= cj_ind_end - cj_ind_start;
1663 /* Return the index of atom a within a cluster */
1664 static inline int cj_mod_cj4(int cj)
1666 return cj & (c_nbnxnGpuJgroupSize - 1);
1669 /* Convert a j-cluster to a cj4 group */
1670 static inline int cj_to_cj4(int cj)
1672 return cj / c_nbnxnGpuJgroupSize;
1675 /* Return the index of an j-atom within a warp */
1676 static inline int a_mod_wj(int a)
1678 return a & (c_nbnxnGpuClusterSize / c_nbnxnGpuClusterpairSplit - 1);
1681 /* As make_fep_list above, but for super/sub lists. */
1682 static void make_fep_list(gmx::ArrayRef<const int> atomIndices,
1683 const nbnxn_atomdata_t* nbat,
1684 NbnxnPairlistGpu* nbl,
1685 gmx_bool bDiagRemoved,
1686 const nbnxn_sci_t* nbl_sci,
1695 const int numJClusterGroups = nbl_sci->numJClusterGroups();
1696 if (numJClusterGroups == 0)
1702 const int sci = nbl_sci->sci;
1704 const int cj4_ind_start = nbl_sci->cj4_ind_start;
1705 const int cj4_ind_end = nbl_sci->cj4_ind_end;
1707 /* Here we process one super-cell, max #atoms na_sc, versus a list
1708 * cj4 entries, each with max c_nbnxnGpuJgroupSize cj's, each
1709 * of size na_cj atoms.
1710 * On the GPU we don't support energy groups (yet).
1711 * So for each of the na_sc i-atoms, we need max one FEP list
1712 * for each max_nrj_fep j-atoms.
1715 nbl->na_sc * nbl->na_cj * (1 + (numJClusterGroups * c_nbnxnGpuJgroupSize) / max_nrj_fep);
1716 if (nlist->nri + nri_max > nlist->maxnri)
1718 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1719 reallocate_nblist(nlist);
1722 /* Loop over the atoms in the i super-cluster */
1723 for (int c = 0; c < c_gpuNumClusterPerCell; c++)
1725 const int c_abs = sci * c_gpuNumClusterPerCell + c;
1727 for (int i = 0; i < nbl->na_ci; i++)
1729 const int ind_i = c_abs * nbl->na_ci + i;
1730 const int ai = atomIndices[ind_i];
1733 int nri = nlist->nri;
1734 nlist->jindex[nri + 1] = nlist->jindex[nri];
1735 nlist->iinr[nri] = ai;
1736 /* With GPUs, energy groups are not supported */
1737 nlist->gid[nri] = 0;
1738 nlist->shift[nri] = nbl_sci->shift & NBNXN_CI_SHIFT;
1741 iGrid.atomIsPerturbed(c_abs - iGrid.cellOffset() * c_gpuNumClusterPerCell, i);
1743 real xi = nbat->x()[ind_i * nbat->xstride + XX] + shx;
1744 real yi = nbat->x()[ind_i * nbat->xstride + YY] + shy;
1745 real zi = nbat->x()[ind_i * nbat->xstride + ZZ] + shz;
1747 const int nrjMax = nlist->nrj + numJClusterGroups * c_nbnxnGpuJgroupSize * nbl->na_cj;
1748 if (nrjMax > nlist->maxnrj)
1750 nlist->maxnrj = over_alloc_small(nrjMax);
1751 nlist->jjnr.resize(nlist->maxnrj);
1752 nlist->excl_fep.resize(nlist->maxnrj);
1755 for (int cj4_ind = cj4_ind_start; cj4_ind < cj4_ind_end; cj4_ind++)
1757 const nbnxn_cj4_t* cj4 = &nbl->cj4[cj4_ind];
1759 for (int gcj = 0; gcj < c_nbnxnGpuJgroupSize; gcj++)
1761 if ((cj4->imei[0].imask & (1U << (gcj * c_gpuNumClusterPerCell + c))) == 0)
1763 /* Skip this ci for this cj */
1767 const int cjr = cj4->cj[gcj] - jGrid.cellOffset() * c_gpuNumClusterPerCell;
1769 if (bFEP_i || jGrid.clusterIsPerturbed(cjr))
1771 for (int j = 0; j < nbl->na_cj; j++)
1773 /* Is this interaction perturbed and not excluded? */
1775 (jGrid.cellOffset() * c_gpuNumClusterPerCell + cjr) * nbl->na_cj + j;
1776 const int aj = atomIndices[ind_j];
1777 if (aj >= 0 && (bFEP_i || jGrid.atomIsPerturbed(cjr, j))
1778 && (!bDiagRemoved || ind_j >= ind_i))
1781 j / (c_nbnxnGpuClusterSize / c_nbnxnGpuClusterpairSplit);
1782 nbnxn_excl_t& excl = get_exclusion_mask(nbl, cj4_ind, jHalf);
1784 int excl_pair = a_mod_wj(j) * nbl->na_ci + i;
1785 unsigned int excl_bit = (1U << (gcj * c_gpuNumClusterPerCell + c));
1787 real dx = nbat->x()[ind_j * nbat->xstride + XX] - xi;
1788 real dy = nbat->x()[ind_j * nbat->xstride + YY] - yi;
1789 real dz = nbat->x()[ind_j * nbat->xstride + ZZ] - zi;
1791 /* The unpruned GPU list has more than 2/3
1792 * of the atom pairs beyond rlist. Using
1793 * this list will cause a lot of overhead
1794 * in the CPU FEP kernels, especially
1795 * relative to the fast GPU kernels.
1796 * So we prune the FEP list here.
1798 if (dx * dx + dy * dy + dz * dz < rlist_fep2)
1800 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1802 fep_list_new_nri_copy(nlist);
1806 /* Add it to the FEP list */
1807 nlist->jjnr[nlist->nrj] = aj;
1808 nlist->excl_fep[nlist->nrj] =
1809 (excl.pair[excl_pair] & excl_bit) ? 1 : 0;
1813 /* Exclude it from the normal list.
1814 * Note that the charge and LJ parameters have
1815 * been set to zero, but we need to avoid 0/0,
1816 * as perturbed atoms can be on top of each other.
1818 excl.pair[excl_pair] &= ~excl_bit;
1822 /* Note that we could mask out this pair in imask
1823 * if all i- and/or all j-particles are perturbed.
1824 * But since the perturbed pairs on the CPU will
1825 * take an order of magnitude more time, the GPU
1826 * will finish before the CPU and there is no gain.
1832 if (nlist->nrj > nlist->jindex[nri])
1834 /* Actually add this new, non-empty, list */
1836 nlist->jindex[nlist->nri] = nlist->nrj;
1843 /* Set all atom-pair exclusions for a GPU type list i-entry
1845 * Sets all atom-pair exclusions from the topology stored in exclusions
1846 * as masks in the pair-list for i-super-cluster list entry iEntry.
1848 static void setExclusionsForIEntry(const Nbnxm::GridSet& gridSet,
1849 NbnxnPairlistGpu* nbl,
1850 gmx_bool diagRemoved,
1851 int gmx_unused na_cj_2log,
1852 const nbnxn_sci_t& iEntry,
1853 const ListOfLists<int>& exclusions)
1855 if (iEntry.numJClusterGroups() == 0)
1861 /* Set the search ranges using start and end j-cluster indices.
1862 * Note that here we can not use cj4_ind_end, since the last cj4
1863 * can be only partially filled, so we use cj_ind.
1865 const JListRanges ranges(iEntry.cj4_ind_start * c_nbnxnGpuJgroupSize,
1867 gmx::makeConstArrayRef(nbl->cj4));
1869 GMX_ASSERT(nbl->na_ci == c_nbnxnGpuClusterSize, "na_ci should match the GPU cluster size");
1870 constexpr int c_clusterSize = c_nbnxnGpuClusterSize;
1871 constexpr int c_superClusterSize = c_nbnxnGpuNumClusterPerSupercluster * c_nbnxnGpuClusterSize;
1873 const int iSuperCluster = iEntry.sci;
1875 gmx::ArrayRef<const int> atomIndices = gridSet.atomIndices();
1876 gmx::ArrayRef<const int> cell = gridSet.cells();
1878 /* Loop over the atoms in the i super-cluster */
1879 for (int i = 0; i < c_superClusterSize; i++)
1881 const int iIndex = iSuperCluster * c_superClusterSize + i;
1882 const int iAtom = atomIndices[iIndex];
1885 const int iCluster = i / c_clusterSize;
1887 /* Loop over the topology-based exclusions for this i-atom */
1888 for (const int jAtom : exclusions[iAtom])
1892 /* The self exclusions are already set, save some time */
1896 /* Get the index of the j-atom in the nbnxn atom data */
1897 const int jIndex = cell[jAtom];
1899 /* Without shifts we only calculate interactions j>i
1900 * for one-way pair-lists.
1902 /* NOTE: We would like to use iIndex on the right hand side,
1903 * but that makes this routine 25% slower with gcc6/7.
1904 * Even using c_superClusterSize makes it slower.
1905 * Either of these changes triggers peeling of the exclIndex
1906 * loop, which apparently leads to far less efficient code.
1908 if (diagRemoved && jIndex <= iSuperCluster * nbl->na_sc + i)
1913 const int jCluster = jIndex / c_clusterSize;
1915 /* Check whether the cluster is in our list? */
1916 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
1919 findJClusterInJList(jCluster, ranges, gmx::makeConstArrayRef(nbl->cj4));
1923 /* We found an exclusion, clear the corresponding
1926 const unsigned int pairMask =
1927 (1U << (cj_mod_cj4(index) * c_gpuNumClusterPerCell + iCluster));
1928 /* Check if the i-cluster interacts with the j-cluster */
1929 if (nbl_imask0(nbl, index) & pairMask)
1931 const int innerI = (i & (c_clusterSize - 1));
1932 const int innerJ = (jIndex & (c_clusterSize - 1));
1934 /* Determine which j-half (CUDA warp) we are in */
1935 const int jHalf = innerJ / (c_clusterSize / c_nbnxnGpuClusterpairSplit);
1937 nbnxn_excl_t& interactionMask =
1938 get_exclusion_mask(nbl, cj_to_cj4(index), jHalf);
1940 interactionMask.pair[a_mod_wj(innerJ) * c_clusterSize + innerI] &= ~pairMask;
1949 /* Make a new ci entry at the back of nbl->ci */
1950 static void addNewIEntry(NbnxnPairlistCpu* nbl, int ci, int shift, int flags)
1954 ciEntry.shift = shift;
1955 /* Store the interaction flags along with the shift */
1956 ciEntry.shift |= flags;
1957 ciEntry.cj_ind_start = nbl->cj.size();
1958 ciEntry.cj_ind_end = nbl->cj.size();
1959 nbl->ci.push_back(ciEntry);
1962 /* Make a new sci entry at index nbl->nsci */
1963 static void addNewIEntry(NbnxnPairlistGpu* nbl, int sci, int shift, int gmx_unused flags)
1965 nbnxn_sci_t sciEntry;
1967 sciEntry.shift = shift;
1968 sciEntry.cj4_ind_start = nbl->cj4.size();
1969 sciEntry.cj4_ind_end = nbl->cj4.size();
1971 nbl->sci.push_back(sciEntry);
1974 /* Sort the simple j-list cj on exclusions.
1975 * Entries with exclusions will all be sorted to the beginning of the list.
1977 static void sort_cj_excl(nbnxn_cj_t* cj, int ncj, NbnxnPairlistCpuWork* work)
1979 work->cj.resize(ncj);
1981 /* Make a list of the j-cells involving exclusions */
1983 for (int j = 0; j < ncj; j++)
1985 if (cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
1987 work->cj[jnew++] = cj[j];
1990 /* Check if there are exclusions at all or not just the first entry */
1991 if (!((jnew == 0) || (jnew == 1 && cj[0].excl != NBNXN_INTERACTION_MASK_ALL)))
1993 for (int j = 0; j < ncj; j++)
1995 if (cj[j].excl == NBNXN_INTERACTION_MASK_ALL)
1997 work->cj[jnew++] = cj[j];
2000 for (int j = 0; j < ncj; j++)
2002 cj[j] = work->cj[j];
2007 /* Close this simple list i entry */
2008 static void closeIEntry(NbnxnPairlistCpu* nbl,
2009 int gmx_unused sp_max_av,
2010 gmx_bool gmx_unused progBal,
2011 float gmx_unused nsp_tot_est,
2012 int gmx_unused thread,
2013 int gmx_unused nthread)
2015 nbnxn_ci_t& ciEntry = nbl->ci.back();
2017 /* All content of the new ci entry have already been filled correctly,
2018 * we only need to sort and increase counts or remove the entry when empty.
2020 const int jlen = ciEntry.cj_ind_end - ciEntry.cj_ind_start;
2023 sort_cj_excl(nbl->cj.data() + ciEntry.cj_ind_start, jlen, nbl->work.get());
2025 /* The counts below are used for non-bonded pair/flop counts
2026 * and should therefore match the available kernel setups.
2028 if (!(ciEntry.shift & NBNXN_CI_DO_COUL(0)))
2030 nbl->work->ncj_noq += jlen;
2032 else if ((ciEntry.shift & NBNXN_CI_HALF_LJ(0)) || !(ciEntry.shift & NBNXN_CI_DO_LJ(0)))
2034 nbl->work->ncj_hlj += jlen;
2039 /* Entry is empty: remove it */
2044 /* Split sci entry for load balancing on the GPU.
2045 * Splitting ensures we have enough lists to fully utilize the whole GPU.
2046 * With progBal we generate progressively smaller lists, which improves
2047 * load balancing. As we only know the current count on our own thread,
2048 * we will need to estimate the current total amount of i-entries.
2049 * As the lists get concatenated later, this estimate depends
2050 * both on nthread and our own thread index.
2052 static void split_sci_entry(NbnxnPairlistGpu* nbl,
2060 int nsp_max = nsp_target_av;
2064 /* Estimate the total numbers of ci's of the nblist combined
2065 * over all threads using the target number of ci's.
2067 float nsp_est = (nsp_tot_est * thread) / nthread + nbl->nci_tot;
2069 /* The first ci blocks should be larger, to avoid overhead.
2070 * The last ci blocks should be smaller, to improve load balancing.
2071 * The factor 3/2 makes the first block 3/2 times the target average
2072 * and ensures that the total number of blocks end up equal to
2073 * that of equally sized blocks of size nsp_target_av.
2075 nsp_max = static_cast<int>(nsp_target_av * (nsp_tot_est * 1.5 / (nsp_est + nsp_tot_est)));
2078 const int cj4_start = nbl->sci.back().cj4_ind_start;
2079 const int cj4_end = nbl->sci.back().cj4_ind_end;
2080 const int j4len = cj4_end - cj4_start;
2082 if (j4len > 1 && j4len * c_gpuNumClusterPerCell * c_nbnxnGpuJgroupSize > nsp_max)
2084 /* Modify the last ci entry and process the cj4's again */
2090 for (int cj4 = cj4_start; cj4 < cj4_end; cj4++)
2092 int nsp_cj4_p = nsp_cj4;
2093 /* Count the number of cluster pairs in this cj4 group */
2095 for (int p = 0; p < c_gpuNumClusterPerCell * c_nbnxnGpuJgroupSize; p++)
2097 nsp_cj4 += (nbl->cj4[cj4].imei[0].imask >> p) & 1;
2100 /* If adding the current cj4 with nsp_cj4 pairs get us further
2101 * away from our target nsp_max, split the list before this cj4.
2103 if (nsp > 0 && nsp_max - nsp < nsp + nsp_cj4 - nsp_max)
2105 /* Split the list at cj4 */
2106 nbl->sci.back().cj4_ind_end = cj4;
2107 /* Create a new sci entry */
2109 sciNew.sci = nbl->sci.back().sci;
2110 sciNew.shift = nbl->sci.back().shift;
2111 sciNew.cj4_ind_start = cj4;
2112 nbl->sci.push_back(sciNew);
2115 nsp_cj4_e = nsp_cj4_p;
2121 /* Put the remaining cj4's in the last sci entry */
2122 nbl->sci.back().cj4_ind_end = cj4_end;
2124 /* Possibly balance out the last two sci's
2125 * by moving the last cj4 of the second last sci.
2127 if (nsp_sci - nsp_cj4_e >= nsp + nsp_cj4_e)
2129 GMX_ASSERT(nbl->sci.size() >= 2, "We expect at least two elements");
2130 nbl->sci[nbl->sci.size() - 2].cj4_ind_end--;
2131 nbl->sci[nbl->sci.size() - 1].cj4_ind_start--;
2136 /* Clost this super/sub list i entry */
2137 static void closeIEntry(NbnxnPairlistGpu* nbl, int nsp_max_av, gmx_bool progBal, float nsp_tot_est, int thread, int nthread)
2139 nbnxn_sci_t& sciEntry = *getOpenIEntry(nbl);
2141 /* All content of the new ci entry have already been filled correctly,
2142 * we only need to, potentially, split or remove the entry when empty.
2144 int j4len = sciEntry.numJClusterGroups();
2147 /* We can only have complete blocks of 4 j-entries in a list,
2148 * so round the count up before closing.
2150 int ncj4 = (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1) / c_nbnxnGpuJgroupSize;
2151 nbl->work->cj_ind = ncj4 * c_nbnxnGpuJgroupSize;
2155 /* Measure the size of the new entry and potentially split it */
2156 split_sci_entry(nbl, nsp_max_av, progBal, nsp_tot_est, thread, nthread);
2161 /* Entry is empty: remove it */
2162 nbl->sci.pop_back();
2166 /* Syncs the working array before adding another grid pair to the GPU list */
2167 static void sync_work(NbnxnPairlistCpu gmx_unused* nbl) {}
2169 /* Syncs the working array before adding another grid pair to the GPU list */
2170 static void sync_work(NbnxnPairlistGpu* nbl)
2172 nbl->work->cj_ind = nbl->cj4.size() * c_nbnxnGpuJgroupSize;
2175 /* Clears an NbnxnPairlistCpu data structure */
2176 static void clear_pairlist(NbnxnPairlistCpu* nbl)
2182 nbl->ciOuter.clear();
2183 nbl->cjOuter.clear();
2185 nbl->work->ncj_noq = 0;
2186 nbl->work->ncj_hlj = 0;
2189 /* Clears an NbnxnPairlistGpu data structure */
2190 static void clear_pairlist(NbnxnPairlistGpu* nbl)
2194 nbl->excl.resize(1);
2198 /* Clears an atom-atom-style pair list */
2199 static void clear_pairlist_fep(t_nblist* nl)
2203 if (nl->jindex.empty())
2205 nl->jindex.resize(1);
2210 /* Sets a simple list i-cell bounding box, including PBC shift */
2212 set_icell_bb_simple(gmx::ArrayRef<const BoundingBox> bb, int ci, real shx, real shy, real shz, BoundingBox* bb_ci)
2214 bb_ci->lower.x = bb[ci].lower.x + shx;
2215 bb_ci->lower.y = bb[ci].lower.y + shy;
2216 bb_ci->lower.z = bb[ci].lower.z + shz;
2217 bb_ci->upper.x = bb[ci].upper.x + shx;
2218 bb_ci->upper.y = bb[ci].upper.y + shy;
2219 bb_ci->upper.z = bb[ci].upper.z + shz;
2222 /* Sets a simple list i-cell bounding box, including PBC shift */
2223 static inline void set_icell_bb(const Grid& iGrid, int ci, real shx, real shy, real shz, NbnxnPairlistCpuWork* work)
2225 set_icell_bb_simple(iGrid.iBoundingBoxes(), ci, shx, shy, shz, &work->iClusterData.bb[0]);
2229 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2230 static void set_icell_bbxxxx_supersub(gmx::ArrayRef<const float> bb, int ci, real shx, real shy, real shz, float* bb_ci)
2232 constexpr int cellBBStride = packedBoundingBoxesIndex(c_gpuNumClusterPerCell);
2233 constexpr int pbbStride = c_packedBoundingBoxesDimSize;
2234 const int ia = ci * cellBBStride;
2235 for (int m = 0; m < cellBBStride; m += c_packedBoundingBoxesSize)
2237 for (int i = 0; i < pbbStride; i++)
2239 bb_ci[m + 0 * pbbStride + i] = bb[ia + m + 0 * pbbStride + i] + shx;
2240 bb_ci[m + 1 * pbbStride + i] = bb[ia + m + 1 * pbbStride + i] + shy;
2241 bb_ci[m + 2 * pbbStride + i] = bb[ia + m + 2 * pbbStride + i] + shz;
2242 bb_ci[m + 3 * pbbStride + i] = bb[ia + m + 3 * pbbStride + i] + shx;
2243 bb_ci[m + 4 * pbbStride + i] = bb[ia + m + 4 * pbbStride + i] + shy;
2244 bb_ci[m + 5 * pbbStride + i] = bb[ia + m + 5 * pbbStride + i] + shz;
2250 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2251 gmx_unused static void set_icell_bb_supersub(gmx::ArrayRef<const BoundingBox> bb,
2258 for (int i = 0; i < c_gpuNumClusterPerCell; i++)
2260 set_icell_bb_simple(bb, ci * c_gpuNumClusterPerCell + i, shx, shy, shz, &bb_ci[i]);
2264 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2265 gmx_unused static void set_icell_bb(const Grid& iGrid, int ci, real shx, real shy, real shz, NbnxnPairlistGpuWork* work)
2268 set_icell_bbxxxx_supersub(
2269 iGrid.packedBoundingBoxes(), ci, shx, shy, shz, work->iSuperClusterData.bbPacked.data());
2271 set_icell_bb_supersub(iGrid.iBoundingBoxes(), ci, shx, shy, shz, work->iSuperClusterData.bb.data());
2275 /* Copies PBC shifted i-cell atom coordinates x,y,z to working array */
2276 static void icell_set_x_simple(int ci,
2282 NbnxnPairlistCpuWork::IClusterData* iClusterData)
2284 const int ia = ci * c_nbnxnCpuIClusterSize;
2286 for (int i = 0; i < c_nbnxnCpuIClusterSize; i++)
2288 iClusterData->x[i * STRIDE_XYZ + XX] = x[(ia + i) * stride + XX] + shx;
2289 iClusterData->x[i * STRIDE_XYZ + YY] = x[(ia + i) * stride + YY] + shy;
2290 iClusterData->x[i * STRIDE_XYZ + ZZ] = x[(ia + i) * stride + ZZ] + shz;
2294 static void icell_set_x(int ci,
2300 const ClusterDistanceKernelType kernelType,
2301 NbnxnPairlistCpuWork* work)
2306 # ifdef GMX_NBNXN_SIMD_4XN
2307 case ClusterDistanceKernelType::CpuSimd_4xM:
2308 icell_set_x_simd_4xn(ci, shx, shy, shz, stride, x, work);
2311 # ifdef GMX_NBNXN_SIMD_2XNN
2312 case ClusterDistanceKernelType::CpuSimd_2xMM:
2313 icell_set_x_simd_2xnn(ci, shx, shy, shz, stride, x, work);
2317 case ClusterDistanceKernelType::CpuPlainC:
2318 icell_set_x_simple(ci, shx, shy, shz, stride, x, &work->iClusterData);
2320 default: GMX_ASSERT(false, "Unhandled case"); break;
2324 /* Copies PBC shifted super-cell atom coordinates x,y,z to working array */
2325 static void icell_set_x(int ci,
2331 ClusterDistanceKernelType gmx_unused kernelType,
2332 NbnxnPairlistGpuWork* work)
2334 #if !GMX_SIMD4_HAVE_REAL
2336 real* x_ci = work->iSuperClusterData.x.data();
2338 int ia = ci * c_gpuNumClusterPerCell * c_nbnxnGpuClusterSize;
2339 for (int i = 0; i < c_gpuNumClusterPerCell * c_nbnxnGpuClusterSize; i++)
2341 x_ci[i * DIM + XX] = x[(ia + i) * stride + XX] + shx;
2342 x_ci[i * DIM + YY] = x[(ia + i) * stride + YY] + shy;
2343 x_ci[i * DIM + ZZ] = x[(ia + i) * stride + ZZ] + shz;
2346 #else /* !GMX_SIMD4_HAVE_REAL */
2348 real* x_ci = work->iSuperClusterData.xSimd.data();
2350 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2352 for (int i = 0; i < c_nbnxnGpuClusterSize; i += GMX_SIMD4_WIDTH)
2354 int io = si * c_nbnxnGpuClusterSize + i;
2355 int ia = ci * c_gpuNumClusterPerCell * c_nbnxnGpuClusterSize + io;
2356 for (int j = 0; j < GMX_SIMD4_WIDTH; j++)
2358 x_ci[io * DIM + j + XX * GMX_SIMD4_WIDTH] = x[(ia + j) * stride + XX] + shx;
2359 x_ci[io * DIM + j + YY * GMX_SIMD4_WIDTH] = x[(ia + j) * stride + YY] + shy;
2360 x_ci[io * DIM + j + ZZ * GMX_SIMD4_WIDTH] = x[(ia + j) * stride + ZZ] + shz;
2365 #endif /* !GMX_SIMD4_HAVE_REAL */
2368 static real minimum_subgrid_size_xy(const Grid& grid)
2370 const Grid::Dimensions& dims = grid.dimensions();
2372 if (grid.geometry().isSimple)
2374 return std::min(dims.cellSize[XX], dims.cellSize[YY]);
2378 return std::min(dims.cellSize[XX] / c_gpuNumClusterPerCellX,
2379 dims.cellSize[YY] / c_gpuNumClusterPerCellY);
2383 static real effective_buffer_1x1_vs_MxN(const Grid& iGrid, const Grid& jGrid)
2385 const real eff_1x1_buffer_fac_overest = 0.1;
2387 /* Determine an atom-pair list cut-off buffer size for atom pairs,
2388 * to be added to rlist (including buffer) used for MxN.
2389 * This is for converting an MxN list to a 1x1 list. This means we can't
2390 * use the normal buffer estimate, as we have an MxN list in which
2391 * some atom pairs beyond rlist are missing. We want to capture
2392 * the beneficial effect of buffering by extra pairs just outside rlist,
2393 * while removing the useless pairs that are further away from rlist.
2394 * (Also the buffer could have been set manually not using the estimate.)
2395 * This buffer size is an overestimate.
2396 * We add 10% of the smallest grid sub-cell dimensions.
2397 * Note that the z-size differs per cell and we don't use this,
2398 * so we overestimate.
2399 * With PME, the 10% value gives a buffer that is somewhat larger
2400 * than the effective buffer with a tolerance of 0.005 kJ/mol/ps.
2401 * Smaller tolerances or using RF lead to a smaller effective buffer,
2402 * so 10% gives a safe overestimate.
2404 return eff_1x1_buffer_fac_overest * (minimum_subgrid_size_xy(iGrid) + minimum_subgrid_size_xy(jGrid));
2407 /* Estimates the interaction volume^2 for non-local interactions */
2408 static real nonlocal_vol2(const struct gmx_domdec_zones_t* zones, const rvec ls, real r)
2410 real vol2_est_tot = 0;
2412 /* Here we simply add up the volumes of 1, 2 or 3 1D decomposition
2413 * not home interaction volume^2. As these volumes are not additive,
2414 * this is an overestimate, but it would only be significant in the limit
2415 * of small cells, where we anyhow need to split the lists into
2416 * as small parts as possible.
2419 for (int z = 0; z < zones->n; z++)
2421 if (zones->shift[z][XX] + zones->shift[z][YY] + zones->shift[z][ZZ] == 1)
2426 for (int d = 0; d < DIM; d++)
2428 if (zones->shift[z][d] == 0)
2432 za *= zones->size[z].x1[d] - zones->size[z].x0[d];
2436 /* 4 octants of a sphere */
2437 real vold_est = 0.25 * M_PI * r * r * r * r;
2438 /* 4 quarter pie slices on the edges */
2439 vold_est += 4 * cl * M_PI / 6.0 * r * r * r;
2440 /* One rectangular volume on a face */
2441 vold_est += ca * 0.5 * r * r;
2443 vol2_est_tot += vold_est * za;
2447 return vol2_est_tot;
2450 /* Estimates the average size of a full j-list for super/sub setup */
2451 static void get_nsubpair_target(const Nbnxm::GridSet& gridSet,
2452 const InteractionLocality iloc,
2454 const int min_ci_balanced,
2455 int* nsubpair_target,
2456 float* nsubpair_tot_est)
2458 /* The target value of 36 seems to be the optimum for Kepler.
2459 * Maxwell is less sensitive to the exact value.
2461 const int nsubpair_target_min = 36;
2463 const Grid& grid = gridSet.grids()[0];
2465 /* We don't need to balance list sizes if:
2466 * - We didn't request balancing.
2467 * - The number of grid cells >= the number of lists requested,
2468 * since we will always generate at least #cells lists.
2469 * - We don't have any cells, since then there won't be any lists.
2471 if (min_ci_balanced <= 0 || grid.numCells() >= min_ci_balanced || grid.numCells() == 0)
2473 /* nsubpair_target==0 signals no balancing */
2474 *nsubpair_target = 0;
2475 *nsubpair_tot_est = 0;
2481 const int numAtomsCluster = grid.geometry().numAtomsICluster;
2482 const Grid::Dimensions& dims = grid.dimensions();
2484 ls[XX] = dims.cellSize[XX] / c_gpuNumClusterPerCellX;
2485 ls[YY] = dims.cellSize[YY] / c_gpuNumClusterPerCellY;
2486 ls[ZZ] = numAtomsCluster / (dims.atomDensity * ls[XX] * ls[YY]);
2488 /* The formulas below are a heuristic estimate of the average nsj per si*/
2489 const real r_eff_sup = rlist + nbnxn_get_rlist_effective_inc(numAtomsCluster, ls);
2491 real nsp_est_nl = 0;
2492 if (gridSet.domainSetup().haveMultipleDomains && gridSet.domainSetup().zones->n != 1)
2494 nsp_est_nl = gmx::square(dims.atomDensity / numAtomsCluster)
2495 * nonlocal_vol2(gridSet.domainSetup().zones, ls, r_eff_sup);
2498 real nsp_est = nsp_est_nl;
2499 if (iloc == InteractionLocality::Local)
2501 /* Sub-cell interacts with itself */
2502 real vol_est = ls[XX] * ls[YY] * ls[ZZ];
2503 /* 6/2 rectangular volume on the faces */
2504 vol_est += (ls[XX] * ls[YY] + ls[XX] * ls[ZZ] + ls[YY] * ls[ZZ]) * r_eff_sup;
2505 /* 12/2 quarter pie slices on the edges */
2506 vol_est += 2 * (ls[XX] + ls[YY] + ls[ZZ]) * 0.25 * M_PI * gmx::square(r_eff_sup);
2507 /* 4 octants of a sphere */
2508 vol_est += 0.5 * 4.0 / 3.0 * M_PI * gmx::power3(r_eff_sup);
2510 /* Estimate the number of cluster pairs as the local number of
2511 * clusters times the volume they interact with times the density.
2513 nsp_est = grid.numClusters() * vol_est * dims.atomDensity / numAtomsCluster;
2515 /* Subtract the non-local pair count */
2516 nsp_est -= nsp_est_nl;
2518 /* For small cut-offs nsp_est will be an underestimate.
2519 * With DD nsp_est_nl is an overestimate so nsp_est can get negative.
2520 * So to avoid too small or negative nsp_est we set a minimum of
2521 * all cells interacting with all 3^3 direct neighbors (3^3-1)/2+1=14.
2522 * This might be a slight overestimate for small non-periodic groups of
2523 * atoms as will occur for a local domain with DD, but for small
2524 * groups of atoms we'll anyhow be limited by nsubpair_target_min,
2525 * so this overestimation will not matter.
2527 nsp_est = std::max(nsp_est, grid.numClusters() * 14._real);
2531 fprintf(debug, "nsp_est local %5.1f non-local %5.1f\n", nsp_est, nsp_est_nl);
2535 /* Thus the (average) maximum j-list size should be as follows.
2536 * Since there is overhead, we shouldn't make the lists too small
2537 * (and we can't chop up j-groups) so we use a minimum target size of 36.
2539 *nsubpair_target = std::max(nsubpair_target_min, roundToInt(nsp_est / min_ci_balanced));
2540 *nsubpair_tot_est = nsp_est;
2544 fprintf(debug, "nbl nsp estimate %.1f, nsubpair_target %d\n", nsp_est, *nsubpair_target);
2548 /* Debug list print function */
2549 static void print_nblist_ci_cj(FILE* fp, const NbnxnPairlistCpu& nbl)
2551 for (const nbnxn_ci_t& ciEntry : nbl.ci)
2553 fprintf(fp, "ci %4d shift %2d ncj %3d\n", ciEntry.ci, ciEntry.shift, ciEntry.cj_ind_end - ciEntry.cj_ind_start);
2555 for (int j = ciEntry.cj_ind_start; j < ciEntry.cj_ind_end; j++)
2557 fprintf(fp, " cj %5d imask %x\n", nbl.cj[j].cj, nbl.cj[j].excl);
2562 /* Debug list print function */
2563 static void print_nblist_sci_cj(FILE* fp, const NbnxnPairlistGpu& nbl)
2565 for (const nbnxn_sci_t& sci : nbl.sci)
2567 fprintf(fp, "ci %4d shift %2d ncj4 %2d\n", sci.sci, sci.shift, sci.numJClusterGroups());
2570 for (int j4 = sci.cj4_ind_start; j4 < sci.cj4_ind_end; j4++)
2572 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
2574 fprintf(fp, " sj %5d imask %x\n", nbl.cj4[j4].cj[j], nbl.cj4[j4].imei[0].imask);
2575 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2577 if (nbl.cj4[j4].imei[0].imask & (1U << (j * c_gpuNumClusterPerCell + si)))
2584 fprintf(fp, "ci %4d shift %2d ncj4 %2d ncp %3d\n", sci.sci, sci.shift, sci.numJClusterGroups(), ncp);
2588 /* Combine pair lists *nbl generated on multiple threads nblc */
2589 static void combine_nblists(gmx::ArrayRef<const NbnxnPairlistGpu> nbls, NbnxnPairlistGpu* nblc)
2591 int nsci = nblc->sci.size();
2592 int ncj4 = nblc->cj4.size();
2593 int nexcl = nblc->excl.size();
2594 for (const auto& nbl : nbls)
2596 nsci += nbl.sci.size();
2597 ncj4 += nbl.cj4.size();
2598 nexcl += nbl.excl.size();
2601 /* Resize with the final, combined size, so we can fill in parallel */
2602 /* NOTE: For better performance we should use default initialization */
2603 nblc->sci.resize(nsci);
2604 nblc->cj4.resize(ncj4);
2605 nblc->excl.resize(nexcl);
2607 /* Each thread should copy its own data to the combined arrays,
2608 * as otherwise data will go back and forth between different caches.
2610 const int gmx_unused nthreads = gmx_omp_nthreads_get(ModuleMultiThread::Pairsearch);
2612 #pragma omp parallel for num_threads(nthreads) schedule(static)
2613 for (gmx::index n = 0; n < nbls.ssize(); n++)
2617 /* Determine the offset in the combined data for our thread.
2618 * Note that the original sizes in nblc are lost.
2620 int sci_offset = nsci;
2621 int cj4_offset = ncj4;
2622 int excl_offset = nexcl;
2624 for (gmx::index i = n; i < nbls.ssize(); i++)
2626 sci_offset -= nbls[i].sci.size();
2627 cj4_offset -= nbls[i].cj4.size();
2628 excl_offset -= nbls[i].excl.size();
2631 const NbnxnPairlistGpu& nbli = nbls[n];
2633 for (size_t i = 0; i < nbli.sci.size(); i++)
2635 nblc->sci[sci_offset + i] = nbli.sci[i];
2636 nblc->sci[sci_offset + i].cj4_ind_start += cj4_offset;
2637 nblc->sci[sci_offset + i].cj4_ind_end += cj4_offset;
2640 for (size_t j4 = 0; j4 < nbli.cj4.size(); j4++)
2642 nblc->cj4[cj4_offset + j4] = nbli.cj4[j4];
2643 nblc->cj4[cj4_offset + j4].imei[0].excl_ind += excl_offset;
2644 nblc->cj4[cj4_offset + j4].imei[1].excl_ind += excl_offset;
2647 for (size_t j4 = 0; j4 < nbli.excl.size(); j4++)
2649 nblc->excl[excl_offset + j4] = nbli.excl[j4];
2652 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2655 for (const auto& nbl : nbls)
2657 nblc->nci_tot += nbl.nci_tot;
2661 static void balance_fep_lists(gmx::ArrayRef<std::unique_ptr<t_nblist>> fepLists,
2662 gmx::ArrayRef<PairsearchWork> work)
2664 const int numLists = fepLists.ssize();
2668 /* Nothing to balance */
2672 /* Count the total i-lists and pairs */
2675 for (const auto& list : fepLists)
2677 nri_tot += list->nri;
2678 nrj_tot += list->nrj;
2681 const int nrj_target = (nrj_tot + numLists - 1) / numLists;
2683 GMX_ASSERT(gmx_omp_nthreads_get(ModuleMultiThread::Nonbonded) == numLists,
2684 "We should have as many work objects as FEP lists");
2686 #pragma omp parallel for schedule(static) num_threads(numLists)
2687 for (int th = 0; th < numLists; th++)
2691 t_nblist* nbl = work[th].nbl_fep.get();
2693 /* Note that here we allocate for the total size, instead of
2694 * a per-thread esimate (which is hard to obtain).
2696 if (nri_tot > nbl->maxnri)
2698 nbl->maxnri = over_alloc_large(nri_tot);
2699 reallocate_nblist(nbl);
2701 if (nri_tot > nbl->maxnri || nrj_tot > nbl->maxnrj)
2703 nbl->maxnrj = over_alloc_small(nrj_tot);
2704 nbl->jjnr.resize(nbl->maxnrj);
2705 nbl->excl_fep.resize(nbl->maxnrj);
2708 clear_pairlist_fep(nbl);
2710 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2713 /* Loop over the source lists and assign and copy i-entries */
2715 t_nblist* nbld = work[th_dest].nbl_fep.get();
2716 for (int th = 0; th < numLists; th++)
2718 const t_nblist* nbls = fepLists[th].get();
2720 for (int i = 0; i < nbls->nri; i++)
2722 /* The number of pairs in this i-entry */
2723 const int nrj = nbls->jindex[i + 1] - nbls->jindex[i];
2725 /* Decide if list th_dest is too large and we should procede
2726 * to the next destination list.
2728 if (th_dest + 1 < numLists && nbld->nrj > 0
2729 && nbld->nrj + nrj - nrj_target > nrj_target - nbld->nrj)
2732 nbld = work[th_dest].nbl_fep.get();
2735 nbld->iinr[nbld->nri] = nbls->iinr[i];
2736 nbld->gid[nbld->nri] = nbls->gid[i];
2737 nbld->shift[nbld->nri] = nbls->shift[i];
2739 for (int j = nbls->jindex[i]; j < nbls->jindex[i + 1]; j++)
2741 nbld->jjnr[nbld->nrj] = nbls->jjnr[j];
2742 nbld->excl_fep[nbld->nrj] = nbls->excl_fep[j];
2746 nbld->jindex[nbld->nri] = nbld->nrj;
2750 /* Swap the list pointers */
2751 for (int th = 0; th < numLists; th++)
2753 fepLists[th].swap(work[th].nbl_fep);
2757 fprintf(debug, "nbl_fep[%d] nri %4d nrj %4d\n", th, fepLists[th]->nri, fepLists[th]->nrj);
2762 /* Returns the next ci to be processes by our thread */
2763 static gmx_bool next_ci(const Grid& grid, int nth, int ci_block, int* ci_x, int* ci_y, int* ci_b, int* ci)
2768 if (*ci_b == ci_block)
2770 /* Jump to the next block assigned to this task */
2771 *ci += (nth - 1) * ci_block;
2775 if (*ci >= grid.numCells())
2780 while (*ci >= grid.firstCellInColumn(*ci_x * grid.dimensions().numCells[YY] + *ci_y + 1))
2783 if (*ci_y == grid.dimensions().numCells[YY])
2793 /* Returns the distance^2 for which we put cell pairs in the list
2794 * without checking atom pair distances. This is usually < rlist^2.
2796 static float boundingbox_only_distance2(const Grid::Dimensions& iGridDims,
2797 const Grid::Dimensions& jGridDims,
2801 /* If the distance between two sub-cell bounding boxes is less
2802 * than this distance, do not check the distance between
2803 * all particle pairs in the sub-cell, since then it is likely
2804 * that the box pair has atom pairs within the cut-off.
2805 * We use the nblist cut-off minus 0.5 times the average x/y diagonal
2806 * spacing of the sub-cells. Around 40% of the checked pairs are pruned.
2807 * Using more than 0.5 gains at most 0.5%.
2808 * If forces are calculated more than twice, the performance gain
2809 * in the force calculation outweighs the cost of checking.
2810 * Note that with subcell lists, the atom-pair distance check
2811 * is only performed when only 1 out of 8 sub-cells in within range,
2812 * this is because the GPU is much faster than the cpu.
2815 real bbx = 0.5 * (iGridDims.cellSize[XX] + jGridDims.cellSize[XX]);
2816 real bby = 0.5 * (iGridDims.cellSize[YY] + jGridDims.cellSize[YY]);
2819 bbx /= c_gpuNumClusterPerCellX;
2820 bby /= c_gpuNumClusterPerCellY;
2823 real rbb2 = std::max(0.0, rlist - 0.5 * std::sqrt(bbx * bbx + bby * bby));
2829 return static_cast<float>((1 + GMX_FLOAT_EPS) * rbb2);
2833 static int get_ci_block_size(const Grid& iGrid, const bool haveMultipleDomains, const int numLists)
2835 const int ci_block_enum = 5;
2836 const int ci_block_denom = 11;
2837 const int ci_block_min_atoms = 16;
2839 /* Here we decide how to distribute the blocks over the threads.
2840 * We use prime numbers to try to avoid that the grid size becomes
2841 * a multiple of the number of threads, which would lead to some
2842 * threads getting "inner" pairs and others getting boundary pairs,
2843 * which in turns will lead to load imbalance between threads.
2844 * Set the block size as 5/11/ntask times the average number of cells
2845 * in a y,z slab. This should ensure a quite uniform distribution
2846 * of the grid parts of the different thread along all three grid
2847 * zone boundaries with 3D domain decomposition. At the same time
2848 * the blocks will not become too small.
2850 GMX_ASSERT(iGrid.dimensions().numCells[XX] > 0, "Grid can't be empty");
2851 GMX_ASSERT(numLists > 0, "We need at least one list");
2852 int ci_block = (iGrid.numCells() * ci_block_enum)
2853 / (ci_block_denom * iGrid.dimensions().numCells[XX] * numLists);
2855 const int numAtomsPerCell = iGrid.geometry().numAtomsPerCell;
2857 /* Ensure the blocks are not too small: avoids cache invalidation */
2858 if (ci_block * numAtomsPerCell < ci_block_min_atoms)
2860 ci_block = (ci_block_min_atoms + numAtomsPerCell - 1) / numAtomsPerCell;
2863 /* Without domain decomposition
2864 * or with less than 3 blocks per task, divide in nth blocks.
2866 if (!haveMultipleDomains || numLists * 3 * ci_block > iGrid.numCells())
2868 ci_block = (iGrid.numCells() + numLists - 1) / numLists;
2871 if (ci_block > 1 && (numLists - 1) * ci_block >= iGrid.numCells())
2873 /* Some threads have no work. Although reducing the block size
2874 * does not decrease the block count on the first few threads,
2875 * with GPUs better mixing of "upper" cells that have more empty
2876 * clusters results in a somewhat lower max load over all threads.
2877 * Without GPUs the regime of so few atoms per thread is less
2878 * performance relevant, but with 8-wide SIMD the same reasoning
2879 * applies, since the pair list uses 4 i-atom "sub-clusters".
2887 /* Returns the number of bits to right-shift a cluster index to obtain
2888 * the corresponding force buffer flag index.
2890 static int getBufferFlagShift(int numAtomsPerCluster)
2892 int bufferFlagShift = 0;
2893 while ((numAtomsPerCluster << bufferFlagShift) < NBNXN_BUFFERFLAG_SIZE)
2898 return bufferFlagShift;
2901 static bool pairlistIsSimple(const NbnxnPairlistCpu gmx_unused& pairlist)
2906 static bool pairlistIsSimple(const NbnxnPairlistGpu gmx_unused& pairlist)
2911 static void makeClusterListWrapper(NbnxnPairlistCpu* nbl,
2912 const Grid gmx_unused& iGrid,
2915 const int firstCell,
2917 const bool excludeSubDiagonal,
2918 const nbnxn_atomdata_t* nbat,
2921 const ClusterDistanceKernelType kernelType,
2922 int* numDistanceChecks)
2926 case ClusterDistanceKernelType::CpuPlainC:
2927 makeClusterListSimple(
2928 jGrid, nbl, ci, firstCell, lastCell, excludeSubDiagonal, nbat->x().data(), rlist2, rbb2, numDistanceChecks);
2930 #ifdef GMX_NBNXN_SIMD_4XN
2931 case ClusterDistanceKernelType::CpuSimd_4xM:
2932 makeClusterListSimd4xn(
2933 jGrid, nbl, ci, firstCell, lastCell, excludeSubDiagonal, nbat->x().data(), rlist2, rbb2, numDistanceChecks);
2936 #ifdef GMX_NBNXN_SIMD_2XNN
2937 case ClusterDistanceKernelType::CpuSimd_2xMM:
2938 makeClusterListSimd2xnn(
2939 jGrid, nbl, ci, firstCell, lastCell, excludeSubDiagonal, nbat->x().data(), rlist2, rbb2, numDistanceChecks);
2942 default: GMX_ASSERT(false, "Unhandled kernel type");
2946 static void makeClusterListWrapper(NbnxnPairlistGpu* nbl,
2947 const Grid& gmx_unused iGrid,
2950 const int firstCell,
2952 const bool excludeSubDiagonal,
2953 const nbnxn_atomdata_t* nbat,
2956 ClusterDistanceKernelType gmx_unused kernelType,
2957 int* numDistanceChecks)
2959 for (int cj = firstCell; cj <= lastCell; cj++)
2961 make_cluster_list_supersub(
2962 iGrid, jGrid, nbl, ci, cj, excludeSubDiagonal, nbat->xstride, nbat->x().data(), rlist2, rbb2, numDistanceChecks);
2966 static int getNumSimpleJClustersInList(const NbnxnPairlistCpu& nbl)
2968 return nbl.cj.size();
2971 static int getNumSimpleJClustersInList(const gmx_unused NbnxnPairlistGpu& nbl)
2976 static void incrementNumSimpleJClustersInList(NbnxnPairlistCpu* nbl, int ncj_old_j)
2978 nbl->ncjInUse += nbl->cj.size();
2979 nbl->ncjInUse -= ncj_old_j;
2982 static void incrementNumSimpleJClustersInList(NbnxnPairlistGpu gmx_unused* nbl, int gmx_unused ncj_old_j)
2986 static void checkListSizeConsistency(const NbnxnPairlistCpu& nbl, const bool haveFreeEnergy)
2988 GMX_RELEASE_ASSERT(static_cast<size_t>(nbl.ncjInUse) == nbl.cj.size() || haveFreeEnergy,
2989 "Without free-energy all cj pair-list entries should be in use. "
2990 "Note that subsequent code does not make use of the equality, "
2991 "this check is only here to catch bugs");
2994 static void checkListSizeConsistency(const NbnxnPairlistGpu gmx_unused& nbl, bool gmx_unused haveFreeEnergy)
2996 /* We currently can not check consistency here */
2999 /* Set the buffer flags for newly added entries in the list */
3000 static void setBufferFlags(const NbnxnPairlistCpu& nbl,
3001 const int ncj_old_j,
3002 const int gridj_flag_shift,
3003 gmx_bitmask_t* gridj_flag,
3006 if (gmx::ssize(nbl.cj) > ncj_old_j)
3008 int cbFirst = nbl.cj[ncj_old_j].cj >> gridj_flag_shift;
3009 int cbLast = nbl.cj.back().cj >> gridj_flag_shift;
3010 for (int cb = cbFirst; cb <= cbLast; cb++)
3012 bitmask_init_bit(&gridj_flag[cb], th);
3017 static void setBufferFlags(const NbnxnPairlistGpu gmx_unused& nbl,
3018 int gmx_unused ncj_old_j,
3019 int gmx_unused gridj_flag_shift,
3020 gmx_bitmask_t gmx_unused* gridj_flag,
3023 GMX_ASSERT(false, "This function should never be called");
3026 /* Generates the part of pair-list nbl assigned to our thread */
3027 template<typename T>
3028 static void nbnxn_make_pairlist_part(const Nbnxm::GridSet& gridSet,
3031 PairsearchWork* work,
3032 const nbnxn_atomdata_t* nbat,
3033 const ListOfLists<int>& exclusions,
3035 const PairlistType pairlistType,
3037 gmx_bool bFBufferFlag,
3040 float nsubpair_tot_est,
3049 int gridi_flag_shift = 0, gridj_flag_shift = 0;
3050 gmx_bitmask_t* gridj_flag = nullptr;
3052 if (jGrid.geometry().isSimple != pairlistIsSimple(*nbl)
3053 || iGrid.geometry().isSimple != pairlistIsSimple(*nbl))
3055 gmx_incons("Grid incompatible with pair-list");
3059 GMX_ASSERT(nbl->na_ci == jGrid.geometry().numAtomsICluster,
3060 "The cluster sizes in the list and grid should match");
3061 nbl->na_cj = JClusterSizePerListType[pairlistType];
3062 const int na_cj_2log = get_2log(nbl->na_cj);
3068 /* Determine conversion of clusters to flag blocks */
3069 gridi_flag_shift = getBufferFlagShift(nbl->na_ci);
3070 gridj_flag_shift = getBufferFlagShift(nbl->na_cj);
3072 gridj_flag = work->buffer_flags.data();
3075 gridSet.getBox(box);
3077 const bool haveFep = gridSet.haveFep();
3079 const real rlist2 = nbl->rlist * nbl->rlist;
3081 // Select the cluster pair distance kernel type
3082 const ClusterDistanceKernelType kernelType = getClusterDistanceKernelType(pairlistType, *nbat);
3084 if (haveFep && !pairlistIsSimple(*nbl))
3086 /* Determine an atom-pair list cut-off distance for FEP atom pairs.
3087 * We should not simply use rlist, since then we would not have
3088 * the small, effective buffering of the NxN lists.
3089 * The buffer is on overestimate, but the resulting cost for pairs
3090 * beyond rlist is negligible compared to the FEP pairs within rlist.
3092 rl_fep2 = nbl->rlist + effective_buffer_1x1_vs_MxN(iGrid, jGrid);
3096 fprintf(debug, "nbl_fep atom-pair rlist %f\n", rl_fep2);
3098 rl_fep2 = rl_fep2 * rl_fep2;
3101 const Grid::Dimensions& iGridDims = iGrid.dimensions();
3102 const Grid::Dimensions& jGridDims = jGrid.dimensions();
3105 boundingbox_only_distance2(iGridDims, jGridDims, nbl->rlist, pairlistIsSimple(*nbl));
3109 fprintf(debug, "nbl bounding box only distance %f\n", std::sqrt(rbb2));
3112 const bool isIntraGridList = (&iGrid == &jGrid);
3114 /* Set the shift range */
3115 for (int d = 0; d < DIM; d++)
3117 /* Check if we need periodicity shifts.
3118 * Without PBC or with domain decomposition we don't need them.
3120 if (d >= numPbcDimensions(gridSet.domainSetup().pbcType)
3121 || gridSet.domainSetup().haveMultipleDomainsPerDim[d])
3127 const real listRangeCellToCell =
3128 listRangeForGridCellToGridCell(rlist, iGrid.dimensions(), jGrid.dimensions());
3129 if (d == XX && box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < listRangeCellToCell)
3139 const bool bSimple = pairlistIsSimple(*nbl);
3140 gmx::ArrayRef<const BoundingBox> bb_i;
3142 gmx::ArrayRef<const float> pbb_i;
3145 bb_i = iGrid.iBoundingBoxes();
3149 pbb_i = iGrid.packedBoundingBoxes();
3152 /* We use the normal bounding box format for both grid types */
3153 bb_i = iGrid.iBoundingBoxes();
3155 gmx::ArrayRef<const BoundingBox1D> bbcz_i = iGrid.zBoundingBoxes();
3156 gmx::ArrayRef<const int> flags_i = iGrid.clusterFlags();
3157 gmx::ArrayRef<const BoundingBox1D> bbcz_j = jGrid.zBoundingBoxes();
3158 int cell0_i = iGrid.cellOffset();
3163 "nbl nc_i %d col.av. %.1f ci_block %d\n",
3165 iGrid.numCells() / static_cast<double>(iGrid.numColumns()),
3169 int numDistanceChecks = 0;
3171 const real listRangeBBToJCell2 =
3172 gmx::square(listRangeForBoundingBoxToGridCell(rlist, jGrid.dimensions()));
3174 /* Initially ci_b and ci to 1 before where we want them to start,
3175 * as they will both be incremented in next_ci.
3178 int ci = th * ci_block - 1;
3181 while (next_ci(iGrid, nth, ci_block, &ci_x, &ci_y, &ci_b, &ci))
3183 if (bSimple && flags_i[ci] == 0)
3187 const int ncj_old_i = getNumSimpleJClustersInList(*nbl);
3190 if (!isIntraGridList && shp[XX] == 0)
3193 bSimple ? bb_i[ci].upper.x
3194 : iGridDims.lowerCorner[XX] + (real(ci_x) + 1) * iGridDims.cellSize[XX];
3195 if (bx1 < jGridDims.lowerCorner[XX])
3197 d2cx = gmx::square(jGridDims.lowerCorner[XX] - bx1);
3199 if (d2cx >= listRangeBBToJCell2)
3206 int ci_xy = ci_x * iGridDims.numCells[YY] + ci_y;
3208 /* Loop over shift vectors in three dimensions */
3209 for (int tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
3211 const real shz = real(tz) * box[ZZ][ZZ];
3213 real bz0 = bbcz_i[ci].lower + shz;
3214 real bz1 = bbcz_i[ci].upper + shz;
3219 d2z = gmx::square(bz1);
3223 d2z = gmx::square(bz0 - box[ZZ][ZZ]);
3226 const real d2z_cx = d2z + d2cx;
3228 if (d2z_cx >= rlist2)
3233 real bz1_frac = bz1 / real(iGrid.numCellsInColumn(ci_xy));
3238 /* The check with bz1_frac close to or larger than 1 comes later */
3240 for (int ty = -shp[YY]; ty <= shp[YY]; ty++)
3242 const real shy = real(ty) * box[YY][YY] + real(tz) * box[ZZ][YY];
3244 const real by0 = bSimple ? bb_i[ci].lower.y + shy
3245 : iGridDims.lowerCorner[YY]
3246 + (real(ci_y)) * iGridDims.cellSize[YY] + shy;
3247 const real by1 = bSimple ? bb_i[ci].upper.y + shy
3248 : iGridDims.lowerCorner[YY]
3249 + (real(ci_y) + 1) * iGridDims.cellSize[YY] + shy;
3251 int cyf, cyl; //NOLINT(cppcoreguidelines-init-variables)
3252 get_cell_range<YY>(by0, by1, jGridDims, d2z_cx, rlist, &cyf, &cyl);
3260 if (by1 < jGridDims.lowerCorner[YY])
3262 d2z_cy += gmx::square(jGridDims.lowerCorner[YY] - by1);
3264 else if (by0 > jGridDims.upperCorner[YY])
3266 d2z_cy += gmx::square(by0 - jGridDims.upperCorner[YY]);
3269 for (int tx = -shp[XX]; tx <= shp[XX]; tx++)
3271 const int shift = xyzToShiftIndex(tx, ty, tz);
3273 const bool excludeSubDiagonal = (isIntraGridList && shift == gmx::c_centralShiftIndex);
3275 if (c_pbcShiftBackward && isIntraGridList && shift > gmx::c_centralShiftIndex)
3281 real(tx) * box[XX][XX] + real(ty) * box[YY][XX] + real(tz) * box[ZZ][XX];
3283 const real bx0 = bSimple ? bb_i[ci].lower.x + shx
3284 : iGridDims.lowerCorner[XX]
3285 + (real(ci_x)) * iGridDims.cellSize[XX] + shx;
3286 const real bx1 = bSimple ? bb_i[ci].upper.x + shx
3287 : iGridDims.lowerCorner[XX]
3288 + (real(ci_x) + 1) * iGridDims.cellSize[XX] + shx;
3290 int cxf, cxl; //NOLINT(cppcoreguidelines-init-variables)
3291 get_cell_range<XX>(bx0, bx1, jGridDims, d2z_cy, rlist, &cxf, &cxl);
3298 addNewIEntry(nbl, cell0_i + ci, shift, flags_i[ci]);
3300 if ((!c_pbcShiftBackward || excludeSubDiagonal) && cxf < ci_x)
3302 /* Leave the pairs with i > j.
3303 * x is the major index, so skip half of it.
3308 set_icell_bb(iGrid, ci, shx, shy, shz, nbl->work.get());
3310 icell_set_x(cell0_i + ci,
3319 for (int cx = cxf; cx <= cxl; cx++)
3321 const real cx_real = cx;
3323 if (jGridDims.lowerCorner[XX] + cx_real * jGridDims.cellSize[XX] > bx1)
3325 d2zx += gmx::square(jGridDims.lowerCorner[XX]
3326 + cx_real * jGridDims.cellSize[XX] - bx1);
3328 else if (jGridDims.lowerCorner[XX] + (cx_real + 1) * jGridDims.cellSize[XX] < bx0)
3330 d2zx += gmx::square(jGridDims.lowerCorner[XX]
3331 + (cx_real + 1) * jGridDims.cellSize[XX] - bx0);
3334 /* When true, leave the pairs with i > j.
3335 * Skip half of y when i and j have the same x.
3337 const bool skipHalfY = (isIntraGridList && cx == 0
3338 && (!c_pbcShiftBackward || shift == gmx::c_centralShiftIndex)
3340 const int cyf_x = skipHalfY ? ci_y : cyf;
3342 for (int cy = cyf_x; cy <= cyl; cy++)
3344 const int columnStart =
3345 jGrid.firstCellInColumn(cx * jGridDims.numCells[YY] + cy);
3346 const int columnEnd =
3347 jGrid.firstCellInColumn(cx * jGridDims.numCells[YY] + cy + 1);
3349 const real cy_real = cy;
3351 if (jGridDims.lowerCorner[YY] + cy_real * jGridDims.cellSize[YY] > by1)
3353 d2zxy += gmx::square(jGridDims.lowerCorner[YY]
3354 + cy_real * jGridDims.cellSize[YY] - by1);
3356 else if (jGridDims.lowerCorner[YY] + (cy_real + 1) * jGridDims.cellSize[YY] < by0)
3358 d2zxy += gmx::square(jGridDims.lowerCorner[YY]
3359 + (cy_real + 1) * jGridDims.cellSize[YY] - by0);
3361 if (columnStart < columnEnd && d2zxy < listRangeBBToJCell2)
3363 /* To improve efficiency in the common case
3364 * of a homogeneous particle distribution,
3365 * we estimate the index of the middle cell
3366 * in range (midCell). We search down and up
3367 * starting from this index.
3369 * Note that the bbcz_j array contains bounds
3370 * for i-clusters, thus for clusters of 4 atoms.
3371 * For the common case where the j-cluster size
3372 * is 8, we could step with a stride of 2,
3373 * but we do not do this because it would
3374 * complicate this code even more.
3379 bz1_frac * static_cast<real>(columnEnd - columnStart));
3380 if (midCell >= columnEnd)
3382 midCell = columnEnd - 1;
3385 const real d2xy = d2zxy - d2z;
3387 /* Find the lowest cell that can possibly
3389 * Check if we hit the bottom of the grid,
3390 * if the j-cell is below the i-cell and if so,
3391 * if it is within range.
3393 int downTestCell = midCell;
3394 while (downTestCell >= columnStart
3395 && (bbcz_j[downTestCell].upper >= bz0
3396 || d2xy + gmx::square(bbcz_j[downTestCell].upper - bz0) < rlist2))
3400 int firstCell = downTestCell + 1;
3402 /* Find the highest cell that can possibly
3404 * Check if we hit the top of the grid,
3405 * if the j-cell is above the i-cell and if so,
3406 * if it is within range.
3408 int upTestCell = midCell + 1;
3409 while (upTestCell < columnEnd
3410 && (bbcz_j[upTestCell].lower <= bz1
3411 || d2xy + gmx::square(bbcz_j[upTestCell].lower - bz1) < rlist2))
3415 int lastCell = upTestCell - 1;
3417 #define NBNXN_REFCODE 0
3420 /* Simple reference code, for debugging,
3421 * overrides the more complex code above.
3423 firstCell = columnEnd;
3425 for (int k = columnStart; k < columnEnd; k++)
3427 if (d2xy + gmx::square(bbcz_j[k * NNBSBB_D + 1] - bz0) < rlist2
3432 if (d2xy + gmx::square(bbcz_j[k * NNBSBB_D] - bz1) < rlist2
3441 if (isIntraGridList)
3443 /* We want each atom/cell pair only once,
3444 * only use cj >= ci.
3446 if (!c_pbcShiftBackward || shift == gmx::c_centralShiftIndex)
3448 firstCell = std::max(firstCell, ci);
3452 if (firstCell <= lastCell)
3454 GMX_ASSERT(firstCell >= columnStart && lastCell < columnEnd,
3455 "The range should reside within the current grid "
3458 /* For f buffer flags with simple lists */
3459 const int ncj_old_j = getNumSimpleJClustersInList(*nbl);
3461 makeClusterListWrapper(nbl,
3472 &numDistanceChecks);
3476 setBufferFlags(*nbl, ncj_old_j, gridj_flag_shift, gridj_flag, th);
3479 incrementNumSimpleJClustersInList(nbl, ncj_old_j);
3485 if (!exclusions.empty())
3487 /* Set the exclusions for this ci list */
3488 setExclusionsForIEntry(
3489 gridSet, nbl, excludeSubDiagonal, na_cj_2log, *getOpenIEntry(nbl), exclusions);
3494 make_fep_list(gridSet.atomIndices(),
3508 /* Close this ci list */
3509 closeIEntry(nbl, nsubpair_max, progBal, nsubpair_tot_est, th, nth);
3514 if (bFBufferFlag && getNumSimpleJClustersInList(*nbl) > ncj_old_i)
3516 bitmask_init_bit(&(work->buffer_flags[(iGrid.cellOffset() + ci) >> gridi_flag_shift]), th);
3520 work->ndistc = numDistanceChecks;
3522 checkListSizeConsistency(*nbl, haveFep);
3526 fprintf(debug, "number of distance checks %d\n", numDistanceChecks);
3528 print_nblist_statistics(debug, *nbl, gridSet, rlist);
3532 fprintf(debug, "nbl FEP list pairs: %d\n", nbl_fep->nrj);
3537 static void reduce_buffer_flags(gmx::ArrayRef<PairsearchWork> searchWork,
3539 gmx::ArrayRef<gmx_bitmask_t> dest)
3541 for (int s = 0; s < nsrc; s++)
3543 gmx::ArrayRef<gmx_bitmask_t> flags(searchWork[s].buffer_flags);
3545 for (size_t b = 0; b < dest.size(); b++)
3547 gmx_bitmask_t& flag = dest[b];
3548 bitmask_union(&flag, flags[b]);
3553 static void print_reduction_cost(gmx::ArrayRef<const gmx_bitmask_t> flags, int nout)
3560 gmx_bitmask_t mask_0; // NOLINT(cppcoreguidelines-init-variables)
3561 bitmask_init_bit(&mask_0, 0);
3562 for (const gmx_bitmask_t& flag_mask : flags)
3564 if (bitmask_is_equal(flag_mask, mask_0))
3566 /* Only flag 0 is set, no copy of reduction required */
3570 else if (!bitmask_is_zero(flag_mask))
3573 for (int out = 0; out < nout; out++)
3575 if (bitmask_is_set(flag_mask, out))
3591 const auto numFlags = static_cast<double>(flags.size());
3593 "nbnxn reduction: #flag %zu #list %d elem %4.2f, keep %4.2f copy %4.2f red %4.2f\n",
3602 /* Copies the list entries from src to dest when cjStart <= *cjGlobal < cjEnd.
3603 * *cjGlobal is updated with the cj count in src.
3604 * When setFlags==true, flag bit t is set in flag for all i and j clusters.
3606 template<bool setFlags>
3607 static void copySelectedListRange(const nbnxn_ci_t* gmx_restrict srcCi,
3608 const NbnxnPairlistCpu* gmx_restrict src,
3609 NbnxnPairlistCpu* gmx_restrict dest,
3610 gmx_bitmask_t* flag,
3615 const int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3617 dest->ci.push_back(*srcCi);
3618 dest->ci.back().cj_ind_start = dest->cj.size();
3619 dest->ci.back().cj_ind_end = dest->ci.back().cj_ind_start + ncj;
3623 bitmask_init_bit(&flag[srcCi->ci >> iFlagShift], t);
3626 for (int j = srcCi->cj_ind_start; j < srcCi->cj_ind_end; j++)
3628 dest->cj.push_back(src->cj[j]);
3632 /* NOTE: This is relatively expensive, since this
3633 * operation is done for all elements in the list,
3634 * whereas at list generation this is done only
3635 * once for each flag entry.
3637 bitmask_init_bit(&flag[src->cj[j].cj >> jFlagShift], t);
3642 #if defined(__GNUC__) && !defined(__clang__) && __GNUC__ == 7
3643 /* Avoid gcc 7 avx512 loop vectorization bug (actually only needed with -mavx512f) */
3644 # pragma GCC push_options
3645 # pragma GCC optimize("no-tree-vectorize")
3648 /* Returns the number of cluster pairs that are in use summed over all lists */
3649 static int countClusterpairs(gmx::ArrayRef<const NbnxnPairlistCpu> pairlists)
3651 /* gcc 7 with -mavx512f can miss the contributions of 16 consecutive
3652 * elements to the sum calculated in this loop. Above we have disabled
3653 * loop vectorization to avoid this bug.
3656 for (const auto& pairlist : pairlists)
3658 ncjTotal += pairlist.ncjInUse;
3663 #if defined(__GNUC__) && !defined(__clang__) && __GNUC__ == 7
3664 # pragma GCC pop_options
3667 /* This routine re-balances the pairlists such that all are nearly equally
3668 * sized. Only whole i-entries are moved between lists. These are moved
3669 * between the ends of the lists, such that the buffer reduction cost should
3670 * not change significantly.
3671 * Note that all original reduction flags are currently kept. This can lead
3672 * to reduction of parts of the force buffer that could be avoided. But since
3673 * the original lists are quite balanced, this will only give minor overhead.
3675 static void rebalanceSimpleLists(gmx::ArrayRef<const NbnxnPairlistCpu> srcSet,
3676 gmx::ArrayRef<NbnxnPairlistCpu> destSet,
3677 gmx::ArrayRef<PairsearchWork> searchWork)
3679 const int ncjTotal = countClusterpairs(srcSet);
3680 const int numLists = srcSet.ssize();
3681 const int ncjTarget = (ncjTotal + numLists - 1) / numLists;
3683 #pragma omp parallel num_threads(numLists)
3685 int t = gmx_omp_get_thread_num();
3687 int cjStart = ncjTarget * t;
3688 int cjEnd = ncjTarget * (t + 1);
3690 /* The destination pair-list for task/thread t */
3691 NbnxnPairlistCpu& dest = destSet[t];
3693 clear_pairlist(&dest);
3694 dest.na_cj = srcSet[0].na_cj;
3696 /* Note that the flags in the work struct (still) contain flags
3697 * for all entries that are present in srcSet->nbl[t].
3699 gmx_bitmask_t* flag = &searchWork[t].buffer_flags[0];
3701 int iFlagShift = getBufferFlagShift(dest.na_ci);
3702 int jFlagShift = getBufferFlagShift(dest.na_cj);
3705 for (int s = 0; s < numLists && cjGlobal < cjEnd; s++)
3707 const NbnxnPairlistCpu* src = &srcSet[s];
3709 if (cjGlobal + src->ncjInUse > cjStart)
3711 for (gmx::index i = 0; i < gmx::ssize(src->ci) && cjGlobal < cjEnd; i++)
3713 const nbnxn_ci_t* srcCi = &src->ci[i];
3714 int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3715 if (cjGlobal >= cjStart)
3717 /* If the source list is not our own, we need to set
3718 * extra flags (the template bool parameter).
3722 copySelectedListRange<true>(srcCi, src, &dest, flag, iFlagShift, jFlagShift, t);
3726 copySelectedListRange<false>(
3727 srcCi, src, &dest, flag, iFlagShift, jFlagShift, t);
3735 cjGlobal += src->ncjInUse;
3739 dest.ncjInUse = dest.cj.size();
3743 const int ncjTotalNew = countClusterpairs(destSet);
3744 GMX_RELEASE_ASSERT(ncjTotalNew == ncjTotal,
3745 "The total size of the lists before and after rebalancing should match");
3749 /* Returns if the pairlists are so imbalanced that it is worth rebalancing. */
3750 static bool checkRebalanceSimpleLists(gmx::ArrayRef<const NbnxnPairlistCpu> lists)
3752 int numLists = lists.ssize();
3755 for (int s = 0; s < numLists; s++)
3757 ncjMax = std::max(ncjMax, lists[s].ncjInUse);
3758 ncjTotal += lists[s].ncjInUse;
3762 fprintf(debug, "Pair-list ncjMax %d ncjTotal %d\n", ncjMax, ncjTotal);
3764 /* The rebalancing adds 3% extra time to the search. Heuristically we
3765 * determined that under common conditions the non-bonded kernel balance
3766 * improvement will outweigh this when the imbalance is more than 3%.
3767 * But this will, obviously, depend on search vs kernel time and nstlist.
3769 const real rebalanceTolerance = 1.03;
3771 return real(numLists * ncjMax) > real(ncjTotal) * rebalanceTolerance;
3774 /* Perform a count (linear) sort to sort the smaller lists to the end.
3775 * This avoids load imbalance on the GPU, as large lists will be
3776 * scheduled and executed first and the smaller lists later.
3777 * Load balancing between multi-processors only happens at the end
3778 * and there smaller lists lead to more effective load balancing.
3779 * The sorting is done on the cj4 count, not on the actual pair counts.
3780 * Not only does this make the sort faster, but it also results in
3781 * better load balancing than using a list sorted on exact load.
3782 * This function swaps the pointer in the pair list to avoid a copy operation.
3784 static void sort_sci(NbnxnPairlistGpu* nbl)
3786 if (nbl->cj4.size() <= nbl->sci.size())
3788 /* nsci = 0 or all sci have size 1, sorting won't change the order */
3792 NbnxnPairlistGpuWork& work = *nbl->work;
3794 /* We will distinguish differences up to double the average */
3795 const int m = static_cast<int>((2 * ssize(nbl->cj4)) / ssize(nbl->sci));
3797 /* Resize work.sci_sort so we can sort into it */
3798 work.sci_sort.resize(nbl->sci.size());
3800 std::vector<int>& sort = work.sortBuffer;
3801 /* Set up m + 1 entries in sort, initialized at 0 */
3803 sort.resize(m + 1, 0);
3804 /* Count the entries of each size */
3805 for (const nbnxn_sci_t& sci : nbl->sci)
3807 int i = std::min(m, sci.numJClusterGroups());
3810 /* Calculate the offset for each count */
3813 for (gmx::index i = m - 1; i >= 0; i--)
3816 sort[i] = sort[i + 1] + s0;
3820 /* Sort entries directly into place */
3821 gmx::ArrayRef<nbnxn_sci_t> sci_sort = work.sci_sort;
3822 for (const nbnxn_sci_t& sci : nbl->sci)
3824 int i = std::min(m, sci.numJClusterGroups());
3825 sci_sort[sort[i]++] = sci;
3828 /* Swap the sci pointers so we use the new, sorted list */
3829 std::swap(nbl->sci, work.sci_sort);
3832 /* Returns the i-zone range for pairlist construction for the give locality */
3833 static Range<int> getIZoneRange(const Nbnxm::GridSet::DomainSetup& domainSetup,
3834 const InteractionLocality locality)
3836 if (domainSetup.doTestParticleInsertion)
3838 /* With TPI we do grid 1, the inserted molecule, versus grid 0, the rest */
3841 else if (locality == InteractionLocality::Local)
3843 /* Local: only zone (grid) 0 vs 0 */
3848 /* Non-local: we need all i-zones */
3849 return { 0, int(domainSetup.zones->iZones.size()) };
3853 /* Returns the j-zone range for pairlist construction for the give locality and i-zone */
3854 static Range<int> getJZoneRange(const gmx_domdec_zones_t* ddZones,
3855 const InteractionLocality locality,
3858 if (locality == InteractionLocality::Local)
3860 /* Local: zone 0 vs 0 or with TPI 1 vs 0 */
3863 else if (iZone == 0)
3865 /* Non-local: we need to avoid the local (zone 0 vs 0) interactions */
3866 return { 1, *ddZones->iZones[iZone].jZoneRange.end() };
3870 /* Non-local with non-local i-zone: use all j-zones */
3871 return ddZones->iZones[iZone].jZoneRange;
3875 //! Prepares CPU lists produced by the search for dynamic pruning
3876 static void prepareListsForDynamicPruning(gmx::ArrayRef<NbnxnPairlistCpu> lists);
3878 void PairlistSet::constructPairlists(gmx::InteractionLocality locality,
3879 const Nbnxm::GridSet& gridSet,
3880 gmx::ArrayRef<PairsearchWork> searchWork,
3881 nbnxn_atomdata_t* nbat,
3882 const ListOfLists<int>& exclusions,
3883 const int minimumIlistCountForGpuBalancing,
3885 SearchCycleCounting* searchCycleCounting)
3887 const real rlist = params_.rlistOuter;
3889 const int numLists = (isCpuType_ ? cpuLists_.size() : gpuLists_.size());
3893 fprintf(debug, "ns making %d nblists\n", numLists);
3896 nbat->bUseBufferFlags = (nbat->out.size() > 1);
3897 /* We should re-init the flags before making the first list */
3898 if (nbat->bUseBufferFlags && locality == InteractionLocality::Local)
3900 resizeAndZeroBufferFlags(&nbat->buffer_flags, nbat->numAtoms());
3903 int nsubpair_target = 0;
3904 float nsubpair_tot_est = 0.0F;
3905 if (!isCpuType_ && minimumIlistCountForGpuBalancing > 0)
3907 get_nsubpair_target(
3908 gridSet, locality, rlist, minimumIlistCountForGpuBalancing, &nsubpair_target, &nsubpair_tot_est);
3911 /* Clear all pair-lists */
3912 for (int th = 0; th < numLists; th++)
3916 clear_pairlist(&cpuLists_[th]);
3920 clear_pairlist(&gpuLists_[th]);
3923 if (params_.haveFep)
3925 clear_pairlist_fep(fepLists_[th].get());
3929 const gmx_domdec_zones_t* ddZones = gridSet.domainSetup().zones;
3930 GMX_ASSERT(locality == InteractionLocality::Local || ddZones != nullptr,
3931 "Nonlocal interaction locality with null ddZones.");
3933 const auto iZoneRange = getIZoneRange(gridSet.domainSetup(), locality);
3935 for (const int iZone : iZoneRange)
3937 const Grid& iGrid = gridSet.grids()[iZone];
3939 const auto jZoneRange = getJZoneRange(ddZones, locality, iZone);
3941 for (int jZone : jZoneRange)
3943 const Grid& jGrid = gridSet.grids()[jZone];
3947 fprintf(debug, "ns search grid %d vs %d\n", iZone, jZone);
3950 searchCycleCounting->start(enbsCCsearch);
3952 const int ci_block =
3953 get_ci_block_size(iGrid, gridSet.domainSetup().haveMultipleDomains, numLists);
3955 /* With GPU: generate progressively smaller lists for
3956 * load balancing for local only or non-local with 2 zones.
3958 const bool progBal = (locality == InteractionLocality::Local || ddZones->n <= 2);
3960 #pragma omp parallel for num_threads(numLists) schedule(static)
3961 for (int th = 0; th < numLists; th++)
3965 /* Re-init the thread-local work flag data before making
3966 * the first list (not an elegant conditional).
3968 if (nbat->bUseBufferFlags && (iZone == 0 && jZone == 0))
3970 resizeAndZeroBufferFlags(&searchWork[th].buffer_flags, nbat->numAtoms());
3973 if (combineLists_ && th > 0)
3975 GMX_ASSERT(!isCpuType_, "Can only combine GPU lists");
3977 clear_pairlist(&gpuLists_[th]);
3980 PairsearchWork& work = searchWork[th];
3982 work.cycleCounter.start();
3984 t_nblist* fepListPtr = (fepLists_.empty() ? nullptr : fepLists_[th].get());
3986 /* Divide the i cells equally over the pairlists */
3989 nbnxn_make_pairlist_part(gridSet,
3996 params_.pairlistType,
3998 nbat->bUseBufferFlags,
4009 nbnxn_make_pairlist_part(gridSet,
4016 params_.pairlistType,
4018 nbat->bUseBufferFlags,
4028 work.cycleCounter.stop();
4030 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
4032 searchCycleCounting->stop(enbsCCsearch);
4037 for (int th = 0; th < numLists; th++)
4039 inc_nrnb(nrnb, eNR_NBNXN_DIST2, searchWork[th].ndistc);
4043 const NbnxnPairlistCpu& nbl = cpuLists_[th];
4044 np_tot += nbl.cj.size();
4045 np_noq += nbl.work->ncj_noq;
4046 np_hlj += nbl.work->ncj_hlj;
4050 const NbnxnPairlistGpu& nbl = gpuLists_[th];
4051 /* This count ignores potential subsequent pair pruning */
4052 np_tot += nbl.nci_tot;
4055 const int nap = isCpuType_ ? cpuLists_[0].na_ci * cpuLists_[0].na_cj
4056 : gmx::square(gpuLists_[0].na_ci);
4058 natpair_ljq_ = (np_tot - np_noq) * nap - np_hlj * nap / 2;
4059 natpair_lj_ = np_noq * nap;
4060 natpair_q_ = np_hlj * nap / 2;
4062 if (combineLists_ && numLists > 1)
4064 GMX_ASSERT(!isCpuType_, "Can only combine GPU lists");
4066 searchCycleCounting->start(enbsCCcombine);
4068 combine_nblists(gmx::constArrayRefFromArray(&gpuLists_[1], numLists - 1), &gpuLists_[0]);
4070 searchCycleCounting->stop(enbsCCcombine);
4077 if (numLists > 1 && checkRebalanceSimpleLists(cpuLists_))
4079 rebalanceSimpleLists(cpuLists_, cpuListsWork_, searchWork);
4081 /* Swap the sets of pair lists */
4082 cpuLists_.swap(cpuListsWork_);
4087 /* Sort the entries on size, large ones first */
4088 if (combineLists_ || gpuLists_.size() == 1)
4090 sort_sci(&gpuLists_[0]);
4094 #pragma omp parallel for num_threads(numLists) schedule(static)
4095 for (int th = 0; th < numLists; th++)
4099 sort_sci(&gpuLists_[th]);
4101 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
4106 if (nbat->bUseBufferFlags)
4108 reduce_buffer_flags(searchWork, numLists, nbat->buffer_flags);
4111 if (gridSet.haveFep())
4113 /* Balance the free-energy lists over all the threads */
4114 balance_fep_lists(fepLists_, searchWork);
4119 /* This is a fresh list, so not pruned, stored using ci.
4120 * ciOuter is invalid at this point.
4122 GMX_ASSERT(cpuLists_[0].ciOuter.empty(), "ciOuter is invalid so it should be empty");
4125 /* If we have more than one list, they either got rebalancing (CPU)
4126 * or combined (GPU), so we should dump the final result to debug.
4130 if (isCpuType_ && cpuLists_.size() > 1)
4132 for (auto& cpuList : cpuLists_)
4134 print_nblist_statistics(debug, cpuList, gridSet, rlist);
4137 else if (!isCpuType_ && gpuLists_.size() > 1)
4139 print_nblist_statistics(debug, gpuLists_[0], gridSet, rlist);
4149 for (auto& cpuList : cpuLists_)
4151 print_nblist_ci_cj(debug, cpuList);
4156 print_nblist_sci_cj(debug, gpuLists_[0]);
4160 if (nbat->bUseBufferFlags)
4162 print_reduction_cost(nbat->buffer_flags, numLists);
4166 if (params_.useDynamicPruning && isCpuType_)
4168 prepareListsForDynamicPruning(cpuLists_);
4172 void PairlistSets::construct(const InteractionLocality iLocality,
4173 PairSearch* pairSearch,
4174 nbnxn_atomdata_t* nbat,
4175 const ListOfLists<int>& exclusions,
4179 const auto& gridSet = pairSearch->gridSet();
4180 const auto* ddZones = gridSet.domainSetup().zones;
4182 /* The Nbnxm code can also work with more exclusions than those in i-zones only
4183 * when using DD, but the equality check can catch more issues.
4186 exclusions.empty() || (!ddZones && exclusions.ssize() == gridSet.numRealAtomsTotal())
4187 || (ddZones && exclusions.ssize() == ddZones->cg_range[ddZones->iZones.size()]),
4188 "exclusions should either be empty or the number of lists should match the number of "
4191 pairlistSet(iLocality).constructPairlists(iLocality,
4196 minimumIlistCountForGpuBalancing_,
4198 &pairSearch->cycleCounting_);
4200 if (iLocality == InteractionLocality::Local)
4202 outerListCreationStep_ = step;
4206 GMX_RELEASE_ASSERT(outerListCreationStep_ == step,
4207 "Outer list should be created at the same step as the inner list");
4210 /* Special performance logging stuff (env.var. GMX_NBNXN_CYCLE) */
4211 if (iLocality == InteractionLocality::Local)
4213 pairSearch->cycleCounting_.searchCount_++;
4215 if (pairSearch->cycleCounting_.recordCycles_
4216 && (!pairSearch->gridSet().domainSetup().haveMultipleDomains || iLocality == InteractionLocality::NonLocal)
4217 && pairSearch->cycleCounting_.searchCount_ % 100 == 0)
4219 pairSearch->cycleCounting_.printCycles(stderr, pairSearch->work());
4223 void nonbonded_verlet_t::constructPairlist(const InteractionLocality iLocality,
4224 const ListOfLists<int>& exclusions,
4228 pairlistSets_->construct(iLocality, pairSearch_.get(), nbat.get(), exclusions, step, nrnb);
4232 /* Launch the transfer of the pairlist to the GPU.
4234 * NOTE: The launch overhead is currently not timed separately
4236 Nbnxm::gpu_init_pairlist(gpu_nbv, pairlistSets().pairlistSet(iLocality).gpuList(), iLocality);
4240 static void prepareListsForDynamicPruning(gmx::ArrayRef<NbnxnPairlistCpu> lists)
4242 /* TODO: Restructure the lists so we have actual outer and inner
4243 * list objects so we can set a single pointer instead of
4244 * swapping several pointers.
4247 for (auto& list : lists)
4249 /* The search produced a list in ci/cj.
4250 * Swap the list pointers so we get the outer list is ciOuter,cjOuter
4251 * and we can prune that to get an inner list in ci/cj.
4253 GMX_RELEASE_ASSERT(list.ciOuter.empty() && list.cjOuter.empty(),
4254 "The outer lists should be empty before preparation");
4256 std::swap(list.ci, list.ciOuter);
4257 std::swap(list.cj, list.cjOuter);