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/utilities.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/mdlib/gmx_omp_nthreads.h"
55 #include "gromacs/mdtypes/group.h"
56 #include "gromacs/mdtypes/md_enums.h"
57 #include "gromacs/mdtypes/nblist.h"
58 #include "gromacs/nbnxm/atomdata.h"
59 #include "gromacs/nbnxm/gpu_data_mgmt.h"
60 #include "gromacs/pbcutil/ishift.h"
61 #include "gromacs/pbcutil/pbc.h"
62 #include "gromacs/simd/simd.h"
63 #include "gromacs/simd/vector_operations.h"
64 #include "gromacs/utility/exceptions.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/gmxomp.h"
67 #include "gromacs/utility/listoflists.h"
68 #include "gromacs/utility/smalloc.h"
70 #include "boundingboxes.h"
71 #include "clusterdistancekerneltype.h"
73 #include "nbnxm_geometry.h"
74 #include "nbnxm_simd.h"
75 #include "pairlistset.h"
76 #include "pairlistsets.h"
77 #include "pairlistwork.h"
78 #include "pairsearch.h"
80 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
82 using BoundingBox = Nbnxm::BoundingBox; // TODO: Remove when refactoring this file
83 using BoundingBox1D = Nbnxm::BoundingBox1D; // TODO: Remove when refactoring this file
85 using Grid = Nbnxm::Grid; // TODO: Remove when refactoring this file
87 // Convience alias for partial Nbnxn namespace usage
88 using InteractionLocality = gmx::InteractionLocality;
90 /* We shift the i-particles backward for PBC.
91 * This leads to more conditionals than shifting forward.
92 * We do this to get more balanced pair lists.
94 constexpr bool c_pbcShiftBackward = true;
96 /* Layout for the nonbonded NxN pair lists */
97 enum class NbnxnLayout
99 NoSimd4x4, // i-cluster size 4, j-cluster size 4
100 Simd4xN, // i-cluster size 4, j-cluster size SIMD width
101 Simd2xNN, // i-cluster size 4, j-cluster size half SIMD width
102 Gpu8x8x8 // i-cluster size 8, j-cluster size 8 + super-clustering
105 #if defined(GMX_NBNXN_SIMD_4XN) || defined(GMX_NBNXN_SIMD_2XNN)
106 /* Returns the j-cluster size */
107 template<NbnxnLayout layout>
108 static constexpr int jClusterSize()
110 static_assert(layout == NbnxnLayout::NoSimd4x4 || layout == NbnxnLayout::Simd4xN
111 || layout == NbnxnLayout::Simd2xNN,
112 "Currently jClusterSize only supports CPU layouts");
114 return layout == NbnxnLayout::Simd4xN
115 ? GMX_SIMD_REAL_WIDTH
116 : (layout == NbnxnLayout::Simd2xNN ? GMX_SIMD_REAL_WIDTH / 2 : c_nbnxnCpuIClusterSize);
119 /*! \brief Returns the j-cluster index given the i-cluster index.
121 * \tparam jClusterSize The number of atoms in a j-cluster
122 * \tparam jSubClusterIndex The j-sub-cluster index (0/1), used when size(j-cluster) <
123 * size(i-cluster) \param[in] ci The i-cluster index
125 template<int jClusterSize, int jSubClusterIndex>
126 static inline int cjFromCi(int ci)
128 static_assert(jClusterSize == c_nbnxnCpuIClusterSize / 2 || jClusterSize == c_nbnxnCpuIClusterSize
129 || jClusterSize == c_nbnxnCpuIClusterSize * 2,
130 "Only j-cluster sizes 2, 4 and 8 are currently implemented");
132 static_assert(jSubClusterIndex == 0 || jSubClusterIndex == 1,
133 "Only sub-cluster indices 0 and 1 are supported");
135 if (jClusterSize == c_nbnxnCpuIClusterSize / 2)
137 if (jSubClusterIndex == 0)
143 return ((ci + 1) << 1) - 1;
146 else if (jClusterSize == c_nbnxnCpuIClusterSize)
156 /*! \brief Returns the j-cluster index given the i-cluster index.
158 * \tparam layout The pair-list layout
159 * \tparam jSubClusterIndex The j-sub-cluster index (0/1), used when size(j-cluster) <
160 * size(i-cluster) \param[in] ci The i-cluster index
162 template<NbnxnLayout layout, int jSubClusterIndex>
163 static inline int cjFromCi(int ci)
165 constexpr int clusterSize = jClusterSize<layout>();
167 return cjFromCi<clusterSize, jSubClusterIndex>(ci);
170 /* Returns the nbnxn coordinate data index given the i-cluster index */
171 template<NbnxnLayout layout>
172 static inline int xIndexFromCi(int ci)
174 constexpr int clusterSize = jClusterSize<layout>();
176 static_assert(clusterSize == c_nbnxnCpuIClusterSize / 2 || clusterSize == c_nbnxnCpuIClusterSize
177 || clusterSize == c_nbnxnCpuIClusterSize * 2,
178 "Only j-cluster sizes 2, 4 and 8 are currently implemented");
180 if (clusterSize <= c_nbnxnCpuIClusterSize)
182 /* Coordinates are stored packed in groups of 4 */
183 return ci * STRIDE_P4;
187 /* Coordinates packed in 8, i-cluster size is half the packing width */
188 return (ci >> 1) * STRIDE_P8 + (ci & 1) * (c_packX8 >> 1);
192 /* Returns the nbnxn coordinate data index given the j-cluster index */
193 template<NbnxnLayout layout>
194 static inline int xIndexFromCj(int cj)
196 constexpr int clusterSize = jClusterSize<layout>();
198 static_assert(clusterSize == c_nbnxnCpuIClusterSize / 2 || clusterSize == c_nbnxnCpuIClusterSize
199 || clusterSize == c_nbnxnCpuIClusterSize * 2,
200 "Only j-cluster sizes 2, 4 and 8 are currently implemented");
202 if (clusterSize == c_nbnxnCpuIClusterSize / 2)
204 /* Coordinates are stored packed in groups of 4 */
205 return (cj >> 1) * STRIDE_P4 + (cj & 1) * (c_packX4 >> 1);
207 else if (clusterSize == c_nbnxnCpuIClusterSize)
209 /* Coordinates are stored packed in groups of 4 */
210 return cj * STRIDE_P4;
214 /* Coordinates are stored packed in groups of 8 */
215 return cj * STRIDE_P8;
218 #endif // defined(GMX_NBNXN_SIMD_4XN) || defined(GMX_NBNXN_SIMD_2XNN)
221 void nbnxn_init_pairlist_fep(t_nblist* nl)
223 nl->type = GMX_NBLIST_INTERACTION_FREE_ENERGY;
224 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
225 /* The interaction functions are set in the free energy kernel fuction */
238 nl->jindex = nullptr;
240 nl->excl_fep = nullptr;
243 static constexpr int sizeNeededForBufferFlags(const int numAtoms)
245 return (numAtoms + NBNXN_BUFFERFLAG_SIZE - 1) / NBNXN_BUFFERFLAG_SIZE;
248 // Resets current flags to 0 and adds more flags if needed.
249 static void resizeAndZeroBufferFlags(std::vector<gmx_bitmask_t>* flags, const int numAtoms)
252 flags->resize(sizeNeededForBufferFlags(numAtoms), gmx_bitmask_t{ 0 });
256 /* Returns the pair-list cutoff between a bounding box and a grid cell given an atom-to-atom pair-list cutoff
258 * Given a cutoff distance between atoms, this functions returns the cutoff
259 * distance2 between a bounding box of a group of atoms and a grid cell.
260 * Since atoms can be geometrically outside of the cell they have been
261 * assigned to (when atom groups instead of individual atoms are assigned
262 * to cells), this distance returned can be larger than the input.
264 static real listRangeForBoundingBoxToGridCell(real rlist, const Grid::Dimensions& gridDims)
266 return rlist + gridDims.maxAtomGroupRadius;
268 /* Returns the pair-list cutoff between a grid cells given an atom-to-atom pair-list cutoff
270 * Given a cutoff distance between atoms, this functions returns the cutoff
271 * distance2 between two grid cells.
272 * Since atoms can be geometrically outside of the cell they have been
273 * assigned to (when atom groups instead of individual atoms are assigned
274 * to cells), this distance returned can be larger than the input.
276 static real listRangeForGridCellToGridCell(real rlist,
277 const Grid::Dimensions& iGridDims,
278 const Grid::Dimensions& jGridDims)
280 return rlist + iGridDims.maxAtomGroupRadius + jGridDims.maxAtomGroupRadius;
283 /* Determines the cell range along one dimension that
284 * the bounding box b0 - b1 sees.
288 get_cell_range(real b0, real b1, const Grid::Dimensions& jGridDims, real d2, real rlist, int* cf, int* cl)
290 real listRangeBBToCell2 = gmx::square(listRangeForBoundingBoxToGridCell(rlist, jGridDims));
291 real distanceInCells = (b0 - jGridDims.lowerCorner[dim]) * jGridDims.invCellSize[dim];
292 *cf = std::max(static_cast<int>(distanceInCells), 0);
295 && d2 + gmx::square((b0 - jGridDims.lowerCorner[dim]) - (*cf - 1 + 1) * jGridDims.cellSize[dim])
296 < listRangeBBToCell2)
301 *cl = std::min(static_cast<int>((b1 - jGridDims.lowerCorner[dim]) * jGridDims.invCellSize[dim]),
302 jGridDims.numCells[dim] - 1);
303 while (*cl < jGridDims.numCells[dim] - 1
304 && d2 + gmx::square((*cl + 1) * jGridDims.cellSize[dim] - (b1 - jGridDims.lowerCorner[dim]))
305 < listRangeBBToCell2)
311 /* Reference code calculating the distance^2 between two bounding boxes */
313 static float box_dist2(float bx0, float bx1, float by0,
314 float by1, float bz0, float bz1,
315 const BoundingBox *bb)
318 float dl, dh, dm, dm0;
322 dl = bx0 - bb->upper.x;
323 dh = bb->lower.x - bx1;
324 dm = std::max(dl, dh);
325 dm0 = std::max(dm, 0.0f);
328 dl = by0 - bb->upper.y;
329 dh = bb->lower.y - by1;
330 dm = std::max(dl, dh);
331 dm0 = std::max(dm, 0.0f);
334 dl = bz0 - bb->upper.z;
335 dh = bb->lower.z - bz1;
336 dm = std::max(dl, dh);
337 dm0 = std::max(dm, 0.0f);
344 #if !NBNXN_SEARCH_BB_SIMD4
346 /*! \brief Plain C code calculating the distance^2 between two bounding boxes in xyz0 format
348 * \param[in] bb_i First bounding box
349 * \param[in] bb_j Second bounding box
351 static float clusterBoundingBoxDistance2(const BoundingBox& bb_i, const BoundingBox& bb_j)
353 float dl = bb_i.lower.x - bb_j.upper.x;
354 float dh = bb_j.lower.x - bb_i.upper.x;
355 float dm = std::max(dl, dh);
356 float dm0 = std::max(dm, 0.0F);
357 float d2 = dm0 * dm0;
359 dl = bb_i.lower.y - bb_j.upper.y;
360 dh = bb_j.lower.y - bb_i.upper.y;
361 dm = std::max(dl, dh);
362 dm0 = std::max(dm, 0.0F);
365 dl = bb_i.lower.z - bb_j.upper.z;
366 dh = bb_j.lower.z - bb_i.upper.z;
367 dm = std::max(dl, dh);
368 dm0 = std::max(dm, 0.0F);
374 #else /* NBNXN_SEARCH_BB_SIMD4 */
376 /*! \brief 4-wide SIMD code calculating the distance^2 between two bounding boxes in xyz0 format
378 * \param[in] bb_i First bounding box, should be aligned for 4-wide SIMD
379 * \param[in] bb_j Second bounding box, should be aligned for 4-wide SIMD
381 static float clusterBoundingBoxDistance2(const BoundingBox& bb_i, const BoundingBox& bb_j)
383 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
386 const Simd4Float bb_i_S0 = load4(bb_i.lower.ptr());
387 const Simd4Float bb_i_S1 = load4(bb_i.upper.ptr());
388 const Simd4Float bb_j_S0 = load4(bb_j.lower.ptr());
389 const Simd4Float bb_j_S1 = load4(bb_j.upper.ptr());
391 const Simd4Float dl_S = bb_i_S0 - bb_j_S1;
392 const Simd4Float dh_S = bb_j_S0 - bb_i_S1;
394 const Simd4Float dm_S = max(dl_S, dh_S);
395 const Simd4Float dm0_S = max(dm_S, simd4SetZeroF());
397 return dotProduct(dm0_S, dm0_S);
400 /* Calculate bb bounding distances of bb_i[si,...,si+3] and store them in d2 */
401 template<int boundingBoxStart>
402 static inline void gmx_simdcall clusterBoundingBoxDistance2_xxxx_simd4_inner(const float* bb_i,
404 const Simd4Float xj_l,
405 const Simd4Float yj_l,
406 const Simd4Float zj_l,
407 const Simd4Float xj_h,
408 const Simd4Float yj_h,
409 const Simd4Float zj_h)
411 constexpr int stride = c_packedBoundingBoxesDimSize;
413 const int shi = boundingBoxStart * Nbnxm::c_numBoundingBoxBounds1D * DIM;
415 const Simd4Float zero = setZero();
417 const Simd4Float xi_l = load4(bb_i + shi + 0 * stride);
418 const Simd4Float yi_l = load4(bb_i + shi + 1 * stride);
419 const Simd4Float zi_l = load4(bb_i + shi + 2 * stride);
420 const Simd4Float xi_h = load4(bb_i + shi + 3 * stride);
421 const Simd4Float yi_h = load4(bb_i + shi + 4 * stride);
422 const Simd4Float zi_h = load4(bb_i + shi + 5 * stride);
424 const Simd4Float dx_0 = xi_l - xj_h;
425 const Simd4Float dy_0 = yi_l - yj_h;
426 const Simd4Float dz_0 = zi_l - zj_h;
428 const Simd4Float dx_1 = xj_l - xi_h;
429 const Simd4Float dy_1 = yj_l - yi_h;
430 const Simd4Float dz_1 = zj_l - zi_h;
432 const Simd4Float mx = max(dx_0, dx_1);
433 const Simd4Float my = max(dy_0, dy_1);
434 const Simd4Float mz = max(dz_0, dz_1);
436 const Simd4Float m0x = max(mx, zero);
437 const Simd4Float m0y = max(my, zero);
438 const Simd4Float m0z = max(mz, zero);
440 const Simd4Float d2x = m0x * m0x;
441 const Simd4Float d2y = m0y * m0y;
442 const Simd4Float d2z = m0z * m0z;
444 const Simd4Float d2s = d2x + d2y;
445 const Simd4Float d2t = d2s + d2z;
447 store4(d2 + boundingBoxStart, d2t);
450 /* 4-wide SIMD code for nsi bb distances for bb format xxxxyyyyzzzz */
451 static void clusterBoundingBoxDistance2_xxxx_simd4(const float* bb_j, const int nsi, const float* bb_i, float* d2)
453 constexpr int stride = c_packedBoundingBoxesDimSize;
455 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
458 const Simd4Float xj_l = Simd4Float(bb_j[0 * stride]);
459 const Simd4Float yj_l = Simd4Float(bb_j[1 * stride]);
460 const Simd4Float zj_l = Simd4Float(bb_j[2 * stride]);
461 const Simd4Float xj_h = Simd4Float(bb_j[3 * stride]);
462 const Simd4Float yj_h = Simd4Float(bb_j[4 * stride]);
463 const Simd4Float zj_h = Simd4Float(bb_j[5 * stride]);
465 /* Here we "loop" over si (0,stride) from 0 to nsi with step stride.
466 * But as we know the number of iterations is 1 or 2, we unroll manually.
468 clusterBoundingBoxDistance2_xxxx_simd4_inner<0>(bb_i, d2, xj_l, yj_l, zj_l, xj_h, yj_h, zj_h);
471 clusterBoundingBoxDistance2_xxxx_simd4_inner<stride>(bb_i, d2, xj_l, yj_l, zj_l, xj_h, yj_h, zj_h);
475 #endif /* NBNXN_SEARCH_BB_SIMD4 */
478 /* Returns if any atom pair from two clusters is within distance sqrt(rlist2) */
479 static inline gmx_bool
480 clusterpair_in_range(const NbnxnPairlistGpuWork& work, int si, int csj, int stride, const real* x_j, real rlist2)
482 #if !GMX_SIMD4_HAVE_REAL
485 * All coordinates are stored as xyzxyz...
488 const real* x_i = work.iSuperClusterData.x.data();
490 for (int i = 0; i < c_nbnxnGpuClusterSize; i++)
492 int i0 = (si * c_nbnxnGpuClusterSize + i) * DIM;
493 for (int j = 0; j < c_nbnxnGpuClusterSize; j++)
495 int j0 = (csj * c_nbnxnGpuClusterSize + j) * stride;
497 real d2 = gmx::square(x_i[i0] - x_j[j0]) + gmx::square(x_i[i0 + 1] - x_j[j0 + 1])
498 + gmx::square(x_i[i0 + 2] - x_j[j0 + 2]);
509 #else /* !GMX_SIMD4_HAVE_REAL */
511 /* 4-wide SIMD version.
512 * The coordinates x_i are stored as xxxxyyyy..., x_j is stored xyzxyz...
513 * Using 8-wide AVX(2) is not faster on Intel Sandy Bridge and Haswell.
515 static_assert(c_nbnxnGpuClusterSize == 8 || c_nbnxnGpuClusterSize == 4,
516 "A cluster is hard-coded to 4/8 atoms.");
518 Simd4Real rc2_S = Simd4Real(rlist2);
520 const real* x_i = work.iSuperClusterData.xSimd.data();
522 int dim_stride = c_nbnxnGpuClusterSize * DIM;
523 Simd4Real ix_S0 = load4(x_i + si * dim_stride + 0 * GMX_SIMD4_WIDTH);
524 Simd4Real iy_S0 = load4(x_i + si * dim_stride + 1 * GMX_SIMD4_WIDTH);
525 Simd4Real iz_S0 = load4(x_i + si * dim_stride + 2 * GMX_SIMD4_WIDTH);
527 Simd4Real ix_S1, iy_S1, iz_S1;
528 if (c_nbnxnGpuClusterSize == 8)
530 ix_S1 = load4(x_i + si * dim_stride + 3 * GMX_SIMD4_WIDTH);
531 iy_S1 = load4(x_i + si * dim_stride + 4 * GMX_SIMD4_WIDTH);
532 iz_S1 = load4(x_i + si * dim_stride + 5 * GMX_SIMD4_WIDTH);
534 /* We loop from the outer to the inner particles to maximize
535 * the chance that we find a pair in range quickly and return.
537 int j0 = csj * c_nbnxnGpuClusterSize;
538 int j1 = j0 + c_nbnxnGpuClusterSize - 1;
541 Simd4Real jx0_S, jy0_S, jz0_S;
542 Simd4Real jx1_S, jy1_S, jz1_S;
544 Simd4Real dx_S0, dy_S0, dz_S0;
545 Simd4Real dx_S1, dy_S1, dz_S1;
546 Simd4Real dx_S2, dy_S2, dz_S2;
547 Simd4Real dx_S3, dy_S3, dz_S3;
558 Simd4Bool wco_any_S01, wco_any_S23, wco_any_S;
560 jx0_S = Simd4Real(x_j[j0 * stride + 0]);
561 jy0_S = Simd4Real(x_j[j0 * stride + 1]);
562 jz0_S = Simd4Real(x_j[j0 * stride + 2]);
564 jx1_S = Simd4Real(x_j[j1 * stride + 0]);
565 jy1_S = Simd4Real(x_j[j1 * stride + 1]);
566 jz1_S = Simd4Real(x_j[j1 * stride + 2]);
568 /* Calculate distance */
569 dx_S0 = ix_S0 - jx0_S;
570 dy_S0 = iy_S0 - jy0_S;
571 dz_S0 = iz_S0 - jz0_S;
572 dx_S2 = ix_S0 - jx1_S;
573 dy_S2 = iy_S0 - jy1_S;
574 dz_S2 = iz_S0 - jz1_S;
575 if (c_nbnxnGpuClusterSize == 8)
577 dx_S1 = ix_S1 - jx0_S;
578 dy_S1 = iy_S1 - jy0_S;
579 dz_S1 = iz_S1 - jz0_S;
580 dx_S3 = ix_S1 - jx1_S;
581 dy_S3 = iy_S1 - jy1_S;
582 dz_S3 = iz_S1 - jz1_S;
585 /* rsq = dx*dx+dy*dy+dz*dz */
586 rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
587 rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
588 if (c_nbnxnGpuClusterSize == 8)
590 rsq_S1 = norm2(dx_S1, dy_S1, dz_S1);
591 rsq_S3 = norm2(dx_S3, dy_S3, dz_S3);
594 wco_S0 = (rsq_S0 < rc2_S);
595 wco_S2 = (rsq_S2 < rc2_S);
596 if (c_nbnxnGpuClusterSize == 8)
598 wco_S1 = (rsq_S1 < rc2_S);
599 wco_S3 = (rsq_S3 < rc2_S);
601 if (c_nbnxnGpuClusterSize == 8)
603 wco_any_S01 = wco_S0 || wco_S1;
604 wco_any_S23 = wco_S2 || wco_S3;
605 wco_any_S = wco_any_S01 || wco_any_S23;
609 wco_any_S = wco_S0 || wco_S2;
612 if (anyTrue(wco_any_S))
623 #endif /* !GMX_SIMD4_HAVE_REAL */
626 /* Returns the j-cluster index for index cjIndex in a cj list */
627 static inline int nblCj(gmx::ArrayRef<const nbnxn_cj_t> cjList, int cjIndex)
629 return cjList[cjIndex].cj;
632 /* Returns the j-cluster index for index cjIndex in a cj4 list */
633 static inline int nblCj(gmx::ArrayRef<const nbnxn_cj4_t> cj4List, int cjIndex)
635 return cj4List[cjIndex / c_nbnxnGpuJgroupSize].cj[cjIndex & (c_nbnxnGpuJgroupSize - 1)];
638 /* Returns the i-interaction mask of the j sub-cell for index cj_ind */
639 static unsigned int nbl_imask0(const NbnxnPairlistGpu* nbl, int cj_ind)
641 return nbl->cj4[cj_ind / c_nbnxnGpuJgroupSize].imei[0].imask;
644 NbnxnPairlistCpu::NbnxnPairlistCpu() :
645 na_ci(c_nbnxnCpuIClusterSize),
650 work(std::make_unique<NbnxnPairlistCpuWork>())
654 NbnxnPairlistGpu::NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy) :
655 na_ci(c_nbnxnGpuClusterSize),
656 na_cj(c_nbnxnGpuClusterSize),
657 na_sc(c_gpuNumClusterPerCell * c_nbnxnGpuClusterSize),
659 sci({}, { pinningPolicy }),
660 cj4({}, { pinningPolicy }),
661 excl({}, { pinningPolicy }),
663 work(std::make_unique<NbnxnPairlistGpuWork>())
665 static_assert(c_nbnxnGpuNumClusterPerSupercluster == c_gpuNumClusterPerCell,
666 "The search code assumes that the a super-cluster matches a search grid cell");
668 static_assert(sizeof(cj4[0].imei[0].imask) * 8 >= c_nbnxnGpuJgroupSize * c_gpuNumClusterPerCell,
669 "The i super-cluster cluster interaction mask does not contain a sufficient "
672 static_assert(sizeof(excl[0]) * 8 >= c_nbnxnGpuJgroupSize * c_gpuNumClusterPerCell,
673 "The GPU exclusion mask does not contain a sufficient number of bits");
675 // We always want a first entry without any exclusions
679 // TODO: Move to pairlistset.cpp
680 PairlistSet::PairlistSet(const InteractionLocality locality, const PairlistParams& pairlistParams) :
682 params_(pairlistParams)
684 isCpuType_ = (params_.pairlistType == PairlistType::Simple4x2
685 || params_.pairlistType == PairlistType::Simple4x4
686 || params_.pairlistType == PairlistType::Simple4x8);
687 // Currently GPU lists are always combined
688 combineLists_ = !isCpuType_;
690 const int numLists = gmx_omp_nthreads_get(emntNonbonded);
692 if (!combineLists_ && numLists > NBNXN_BUFFERFLAG_MAX_THREADS)
695 "%d OpenMP threads were requested. Since the non-bonded force buffer reduction "
696 "is prohibitively slow with more than %d threads, we do not allow this. Use %d "
697 "or less OpenMP threads.",
699 NBNXN_BUFFERFLAG_MAX_THREADS,
700 NBNXN_BUFFERFLAG_MAX_THREADS);
705 cpuLists_.resize(numLists);
708 cpuListsWork_.resize(numLists);
713 /* Only list 0 is used on the GPU, use normal allocation for i>0 */
714 gpuLists_.emplace_back(gmx::PinningPolicy::PinnedIfSupported);
715 /* Lists 0 to numLists are use for constructing lists in parallel
716 * on the CPU using numLists threads (and then merged into list 0).
718 for (int i = 1; i < numLists; i++)
720 gpuLists_.emplace_back(gmx::PinningPolicy::CannotBePinned);
725 fepLists_.resize(numLists);
727 /* Execute in order to avoid memory interleaving between threads */
728 #pragma omp parallel for num_threads(numLists) schedule(static)
729 for (int i = 0; i < numLists; i++)
733 /* We used to allocate all normal lists locally on each thread
734 * as well. The question is if allocating the object on the
735 * master thread (but all contained list memory thread local)
736 * impacts performance.
738 fepLists_[i] = std::make_unique<t_nblist>();
739 nbnxn_init_pairlist_fep(fepLists_[i].get());
741 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
746 /* Print statistics of a pair list, used for debug output */
747 static void print_nblist_statistics(FILE* fp,
748 const NbnxnPairlistCpu& nbl,
749 const Nbnxm::GridSet& gridSet,
752 const Grid& grid = gridSet.grids()[0];
753 const Grid::Dimensions& dims = grid.dimensions();
755 fprintf(fp, "nbl nci %zu ncj %d\n", nbl.ci.size(), nbl.ncjInUse);
756 const int numAtomsJCluster = grid.geometry().numAtomsJCluster;
757 const double numAtomsPerCell = nbl.ncjInUse / static_cast<double>(grid.numCells()) * numAtomsJCluster;
759 "nbl na_cj %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
763 nbl.ncjInUse / static_cast<double>(grid.numCells()),
766 / (0.5 * 4.0 / 3.0 * M_PI * rl * rl * rl * grid.numCells() * numAtomsJCluster
767 / (dims.gridSize[XX] * dims.gridSize[YY] * dims.gridSize[ZZ])));
770 "nbl average j cell list length %.1f\n",
771 0.25 * nbl.ncjInUse / std::max(static_cast<double>(nbl.ci.size()), 1.0));
773 int cs[SHIFTS] = { 0 };
775 for (const nbnxn_ci_t& ciEntry : nbl.ci)
777 cs[ciEntry.shift & NBNXN_CI_SHIFT] += ciEntry.cj_ind_end - ciEntry.cj_ind_start;
779 int j = ciEntry.cj_ind_start;
780 while (j < ciEntry.cj_ind_end && nbl.cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
787 "nbl cell pairs, total: %zu excl: %d %.1f%%\n",
790 100 * npexcl / std::max(static_cast<double>(nbl.cj.size()), 1.0));
791 for (int s = 0; s < SHIFTS; s++)
795 fprintf(fp, "nbl shift %2d ncj %3d\n", s, cs[s]);
800 /* Print statistics of a pair lists, used for debug output */
801 static void print_nblist_statistics(FILE* fp,
802 const NbnxnPairlistGpu& nbl,
803 const Nbnxm::GridSet& gridSet,
806 const Grid& grid = gridSet.grids()[0];
807 const Grid::Dimensions& dims = grid.dimensions();
810 "nbl nsci %zu ncj4 %zu nsi %d excl4 %zu\n",
815 const int numAtomsCluster = grid.geometry().numAtomsICluster;
816 const double numAtomsPerCell = nbl.nci_tot / static_cast<double>(grid.numClusters()) * numAtomsCluster;
818 "nbl na_c %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
822 nbl.nci_tot / static_cast<double>(grid.numClusters()),
825 / (0.5 * 4.0 / 3.0 * M_PI * rl * rl * rl * grid.numClusters() * numAtomsCluster
826 / (dims.gridSize[XX] * dims.gridSize[YY] * dims.gridSize[ZZ])));
831 int c[c_gpuNumClusterPerCell + 1] = { 0 };
832 for (const nbnxn_sci_t& sci : nbl.sci)
835 for (int j4 = sci.cj4_ind_start; j4 < sci.cj4_ind_end; j4++)
837 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
840 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
842 if (nbl.cj4[j4].imei[0].imask & (1U << (j * c_gpuNumClusterPerCell + si)))
852 sum_nsp2 += nsp * nsp;
853 nsp_max = std::max(nsp_max, nsp);
855 if (!nbl.sci.empty())
857 sum_nsp /= nbl.sci.size();
858 sum_nsp2 /= nbl.sci.size();
861 "nbl #cluster-pairs: av %.1f stddev %.1f max %d\n",
863 std::sqrt(sum_nsp2 - sum_nsp * sum_nsp),
866 if (!nbl.cj4.empty())
868 for (int b = 0; b <= c_gpuNumClusterPerCell; b++)
871 "nbl j-list #i-subcell %d %7d %4.1f\n",
874 100.0 * c[b] / size_t{ nbl.cj4.size() * c_nbnxnGpuJgroupSize });
879 /* Returns a reference to the exclusion mask for j-cluster-group \p cj4 and warp \p warp
880 * Generates a new exclusion entry when the j-cluster-group uses
881 * the default all-interaction mask at call time, so the returned mask
882 * can be modified when needed.
884 static nbnxn_excl_t& get_exclusion_mask(NbnxnPairlistGpu* nbl, int cj4, int warp)
886 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
888 /* No exclusions set, make a new list entry */
889 const size_t oldSize = nbl->excl.size();
890 GMX_ASSERT(oldSize >= 1, "We should always have entry [0]");
891 /* Add entry with default values: no exclusions */
892 nbl->excl.resize(oldSize + 1);
893 nbl->cj4[cj4].imei[warp].excl_ind = oldSize;
896 return nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
899 /* Sets self exclusions and excludes half of the double pairs in the self cluster-pair \p nbl->cj4[cj4Index].cj[jOffsetInGroup]
901 * \param[in,out] nbl The cluster pair list
902 * \param[in] cj4Index The j-cluster group index into \p nbl->cj4
903 * \param[in] jOffsetInGroup The j-entry offset in \p nbl->cj4[cj4Index]
904 * \param[in] iClusterInCell The i-cluster index in the cell
906 static void setSelfAndNewtonExclusionsGpu(NbnxnPairlistGpu* nbl,
908 const int jOffsetInGroup,
909 const int iClusterInCell)
911 constexpr int numJatomsPerPart = c_nbnxnGpuClusterSize / c_nbnxnGpuClusterpairSplit;
913 /* The exclusions are stored separately for each part of the split */
914 for (int part = 0; part < c_nbnxnGpuClusterpairSplit; part++)
916 const int jOffset = part * numJatomsPerPart;
917 /* Make a new exclusion mask entry for each part, if we don't already have one yet */
918 nbnxn_excl_t& excl = get_exclusion_mask(nbl, cj4Index, part);
920 /* Set all bits with j-index <= i-index */
921 for (int jIndexInPart = 0; jIndexInPart < numJatomsPerPart; jIndexInPart++)
923 for (int i = jOffset + jIndexInPart; i < c_nbnxnGpuClusterSize; i++)
925 excl.pair[jIndexInPart * c_nbnxnGpuClusterSize + i] &=
926 ~(1U << (jOffsetInGroup * c_gpuNumClusterPerCell + iClusterInCell));
932 /* Returns a diagonal or off-diagonal interaction mask for plain C lists */
933 static unsigned int get_imask(gmx_bool rdiag, int ci, int cj)
935 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
938 /* Returns a diagonal or off-diagonal interaction mask for cj-size=2 */
939 gmx_unused static unsigned int get_imask_simd_j2(gmx_bool rdiag, int ci, int cj)
941 return (rdiag && ci * 2 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_0
942 : (rdiag && ci * 2 + 1 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_1
943 : NBNXN_INTERACTION_MASK_ALL));
946 /* Returns a diagonal or off-diagonal interaction mask for cj-size=4 */
947 gmx_unused static unsigned int get_imask_simd_j4(gmx_bool rdiag, int ci, int cj)
949 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
952 /* Returns a diagonal or off-diagonal interaction mask for cj-size=8 */
953 gmx_unused static unsigned int get_imask_simd_j8(gmx_bool rdiag, int ci, int cj)
955 return (rdiag && ci == cj * 2 ? NBNXN_INTERACTION_MASK_DIAG_J8_0
956 : (rdiag && ci == cj * 2 + 1 ? NBNXN_INTERACTION_MASK_DIAG_J8_1
957 : NBNXN_INTERACTION_MASK_ALL));
961 # if GMX_SIMD_REAL_WIDTH == 2
962 # define get_imask_simd_4xn get_imask_simd_j2
964 # if GMX_SIMD_REAL_WIDTH == 4
965 # define get_imask_simd_4xn get_imask_simd_j4
967 # if GMX_SIMD_REAL_WIDTH == 8
968 # define get_imask_simd_4xn get_imask_simd_j8
969 # define get_imask_simd_2xnn get_imask_simd_j4
971 # if GMX_SIMD_REAL_WIDTH == 16
972 # define get_imask_simd_2xnn get_imask_simd_j8
976 /* Plain C code for checking and adding cluster-pairs to the list.
978 * \param[in] gridj The j-grid
979 * \param[in,out] nbl The pair-list to store the cluster pairs in
980 * \param[in] icluster The index of the i-cluster
981 * \param[in] jclusterFirst The first cluster in the j-range
982 * \param[in] jclusterLast The last cluster in the j-range
983 * \param[in] excludeSubDiagonal Exclude atom pairs with i-index > j-index
984 * \param[in] x_j Coordinates for the j-atom, in xyz format
985 * \param[in] rlist2 The squared list cut-off
986 * \param[in] rbb2 The squared cut-off for putting cluster-pairs in the list based on bounding box distance only
987 * \param[in,out] numDistanceChecks The number of distance checks performed
989 static void makeClusterListSimple(const Grid& jGrid,
990 NbnxnPairlistCpu* nbl,
994 bool excludeSubDiagonal,
995 const real* gmx_restrict x_j,
998 int* gmx_restrict numDistanceChecks)
1000 const BoundingBox* gmx_restrict bb_ci = nbl->work->iClusterData.bb.data();
1001 const real* gmx_restrict x_ci = nbl->work->iClusterData.x.data();
1006 while (!InRange && jclusterFirst <= jclusterLast)
1008 real d2 = clusterBoundingBoxDistance2(bb_ci[0], jGrid.jBoundingBoxes()[jclusterFirst]);
1009 *numDistanceChecks += 2;
1011 /* Check if the distance is within the distance where
1012 * we use only the bounding box distance rbb,
1013 * or within the cut-off and there is at least one atom pair
1014 * within the cut-off.
1020 else if (d2 < rlist2)
1022 int cjf_gl = jGrid.cellOffset() + jclusterFirst;
1023 for (int i = 0; i < c_nbnxnCpuIClusterSize && !InRange; i++)
1025 for (int j = 0; j < c_nbnxnCpuIClusterSize; j++)
1029 || (gmx::square(x_ci[i * STRIDE_XYZ + XX]
1030 - x_j[(cjf_gl * c_nbnxnCpuIClusterSize + j) * STRIDE_XYZ + XX])
1031 + gmx::square(x_ci[i * STRIDE_XYZ + YY]
1032 - x_j[(cjf_gl * c_nbnxnCpuIClusterSize + j) * STRIDE_XYZ + YY])
1033 + gmx::square(x_ci[i * STRIDE_XYZ + ZZ]
1034 - x_j[(cjf_gl * c_nbnxnCpuIClusterSize + j) * STRIDE_XYZ + ZZ])
1038 *numDistanceChecks += c_nbnxnCpuIClusterSize * c_nbnxnCpuIClusterSize;
1051 while (!InRange && jclusterLast > jclusterFirst)
1053 real d2 = clusterBoundingBoxDistance2(bb_ci[0], jGrid.jBoundingBoxes()[jclusterLast]);
1054 *numDistanceChecks += 2;
1056 /* Check if the distance is within the distance where
1057 * we use only the bounding box distance rbb,
1058 * or within the cut-off and there is at least one atom pair
1059 * within the cut-off.
1065 else if (d2 < rlist2)
1067 int cjl_gl = jGrid.cellOffset() + jclusterLast;
1068 for (int i = 0; i < c_nbnxnCpuIClusterSize && !InRange; i++)
1070 for (int j = 0; j < c_nbnxnCpuIClusterSize; j++)
1074 || (gmx::square(x_ci[i * STRIDE_XYZ + XX]
1075 - x_j[(cjl_gl * c_nbnxnCpuIClusterSize + j) * STRIDE_XYZ + XX])
1076 + gmx::square(x_ci[i * STRIDE_XYZ + YY]
1077 - x_j[(cjl_gl * c_nbnxnCpuIClusterSize + j) * STRIDE_XYZ + YY])
1078 + gmx::square(x_ci[i * STRIDE_XYZ + ZZ]
1079 - x_j[(cjl_gl * c_nbnxnCpuIClusterSize + j) * STRIDE_XYZ + ZZ])
1083 *numDistanceChecks += c_nbnxnCpuIClusterSize * c_nbnxnCpuIClusterSize;
1091 if (jclusterFirst <= jclusterLast)
1093 for (int jcluster = jclusterFirst; jcluster <= jclusterLast; jcluster++)
1095 /* Store cj and the interaction mask */
1097 cjEntry.cj = jGrid.cellOffset() + jcluster;
1098 cjEntry.excl = get_imask(excludeSubDiagonal, icluster, jcluster);
1099 nbl->cj.push_back(cjEntry);
1101 /* Increase the closing index in the i list */
1102 nbl->ci.back().cj_ind_end = nbl->cj.size();
1106 #ifdef GMX_NBNXN_SIMD_4XN
1107 # include "pairlist_simd_4xm.h"
1109 #ifdef GMX_NBNXN_SIMD_2XNN
1110 # include "pairlist_simd_2xmm.h"
1113 /* Plain C or SIMD4 code for making a pair list of super-cell sci vs scj.
1114 * Checks bounding box distances and possibly atom pair distances.
1116 static void make_cluster_list_supersub(const Grid& iGrid,
1118 NbnxnPairlistGpu* nbl,
1121 const bool excludeSubDiagonal,
1126 int* numDistanceChecks)
1128 NbnxnPairlistGpuWork& work = *nbl->work;
1131 const float* pbb_ci = work.iSuperClusterData.bbPacked.data();
1133 const BoundingBox* bb_ci = work.iSuperClusterData.bb.data();
1136 assert(c_nbnxnGpuClusterSize == iGrid.geometry().numAtomsICluster);
1137 assert(c_nbnxnGpuClusterSize == jGrid.geometry().numAtomsICluster);
1139 /* We generate the pairlist mainly based on bounding-box distances
1140 * and do atom pair distance based pruning on the GPU.
1141 * Only if a j-group contains a single cluster-pair, we try to prune
1142 * that pair based on atom distances on the CPU to avoid empty j-groups.
1144 #define PRUNE_LIST_CPU_ONE 1
1145 #define PRUNE_LIST_CPU_ALL 0
1147 #if PRUNE_LIST_CPU_ONE
1151 float* d2l = work.distanceBuffer.data();
1153 for (int subc = 0; subc < jGrid.numClustersPerCell()[scj]; subc++)
1155 const int cj4_ind = work.cj_ind / c_nbnxnGpuJgroupSize;
1156 const int cj_offset = work.cj_ind - cj4_ind * c_nbnxnGpuJgroupSize;
1157 const int cj = scj * c_gpuNumClusterPerCell + subc;
1159 const int cj_gl = jGrid.cellOffset() * c_gpuNumClusterPerCell + cj;
1162 if (excludeSubDiagonal && sci == scj)
1168 ci1 = iGrid.numClustersPerCell()[sci];
1172 /* Determine all ci1 bb distances in one call with SIMD4 */
1173 const int offset = packedBoundingBoxesIndex(cj) + (cj & (c_packedBoundingBoxesDimSize - 1));
1174 clusterBoundingBoxDistance2_xxxx_simd4(
1175 jGrid.packedBoundingBoxes().data() + offset, ci1, pbb_ci, d2l);
1176 *numDistanceChecks += c_nbnxnGpuClusterSize * 2;
1180 unsigned int imask = 0;
1181 /* We use a fixed upper-bound instead of ci1 to help optimization */
1182 for (int ci = 0; ci < c_gpuNumClusterPerCell; ci++)
1190 /* Determine the bb distance between ci and cj */
1191 d2l[ci] = clusterBoundingBoxDistance2(bb_ci[ci], jGrid.jBoundingBoxes()[cj]);
1192 *numDistanceChecks += 2;
1196 #if PRUNE_LIST_CPU_ALL
1197 /* Check if the distance is within the distance where
1198 * we use only the bounding box distance rbb,
1199 * or within the cut-off and there is at least one atom pair
1200 * within the cut-off. This check is very costly.
1202 *numDistanceChecks += c_nbnxnGpuClusterSize * c_nbnxnGpuClusterSize;
1203 if (d2 < rbb2 || (d2 < rlist2 && clusterpair_in_range(work, ci, cj_gl, stride, x, rlist2)))
1205 /* Check if the distance between the two bounding boxes
1206 * in within the pair-list cut-off.
1211 /* Flag this i-subcell to be taken into account */
1212 imask |= (1U << (cj_offset * c_gpuNumClusterPerCell + ci));
1214 #if PRUNE_LIST_CPU_ONE
1222 #if PRUNE_LIST_CPU_ONE
1223 /* If we only found 1 pair, check if any atoms are actually
1224 * within the cut-off, so we could get rid of it.
1226 if (npair == 1 && d2l[ci_last] >= rbb2
1227 && !clusterpair_in_range(work, ci_last, cj_gl, stride, x, rlist2))
1229 imask &= ~(1U << (cj_offset * c_gpuNumClusterPerCell + ci_last));
1236 /* We have at least one cluster pair: add a j-entry */
1237 if (static_cast<size_t>(cj4_ind) == nbl->cj4.size())
1239 nbl->cj4.resize(nbl->cj4.size() + 1);
1241 nbnxn_cj4_t* cj4 = &nbl->cj4[cj4_ind];
1243 cj4->cj[cj_offset] = cj_gl;
1245 /* Set the exclusions for the ci==sj entry.
1246 * Here we don't bother to check if this entry is actually flagged,
1247 * as it will nearly always be in the list.
1249 if (excludeSubDiagonal && sci == scj)
1251 setSelfAndNewtonExclusionsGpu(nbl, cj4_ind, cj_offset, subc);
1254 /* Copy the cluster interaction mask to the list */
1255 for (int w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
1257 cj4->imei[w].imask |= imask;
1260 nbl->work->cj_ind++;
1262 /* Keep the count */
1263 nbl->nci_tot += npair;
1265 /* Increase the closing index in i super-cell list */
1266 nbl->sci.back().cj4_ind_end =
1267 (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1) / c_nbnxnGpuJgroupSize;
1272 /* Returns how many contiguous j-clusters we have starting in the i-list */
1273 template<typename CjListType>
1274 static int numContiguousJClusters(const int cjIndexStart,
1275 const int cjIndexEnd,
1276 gmx::ArrayRef<const CjListType> cjList)
1278 const int firstJCluster = nblCj(cjList, cjIndexStart);
1280 int numContiguous = 0;
1282 while (cjIndexStart + numContiguous < cjIndexEnd
1283 && nblCj(cjList, cjIndexStart + numContiguous) == firstJCluster + numContiguous)
1288 return numContiguous;
1292 * \brief Helper struct for efficient searching for excluded atoms in a j-list
1296 /*! \brief Constructs a j-list range from \p cjList with the given index range */
1297 template<typename CjListType>
1298 JListRanges(int cjIndexStart, int cjIndexEnd, gmx::ArrayRef<const CjListType> cjList);
1300 int cjIndexStart; //!< The start index in the j-list
1301 int cjIndexEnd; //!< The end index in the j-list
1302 int cjFirst; //!< The j-cluster with index cjIndexStart
1303 int cjLast; //!< The j-cluster with index cjIndexEnd-1
1304 int numDirect; //!< Up to cjIndexStart+numDirect the j-clusters are cjFirst + the index offset
1308 template<typename CjListType>
1309 JListRanges::JListRanges(int cjIndexStart, int cjIndexEnd, gmx::ArrayRef<const CjListType> cjList) :
1310 cjIndexStart(cjIndexStart),
1311 cjIndexEnd(cjIndexEnd)
1313 GMX_ASSERT(cjIndexEnd > cjIndexStart, "JListRanges should only be called with non-empty lists");
1315 cjFirst = nblCj(cjList, cjIndexStart);
1316 cjLast = nblCj(cjList, cjIndexEnd - 1);
1318 /* Determine how many contiguous j-cells we have starting
1319 * from the first i-cell. This number can be used to directly
1320 * calculate j-cell indices for excluded atoms.
1322 numDirect = numContiguousJClusters(cjIndexStart, cjIndexEnd, cjList);
1326 /* Return the index of \p jCluster in the given range or -1 when not present
1328 * Note: This code is executed very often and therefore performance is
1329 * important. It should be inlined and fully optimized.
1331 template<typename CjListType>
1332 static inline int findJClusterInJList(int jCluster,
1333 const JListRanges& ranges,
1334 gmx::ArrayRef<const CjListType> cjList)
1338 if (jCluster < ranges.cjFirst + ranges.numDirect)
1340 /* We can calculate the index directly using the offset */
1341 index = ranges.cjIndexStart + jCluster - ranges.cjFirst;
1345 /* Search for jCluster using bisection */
1347 int rangeStart = ranges.cjIndexStart + ranges.numDirect;
1348 int rangeEnd = ranges.cjIndexEnd;
1350 while (index == -1 && rangeStart < rangeEnd)
1352 rangeMiddle = (rangeStart + rangeEnd) >> 1;
1354 const int clusterMiddle = nblCj(cjList, rangeMiddle);
1356 if (jCluster == clusterMiddle)
1358 index = rangeMiddle;
1360 else if (jCluster < clusterMiddle)
1362 rangeEnd = rangeMiddle;
1366 rangeStart = rangeMiddle + 1;
1374 // TODO: Get rid of the two functions below by renaming sci to ci (or something better)
1376 /* Return the i-entry in the list we are currently operating on */
1377 static nbnxn_ci_t* getOpenIEntry(NbnxnPairlistCpu* nbl)
1379 return &nbl->ci.back();
1382 /* Return the i-entry in the list we are currently operating on */
1383 static nbnxn_sci_t* getOpenIEntry(NbnxnPairlistGpu* nbl)
1385 return &nbl->sci.back();
1388 /* Set all atom-pair exclusions for a simple type list i-entry
1390 * Set all atom-pair exclusions from the topology stored in exclusions
1391 * as masks in the pair-list for simple list entry iEntry.
1393 static void setExclusionsForIEntry(const Nbnxm::GridSet& gridSet,
1394 NbnxnPairlistCpu* nbl,
1395 gmx_bool diagRemoved,
1397 const nbnxn_ci_t& iEntry,
1398 const ListOfLists<int>& exclusions)
1400 if (iEntry.cj_ind_end == iEntry.cj_ind_start)
1402 /* Empty list: no exclusions */
1406 const JListRanges ranges(iEntry.cj_ind_start, iEntry.cj_ind_end, gmx::makeConstArrayRef(nbl->cj));
1408 const int iCluster = iEntry.ci;
1410 gmx::ArrayRef<const int> cell = gridSet.cells();
1411 gmx::ArrayRef<const int> atomIndices = gridSet.atomIndices();
1413 /* Loop over the atoms in the i-cluster */
1414 for (int i = 0; i < nbl->na_ci; i++)
1416 const int iIndex = iCluster * nbl->na_ci + i;
1417 const int iAtom = atomIndices[iIndex];
1420 /* Loop over the topology-based exclusions for this i-atom */
1421 for (const int jAtom : exclusions[iAtom])
1425 /* The self exclusion are already set, save some time */
1429 /* Get the index of the j-atom in the nbnxn atom data */
1430 const int jIndex = cell[jAtom];
1432 /* Without shifts we only calculate interactions j>i
1433 * for one-way pair-lists.
1435 if (diagRemoved && jIndex <= iIndex)
1440 const int jCluster = (jIndex >> na_cj_2log);
1442 /* Could the cluster se be in our list? */
1443 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
1446 findJClusterInJList(jCluster, ranges, gmx::makeConstArrayRef(nbl->cj));
1450 /* We found an exclusion, clear the corresponding
1453 const int innerJ = jIndex - (jCluster << na_cj_2log);
1455 nbl->cj[index].excl &= ~(1U << ((i << na_cj_2log) + innerJ));
1463 /* Add a new i-entry to the FEP list and copy the i-properties */
1464 static inline void fep_list_new_nri_copy(t_nblist* nlist)
1466 /* Add a new i-entry */
1469 assert(nlist->nri < nlist->maxnri);
1471 /* Duplicate the last i-entry, except for jindex, which continues */
1472 nlist->iinr[nlist->nri] = nlist->iinr[nlist->nri - 1];
1473 nlist->shift[nlist->nri] = nlist->shift[nlist->nri - 1];
1474 nlist->gid[nlist->nri] = nlist->gid[nlist->nri - 1];
1475 nlist->jindex[nlist->nri] = nlist->nrj;
1478 /* Rellocate FEP list for size nl->maxnri, TODO: replace by C++ */
1479 static void reallocate_nblist(t_nblist* nl)
1484 "reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, type=%d), maxnri=%d\n",
1491 srenew(nl->iinr, nl->maxnri);
1492 srenew(nl->gid, nl->maxnri);
1493 srenew(nl->shift, nl->maxnri);
1494 srenew(nl->jindex, nl->maxnri + 1);
1497 /* For load balancing of the free-energy lists over threads, we set
1498 * the maximum nrj size of an i-entry to 40. This leads to good
1499 * load balancing in the worst case scenario of a single perturbed
1500 * particle on 16 threads, while not introducing significant overhead.
1501 * Note that half of the perturbed pairs will anyhow end up in very small lists,
1502 * since non perturbed i-particles will see few perturbed j-particles).
1504 const int max_nrj_fep = 40;
1506 /* Exclude the perturbed pairs from the Verlet list. This is only done to avoid
1507 * singularities for overlapping particles (0/0), since the charges and
1508 * LJ parameters have been zeroed in the nbnxn data structure.
1509 * Simultaneously make a group pair list for the perturbed pairs.
1511 static void make_fep_list(gmx::ArrayRef<const int> atomIndices,
1512 const nbnxn_atomdata_t* nbat,
1513 NbnxnPairlistCpu* nbl,
1514 gmx_bool bDiagRemoved,
1516 real gmx_unused shx,
1517 real gmx_unused shy,
1518 real gmx_unused shz,
1519 real gmx_unused rlist_fep2,
1524 int ci, cj_ind_start, cj_ind_end, cja, cjr;
1526 int gid_i = 0, gid_j, gid;
1527 int egp_shift, egp_mask;
1529 int ind_i, ind_j, ai, aj;
1531 gmx_bool bFEP_i, bFEP_i_all;
1533 if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
1541 cj_ind_start = nbl_ci->cj_ind_start;
1542 cj_ind_end = nbl_ci->cj_ind_end;
1544 /* In worst case we have alternating energy groups
1545 * and create #atom-pair lists, which means we need the size
1546 * of a cluster pair (na_ci*na_cj) times the number of cj's.
1548 nri_max = nbl->na_ci * nbl->na_cj * (cj_ind_end - cj_ind_start);
1549 if (nlist->nri + nri_max > nlist->maxnri)
1551 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1552 reallocate_nblist(nlist);
1555 const int numAtomsJCluster = jGrid.geometry().numAtomsJCluster;
1557 const nbnxn_atomdata_t::Params& nbatParams = nbat->params();
1559 const int ngid = nbatParams.nenergrp;
1561 /* TODO: Consider adding a check in grompp and changing this to an assert */
1562 const int numBitsInEnergyGroupIdsForAtomsInJCluster = sizeof(gid_cj) * 8;
1563 if (ngid * numAtomsJCluster > numBitsInEnergyGroupIdsForAtomsInJCluster)
1566 "The Verlet scheme with %dx%d kernels and free-energy only supports up to %zu "
1568 iGrid.geometry().numAtomsICluster,
1570 (sizeof(gid_cj) * 8) / numAtomsJCluster);
1573 egp_shift = nbatParams.neg_2log;
1574 egp_mask = (1 << egp_shift) - 1;
1576 /* Loop over the atoms in the i sub-cell */
1578 for (int i = 0; i < nbl->na_ci; i++)
1580 ind_i = ci * nbl->na_ci + i;
1581 ai = atomIndices[ind_i];
1585 nlist->jindex[nri + 1] = nlist->jindex[nri];
1586 nlist->iinr[nri] = ai;
1587 /* The actual energy group pair index is set later */
1588 nlist->gid[nri] = 0;
1589 nlist->shift[nri] = nbl_ci->shift & NBNXN_CI_SHIFT;
1591 bFEP_i = iGrid.atomIsPerturbed(ci - iGrid.cellOffset(), i);
1593 bFEP_i_all = bFEP_i_all && bFEP_i;
1595 if (nlist->nrj + (cj_ind_end - cj_ind_start) * nbl->na_cj > nlist->maxnrj)
1597 nlist->maxnrj = over_alloc_small(nlist->nrj + (cj_ind_end - cj_ind_start) * nbl->na_cj);
1598 srenew(nlist->jjnr, nlist->maxnrj);
1599 srenew(nlist->excl_fep, nlist->maxnrj);
1604 gid_i = (nbatParams.energrp[ci] >> (egp_shift * i)) & egp_mask;
1607 for (int cj_ind = cj_ind_start; cj_ind < cj_ind_end; cj_ind++)
1609 unsigned int fep_cj;
1611 cja = nbl->cj[cj_ind].cj;
1613 if (numAtomsJCluster == jGrid.geometry().numAtomsICluster)
1615 cjr = cja - jGrid.cellOffset();
1616 fep_cj = jGrid.fepBits(cjr);
1619 gid_cj = nbatParams.energrp[cja];
1622 else if (2 * numAtomsJCluster == jGrid.geometry().numAtomsICluster)
1624 cjr = cja - jGrid.cellOffset() * 2;
1625 /* Extract half of the ci fep/energrp mask */
1626 fep_cj = (jGrid.fepBits(cjr >> 1) >> ((cjr & 1) * numAtomsJCluster))
1627 & ((1 << numAtomsJCluster) - 1);
1630 gid_cj = nbatParams.energrp[cja >> 1] >> ((cja & 1) * numAtomsJCluster * egp_shift)
1631 & ((1 << (numAtomsJCluster * egp_shift)) - 1);
1636 cjr = cja - (jGrid.cellOffset() >> 1);
1637 /* Combine two ci fep masks/energrp */
1638 fep_cj = jGrid.fepBits(cjr * 2)
1639 + (jGrid.fepBits(cjr * 2 + 1) << jGrid.geometry().numAtomsICluster);
1642 gid_cj = nbatParams.energrp[cja * 2]
1643 + (nbatParams.energrp[cja * 2 + 1]
1644 << (jGrid.geometry().numAtomsICluster * egp_shift));
1648 if (bFEP_i || fep_cj != 0)
1650 for (int j = 0; j < nbl->na_cj; j++)
1652 /* Is this interaction perturbed and not excluded? */
1653 ind_j = cja * nbl->na_cj + j;
1654 aj = atomIndices[ind_j];
1655 if (aj >= 0 && (bFEP_i || (fep_cj & (1 << j))) && (!bDiagRemoved || ind_j >= ind_i))
1659 gid_j = (gid_cj >> (j * egp_shift)) & egp_mask;
1660 gid = GID(gid_i, gid_j, ngid);
1662 if (nlist->nrj > nlist->jindex[nri] && nlist->gid[nri] != gid)
1664 /* Energy group pair changed: new list */
1665 fep_list_new_nri_copy(nlist);
1668 nlist->gid[nri] = gid;
1671 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1673 fep_list_new_nri_copy(nlist);
1677 /* Add it to the FEP list */
1678 nlist->jjnr[nlist->nrj] = aj;
1679 nlist->excl_fep[nlist->nrj] =
1680 (nbl->cj[cj_ind].excl >> (i * nbl->na_cj + j)) & 1;
1683 /* Exclude it from the normal list.
1684 * Note that the charge has been set to zero,
1685 * but we need to avoid 0/0, as perturbed atoms
1686 * can be on top of each other.
1688 nbl->cj[cj_ind].excl &= ~(1U << (i * nbl->na_cj + j));
1694 if (nlist->nrj > nlist->jindex[nri])
1696 /* Actually add this new, non-empty, list */
1698 nlist->jindex[nlist->nri] = nlist->nrj;
1705 /* All interactions are perturbed, we can skip this entry */
1706 nbl_ci->cj_ind_end = cj_ind_start;
1707 nbl->ncjInUse -= cj_ind_end - cj_ind_start;
1711 /* Return the index of atom a within a cluster */
1712 static inline int cj_mod_cj4(int cj)
1714 return cj & (c_nbnxnGpuJgroupSize - 1);
1717 /* Convert a j-cluster to a cj4 group */
1718 static inline int cj_to_cj4(int cj)
1720 return cj / c_nbnxnGpuJgroupSize;
1723 /* Return the index of an j-atom within a warp */
1724 static inline int a_mod_wj(int a)
1726 return a & (c_nbnxnGpuClusterSize / c_nbnxnGpuClusterpairSplit - 1);
1729 /* As make_fep_list above, but for super/sub lists. */
1730 static void make_fep_list(gmx::ArrayRef<const int> atomIndices,
1731 const nbnxn_atomdata_t* nbat,
1732 NbnxnPairlistGpu* nbl,
1733 gmx_bool bDiagRemoved,
1734 const nbnxn_sci_t* nbl_sci,
1745 int ind_i, ind_j, ai, aj;
1749 const nbnxn_cj4_t* cj4;
1751 const int numJClusterGroups = nbl_sci->numJClusterGroups();
1752 if (numJClusterGroups == 0)
1758 const int sci = nbl_sci->sci;
1760 const int cj4_ind_start = nbl_sci->cj4_ind_start;
1761 const int cj4_ind_end = nbl_sci->cj4_ind_end;
1763 /* Here we process one super-cell, max #atoms na_sc, versus a list
1764 * cj4 entries, each with max c_nbnxnGpuJgroupSize cj's, each
1765 * of size na_cj atoms.
1766 * On the GPU we don't support energy groups (yet).
1767 * So for each of the na_sc i-atoms, we need max one FEP list
1768 * for each max_nrj_fep j-atoms.
1770 nri_max = nbl->na_sc * nbl->na_cj * (1 + (numJClusterGroups * c_nbnxnGpuJgroupSize) / max_nrj_fep);
1771 if (nlist->nri + nri_max > nlist->maxnri)
1773 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1774 reallocate_nblist(nlist);
1777 /* Loop over the atoms in the i super-cluster */
1778 for (int c = 0; c < c_gpuNumClusterPerCell; c++)
1780 c_abs = sci * c_gpuNumClusterPerCell + c;
1782 for (int i = 0; i < nbl->na_ci; i++)
1784 ind_i = c_abs * nbl->na_ci + i;
1785 ai = atomIndices[ind_i];
1789 nlist->jindex[nri + 1] = nlist->jindex[nri];
1790 nlist->iinr[nri] = ai;
1791 /* With GPUs, energy groups are not supported */
1792 nlist->gid[nri] = 0;
1793 nlist->shift[nri] = nbl_sci->shift & NBNXN_CI_SHIFT;
1795 bFEP_i = iGrid.atomIsPerturbed(c_abs - iGrid.cellOffset() * c_gpuNumClusterPerCell, i);
1797 xi = nbat->x()[ind_i * nbat->xstride + XX] + shx;
1798 yi = nbat->x()[ind_i * nbat->xstride + YY] + shy;
1799 zi = nbat->x()[ind_i * nbat->xstride + ZZ] + shz;
1801 const int nrjMax = nlist->nrj + numJClusterGroups * c_nbnxnGpuJgroupSize * nbl->na_cj;
1802 if (nrjMax > nlist->maxnrj)
1804 nlist->maxnrj = over_alloc_small(nrjMax);
1805 srenew(nlist->jjnr, nlist->maxnrj);
1806 srenew(nlist->excl_fep, nlist->maxnrj);
1809 for (int cj4_ind = cj4_ind_start; cj4_ind < cj4_ind_end; cj4_ind++)
1811 cj4 = &nbl->cj4[cj4_ind];
1813 for (int gcj = 0; gcj < c_nbnxnGpuJgroupSize; gcj++)
1815 if ((cj4->imei[0].imask & (1U << (gcj * c_gpuNumClusterPerCell + c))) == 0)
1817 /* Skip this ci for this cj */
1821 const int cjr = cj4->cj[gcj] - jGrid.cellOffset() * c_gpuNumClusterPerCell;
1823 if (bFEP_i || jGrid.clusterIsPerturbed(cjr))
1825 for (int j = 0; j < nbl->na_cj; j++)
1827 /* Is this interaction perturbed and not excluded? */
1828 ind_j = (jGrid.cellOffset() * c_gpuNumClusterPerCell + cjr) * nbl->na_cj + j;
1829 aj = atomIndices[ind_j];
1830 if (aj >= 0 && (bFEP_i || jGrid.atomIsPerturbed(cjr, j))
1831 && (!bDiagRemoved || ind_j >= ind_i))
1834 unsigned int excl_bit;
1838 j / (c_nbnxnGpuClusterSize / c_nbnxnGpuClusterpairSplit);
1839 nbnxn_excl_t& excl = get_exclusion_mask(nbl, cj4_ind, jHalf);
1841 excl_pair = a_mod_wj(j) * nbl->na_ci + i;
1842 excl_bit = (1U << (gcj * c_gpuNumClusterPerCell + c));
1844 dx = nbat->x()[ind_j * nbat->xstride + XX] - xi;
1845 dy = nbat->x()[ind_j * nbat->xstride + YY] - yi;
1846 dz = nbat->x()[ind_j * nbat->xstride + ZZ] - zi;
1848 /* The unpruned GPU list has more than 2/3
1849 * of the atom pairs beyond rlist. Using
1850 * this list will cause a lot of overhead
1851 * in the CPU FEP kernels, especially
1852 * relative to the fast GPU kernels.
1853 * So we prune the FEP list here.
1855 if (dx * dx + dy * dy + dz * dz < rlist_fep2)
1857 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1859 fep_list_new_nri_copy(nlist);
1863 /* Add it to the FEP list */
1864 nlist->jjnr[nlist->nrj] = aj;
1865 nlist->excl_fep[nlist->nrj] =
1866 (excl.pair[excl_pair] & excl_bit) ? 1 : 0;
1870 /* Exclude it from the normal list.
1871 * Note that the charge and LJ parameters have
1872 * been set to zero, but we need to avoid 0/0,
1873 * as perturbed atoms can be on top of each other.
1875 excl.pair[excl_pair] &= ~excl_bit;
1879 /* Note that we could mask out this pair in imask
1880 * if all i- and/or all j-particles are perturbed.
1881 * But since the perturbed pairs on the CPU will
1882 * take an order of magnitude more time, the GPU
1883 * will finish before the CPU and there is no gain.
1889 if (nlist->nrj > nlist->jindex[nri])
1891 /* Actually add this new, non-empty, list */
1893 nlist->jindex[nlist->nri] = nlist->nrj;
1900 /* Set all atom-pair exclusions for a GPU type list i-entry
1902 * Sets all atom-pair exclusions from the topology stored in exclusions
1903 * as masks in the pair-list for i-super-cluster list entry iEntry.
1905 static void setExclusionsForIEntry(const Nbnxm::GridSet& gridSet,
1906 NbnxnPairlistGpu* nbl,
1907 gmx_bool diagRemoved,
1908 int gmx_unused na_cj_2log,
1909 const nbnxn_sci_t& iEntry,
1910 const ListOfLists<int>& exclusions)
1912 if (iEntry.numJClusterGroups() == 0)
1918 /* Set the search ranges using start and end j-cluster indices.
1919 * Note that here we can not use cj4_ind_end, since the last cj4
1920 * can be only partially filled, so we use cj_ind.
1922 const JListRanges ranges(iEntry.cj4_ind_start * c_nbnxnGpuJgroupSize,
1924 gmx::makeConstArrayRef(nbl->cj4));
1926 GMX_ASSERT(nbl->na_ci == c_nbnxnGpuClusterSize, "na_ci should match the GPU cluster size");
1927 constexpr int c_clusterSize = c_nbnxnGpuClusterSize;
1928 constexpr int c_superClusterSize = c_nbnxnGpuNumClusterPerSupercluster * c_nbnxnGpuClusterSize;
1930 const int iSuperCluster = iEntry.sci;
1932 gmx::ArrayRef<const int> atomIndices = gridSet.atomIndices();
1933 gmx::ArrayRef<const int> cell = gridSet.cells();
1935 /* Loop over the atoms in the i super-cluster */
1936 for (int i = 0; i < c_superClusterSize; i++)
1938 const int iIndex = iSuperCluster * c_superClusterSize + i;
1939 const int iAtom = atomIndices[iIndex];
1942 const int iCluster = i / c_clusterSize;
1944 /* Loop over the topology-based exclusions for this i-atom */
1945 for (const int jAtom : exclusions[iAtom])
1949 /* The self exclusions are already set, save some time */
1953 /* Get the index of the j-atom in the nbnxn atom data */
1954 const int jIndex = cell[jAtom];
1956 /* Without shifts we only calculate interactions j>i
1957 * for one-way pair-lists.
1959 /* NOTE: We would like to use iIndex on the right hand side,
1960 * but that makes this routine 25% slower with gcc6/7.
1961 * Even using c_superClusterSize makes it slower.
1962 * Either of these changes triggers peeling of the exclIndex
1963 * loop, which apparently leads to far less efficient code.
1965 if (diagRemoved && jIndex <= iSuperCluster * nbl->na_sc + i)
1970 const int jCluster = jIndex / c_clusterSize;
1972 /* Check whether the cluster is in our list? */
1973 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
1976 findJClusterInJList(jCluster, ranges, gmx::makeConstArrayRef(nbl->cj4));
1980 /* We found an exclusion, clear the corresponding
1983 const unsigned int pairMask =
1984 (1U << (cj_mod_cj4(index) * c_gpuNumClusterPerCell + iCluster));
1985 /* Check if the i-cluster interacts with the j-cluster */
1986 if (nbl_imask0(nbl, index) & pairMask)
1988 const int innerI = (i & (c_clusterSize - 1));
1989 const int innerJ = (jIndex & (c_clusterSize - 1));
1991 /* Determine which j-half (CUDA warp) we are in */
1992 const int jHalf = innerJ / (c_clusterSize / c_nbnxnGpuClusterpairSplit);
1994 nbnxn_excl_t& interactionMask =
1995 get_exclusion_mask(nbl, cj_to_cj4(index), jHalf);
1997 interactionMask.pair[a_mod_wj(innerJ) * c_clusterSize + innerI] &= ~pairMask;
2006 /* Make a new ci entry at the back of nbl->ci */
2007 static void addNewIEntry(NbnxnPairlistCpu* nbl, int ci, int shift, int flags)
2011 ciEntry.shift = shift;
2012 /* Store the interaction flags along with the shift */
2013 ciEntry.shift |= flags;
2014 ciEntry.cj_ind_start = nbl->cj.size();
2015 ciEntry.cj_ind_end = nbl->cj.size();
2016 nbl->ci.push_back(ciEntry);
2019 /* Make a new sci entry at index nbl->nsci */
2020 static void addNewIEntry(NbnxnPairlistGpu* nbl, int sci, int shift, int gmx_unused flags)
2022 nbnxn_sci_t sciEntry;
2024 sciEntry.shift = shift;
2025 sciEntry.cj4_ind_start = nbl->cj4.size();
2026 sciEntry.cj4_ind_end = nbl->cj4.size();
2028 nbl->sci.push_back(sciEntry);
2031 /* Sort the simple j-list cj on exclusions.
2032 * Entries with exclusions will all be sorted to the beginning of the list.
2034 static void sort_cj_excl(nbnxn_cj_t* cj, int ncj, NbnxnPairlistCpuWork* work)
2036 work->cj.resize(ncj);
2038 /* Make a list of the j-cells involving exclusions */
2040 for (int j = 0; j < ncj; j++)
2042 if (cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
2044 work->cj[jnew++] = cj[j];
2047 /* Check if there are exclusions at all or not just the first entry */
2048 if (!((jnew == 0) || (jnew == 1 && cj[0].excl != NBNXN_INTERACTION_MASK_ALL)))
2050 for (int j = 0; j < ncj; j++)
2052 if (cj[j].excl == NBNXN_INTERACTION_MASK_ALL)
2054 work->cj[jnew++] = cj[j];
2057 for (int j = 0; j < ncj; j++)
2059 cj[j] = work->cj[j];
2064 /* Close this simple list i entry */
2065 static void closeIEntry(NbnxnPairlistCpu* nbl,
2066 int gmx_unused sp_max_av,
2067 gmx_bool gmx_unused progBal,
2068 float gmx_unused nsp_tot_est,
2069 int gmx_unused thread,
2070 int gmx_unused nthread)
2072 nbnxn_ci_t& ciEntry = nbl->ci.back();
2074 /* All content of the new ci entry have already been filled correctly,
2075 * we only need to sort and increase counts or remove the entry when empty.
2077 const int jlen = ciEntry.cj_ind_end - ciEntry.cj_ind_start;
2080 sort_cj_excl(nbl->cj.data() + ciEntry.cj_ind_start, jlen, nbl->work.get());
2082 /* The counts below are used for non-bonded pair/flop counts
2083 * and should therefore match the available kernel setups.
2085 if (!(ciEntry.shift & NBNXN_CI_DO_COUL(0)))
2087 nbl->work->ncj_noq += jlen;
2089 else if ((ciEntry.shift & NBNXN_CI_HALF_LJ(0)) || !(ciEntry.shift & NBNXN_CI_DO_LJ(0)))
2091 nbl->work->ncj_hlj += jlen;
2096 /* Entry is empty: remove it */
2101 /* Split sci entry for load balancing on the GPU.
2102 * Splitting ensures we have enough lists to fully utilize the whole GPU.
2103 * With progBal we generate progressively smaller lists, which improves
2104 * load balancing. As we only know the current count on our own thread,
2105 * we will need to estimate the current total amount of i-entries.
2106 * As the lists get concatenated later, this estimate depends
2107 * both on nthread and our own thread index.
2109 static void split_sci_entry(NbnxnPairlistGpu* nbl,
2122 /* Estimate the total numbers of ci's of the nblist combined
2123 * over all threads using the target number of ci's.
2125 nsp_est = (nsp_tot_est * thread) / nthread + nbl->nci_tot;
2127 /* The first ci blocks should be larger, to avoid overhead.
2128 * The last ci blocks should be smaller, to improve load balancing.
2129 * The factor 3/2 makes the first block 3/2 times the target average
2130 * and ensures that the total number of blocks end up equal to
2131 * that of equally sized blocks of size nsp_target_av.
2133 nsp_max = static_cast<int>(nsp_target_av * (nsp_tot_est * 1.5 / (nsp_est + nsp_tot_est)));
2137 nsp_max = nsp_target_av;
2140 const int cj4_start = nbl->sci.back().cj4_ind_start;
2141 const int cj4_end = nbl->sci.back().cj4_ind_end;
2142 const int j4len = cj4_end - cj4_start;
2144 if (j4len > 1 && j4len * c_gpuNumClusterPerCell * c_nbnxnGpuJgroupSize > nsp_max)
2146 /* Modify the last ci entry and process the cj4's again */
2152 for (int cj4 = cj4_start; cj4 < cj4_end; cj4++)
2154 int nsp_cj4_p = nsp_cj4;
2155 /* Count the number of cluster pairs in this cj4 group */
2157 for (int p = 0; p < c_gpuNumClusterPerCell * c_nbnxnGpuJgroupSize; p++)
2159 nsp_cj4 += (nbl->cj4[cj4].imei[0].imask >> p) & 1;
2162 /* If adding the current cj4 with nsp_cj4 pairs get us further
2163 * away from our target nsp_max, split the list before this cj4.
2165 if (nsp > 0 && nsp_max - nsp < nsp + nsp_cj4 - nsp_max)
2167 /* Split the list at cj4 */
2168 nbl->sci.back().cj4_ind_end = cj4;
2169 /* Create a new sci entry */
2171 sciNew.sci = nbl->sci.back().sci;
2172 sciNew.shift = nbl->sci.back().shift;
2173 sciNew.cj4_ind_start = cj4;
2174 nbl->sci.push_back(sciNew);
2177 nsp_cj4_e = nsp_cj4_p;
2183 /* Put the remaining cj4's in the last sci entry */
2184 nbl->sci.back().cj4_ind_end = cj4_end;
2186 /* Possibly balance out the last two sci's
2187 * by moving the last cj4 of the second last sci.
2189 if (nsp_sci - nsp_cj4_e >= nsp + nsp_cj4_e)
2191 GMX_ASSERT(nbl->sci.size() >= 2, "We expect at least two elements");
2192 nbl->sci[nbl->sci.size() - 2].cj4_ind_end--;
2193 nbl->sci[nbl->sci.size() - 1].cj4_ind_start--;
2198 /* Clost this super/sub list i entry */
2199 static void closeIEntry(NbnxnPairlistGpu* nbl, int nsp_max_av, gmx_bool progBal, float nsp_tot_est, int thread, int nthread)
2201 nbnxn_sci_t& sciEntry = *getOpenIEntry(nbl);
2203 /* All content of the new ci entry have already been filled correctly,
2204 * we only need to, potentially, split or remove the entry when empty.
2206 int j4len = sciEntry.numJClusterGroups();
2209 /* We can only have complete blocks of 4 j-entries in a list,
2210 * so round the count up before closing.
2212 int ncj4 = (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1) / c_nbnxnGpuJgroupSize;
2213 nbl->work->cj_ind = ncj4 * c_nbnxnGpuJgroupSize;
2217 /* Measure the size of the new entry and potentially split it */
2218 split_sci_entry(nbl, nsp_max_av, progBal, nsp_tot_est, thread, nthread);
2223 /* Entry is empty: remove it */
2224 nbl->sci.pop_back();
2228 /* Syncs the working array before adding another grid pair to the GPU list */
2229 static void sync_work(NbnxnPairlistCpu gmx_unused* nbl) {}
2231 /* Syncs the working array before adding another grid pair to the GPU list */
2232 static void sync_work(NbnxnPairlistGpu* nbl)
2234 nbl->work->cj_ind = nbl->cj4.size() * c_nbnxnGpuJgroupSize;
2237 /* Clears an NbnxnPairlistCpu data structure */
2238 static void clear_pairlist(NbnxnPairlistCpu* nbl)
2244 nbl->ciOuter.clear();
2245 nbl->cjOuter.clear();
2247 nbl->work->ncj_noq = 0;
2248 nbl->work->ncj_hlj = 0;
2251 /* Clears an NbnxnPairlistGpu data structure */
2252 static void clear_pairlist(NbnxnPairlistGpu* nbl)
2256 nbl->excl.resize(1);
2260 /* Clears an atom-atom-style pair list */
2261 static void clear_pairlist_fep(t_nblist* nl)
2265 if (nl->jindex == nullptr)
2267 snew(nl->jindex, 1);
2272 /* Sets a simple list i-cell bounding box, including PBC shift */
2274 set_icell_bb_simple(gmx::ArrayRef<const BoundingBox> bb, int ci, real shx, real shy, real shz, BoundingBox* bb_ci)
2276 bb_ci->lower.x = bb[ci].lower.x + shx;
2277 bb_ci->lower.y = bb[ci].lower.y + shy;
2278 bb_ci->lower.z = bb[ci].lower.z + shz;
2279 bb_ci->upper.x = bb[ci].upper.x + shx;
2280 bb_ci->upper.y = bb[ci].upper.y + shy;
2281 bb_ci->upper.z = bb[ci].upper.z + shz;
2284 /* Sets a simple list i-cell bounding box, including PBC shift */
2285 static inline void set_icell_bb(const Grid& iGrid, int ci, real shx, real shy, real shz, NbnxnPairlistCpuWork* work)
2287 set_icell_bb_simple(iGrid.iBoundingBoxes(), ci, shx, shy, shz, &work->iClusterData.bb[0]);
2291 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2292 static void set_icell_bbxxxx_supersub(gmx::ArrayRef<const float> bb, int ci, real shx, real shy, real shz, float* bb_ci)
2294 constexpr int cellBBStride = packedBoundingBoxesIndex(c_gpuNumClusterPerCell);
2295 constexpr int pbbStride = c_packedBoundingBoxesDimSize;
2296 const int ia = ci * cellBBStride;
2297 for (int m = 0; m < cellBBStride; m += c_packedBoundingBoxesSize)
2299 for (int i = 0; i < pbbStride; i++)
2301 bb_ci[m + 0 * pbbStride + i] = bb[ia + m + 0 * pbbStride + i] + shx;
2302 bb_ci[m + 1 * pbbStride + i] = bb[ia + m + 1 * pbbStride + i] + shy;
2303 bb_ci[m + 2 * pbbStride + i] = bb[ia + m + 2 * pbbStride + i] + shz;
2304 bb_ci[m + 3 * pbbStride + i] = bb[ia + m + 3 * pbbStride + i] + shx;
2305 bb_ci[m + 4 * pbbStride + i] = bb[ia + m + 4 * pbbStride + i] + shy;
2306 bb_ci[m + 5 * pbbStride + i] = bb[ia + m + 5 * pbbStride + i] + shz;
2312 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2313 gmx_unused static void set_icell_bb_supersub(gmx::ArrayRef<const BoundingBox> bb,
2320 for (int i = 0; i < c_gpuNumClusterPerCell; i++)
2322 set_icell_bb_simple(bb, ci * c_gpuNumClusterPerCell + i, shx, shy, shz, &bb_ci[i]);
2326 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2327 gmx_unused static void set_icell_bb(const Grid& iGrid, int ci, real shx, real shy, real shz, NbnxnPairlistGpuWork* work)
2330 set_icell_bbxxxx_supersub(
2331 iGrid.packedBoundingBoxes(), ci, shx, shy, shz, work->iSuperClusterData.bbPacked.data());
2333 set_icell_bb_supersub(iGrid.iBoundingBoxes(), ci, shx, shy, shz, work->iSuperClusterData.bb.data());
2337 /* Copies PBC shifted i-cell atom coordinates x,y,z to working array */
2338 static void icell_set_x_simple(int ci,
2344 NbnxnPairlistCpuWork::IClusterData* iClusterData)
2346 const int ia = ci * c_nbnxnCpuIClusterSize;
2348 for (int i = 0; i < c_nbnxnCpuIClusterSize; i++)
2350 iClusterData->x[i * STRIDE_XYZ + XX] = x[(ia + i) * stride + XX] + shx;
2351 iClusterData->x[i * STRIDE_XYZ + YY] = x[(ia + i) * stride + YY] + shy;
2352 iClusterData->x[i * STRIDE_XYZ + ZZ] = x[(ia + i) * stride + ZZ] + shz;
2356 static void icell_set_x(int ci,
2362 const ClusterDistanceKernelType kernelType,
2363 NbnxnPairlistCpuWork* work)
2368 # ifdef GMX_NBNXN_SIMD_4XN
2369 case ClusterDistanceKernelType::CpuSimd_4xM:
2370 icell_set_x_simd_4xn(ci, shx, shy, shz, stride, x, work);
2373 # ifdef GMX_NBNXN_SIMD_2XNN
2374 case ClusterDistanceKernelType::CpuSimd_2xMM:
2375 icell_set_x_simd_2xnn(ci, shx, shy, shz, stride, x, work);
2379 case ClusterDistanceKernelType::CpuPlainC:
2380 icell_set_x_simple(ci, shx, shy, shz, stride, x, &work->iClusterData);
2382 default: GMX_ASSERT(false, "Unhandled case"); break;
2386 /* Copies PBC shifted super-cell atom coordinates x,y,z to working array */
2387 static void icell_set_x(int ci,
2393 ClusterDistanceKernelType gmx_unused kernelType,
2394 NbnxnPairlistGpuWork* work)
2396 #if !GMX_SIMD4_HAVE_REAL
2398 real* x_ci = work->iSuperClusterData.x.data();
2400 int ia = ci * c_gpuNumClusterPerCell * c_nbnxnGpuClusterSize;
2401 for (int i = 0; i < c_gpuNumClusterPerCell * c_nbnxnGpuClusterSize; i++)
2403 x_ci[i * DIM + XX] = x[(ia + i) * stride + XX] + shx;
2404 x_ci[i * DIM + YY] = x[(ia + i) * stride + YY] + shy;
2405 x_ci[i * DIM + ZZ] = x[(ia + i) * stride + ZZ] + shz;
2408 #else /* !GMX_SIMD4_HAVE_REAL */
2410 real* x_ci = work->iSuperClusterData.xSimd.data();
2412 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2414 for (int i = 0; i < c_nbnxnGpuClusterSize; i += GMX_SIMD4_WIDTH)
2416 int io = si * c_nbnxnGpuClusterSize + i;
2417 int ia = ci * c_gpuNumClusterPerCell * c_nbnxnGpuClusterSize + io;
2418 for (int j = 0; j < GMX_SIMD4_WIDTH; j++)
2420 x_ci[io * DIM + j + XX * GMX_SIMD4_WIDTH] = x[(ia + j) * stride + XX] + shx;
2421 x_ci[io * DIM + j + YY * GMX_SIMD4_WIDTH] = x[(ia + j) * stride + YY] + shy;
2422 x_ci[io * DIM + j + ZZ * GMX_SIMD4_WIDTH] = x[(ia + j) * stride + ZZ] + shz;
2427 #endif /* !GMX_SIMD4_HAVE_REAL */
2430 static real minimum_subgrid_size_xy(const Grid& grid)
2432 const Grid::Dimensions& dims = grid.dimensions();
2434 if (grid.geometry().isSimple)
2436 return std::min(dims.cellSize[XX], dims.cellSize[YY]);
2440 return std::min(dims.cellSize[XX] / c_gpuNumClusterPerCellX,
2441 dims.cellSize[YY] / c_gpuNumClusterPerCellY);
2445 static real effective_buffer_1x1_vs_MxN(const Grid& iGrid, const Grid& jGrid)
2447 const real eff_1x1_buffer_fac_overest = 0.1;
2449 /* Determine an atom-pair list cut-off buffer size for atom pairs,
2450 * to be added to rlist (including buffer) used for MxN.
2451 * This is for converting an MxN list to a 1x1 list. This means we can't
2452 * use the normal buffer estimate, as we have an MxN list in which
2453 * some atom pairs beyond rlist are missing. We want to capture
2454 * the beneficial effect of buffering by extra pairs just outside rlist,
2455 * while removing the useless pairs that are further away from rlist.
2456 * (Also the buffer could have been set manually not using the estimate.)
2457 * This buffer size is an overestimate.
2458 * We add 10% of the smallest grid sub-cell dimensions.
2459 * Note that the z-size differs per cell and we don't use this,
2460 * so we overestimate.
2461 * With PME, the 10% value gives a buffer that is somewhat larger
2462 * than the effective buffer with a tolerance of 0.005 kJ/mol/ps.
2463 * Smaller tolerances or using RF lead to a smaller effective buffer,
2464 * so 10% gives a safe overestimate.
2466 return eff_1x1_buffer_fac_overest * (minimum_subgrid_size_xy(iGrid) + minimum_subgrid_size_xy(jGrid));
2469 /* Estimates the interaction volume^2 for non-local interactions */
2470 static real nonlocal_vol2(const struct gmx_domdec_zones_t* zones, const rvec ls, real r)
2478 /* Here we simply add up the volumes of 1, 2 or 3 1D decomposition
2479 * not home interaction volume^2. As these volumes are not additive,
2480 * this is an overestimate, but it would only be significant in the limit
2481 * of small cells, where we anyhow need to split the lists into
2482 * as small parts as possible.
2485 for (int z = 0; z < zones->n; z++)
2487 if (zones->shift[z][XX] + zones->shift[z][YY] + zones->shift[z][ZZ] == 1)
2492 for (int d = 0; d < DIM; d++)
2494 if (zones->shift[z][d] == 0)
2498 za *= zones->size[z].x1[d] - zones->size[z].x0[d];
2502 /* 4 octants of a sphere */
2503 vold_est = 0.25 * M_PI * r * r * r * r;
2504 /* 4 quarter pie slices on the edges */
2505 vold_est += 4 * cl * M_PI / 6.0 * r * r * r;
2506 /* One rectangular volume on a face */
2507 vold_est += ca * 0.5 * r * r;
2509 vol2_est_tot += vold_est * za;
2513 return vol2_est_tot;
2516 /* Estimates the average size of a full j-list for super/sub setup */
2517 static void get_nsubpair_target(const Nbnxm::GridSet& gridSet,
2518 const InteractionLocality iloc,
2520 const int min_ci_balanced,
2521 int* nsubpair_target,
2522 float* nsubpair_tot_est)
2524 /* The target value of 36 seems to be the optimum for Kepler.
2525 * Maxwell is less sensitive to the exact value.
2527 const int nsubpair_target_min = 36;
2528 real r_eff_sup, vol_est, nsp_est, nsp_est_nl;
2530 const Grid& grid = gridSet.grids()[0];
2532 /* We don't need to balance list sizes if:
2533 * - We didn't request balancing.
2534 * - The number of grid cells >= the number of lists requested,
2535 * since we will always generate at least #cells lists.
2536 * - We don't have any cells, since then there won't be any lists.
2538 if (min_ci_balanced <= 0 || grid.numCells() >= min_ci_balanced || grid.numCells() == 0)
2540 /* nsubpair_target==0 signals no balancing */
2541 *nsubpair_target = 0;
2542 *nsubpair_tot_est = 0;
2548 const int numAtomsCluster = grid.geometry().numAtomsICluster;
2549 const Grid::Dimensions& dims = grid.dimensions();
2551 ls[XX] = dims.cellSize[XX] / c_gpuNumClusterPerCellX;
2552 ls[YY] = dims.cellSize[YY] / c_gpuNumClusterPerCellY;
2553 ls[ZZ] = numAtomsCluster / (dims.atomDensity * ls[XX] * ls[YY]);
2555 /* The formulas below are a heuristic estimate of the average nsj per si*/
2556 r_eff_sup = rlist + nbnxn_get_rlist_effective_inc(numAtomsCluster, ls);
2558 if (!gridSet.domainSetup().haveMultipleDomains || gridSet.domainSetup().zones->n == 1)
2564 nsp_est_nl = gmx::square(dims.atomDensity / numAtomsCluster)
2565 * nonlocal_vol2(gridSet.domainSetup().zones, ls, r_eff_sup);
2568 if (iloc == InteractionLocality::Local)
2570 /* Sub-cell interacts with itself */
2571 vol_est = ls[XX] * ls[YY] * ls[ZZ];
2572 /* 6/2 rectangular volume on the faces */
2573 vol_est += (ls[XX] * ls[YY] + ls[XX] * ls[ZZ] + ls[YY] * ls[ZZ]) * r_eff_sup;
2574 /* 12/2 quarter pie slices on the edges */
2575 vol_est += 2 * (ls[XX] + ls[YY] + ls[ZZ]) * 0.25 * M_PI * gmx::square(r_eff_sup);
2576 /* 4 octants of a sphere */
2577 vol_est += 0.5 * 4.0 / 3.0 * M_PI * gmx::power3(r_eff_sup);
2579 /* Estimate the number of cluster pairs as the local number of
2580 * clusters times the volume they interact with times the density.
2582 nsp_est = grid.numClusters() * vol_est * dims.atomDensity / numAtomsCluster;
2584 /* Subtract the non-local pair count */
2585 nsp_est -= nsp_est_nl;
2587 /* For small cut-offs nsp_est will be an underesimate.
2588 * With DD nsp_est_nl is an overestimate so nsp_est can get negative.
2589 * So to avoid too small or negative nsp_est we set a minimum of
2590 * all cells interacting with all 3^3 direct neighbors (3^3-1)/2+1=14.
2591 * This might be a slight overestimate for small non-periodic groups of
2592 * atoms as will occur for a local domain with DD, but for small
2593 * groups of atoms we'll anyhow be limited by nsubpair_target_min,
2594 * so this overestimation will not matter.
2596 nsp_est = std::max(nsp_est, grid.numClusters() * 14._real);
2600 fprintf(debug, "nsp_est local %5.1f non-local %5.1f\n", nsp_est, nsp_est_nl);
2605 nsp_est = nsp_est_nl;
2608 /* Thus the (average) maximum j-list size should be as follows.
2609 * Since there is overhead, we shouldn't make the lists too small
2610 * (and we can't chop up j-groups) so we use a minimum target size of 36.
2612 *nsubpair_target = std::max(nsubpair_target_min, roundToInt(nsp_est / min_ci_balanced));
2613 *nsubpair_tot_est = nsp_est;
2617 fprintf(debug, "nbl nsp estimate %.1f, nsubpair_target %d\n", nsp_est, *nsubpair_target);
2621 /* Debug list print function */
2622 static void print_nblist_ci_cj(FILE* fp, const NbnxnPairlistCpu& nbl)
2624 for (const nbnxn_ci_t& ciEntry : nbl.ci)
2626 fprintf(fp, "ci %4d shift %2d ncj %3d\n", ciEntry.ci, ciEntry.shift, ciEntry.cj_ind_end - ciEntry.cj_ind_start);
2628 for (int j = ciEntry.cj_ind_start; j < ciEntry.cj_ind_end; j++)
2630 fprintf(fp, " cj %5d imask %x\n", nbl.cj[j].cj, nbl.cj[j].excl);
2635 /* Debug list print function */
2636 static void print_nblist_sci_cj(FILE* fp, const NbnxnPairlistGpu& nbl)
2638 for (const nbnxn_sci_t& sci : nbl.sci)
2640 fprintf(fp, "ci %4d shift %2d ncj4 %2d\n", sci.sci, sci.shift, sci.numJClusterGroups());
2643 for (int j4 = sci.cj4_ind_start; j4 < sci.cj4_ind_end; j4++)
2645 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
2647 fprintf(fp, " sj %5d imask %x\n", nbl.cj4[j4].cj[j], nbl.cj4[j4].imei[0].imask);
2648 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2650 if (nbl.cj4[j4].imei[0].imask & (1U << (j * c_gpuNumClusterPerCell + si)))
2657 fprintf(fp, "ci %4d shift %2d ncj4 %2d ncp %3d\n", sci.sci, sci.shift, sci.numJClusterGroups(), ncp);
2661 /* Combine pair lists *nbl generated on multiple threads nblc */
2662 static void combine_nblists(gmx::ArrayRef<const NbnxnPairlistGpu> nbls, NbnxnPairlistGpu* nblc)
2664 int nsci = nblc->sci.size();
2665 int ncj4 = nblc->cj4.size();
2666 int nexcl = nblc->excl.size();
2667 for (auto& nbl : nbls)
2669 nsci += nbl.sci.size();
2670 ncj4 += nbl.cj4.size();
2671 nexcl += nbl.excl.size();
2674 /* Resize with the final, combined size, so we can fill in parallel */
2675 /* NOTE: For better performance we should use default initialization */
2676 nblc->sci.resize(nsci);
2677 nblc->cj4.resize(ncj4);
2678 nblc->excl.resize(nexcl);
2680 /* Each thread should copy its own data to the combined arrays,
2681 * as otherwise data will go back and forth between different caches.
2683 const int gmx_unused nthreads = gmx_omp_nthreads_get(emntPairsearch);
2685 #pragma omp parallel for num_threads(nthreads) schedule(static)
2686 for (gmx::index n = 0; n < nbls.ssize(); n++)
2690 /* Determine the offset in the combined data for our thread.
2691 * Note that the original sizes in nblc are lost.
2693 int sci_offset = nsci;
2694 int cj4_offset = ncj4;
2695 int excl_offset = nexcl;
2697 for (gmx::index i = n; i < nbls.ssize(); i++)
2699 sci_offset -= nbls[i].sci.size();
2700 cj4_offset -= nbls[i].cj4.size();
2701 excl_offset -= nbls[i].excl.size();
2704 const NbnxnPairlistGpu& nbli = nbls[n];
2706 for (size_t i = 0; i < nbli.sci.size(); i++)
2708 nblc->sci[sci_offset + i] = nbli.sci[i];
2709 nblc->sci[sci_offset + i].cj4_ind_start += cj4_offset;
2710 nblc->sci[sci_offset + i].cj4_ind_end += cj4_offset;
2713 for (size_t j4 = 0; j4 < nbli.cj4.size(); j4++)
2715 nblc->cj4[cj4_offset + j4] = nbli.cj4[j4];
2716 nblc->cj4[cj4_offset + j4].imei[0].excl_ind += excl_offset;
2717 nblc->cj4[cj4_offset + j4].imei[1].excl_ind += excl_offset;
2720 for (size_t j4 = 0; j4 < nbli.excl.size(); j4++)
2722 nblc->excl[excl_offset + j4] = nbli.excl[j4];
2725 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2728 for (auto& nbl : nbls)
2730 nblc->nci_tot += nbl.nci_tot;
2734 static void balance_fep_lists(gmx::ArrayRef<std::unique_ptr<t_nblist>> fepLists,
2735 gmx::ArrayRef<PairsearchWork> work)
2737 const int numLists = fepLists.ssize();
2741 /* Nothing to balance */
2745 /* Count the total i-lists and pairs */
2748 for (const auto& list : fepLists)
2750 nri_tot += list->nri;
2751 nrj_tot += list->nrj;
2754 const int nrj_target = (nrj_tot + numLists - 1) / numLists;
2756 GMX_ASSERT(gmx_omp_nthreads_get(emntNonbonded) == numLists,
2757 "We should have as many work objects as FEP lists");
2759 #pragma omp parallel for schedule(static) num_threads(numLists)
2760 for (int th = 0; th < numLists; th++)
2764 t_nblist* nbl = work[th].nbl_fep.get();
2766 /* Note that here we allocate for the total size, instead of
2767 * a per-thread esimate (which is hard to obtain).
2769 if (nri_tot > nbl->maxnri)
2771 nbl->maxnri = over_alloc_large(nri_tot);
2772 reallocate_nblist(nbl);
2774 if (nri_tot > nbl->maxnri || nrj_tot > nbl->maxnrj)
2776 nbl->maxnrj = over_alloc_small(nrj_tot);
2777 srenew(nbl->jjnr, nbl->maxnrj);
2778 srenew(nbl->excl_fep, nbl->maxnrj);
2781 clear_pairlist_fep(nbl);
2783 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
2786 /* Loop over the source lists and assign and copy i-entries */
2788 t_nblist* nbld = work[th_dest].nbl_fep.get();
2789 for (int th = 0; th < numLists; th++)
2791 const t_nblist* nbls = fepLists[th].get();
2793 for (int i = 0; i < nbls->nri; i++)
2797 /* The number of pairs in this i-entry */
2798 nrj = nbls->jindex[i + 1] - nbls->jindex[i];
2800 /* Decide if list th_dest is too large and we should procede
2801 * to the next destination list.
2803 if (th_dest + 1 < numLists && nbld->nrj > 0
2804 && nbld->nrj + nrj - nrj_target > nrj_target - nbld->nrj)
2807 nbld = work[th_dest].nbl_fep.get();
2810 nbld->iinr[nbld->nri] = nbls->iinr[i];
2811 nbld->gid[nbld->nri] = nbls->gid[i];
2812 nbld->shift[nbld->nri] = nbls->shift[i];
2814 for (int j = nbls->jindex[i]; j < nbls->jindex[i + 1]; j++)
2816 nbld->jjnr[nbld->nrj] = nbls->jjnr[j];
2817 nbld->excl_fep[nbld->nrj] = nbls->excl_fep[j];
2821 nbld->jindex[nbld->nri] = nbld->nrj;
2825 /* Swap the list pointers */
2826 for (int th = 0; th < numLists; th++)
2828 fepLists[th].swap(work[th].nbl_fep);
2832 fprintf(debug, "nbl_fep[%d] nri %4d nrj %4d\n", th, fepLists[th]->nri, fepLists[th]->nrj);
2837 /* Returns the next ci to be processes by our thread */
2838 static gmx_bool next_ci(const Grid& grid, int nth, int ci_block, int* ci_x, int* ci_y, int* ci_b, int* ci)
2843 if (*ci_b == ci_block)
2845 /* Jump to the next block assigned to this task */
2846 *ci += (nth - 1) * ci_block;
2850 if (*ci >= grid.numCells())
2855 while (*ci >= grid.firstCellInColumn(*ci_x * grid.dimensions().numCells[YY] + *ci_y + 1))
2858 if (*ci_y == grid.dimensions().numCells[YY])
2868 /* Returns the distance^2 for which we put cell pairs in the list
2869 * without checking atom pair distances. This is usually < rlist^2.
2871 static float boundingbox_only_distance2(const Grid::Dimensions& iGridDims,
2872 const Grid::Dimensions& jGridDims,
2876 /* If the distance between two sub-cell bounding boxes is less
2877 * than this distance, do not check the distance between
2878 * all particle pairs in the sub-cell, since then it is likely
2879 * that the box pair has atom pairs within the cut-off.
2880 * We use the nblist cut-off minus 0.5 times the average x/y diagonal
2881 * spacing of the sub-cells. Around 40% of the checked pairs are pruned.
2882 * Using more than 0.5 gains at most 0.5%.
2883 * If forces are calculated more than twice, the performance gain
2884 * in the force calculation outweighs the cost of checking.
2885 * Note that with subcell lists, the atom-pair distance check
2886 * is only performed when only 1 out of 8 sub-cells in within range,
2887 * this is because the GPU is much faster than the cpu.
2892 bbx = 0.5 * (iGridDims.cellSize[XX] + jGridDims.cellSize[XX]);
2893 bby = 0.5 * (iGridDims.cellSize[YY] + jGridDims.cellSize[YY]);
2896 bbx /= c_gpuNumClusterPerCellX;
2897 bby /= c_gpuNumClusterPerCellY;
2900 rbb2 = std::max(0.0, rlist - 0.5 * std::sqrt(bbx * bbx + bby * bby));
2906 return static_cast<float>((1 + GMX_FLOAT_EPS) * rbb2);
2910 static int get_ci_block_size(const Grid& iGrid, const bool haveMultipleDomains, const int numLists)
2912 const int ci_block_enum = 5;
2913 const int ci_block_denom = 11;
2914 const int ci_block_min_atoms = 16;
2917 /* Here we decide how to distribute the blocks over the threads.
2918 * We use prime numbers to try to avoid that the grid size becomes
2919 * a multiple of the number of threads, which would lead to some
2920 * threads getting "inner" pairs and others getting boundary pairs,
2921 * which in turns will lead to load imbalance between threads.
2922 * Set the block size as 5/11/ntask times the average number of cells
2923 * in a y,z slab. This should ensure a quite uniform distribution
2924 * of the grid parts of the different thread along all three grid
2925 * zone boundaries with 3D domain decomposition. At the same time
2926 * the blocks will not become too small.
2928 GMX_ASSERT(iGrid.dimensions().numCells[XX] > 0, "Grid can't be empty");
2929 GMX_ASSERT(numLists > 0, "We need at least one list");
2930 ci_block = (iGrid.numCells() * ci_block_enum)
2931 / (ci_block_denom * iGrid.dimensions().numCells[XX] * numLists);
2933 const int numAtomsPerCell = iGrid.geometry().numAtomsPerCell;
2935 /* Ensure the blocks are not too small: avoids cache invalidation */
2936 if (ci_block * numAtomsPerCell < ci_block_min_atoms)
2938 ci_block = (ci_block_min_atoms + numAtomsPerCell - 1) / numAtomsPerCell;
2941 /* Without domain decomposition
2942 * or with less than 3 blocks per task, divide in nth blocks.
2944 if (!haveMultipleDomains || numLists * 3 * ci_block > iGrid.numCells())
2946 ci_block = (iGrid.numCells() + numLists - 1) / numLists;
2949 if (ci_block > 1 && (numLists - 1) * ci_block >= iGrid.numCells())
2951 /* Some threads have no work. Although reducing the block size
2952 * does not decrease the block count on the first few threads,
2953 * with GPUs better mixing of "upper" cells that have more empty
2954 * clusters results in a somewhat lower max load over all threads.
2955 * Without GPUs the regime of so few atoms per thread is less
2956 * performance relevant, but with 8-wide SIMD the same reasoning
2957 * applies, since the pair list uses 4 i-atom "sub-clusters".
2965 /* Returns the number of bits to right-shift a cluster index to obtain
2966 * the corresponding force buffer flag index.
2968 static int getBufferFlagShift(int numAtomsPerCluster)
2970 int bufferFlagShift = 0;
2971 while ((numAtomsPerCluster << bufferFlagShift) < NBNXN_BUFFERFLAG_SIZE)
2976 return bufferFlagShift;
2979 static bool pairlistIsSimple(const NbnxnPairlistCpu gmx_unused& pairlist)
2984 static bool pairlistIsSimple(const NbnxnPairlistGpu gmx_unused& pairlist)
2989 static void makeClusterListWrapper(NbnxnPairlistCpu* nbl,
2990 const Grid gmx_unused& iGrid,
2993 const int firstCell,
2995 const bool excludeSubDiagonal,
2996 const nbnxn_atomdata_t* nbat,
2999 const ClusterDistanceKernelType kernelType,
3000 int* numDistanceChecks)
3004 case ClusterDistanceKernelType::CpuPlainC:
3005 makeClusterListSimple(
3006 jGrid, nbl, ci, firstCell, lastCell, excludeSubDiagonal, nbat->x().data(), rlist2, rbb2, numDistanceChecks);
3008 #ifdef GMX_NBNXN_SIMD_4XN
3009 case ClusterDistanceKernelType::CpuSimd_4xM:
3010 makeClusterListSimd4xn(
3011 jGrid, nbl, ci, firstCell, lastCell, excludeSubDiagonal, nbat->x().data(), rlist2, rbb2, numDistanceChecks);
3014 #ifdef GMX_NBNXN_SIMD_2XNN
3015 case ClusterDistanceKernelType::CpuSimd_2xMM:
3016 makeClusterListSimd2xnn(
3017 jGrid, nbl, ci, firstCell, lastCell, excludeSubDiagonal, nbat->x().data(), rlist2, rbb2, numDistanceChecks);
3020 default: GMX_ASSERT(false, "Unhandled kernel type");
3024 static void makeClusterListWrapper(NbnxnPairlistGpu* nbl,
3025 const Grid& gmx_unused iGrid,
3028 const int firstCell,
3030 const bool excludeSubDiagonal,
3031 const nbnxn_atomdata_t* nbat,
3034 ClusterDistanceKernelType gmx_unused kernelType,
3035 int* numDistanceChecks)
3037 for (int cj = firstCell; cj <= lastCell; cj++)
3039 make_cluster_list_supersub(
3040 iGrid, jGrid, nbl, ci, cj, excludeSubDiagonal, nbat->xstride, nbat->x().data(), rlist2, rbb2, numDistanceChecks);
3044 static int getNumSimpleJClustersInList(const NbnxnPairlistCpu& nbl)
3046 return nbl.cj.size();
3049 static int getNumSimpleJClustersInList(const gmx_unused NbnxnPairlistGpu& nbl)
3054 static void incrementNumSimpleJClustersInList(NbnxnPairlistCpu* nbl, int ncj_old_j)
3056 nbl->ncjInUse += nbl->cj.size();
3057 nbl->ncjInUse -= ncj_old_j;
3060 static void incrementNumSimpleJClustersInList(NbnxnPairlistGpu gmx_unused* nbl, int gmx_unused ncj_old_j)
3064 static void checkListSizeConsistency(const NbnxnPairlistCpu& nbl, const bool haveFreeEnergy)
3066 GMX_RELEASE_ASSERT(static_cast<size_t>(nbl.ncjInUse) == nbl.cj.size() || haveFreeEnergy,
3067 "Without free-energy all cj pair-list entries should be in use. "
3068 "Note that subsequent code does not make use of the equality, "
3069 "this check is only here to catch bugs");
3072 static void checkListSizeConsistency(const NbnxnPairlistGpu gmx_unused& nbl, bool gmx_unused haveFreeEnergy)
3074 /* We currently can not check consistency here */
3077 /* Set the buffer flags for newly added entries in the list */
3078 static void setBufferFlags(const NbnxnPairlistCpu& nbl,
3079 const int ncj_old_j,
3080 const int gridj_flag_shift,
3081 gmx_bitmask_t* gridj_flag,
3084 if (gmx::ssize(nbl.cj) > ncj_old_j)
3086 int cbFirst = nbl.cj[ncj_old_j].cj >> gridj_flag_shift;
3087 int cbLast = nbl.cj.back().cj >> gridj_flag_shift;
3088 for (int cb = cbFirst; cb <= cbLast; cb++)
3090 bitmask_init_bit(&gridj_flag[cb], th);
3095 static void setBufferFlags(const NbnxnPairlistGpu gmx_unused& nbl,
3096 int gmx_unused ncj_old_j,
3097 int gmx_unused gridj_flag_shift,
3098 gmx_bitmask_t gmx_unused* gridj_flag,
3101 GMX_ASSERT(false, "This function should never be called");
3104 /* Generates the part of pair-list nbl assigned to our thread */
3105 template<typename T>
3106 static void nbnxn_make_pairlist_part(const Nbnxm::GridSet& gridSet,
3109 PairsearchWork* work,
3110 const nbnxn_atomdata_t* nbat,
3111 const ListOfLists<int>& exclusions,
3113 const PairlistType pairlistType,
3115 gmx_bool bFBufferFlag,
3118 float nsubpair_tot_est,
3128 int ci_b, ci, ci_x, ci_y, ci_xy;
3130 real bx0, bx1, by0, by1, bz0, bz1;
3132 real d2cx, d2z, d2z_cx, d2z_cy, d2zx, d2zxy, d2xy;
3133 int cxf, cxl, cyf, cyf_x, cyl;
3134 int numDistanceChecks;
3135 int gridi_flag_shift = 0, gridj_flag_shift = 0;
3136 gmx_bitmask_t* gridj_flag = nullptr;
3137 int ncj_old_i, ncj_old_j;
3139 if (jGrid.geometry().isSimple != pairlistIsSimple(*nbl)
3140 || iGrid.geometry().isSimple != pairlistIsSimple(*nbl))
3142 gmx_incons("Grid incompatible with pair-list");
3146 GMX_ASSERT(nbl->na_ci == jGrid.geometry().numAtomsICluster,
3147 "The cluster sizes in the list and grid should match");
3148 nbl->na_cj = JClusterSizePerListType[pairlistType];
3149 na_cj_2log = get_2log(nbl->na_cj);
3155 /* Determine conversion of clusters to flag blocks */
3156 gridi_flag_shift = getBufferFlagShift(nbl->na_ci);
3157 gridj_flag_shift = getBufferFlagShift(nbl->na_cj);
3159 gridj_flag = work->buffer_flags.data();
3162 gridSet.getBox(box);
3164 const bool haveFep = gridSet.haveFep();
3166 const real rlist2 = nbl->rlist * nbl->rlist;
3168 // Select the cluster pair distance kernel type
3169 const ClusterDistanceKernelType kernelType = getClusterDistanceKernelType(pairlistType, *nbat);
3171 if (haveFep && !pairlistIsSimple(*nbl))
3173 /* Determine an atom-pair list cut-off distance for FEP atom pairs.
3174 * We should not simply use rlist, since then we would not have
3175 * the small, effective buffering of the NxN lists.
3176 * The buffer is on overestimate, but the resulting cost for pairs
3177 * beyond rlist is neglible compared to the FEP pairs within rlist.
3179 rl_fep2 = nbl->rlist + effective_buffer_1x1_vs_MxN(iGrid, jGrid);
3183 fprintf(debug, "nbl_fep atom-pair rlist %f\n", rl_fep2);
3185 rl_fep2 = rl_fep2 * rl_fep2;
3188 const Grid::Dimensions& iGridDims = iGrid.dimensions();
3189 const Grid::Dimensions& jGridDims = jGrid.dimensions();
3191 rbb2 = boundingbox_only_distance2(iGridDims, jGridDims, nbl->rlist, pairlistIsSimple(*nbl));
3195 fprintf(debug, "nbl bounding box only distance %f\n", std::sqrt(rbb2));
3198 const bool isIntraGridList = (&iGrid == &jGrid);
3200 /* Set the shift range */
3201 for (int d = 0; d < DIM; d++)
3203 /* Check if we need periodicity shifts.
3204 * Without PBC or with domain decomposition we don't need them.
3206 if (d >= numPbcDimensions(gridSet.domainSetup().pbcType)
3207 || gridSet.domainSetup().haveMultipleDomainsPerDim[d])
3213 const real listRangeCellToCell =
3214 listRangeForGridCellToGridCell(rlist, iGrid.dimensions(), jGrid.dimensions());
3215 if (d == XX && box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < listRangeCellToCell)
3225 const bool bSimple = pairlistIsSimple(*nbl);
3226 gmx::ArrayRef<const BoundingBox> bb_i;
3228 gmx::ArrayRef<const float> pbb_i;
3231 bb_i = iGrid.iBoundingBoxes();
3235 pbb_i = iGrid.packedBoundingBoxes();
3238 /* We use the normal bounding box format for both grid types */
3239 bb_i = iGrid.iBoundingBoxes();
3241 gmx::ArrayRef<const BoundingBox1D> bbcz_i = iGrid.zBoundingBoxes();
3242 gmx::ArrayRef<const int> flags_i = iGrid.clusterFlags();
3243 gmx::ArrayRef<const BoundingBox1D> bbcz_j = jGrid.zBoundingBoxes();
3244 int cell0_i = iGrid.cellOffset();
3249 "nbl nc_i %d col.av. %.1f ci_block %d\n",
3251 iGrid.numCells() / static_cast<double>(iGrid.numColumns()),
3255 numDistanceChecks = 0;
3257 const real listRangeBBToJCell2 =
3258 gmx::square(listRangeForBoundingBoxToGridCell(rlist, jGrid.dimensions()));
3260 /* Initially ci_b and ci to 1 before where we want them to start,
3261 * as they will both be incremented in next_ci.
3264 ci = th * ci_block - 1;
3267 while (next_ci(iGrid, nth, ci_block, &ci_x, &ci_y, &ci_b, &ci))
3269 if (bSimple && flags_i[ci] == 0)
3273 ncj_old_i = getNumSimpleJClustersInList(*nbl);
3276 if (!isIntraGridList && shp[XX] == 0)
3280 bx1 = bb_i[ci].upper.x;
3284 bx1 = iGridDims.lowerCorner[XX] + (real(ci_x) + 1) * iGridDims.cellSize[XX];
3286 if (bx1 < jGridDims.lowerCorner[XX])
3288 d2cx = gmx::square(jGridDims.lowerCorner[XX] - bx1);
3290 if (d2cx >= listRangeBBToJCell2)
3297 ci_xy = ci_x * iGridDims.numCells[YY] + ci_y;
3299 /* Loop over shift vectors in three dimensions */
3300 for (int tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
3302 const real shz = real(tz) * box[ZZ][ZZ];
3304 bz0 = bbcz_i[ci].lower + shz;
3305 bz1 = bbcz_i[ci].upper + shz;
3313 d2z = gmx::square(bz1);
3317 d2z = gmx::square(bz0 - box[ZZ][ZZ]);
3320 d2z_cx = d2z + d2cx;
3322 if (d2z_cx >= rlist2)
3327 bz1_frac = bz1 / real(iGrid.numCellsInColumn(ci_xy));
3332 /* The check with bz1_frac close to or larger than 1 comes later */
3334 for (int ty = -shp[YY]; ty <= shp[YY]; ty++)
3336 const real shy = real(ty) * box[YY][YY] + real(tz) * box[ZZ][YY];
3340 by0 = bb_i[ci].lower.y + shy;
3341 by1 = bb_i[ci].upper.y + shy;
3345 by0 = iGridDims.lowerCorner[YY] + (real(ci_y)) * iGridDims.cellSize[YY] + shy;
3346 by1 = iGridDims.lowerCorner[YY] + (real(ci_y) + 1) * iGridDims.cellSize[YY] + shy;
3349 get_cell_range<YY>(by0, by1, jGridDims, d2z_cx, rlist, &cyf, &cyl);
3357 if (by1 < jGridDims.lowerCorner[YY])
3359 d2z_cy += gmx::square(jGridDims.lowerCorner[YY] - by1);
3361 else if (by0 > jGridDims.upperCorner[YY])
3363 d2z_cy += gmx::square(by0 - jGridDims.upperCorner[YY]);
3366 for (int tx = -shp[XX]; tx <= shp[XX]; tx++)
3368 const int shift = XYZ2IS(tx, ty, tz);
3370 const bool excludeSubDiagonal = (isIntraGridList && shift == CENTRAL);
3372 if (c_pbcShiftBackward && isIntraGridList && shift > CENTRAL)
3378 real(tx) * box[XX][XX] + real(ty) * box[YY][XX] + real(tz) * box[ZZ][XX];
3382 bx0 = bb_i[ci].lower.x + shx;
3383 bx1 = bb_i[ci].upper.x + shx;
3387 bx0 = iGridDims.lowerCorner[XX] + (real(ci_x)) * iGridDims.cellSize[XX] + shx;
3388 bx1 = iGridDims.lowerCorner[XX] + (real(ci_x) + 1) * iGridDims.cellSize[XX] + shx;
3391 get_cell_range<XX>(bx0, bx1, jGridDims, d2z_cy, rlist, &cxf, &cxl);
3398 addNewIEntry(nbl, cell0_i + ci, shift, flags_i[ci]);
3400 if ((!c_pbcShiftBackward || excludeSubDiagonal) && cxf < ci_x)
3402 /* Leave the pairs with i > j.
3403 * x is the major index, so skip half of it.
3408 set_icell_bb(iGrid, ci, shx, shy, shz, nbl->work.get());
3410 icell_set_x(cell0_i + ci,
3419 for (int cx = cxf; cx <= cxl; cx++)
3421 const real cx_real = cx;
3423 if (jGridDims.lowerCorner[XX] + cx_real * jGridDims.cellSize[XX] > bx1)
3425 d2zx += gmx::square(jGridDims.lowerCorner[XX]
3426 + cx_real * jGridDims.cellSize[XX] - bx1);
3428 else if (jGridDims.lowerCorner[XX] + (cx_real + 1) * jGridDims.cellSize[XX] < bx0)
3430 d2zx += gmx::square(jGridDims.lowerCorner[XX]
3431 + (cx_real + 1) * jGridDims.cellSize[XX] - bx0);
3434 if (isIntraGridList && cx == 0 && (!c_pbcShiftBackward || shift == CENTRAL)
3437 /* Leave the pairs with i > j.
3438 * Skip half of y when i and j have the same x.
3447 for (int cy = cyf_x; cy <= cyl; cy++)
3449 const int columnStart =
3450 jGrid.firstCellInColumn(cx * jGridDims.numCells[YY] + cy);
3451 const int columnEnd =
3452 jGrid.firstCellInColumn(cx * jGridDims.numCells[YY] + cy + 1);
3454 const real cy_real = cy;
3456 if (jGridDims.lowerCorner[YY] + cy_real * jGridDims.cellSize[YY] > by1)
3458 d2zxy += gmx::square(jGridDims.lowerCorner[YY]
3459 + cy_real * jGridDims.cellSize[YY] - by1);
3461 else if (jGridDims.lowerCorner[YY] + (cy_real + 1) * jGridDims.cellSize[YY] < by0)
3463 d2zxy += gmx::square(jGridDims.lowerCorner[YY]
3464 + (cy_real + 1) * jGridDims.cellSize[YY] - by0);
3466 if (columnStart < columnEnd && d2zxy < listRangeBBToJCell2)
3468 /* To improve efficiency in the common case
3469 * of a homogeneous particle distribution,
3470 * we estimate the index of the middle cell
3471 * in range (midCell). We search down and up
3472 * starting from this index.
3474 * Note that the bbcz_j array contains bounds
3475 * for i-clusters, thus for clusters of 4 atoms.
3476 * For the common case where the j-cluster size
3477 * is 8, we could step with a stride of 2,
3478 * but we do not do this because it would
3479 * complicate this code even more.
3483 + static_cast<int>(bz1_frac
3484 * static_cast<real>(columnEnd - columnStart));
3485 if (midCell >= columnEnd)
3487 midCell = columnEnd - 1;
3492 /* Find the lowest cell that can possibly
3494 * Check if we hit the bottom of the grid,
3495 * if the j-cell is below the i-cell and if so,
3496 * if it is within range.
3498 int downTestCell = midCell;
3499 while (downTestCell >= columnStart
3500 && (bbcz_j[downTestCell].upper >= bz0
3501 || d2xy + gmx::square(bbcz_j[downTestCell].upper - bz0) < rlist2))
3505 int firstCell = downTestCell + 1;
3507 /* Find the highest cell that can possibly
3509 * Check if we hit the top of the grid,
3510 * if the j-cell is above the i-cell and if so,
3511 * if it is within range.
3513 int upTestCell = midCell + 1;
3514 while (upTestCell < columnEnd
3515 && (bbcz_j[upTestCell].lower <= bz1
3516 || d2xy + gmx::square(bbcz_j[upTestCell].lower - bz1) < rlist2))
3520 int lastCell = upTestCell - 1;
3522 #define NBNXN_REFCODE 0
3525 /* Simple reference code, for debugging,
3526 * overrides the more complex code above.
3528 firstCell = columnEnd;
3530 for (int k = columnStart; k < columnEnd; k++)
3532 if (d2xy + gmx::square(bbcz_j[k * NNBSBB_D + 1] - bz0) < rlist2
3537 if (d2xy + gmx::square(bbcz_j[k * NNBSBB_D] - bz1) < rlist2
3546 if (isIntraGridList)
3548 /* We want each atom/cell pair only once,
3549 * only use cj >= ci.
3551 if (!c_pbcShiftBackward || shift == CENTRAL)
3553 firstCell = std::max(firstCell, ci);
3557 if (firstCell <= lastCell)
3559 GMX_ASSERT(firstCell >= columnStart && lastCell < columnEnd,
3560 "The range should reside within the current grid "
3563 /* For f buffer flags with simple lists */
3564 ncj_old_j = getNumSimpleJClustersInList(*nbl);
3566 makeClusterListWrapper(nbl,
3577 &numDistanceChecks);
3581 setBufferFlags(*nbl, ncj_old_j, gridj_flag_shift, gridj_flag, th);
3584 incrementNumSimpleJClustersInList(nbl, ncj_old_j);
3590 if (!exclusions.empty())
3592 /* Set the exclusions for this ci list */
3593 setExclusionsForIEntry(
3594 gridSet, nbl, excludeSubDiagonal, na_cj_2log, *getOpenIEntry(nbl), exclusions);
3599 make_fep_list(gridSet.atomIndices(),
3613 /* Close this ci list */
3614 closeIEntry(nbl, nsubpair_max, progBal, nsubpair_tot_est, th, nth);
3619 if (bFBufferFlag && getNumSimpleJClustersInList(*nbl) > ncj_old_i)
3621 bitmask_init_bit(&(work->buffer_flags[(iGrid.cellOffset() + ci) >> gridi_flag_shift]), th);
3625 work->ndistc = numDistanceChecks;
3627 checkListSizeConsistency(*nbl, haveFep);
3631 fprintf(debug, "number of distance checks %d\n", numDistanceChecks);
3633 print_nblist_statistics(debug, *nbl, gridSet, rlist);
3637 fprintf(debug, "nbl FEP list pairs: %d\n", nbl_fep->nrj);
3642 static void reduce_buffer_flags(gmx::ArrayRef<PairsearchWork> searchWork,
3644 gmx::ArrayRef<gmx_bitmask_t> dest)
3646 for (int s = 0; s < nsrc; s++)
3648 gmx::ArrayRef<gmx_bitmask_t> flags(searchWork[s].buffer_flags);
3650 for (size_t b = 0; b < dest.size(); b++)
3652 gmx_bitmask_t& flag = dest[b];
3653 bitmask_union(&flag, flags[b]);
3658 static void print_reduction_cost(gmx::ArrayRef<const gmx_bitmask_t> flags, int nout)
3660 int nelem, nkeep, ncopy, nred, out;
3661 gmx_bitmask_t mask_0;
3667 bitmask_init_bit(&mask_0, 0);
3668 for (const gmx_bitmask_t& flag_mask : flags)
3670 if (bitmask_is_equal(flag_mask, mask_0))
3672 /* Only flag 0 is set, no copy of reduction required */
3676 else if (!bitmask_is_zero(flag_mask))
3679 for (out = 0; out < nout; out++)
3681 if (bitmask_is_set(flag_mask, out))
3697 const auto numFlags = static_cast<double>(flags.size());
3699 "nbnxn reduction: #flag %zu #list %d elem %4.2f, keep %4.2f copy %4.2f red %4.2f\n",
3708 /* Copies the list entries from src to dest when cjStart <= *cjGlobal < cjEnd.
3709 * *cjGlobal is updated with the cj count in src.
3710 * When setFlags==true, flag bit t is set in flag for all i and j clusters.
3712 template<bool setFlags>
3713 static void copySelectedListRange(const nbnxn_ci_t* gmx_restrict srcCi,
3714 const NbnxnPairlistCpu* gmx_restrict src,
3715 NbnxnPairlistCpu* gmx_restrict dest,
3716 gmx_bitmask_t* flag,
3721 const int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3723 dest->ci.push_back(*srcCi);
3724 dest->ci.back().cj_ind_start = dest->cj.size();
3725 dest->ci.back().cj_ind_end = dest->ci.back().cj_ind_start + ncj;
3729 bitmask_init_bit(&flag[srcCi->ci >> iFlagShift], t);
3732 for (int j = srcCi->cj_ind_start; j < srcCi->cj_ind_end; j++)
3734 dest->cj.push_back(src->cj[j]);
3738 /* NOTE: This is relatively expensive, since this
3739 * operation is done for all elements in the list,
3740 * whereas at list generation this is done only
3741 * once for each flag entry.
3743 bitmask_init_bit(&flag[src->cj[j].cj >> jFlagShift], t);
3748 #if defined(__GNUC__) && !defined(__clang__) && __GNUC__ == 7
3749 /* Avoid gcc 7 avx512 loop vectorization bug (actually only needed with -mavx512f) */
3750 # pragma GCC push_options
3751 # pragma GCC optimize("no-tree-vectorize")
3754 /* Returns the number of cluster pairs that are in use summed over all lists */
3755 static int countClusterpairs(gmx::ArrayRef<const NbnxnPairlistCpu> pairlists)
3757 /* gcc 7 with -mavx512f can miss the contributions of 16 consecutive
3758 * elements to the sum calculated in this loop. Above we have disabled
3759 * loop vectorization to avoid this bug.
3762 for (const auto& pairlist : pairlists)
3764 ncjTotal += pairlist.ncjInUse;
3769 #if defined(__GNUC__) && !defined(__clang__) && __GNUC__ == 7
3770 # pragma GCC pop_options
3773 /* This routine re-balances the pairlists such that all are nearly equally
3774 * sized. Only whole i-entries are moved between lists. These are moved
3775 * between the ends of the lists, such that the buffer reduction cost should
3776 * not change significantly.
3777 * Note that all original reduction flags are currently kept. This can lead
3778 * to reduction of parts of the force buffer that could be avoided. But since
3779 * the original lists are quite balanced, this will only give minor overhead.
3781 static void rebalanceSimpleLists(gmx::ArrayRef<const NbnxnPairlistCpu> srcSet,
3782 gmx::ArrayRef<NbnxnPairlistCpu> destSet,
3783 gmx::ArrayRef<PairsearchWork> searchWork)
3785 const int ncjTotal = countClusterpairs(srcSet);
3786 const int numLists = srcSet.ssize();
3787 const int ncjTarget = (ncjTotal + numLists - 1) / numLists;
3789 #pragma omp parallel num_threads(numLists)
3791 int t = gmx_omp_get_thread_num();
3793 int cjStart = ncjTarget * t;
3794 int cjEnd = ncjTarget * (t + 1);
3796 /* The destination pair-list for task/thread t */
3797 NbnxnPairlistCpu& dest = destSet[t];
3799 clear_pairlist(&dest);
3800 dest.na_cj = srcSet[0].na_cj;
3802 /* Note that the flags in the work struct (still) contain flags
3803 * for all entries that are present in srcSet->nbl[t].
3805 gmx_bitmask_t* flag = &searchWork[t].buffer_flags[0];
3807 int iFlagShift = getBufferFlagShift(dest.na_ci);
3808 int jFlagShift = getBufferFlagShift(dest.na_cj);
3811 for (int s = 0; s < numLists && cjGlobal < cjEnd; s++)
3813 const NbnxnPairlistCpu* src = &srcSet[s];
3815 if (cjGlobal + src->ncjInUse > cjStart)
3817 for (gmx::index i = 0; i < gmx::ssize(src->ci) && cjGlobal < cjEnd; i++)
3819 const nbnxn_ci_t* srcCi = &src->ci[i];
3820 int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3821 if (cjGlobal >= cjStart)
3823 /* If the source list is not our own, we need to set
3824 * extra flags (the template bool parameter).
3828 copySelectedListRange<true>(srcCi, src, &dest, flag, iFlagShift, jFlagShift, t);
3832 copySelectedListRange<false>(
3833 srcCi, src, &dest, flag, iFlagShift, jFlagShift, t);
3841 cjGlobal += src->ncjInUse;
3845 dest.ncjInUse = dest.cj.size();
3849 const int ncjTotalNew = countClusterpairs(destSet);
3850 GMX_RELEASE_ASSERT(ncjTotalNew == ncjTotal,
3851 "The total size of the lists before and after rebalancing should match");
3855 /* Returns if the pairlists are so imbalanced that it is worth rebalancing. */
3856 static bool checkRebalanceSimpleLists(gmx::ArrayRef<const NbnxnPairlistCpu> lists)
3858 int numLists = lists.ssize();
3861 for (int s = 0; s < numLists; s++)
3863 ncjMax = std::max(ncjMax, lists[s].ncjInUse);
3864 ncjTotal += lists[s].ncjInUse;
3868 fprintf(debug, "Pair-list ncjMax %d ncjTotal %d\n", ncjMax, ncjTotal);
3870 /* The rebalancing adds 3% extra time to the search. Heuristically we
3871 * determined that under common conditions the non-bonded kernel balance
3872 * improvement will outweigh this when the imbalance is more than 3%.
3873 * But this will, obviously, depend on search vs kernel time and nstlist.
3875 const real rebalanceTolerance = 1.03;
3877 return real(numLists * ncjMax) > real(ncjTotal) * rebalanceTolerance;
3880 /* Perform a count (linear) sort to sort the smaller lists to the end.
3881 * This avoids load imbalance on the GPU, as large lists will be
3882 * scheduled and executed first and the smaller lists later.
3883 * Load balancing between multi-processors only happens at the end
3884 * and there smaller lists lead to more effective load balancing.
3885 * The sorting is done on the cj4 count, not on the actual pair counts.
3886 * Not only does this make the sort faster, but it also results in
3887 * better load balancing than using a list sorted on exact load.
3888 * This function swaps the pointer in the pair list to avoid a copy operation.
3890 static void sort_sci(NbnxnPairlistGpu* nbl)
3892 if (nbl->cj4.size() <= nbl->sci.size())
3894 /* nsci = 0 or all sci have size 1, sorting won't change the order */
3898 NbnxnPairlistGpuWork& work = *nbl->work;
3900 /* We will distinguish differences up to double the average */
3901 const int m = static_cast<int>((2 * ssize(nbl->cj4)) / ssize(nbl->sci));
3903 /* Resize work.sci_sort so we can sort into it */
3904 work.sci_sort.resize(nbl->sci.size());
3906 std::vector<int>& sort = work.sortBuffer;
3907 /* Set up m + 1 entries in sort, initialized at 0 */
3909 sort.resize(m + 1, 0);
3910 /* Count the entries of each size */
3911 for (const nbnxn_sci_t& sci : nbl->sci)
3913 int i = std::min(m, sci.numJClusterGroups());
3916 /* Calculate the offset for each count */
3919 for (gmx::index i = m - 1; i >= 0; i--)
3922 sort[i] = sort[i + 1] + s0;
3926 /* Sort entries directly into place */
3927 gmx::ArrayRef<nbnxn_sci_t> sci_sort = work.sci_sort;
3928 for (const nbnxn_sci_t& sci : nbl->sci)
3930 int i = std::min(m, sci.numJClusterGroups());
3931 sci_sort[sort[i]++] = sci;
3934 /* Swap the sci pointers so we use the new, sorted list */
3935 std::swap(nbl->sci, work.sci_sort);
3938 /* Returns the i-zone range for pairlist construction for the give locality */
3939 static Range<int> getIZoneRange(const Nbnxm::GridSet::DomainSetup& domainSetup,
3940 const InteractionLocality locality)
3942 if (domainSetup.doTestParticleInsertion)
3944 /* With TPI we do grid 1, the inserted molecule, versus grid 0, the rest */
3947 else if (locality == InteractionLocality::Local)
3949 /* Local: only zone (grid) 0 vs 0 */
3954 /* Non-local: we need all i-zones */
3955 return { 0, int(domainSetup.zones->iZones.size()) };
3959 /* Returns the j-zone range for pairlist construction for the give locality and i-zone */
3960 static Range<int> getJZoneRange(const gmx_domdec_zones_t* ddZones,
3961 const InteractionLocality locality,
3964 if (locality == InteractionLocality::Local)
3966 /* Local: zone 0 vs 0 or with TPI 1 vs 0 */
3969 else if (iZone == 0)
3971 /* Non-local: we need to avoid the local (zone 0 vs 0) interactions */
3972 return { 1, *ddZones->iZones[iZone].jZoneRange.end() };
3976 /* Non-local with non-local i-zone: use all j-zones */
3977 return ddZones->iZones[iZone].jZoneRange;
3981 //! Prepares CPU lists produced by the search for dynamic pruning
3982 static void prepareListsForDynamicPruning(gmx::ArrayRef<NbnxnPairlistCpu> lists);
3984 void PairlistSet::constructPairlists(const Nbnxm::GridSet& gridSet,
3985 gmx::ArrayRef<PairsearchWork> searchWork,
3986 nbnxn_atomdata_t* nbat,
3987 const ListOfLists<int>& exclusions,
3988 const int minimumIlistCountForGpuBalancing,
3990 SearchCycleCounting* searchCycleCounting)
3992 const real rlist = params_.rlistOuter;
3994 int nsubpair_target;
3995 float nsubpair_tot_est;
3998 int np_tot, np_noq, np_hlj, nap;
4000 const int numLists = (isCpuType_ ? cpuLists_.size() : gpuLists_.size());
4004 fprintf(debug, "ns making %d nblists\n", numLists);
4007 nbat->bUseBufferFlags = (nbat->out.size() > 1);
4008 /* We should re-init the flags before making the first list */
4009 if (nbat->bUseBufferFlags && locality_ == InteractionLocality::Local)
4011 resizeAndZeroBufferFlags(&nbat->buffer_flags, nbat->numAtoms());
4014 if (!isCpuType_ && minimumIlistCountForGpuBalancing > 0)
4016 get_nsubpair_target(
4017 gridSet, locality_, rlist, minimumIlistCountForGpuBalancing, &nsubpair_target, &nsubpair_tot_est);
4021 nsubpair_target = 0;
4022 nsubpair_tot_est = 0;
4025 /* Clear all pair-lists */
4026 for (int th = 0; th < numLists; th++)
4030 clear_pairlist(&cpuLists_[th]);
4034 clear_pairlist(&gpuLists_[th]);
4037 if (params_.haveFep)
4039 clear_pairlist_fep(fepLists_[th].get());
4043 const gmx_domdec_zones_t* ddZones = gridSet.domainSetup().zones;
4044 GMX_ASSERT(locality_ == InteractionLocality::Local || ddZones != nullptr,
4045 "Nonlocal interaction locality with null ddZones.");
4047 const auto iZoneRange = getIZoneRange(gridSet.domainSetup(), locality_);
4049 for (const int iZone : iZoneRange)
4051 const Grid& iGrid = gridSet.grids()[iZone];
4053 const auto jZoneRange = getJZoneRange(ddZones, locality_, iZone);
4055 for (int jZone : jZoneRange)
4057 const Grid& jGrid = gridSet.grids()[jZone];
4061 fprintf(debug, "ns search grid %d vs %d\n", iZone, jZone);
4064 searchCycleCounting->start(enbsCCsearch);
4066 ci_block = get_ci_block_size(iGrid, gridSet.domainSetup().haveMultipleDomains, numLists);
4068 /* With GPU: generate progressively smaller lists for
4069 * load balancing for local only or non-local with 2 zones.
4071 progBal = (locality_ == InteractionLocality::Local || ddZones->n <= 2);
4073 #pragma omp parallel for num_threads(numLists) schedule(static)
4074 for (int th = 0; th < numLists; th++)
4078 /* Re-init the thread-local work flag data before making
4079 * the first list (not an elegant conditional).
4081 if (nbat->bUseBufferFlags && (iZone == 0 && jZone == 0))
4083 resizeAndZeroBufferFlags(&searchWork[th].buffer_flags, nbat->numAtoms());
4086 if (combineLists_ && th > 0)
4088 GMX_ASSERT(!isCpuType_, "Can only combine GPU lists");
4090 clear_pairlist(&gpuLists_[th]);
4093 PairsearchWork& work = searchWork[th];
4095 work.cycleCounter.start();
4097 t_nblist* fepListPtr = (fepLists_.empty() ? nullptr : fepLists_[th].get());
4099 /* Divide the i cells equally over the pairlists */
4102 nbnxn_make_pairlist_part(gridSet,
4109 params_.pairlistType,
4111 nbat->bUseBufferFlags,
4122 nbnxn_make_pairlist_part(gridSet,
4129 params_.pairlistType,
4131 nbat->bUseBufferFlags,
4141 work.cycleCounter.stop();
4143 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
4145 searchCycleCounting->stop(enbsCCsearch);
4150 for (int th = 0; th < numLists; th++)
4152 inc_nrnb(nrnb, eNR_NBNXN_DIST2, searchWork[th].ndistc);
4156 const NbnxnPairlistCpu& nbl = cpuLists_[th];
4157 np_tot += nbl.cj.size();
4158 np_noq += nbl.work->ncj_noq;
4159 np_hlj += nbl.work->ncj_hlj;
4163 const NbnxnPairlistGpu& nbl = gpuLists_[th];
4164 /* This count ignores potential subsequent pair pruning */
4165 np_tot += nbl.nci_tot;
4170 nap = cpuLists_[0].na_ci * cpuLists_[0].na_cj;
4174 nap = gmx::square(gpuLists_[0].na_ci);
4176 natpair_ljq_ = (np_tot - np_noq) * nap - np_hlj * nap / 2;
4177 natpair_lj_ = np_noq * nap;
4178 natpair_q_ = np_hlj * nap / 2;
4180 if (combineLists_ && numLists > 1)
4182 GMX_ASSERT(!isCpuType_, "Can only combine GPU lists");
4184 searchCycleCounting->start(enbsCCcombine);
4186 combine_nblists(gmx::constArrayRefFromArray(&gpuLists_[1], numLists - 1), &gpuLists_[0]);
4188 searchCycleCounting->stop(enbsCCcombine);
4195 if (numLists > 1 && checkRebalanceSimpleLists(cpuLists_))
4197 rebalanceSimpleLists(cpuLists_, cpuListsWork_, searchWork);
4199 /* Swap the sets of pair lists */
4200 cpuLists_.swap(cpuListsWork_);
4205 /* Sort the entries on size, large ones first */
4206 if (combineLists_ || gpuLists_.size() == 1)
4208 sort_sci(&gpuLists_[0]);
4212 #pragma omp parallel for num_threads(numLists) schedule(static)
4213 for (int th = 0; th < numLists; th++)
4217 sort_sci(&gpuLists_[th]);
4219 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
4224 if (nbat->bUseBufferFlags)
4226 reduce_buffer_flags(searchWork, numLists, nbat->buffer_flags);
4229 if (gridSet.haveFep())
4231 /* Balance the free-energy lists over all the threads */
4232 balance_fep_lists(fepLists_, searchWork);
4237 /* This is a fresh list, so not pruned, stored using ci.
4238 * ciOuter is invalid at this point.
4240 GMX_ASSERT(cpuLists_[0].ciOuter.empty(), "ciOuter is invalid so it should be empty");
4243 /* If we have more than one list, they either got rebalancing (CPU)
4244 * or combined (GPU), so we should dump the final result to debug.
4248 if (isCpuType_ && cpuLists_.size() > 1)
4250 for (auto& cpuList : cpuLists_)
4252 print_nblist_statistics(debug, cpuList, gridSet, rlist);
4255 else if (!isCpuType_ && gpuLists_.size() > 1)
4257 print_nblist_statistics(debug, gpuLists_[0], gridSet, rlist);
4267 for (auto& cpuList : cpuLists_)
4269 print_nblist_ci_cj(debug, cpuList);
4274 print_nblist_sci_cj(debug, gpuLists_[0]);
4278 if (nbat->bUseBufferFlags)
4280 print_reduction_cost(nbat->buffer_flags, numLists);
4284 if (params_.useDynamicPruning && isCpuType_)
4286 prepareListsForDynamicPruning(cpuLists_);
4290 void PairlistSets::construct(const InteractionLocality iLocality,
4291 PairSearch* pairSearch,
4292 nbnxn_atomdata_t* nbat,
4293 const ListOfLists<int>& exclusions,
4297 const auto& gridSet = pairSearch->gridSet();
4298 const auto* ddZones = gridSet.domainSetup().zones;
4300 /* The Nbnxm code can also work with more exclusions than those in i-zones only
4301 * when using DD, but the equality check can catch more issues.
4304 exclusions.empty() || (!ddZones && exclusions.ssize() == gridSet.numRealAtomsTotal())
4305 || (ddZones && exclusions.ssize() == ddZones->cg_range[ddZones->iZones.size()]),
4306 "exclusions should either be empty or the number of lists should match the number of "
4309 pairlistSet(iLocality).constructPairlists(gridSet,
4313 minimumIlistCountForGpuBalancing_,
4315 &pairSearch->cycleCounting_);
4317 if (iLocality == InteractionLocality::Local)
4319 outerListCreationStep_ = step;
4323 GMX_RELEASE_ASSERT(outerListCreationStep_ == step,
4324 "Outer list should be created at the same step as the inner list");
4327 /* Special performance logging stuff (env.var. GMX_NBNXN_CYCLE) */
4328 if (iLocality == InteractionLocality::Local)
4330 pairSearch->cycleCounting_.searchCount_++;
4332 if (pairSearch->cycleCounting_.recordCycles_
4333 && (!pairSearch->gridSet().domainSetup().haveMultipleDomains || iLocality == InteractionLocality::NonLocal)
4334 && pairSearch->cycleCounting_.searchCount_ % 100 == 0)
4336 pairSearch->cycleCounting_.printCycles(stderr, pairSearch->work());
4340 void nonbonded_verlet_t::constructPairlist(const InteractionLocality iLocality,
4341 const ListOfLists<int>& exclusions,
4345 pairlistSets_->construct(iLocality, pairSearch_.get(), nbat.get(), exclusions, step, nrnb);
4349 /* Launch the transfer of the pairlist to the GPU.
4351 * NOTE: The launch overhead is currently not timed separately
4353 Nbnxm::gpu_init_pairlist(gpu_nbv, pairlistSets().pairlistSet(iLocality).gpuList(), iLocality);
4357 static void prepareListsForDynamicPruning(gmx::ArrayRef<NbnxnPairlistCpu> lists)
4359 /* TODO: Restructure the lists so we have actual outer and inner
4360 * list objects so we can set a single pointer instead of
4361 * swapping several pointers.
4364 for (auto& list : lists)
4366 /* The search produced a list in ci/cj.
4367 * Swap the list pointers so we get the outer list is ciOuter,cjOuter
4368 * and we can prune that to get an inner list in ci/cj.
4370 GMX_RELEASE_ASSERT(list.ciOuter.empty() && list.cjOuter.empty(),
4371 "The outer lists should be empty before preparation");
4373 std::swap(list.ci, list.ciOuter);
4374 std::swap(list.cj, list.cjOuter);