2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
48 #include "gromacs/domdec/domdec_struct.h"
49 #include "gromacs/gmxlib/nrnb.h"
50 #include "gromacs/math/functions.h"
51 #include "gromacs/math/utilities.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/mdlib/gmx_omp_nthreads.h"
54 #include "gromacs/mdtypes/group.h"
55 #include "gromacs/mdtypes/md_enums.h"
56 #include "gromacs/nbnxm/atomdata.h"
57 #include "gromacs/nbnxm/gpu_data_mgmt.h"
58 #include "gromacs/nbnxm/nbnxm_geometry.h"
59 #include "gromacs/nbnxm/nbnxm_simd.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/topology/block.h"
65 #include "gromacs/utility/exceptions.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/gmxomp.h"
68 #include "gromacs/utility/smalloc.h"
70 #include "clusterdistancekerneltype.h"
72 #include "pairlistset.h"
73 #include "pairlistsets.h"
74 #include "pairlistwork.h"
75 #include "pairsearch.h"
77 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
79 using BoundingBox = Nbnxm::BoundingBox; // TODO: Remove when refactoring this file
80 using BoundingBox1D = Nbnxm::BoundingBox1D; // TODO: Remove when refactoring this file
82 using Grid = Nbnxm::Grid; // TODO: Remove when refactoring this file
84 // Convience alias for partial Nbnxn namespace usage
85 using InteractionLocality = Nbnxm::InteractionLocality;
87 /* We shift the i-particles backward for PBC.
88 * This leads to more conditionals than shifting forward.
89 * We do this to get more balanced pair lists.
91 constexpr bool c_pbcShiftBackward = true;
93 /* Layout for the nonbonded NxN pair lists */
94 enum class NbnxnLayout
96 NoSimd4x4, // i-cluster size 4, j-cluster size 4
97 Simd4xN, // i-cluster size 4, j-cluster size SIMD width
98 Simd2xNN, // i-cluster size 4, j-cluster size half SIMD width
99 Gpu8x8x8 // i-cluster size 8, j-cluster size 8 + super-clustering
103 /* Returns the j-cluster size */
104 template <NbnxnLayout layout>
105 static constexpr int jClusterSize()
107 static_assert(layout == NbnxnLayout::NoSimd4x4 || layout == NbnxnLayout::Simd4xN || layout == NbnxnLayout::Simd2xNN, "Currently jClusterSize only supports CPU layouts");
109 return layout == NbnxnLayout::Simd4xN ? GMX_SIMD_REAL_WIDTH : (layout == NbnxnLayout::Simd2xNN ? GMX_SIMD_REAL_WIDTH/2 : c_nbnxnCpuIClusterSize);
112 /*! \brief Returns the j-cluster index given the i-cluster index.
114 * \tparam jClusterSize The number of atoms in a j-cluster
115 * \tparam jSubClusterIndex The j-sub-cluster index (0/1), used when size(j-cluster) < size(i-cluster)
116 * \param[in] ci The i-cluster index
118 template <int jClusterSize, int jSubClusterIndex>
119 static inline int cjFromCi(int ci)
121 static_assert(jClusterSize == c_nbnxnCpuIClusterSize/2 || jClusterSize == c_nbnxnCpuIClusterSize || jClusterSize == c_nbnxnCpuIClusterSize*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
123 static_assert(jSubClusterIndex == 0 || jSubClusterIndex == 1,
124 "Only sub-cluster indices 0 and 1 are supported");
126 if (jClusterSize == c_nbnxnCpuIClusterSize/2)
128 if (jSubClusterIndex == 0)
134 return ((ci + 1) << 1) - 1;
137 else if (jClusterSize == c_nbnxnCpuIClusterSize)
147 /*! \brief Returns the j-cluster index given the i-cluster index.
149 * \tparam layout The pair-list layout
150 * \tparam jSubClusterIndex The j-sub-cluster index (0/1), used when size(j-cluster) < size(i-cluster)
151 * \param[in] ci The i-cluster index
153 template <NbnxnLayout layout, int jSubClusterIndex>
154 static inline int cjFromCi(int ci)
156 constexpr int clusterSize = jClusterSize<layout>();
158 return cjFromCi<clusterSize, jSubClusterIndex>(ci);
161 /* Returns the nbnxn coordinate data index given the i-cluster index */
162 template <NbnxnLayout layout>
163 static inline int xIndexFromCi(int ci)
165 constexpr int clusterSize = jClusterSize<layout>();
167 static_assert(clusterSize == c_nbnxnCpuIClusterSize/2 || clusterSize == c_nbnxnCpuIClusterSize || clusterSize == c_nbnxnCpuIClusterSize*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
169 if (clusterSize <= c_nbnxnCpuIClusterSize)
171 /* Coordinates are stored packed in groups of 4 */
176 /* Coordinates packed in 8, i-cluster size is half the packing width */
177 return (ci >> 1)*STRIDE_P8 + (ci & 1)*(c_packX8 >> 1);
181 /* Returns the nbnxn coordinate data index given the j-cluster index */
182 template <NbnxnLayout layout>
183 static inline int xIndexFromCj(int cj)
185 constexpr int clusterSize = jClusterSize<layout>();
187 static_assert(clusterSize == c_nbnxnCpuIClusterSize/2 || clusterSize == c_nbnxnCpuIClusterSize || clusterSize == c_nbnxnCpuIClusterSize*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
189 if (clusterSize == c_nbnxnCpuIClusterSize/2)
191 /* Coordinates are stored packed in groups of 4 */
192 return (cj >> 1)*STRIDE_P4 + (cj & 1)*(c_packX4 >> 1);
194 else if (clusterSize == c_nbnxnCpuIClusterSize)
196 /* Coordinates are stored packed in groups of 4 */
201 /* Coordinates are stored packed in groups of 8 */
208 void nbnxn_init_pairlist_fep(t_nblist *nl)
210 nl->type = GMX_NBLIST_INTERACTION_FREE_ENERGY;
211 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
212 /* The interaction functions are set in the free energy kernel fuction */
225 nl->jindex = nullptr;
227 nl->excl_fep = nullptr;
231 static void init_buffer_flags(nbnxn_buffer_flags_t *flags,
234 flags->nflag = (natoms + NBNXN_BUFFERFLAG_SIZE - 1)/NBNXN_BUFFERFLAG_SIZE;
235 if (flags->nflag > flags->flag_nalloc)
237 flags->flag_nalloc = over_alloc_large(flags->nflag);
238 srenew(flags->flag, flags->flag_nalloc);
240 for (int b = 0; b < flags->nflag; b++)
242 bitmask_clear(&(flags->flag[b]));
246 /* Returns the pair-list cutoff between a bounding box and a grid cell given an atom-to-atom pair-list cutoff
248 * Given a cutoff distance between atoms, this functions returns the cutoff
249 * distance2 between a bounding box of a group of atoms and a grid cell.
250 * Since atoms can be geometrically outside of the cell they have been
251 * assigned to (when atom groups instead of individual atoms are assigned
252 * to cells), this distance returned can be larger than the input.
255 listRangeForBoundingBoxToGridCell(real rlist,
256 const Grid::Dimensions &gridDims)
258 return rlist + gridDims.maxAtomGroupRadius;
261 /* Returns the pair-list cutoff between a grid cells given an atom-to-atom pair-list cutoff
263 * Given a cutoff distance between atoms, this functions returns the cutoff
264 * distance2 between two grid cells.
265 * Since atoms can be geometrically outside of the cell they have been
266 * assigned to (when atom groups instead of individual atoms are assigned
267 * to cells), this distance returned can be larger than the input.
270 listRangeForGridCellToGridCell(real rlist,
271 const Grid::Dimensions &iGridDims,
272 const Grid::Dimensions &jGridDims)
274 return rlist + iGridDims.maxAtomGroupRadius + jGridDims.maxAtomGroupRadius;
277 /* Determines the cell range along one dimension that
278 * the bounding box b0 - b1 sees.
281 static void get_cell_range(real b0, real b1,
282 const Grid::Dimensions &jGridDims,
283 real d2, real rlist, int *cf, int *cl)
285 real listRangeBBToCell2 = gmx::square(listRangeForBoundingBoxToGridCell(rlist, jGridDims));
286 real distanceInCells = (b0 - jGridDims.lowerCorner[dim])*jGridDims.invCellSize[dim];
287 *cf = std::max(static_cast<int>(distanceInCells), 0);
290 d2 + gmx::square((b0 - jGridDims.lowerCorner[dim]) - (*cf - 1 + 1)*jGridDims.cellSize[dim]) < listRangeBBToCell2)
295 *cl = std::min(static_cast<int>((b1 - jGridDims.lowerCorner[dim])*jGridDims.invCellSize[dim]), jGridDims.numCells[dim] - 1);
296 while (*cl < jGridDims.numCells[dim] - 1 &&
297 d2 + gmx::square((*cl + 1)*jGridDims.cellSize[dim] - (b1 - jGridDims.lowerCorner[dim])) < listRangeBBToCell2)
303 /* Reference code calculating the distance^2 between two bounding boxes */
305 static float box_dist2(float bx0, float bx1, float by0,
306 float by1, float bz0, float bz1,
307 const BoundingBox *bb)
310 float dl, dh, dm, dm0;
314 dl = bx0 - bb->upper.x;
315 dh = bb->lower.x - bx1;
316 dm = std::max(dl, dh);
317 dm0 = std::max(dm, 0.0f);
320 dl = by0 - bb->upper.y;
321 dh = bb->lower.y - by1;
322 dm = std::max(dl, dh);
323 dm0 = std::max(dm, 0.0f);
326 dl = bz0 - bb->upper.z;
327 dh = bb->lower.z - bz1;
328 dm = std::max(dl, dh);
329 dm0 = std::max(dm, 0.0f);
336 #if !NBNXN_SEARCH_BB_SIMD4
338 /*! \brief Plain C code calculating the distance^2 between two bounding boxes in xyz0 format
340 * \param[in] bb_i First bounding box
341 * \param[in] bb_j Second bounding box
343 static float clusterBoundingBoxDistance2(const BoundingBox &bb_i,
344 const BoundingBox &bb_j)
346 float dl = bb_i.lower.x - bb_j.upper.x;
347 float dh = bb_j.lower.x - bb_i.upper.x;
348 float dm = std::max(dl, dh);
349 float dm0 = std::max(dm, 0.0f);
352 dl = bb_i.lower.y - bb_j.upper.y;
353 dh = bb_j.lower.y - bb_i.upper.y;
354 dm = std::max(dl, dh);
355 dm0 = std::max(dm, 0.0f);
358 dl = bb_i.lower.z - bb_j.upper.z;
359 dh = bb_j.lower.z - bb_i.upper.z;
360 dm = std::max(dl, dh);
361 dm0 = std::max(dm, 0.0f);
367 #else /* NBNXN_SEARCH_BB_SIMD4 */
369 /*! \brief 4-wide SIMD code calculating the distance^2 between two bounding boxes in xyz0 format
371 * \param[in] bb_i First bounding box, should be aligned for 4-wide SIMD
372 * \param[in] bb_j Second bounding box, should be aligned for 4-wide SIMD
374 static float clusterBoundingBoxDistance2(const BoundingBox &bb_i,
375 const BoundingBox &bb_j)
377 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
380 const Simd4Float bb_i_S0 = load4(bb_i.lower.ptr());
381 const Simd4Float bb_i_S1 = load4(bb_i.upper.ptr());
382 const Simd4Float bb_j_S0 = load4(bb_j.lower.ptr());
383 const Simd4Float bb_j_S1 = load4(bb_j.upper.ptr());
385 const Simd4Float dl_S = bb_i_S0 - bb_j_S1;
386 const Simd4Float dh_S = bb_j_S0 - bb_i_S1;
388 const Simd4Float dm_S = max(dl_S, dh_S);
389 const Simd4Float dm0_S = max(dm_S, simd4SetZeroF());
391 return dotProduct(dm0_S, dm0_S);
394 /* Calculate bb bounding distances of bb_i[si,...,si+3] and store them in d2 */
395 template <int boundingBoxStart>
396 static inline void gmx_simdcall
397 clusterBoundingBoxDistance2_xxxx_simd4_inner(const float *bb_i,
399 const Simd4Float xj_l,
400 const Simd4Float yj_l,
401 const Simd4Float zj_l,
402 const Simd4Float xj_h,
403 const Simd4Float yj_h,
404 const Simd4Float zj_h)
406 constexpr int stride = c_packedBoundingBoxesDimSize;
408 const int shi = boundingBoxStart*Nbnxm::c_numBoundingBoxBounds1D*DIM;
410 const Simd4Float zero = setZero();
412 const Simd4Float xi_l = load4(bb_i + shi + 0*stride);
413 const Simd4Float yi_l = load4(bb_i + shi + 1*stride);
414 const Simd4Float zi_l = load4(bb_i + shi + 2*stride);
415 const Simd4Float xi_h = load4(bb_i + shi + 3*stride);
416 const Simd4Float yi_h = load4(bb_i + shi + 4*stride);
417 const Simd4Float zi_h = load4(bb_i + shi + 5*stride);
419 const Simd4Float dx_0 = xi_l - xj_h;
420 const Simd4Float dy_0 = yi_l - yj_h;
421 const Simd4Float dz_0 = zi_l - zj_h;
423 const Simd4Float dx_1 = xj_l - xi_h;
424 const Simd4Float dy_1 = yj_l - yi_h;
425 const Simd4Float dz_1 = zj_l - zi_h;
427 const Simd4Float mx = max(dx_0, dx_1);
428 const Simd4Float my = max(dy_0, dy_1);
429 const Simd4Float mz = max(dz_0, dz_1);
431 const Simd4Float m0x = max(mx, zero);
432 const Simd4Float m0y = max(my, zero);
433 const Simd4Float m0z = max(mz, zero);
435 const Simd4Float d2x = m0x * m0x;
436 const Simd4Float d2y = m0y * m0y;
437 const Simd4Float d2z = m0z * m0z;
439 const Simd4Float d2s = d2x + d2y;
440 const Simd4Float d2t = d2s + d2z;
442 store4(d2 + boundingBoxStart, d2t);
445 /* 4-wide SIMD code for nsi bb distances for bb format xxxxyyyyzzzz */
447 clusterBoundingBoxDistance2_xxxx_simd4(const float *bb_j,
452 constexpr int stride = c_packedBoundingBoxesDimSize;
454 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
457 const Simd4Float xj_l = Simd4Float(bb_j[0*stride]);
458 const Simd4Float yj_l = Simd4Float(bb_j[1*stride]);
459 const Simd4Float zj_l = Simd4Float(bb_j[2*stride]);
460 const Simd4Float xj_h = Simd4Float(bb_j[3*stride]);
461 const Simd4Float yj_h = Simd4Float(bb_j[4*stride]);
462 const Simd4Float zj_h = Simd4Float(bb_j[5*stride]);
464 /* Here we "loop" over si (0,stride) from 0 to nsi with step stride.
465 * But as we know the number of iterations is 1 or 2, we unroll manually.
467 clusterBoundingBoxDistance2_xxxx_simd4_inner<0>(bb_i, d2,
472 clusterBoundingBoxDistance2_xxxx_simd4_inner<stride>(bb_i, d2,
478 #endif /* NBNXN_SEARCH_BB_SIMD4 */
481 /* Returns if any atom pair from two clusters is within distance sqrt(rlist2) */
482 static inline gmx_bool
483 clusterpair_in_range(const NbnxnPairlistGpuWork &work,
485 int csj, int stride, const real *x_j,
488 #if !GMX_SIMD4_HAVE_REAL
491 * All coordinates are stored as xyzxyz...
494 const real *x_i = work.iSuperClusterData.x.data();
496 for (int i = 0; i < c_nbnxnGpuClusterSize; i++)
498 int i0 = (si*c_nbnxnGpuClusterSize + i)*DIM;
499 for (int j = 0; j < c_nbnxnGpuClusterSize; j++)
501 int j0 = (csj*c_nbnxnGpuClusterSize + j)*stride;
503 real d2 = gmx::square(x_i[i0 ] - x_j[j0 ]) + gmx::square(x_i[i0+1] - x_j[j0+1]) + gmx::square(x_i[i0+2] - x_j[j0+2]);
514 #else /* !GMX_SIMD4_HAVE_REAL */
516 /* 4-wide SIMD version.
517 * The coordinates x_i are stored as xxxxyyyy..., x_j is stored xyzxyz...
518 * Using 8-wide AVX(2) is not faster on Intel Sandy Bridge and Haswell.
520 static_assert(c_nbnxnGpuClusterSize == 8 || c_nbnxnGpuClusterSize == 4,
521 "A cluster is hard-coded to 4/8 atoms.");
523 Simd4Real rc2_S = Simd4Real(rlist2);
525 const real *x_i = work.iSuperClusterData.xSimd.data();
527 int dim_stride = c_nbnxnGpuClusterSize*DIM;
528 Simd4Real ix_S0 = load4(x_i + si*dim_stride + 0*GMX_SIMD4_WIDTH);
529 Simd4Real iy_S0 = load4(x_i + si*dim_stride + 1*GMX_SIMD4_WIDTH);
530 Simd4Real iz_S0 = load4(x_i + si*dim_stride + 2*GMX_SIMD4_WIDTH);
532 Simd4Real ix_S1, iy_S1, iz_S1;
533 if (c_nbnxnGpuClusterSize == 8)
535 ix_S1 = load4(x_i + si*dim_stride + 3*GMX_SIMD4_WIDTH);
536 iy_S1 = load4(x_i + si*dim_stride + 4*GMX_SIMD4_WIDTH);
537 iz_S1 = load4(x_i + si*dim_stride + 5*GMX_SIMD4_WIDTH);
539 /* We loop from the outer to the inner particles to maximize
540 * the chance that we find a pair in range quickly and return.
542 int j0 = csj*c_nbnxnGpuClusterSize;
543 int j1 = j0 + c_nbnxnGpuClusterSize - 1;
546 Simd4Real jx0_S, jy0_S, jz0_S;
547 Simd4Real jx1_S, jy1_S, jz1_S;
549 Simd4Real dx_S0, dy_S0, dz_S0;
550 Simd4Real dx_S1, dy_S1, dz_S1;
551 Simd4Real dx_S2, dy_S2, dz_S2;
552 Simd4Real dx_S3, dy_S3, dz_S3;
563 Simd4Bool wco_any_S01, wco_any_S23, wco_any_S;
565 jx0_S = Simd4Real(x_j[j0*stride+0]);
566 jy0_S = Simd4Real(x_j[j0*stride+1]);
567 jz0_S = Simd4Real(x_j[j0*stride+2]);
569 jx1_S = Simd4Real(x_j[j1*stride+0]);
570 jy1_S = Simd4Real(x_j[j1*stride+1]);
571 jz1_S = Simd4Real(x_j[j1*stride+2]);
573 /* Calculate distance */
574 dx_S0 = ix_S0 - jx0_S;
575 dy_S0 = iy_S0 - jy0_S;
576 dz_S0 = iz_S0 - jz0_S;
577 dx_S2 = ix_S0 - jx1_S;
578 dy_S2 = iy_S0 - jy1_S;
579 dz_S2 = iz_S0 - jz1_S;
580 if (c_nbnxnGpuClusterSize == 8)
582 dx_S1 = ix_S1 - jx0_S;
583 dy_S1 = iy_S1 - jy0_S;
584 dz_S1 = iz_S1 - jz0_S;
585 dx_S3 = ix_S1 - jx1_S;
586 dy_S3 = iy_S1 - jy1_S;
587 dz_S3 = iz_S1 - jz1_S;
590 /* rsq = dx*dx+dy*dy+dz*dz */
591 rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
592 rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
593 if (c_nbnxnGpuClusterSize == 8)
595 rsq_S1 = norm2(dx_S1, dy_S1, dz_S1);
596 rsq_S3 = norm2(dx_S3, dy_S3, dz_S3);
599 wco_S0 = (rsq_S0 < rc2_S);
600 wco_S2 = (rsq_S2 < rc2_S);
601 if (c_nbnxnGpuClusterSize == 8)
603 wco_S1 = (rsq_S1 < rc2_S);
604 wco_S3 = (rsq_S3 < rc2_S);
606 if (c_nbnxnGpuClusterSize == 8)
608 wco_any_S01 = wco_S0 || wco_S1;
609 wco_any_S23 = wco_S2 || wco_S3;
610 wco_any_S = wco_any_S01 || wco_any_S23;
614 wco_any_S = wco_S0 || wco_S2;
617 if (anyTrue(wco_any_S))
628 #endif /* !GMX_SIMD4_HAVE_REAL */
631 /* Returns the j-cluster index for index cjIndex in a cj list */
632 static inline int nblCj(gmx::ArrayRef<const nbnxn_cj_t> cjList,
635 return cjList[cjIndex].cj;
638 /* Returns the j-cluster index for index cjIndex in a cj4 list */
639 static inline int nblCj(gmx::ArrayRef<const nbnxn_cj4_t> cj4List,
642 return cj4List[cjIndex/c_nbnxnGpuJgroupSize].cj[cjIndex & (c_nbnxnGpuJgroupSize - 1)];
645 /* Returns the i-interaction mask of the j sub-cell for index cj_ind */
646 static unsigned int nbl_imask0(const NbnxnPairlistGpu *nbl, int cj_ind)
648 return nbl->cj4[cj_ind/c_nbnxnGpuJgroupSize].imei[0].imask;
651 NbnxnPairlistCpu::NbnxnPairlistCpu() :
652 na_ci(c_nbnxnCpuIClusterSize),
657 work(std::make_unique<NbnxnPairlistCpuWork>())
661 NbnxnPairlistGpu::NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy) :
662 na_ci(c_nbnxnGpuClusterSize),
663 na_cj(c_nbnxnGpuClusterSize),
664 na_sc(c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize),
666 sci({}, {pinningPolicy}),
667 cj4({}, {pinningPolicy}),
668 excl({}, {pinningPolicy}),
670 work(std::make_unique<NbnxnPairlistGpuWork>())
672 static_assert(c_nbnxnGpuNumClusterPerSupercluster == c_gpuNumClusterPerCell,
673 "The search code assumes that the a super-cluster matches a search grid cell");
675 static_assert(sizeof(cj4[0].imei[0].imask)*8 >= c_nbnxnGpuJgroupSize*c_gpuNumClusterPerCell,
676 "The i super-cluster cluster interaction mask does not contain a sufficient number of bits");
678 static_assert(sizeof(excl[0])*8 >= c_nbnxnGpuJgroupSize*c_gpuNumClusterPerCell, "The GPU exclusion mask does not contain a sufficient number of bits");
680 // We always want a first entry without any exclusions
684 // TODO: Move to pairlistset.cpp
685 PairlistSet::PairlistSet(const Nbnxm::InteractionLocality locality,
686 const PairlistParams &pairlistParams) :
688 params_(pairlistParams)
691 (params_.pairlistType == PairlistType::Simple4x2 ||
692 params_.pairlistType == PairlistType::Simple4x4 ||
693 params_.pairlistType == PairlistType::Simple4x8);
694 // Currently GPU lists are always combined
695 combineLists_ = !isCpuType_;
697 const int numLists = gmx_omp_nthreads_get(emntNonbonded);
699 if (!combineLists_ &&
700 numLists > NBNXN_BUFFERFLAG_MAX_THREADS)
702 gmx_fatal(FARGS, "%d OpenMP threads were requested. Since the non-bonded force buffer reduction is prohibitively slow with more than %d threads, we do not allow this. Use %d or less OpenMP threads.",
703 numLists, NBNXN_BUFFERFLAG_MAX_THREADS, NBNXN_BUFFERFLAG_MAX_THREADS);
708 cpuLists_.resize(numLists);
711 cpuListsWork_.resize(numLists);
716 /* Only list 0 is used on the GPU, use normal allocation for i>0 */
717 gpuLists_.emplace_back(gmx::PinningPolicy::PinnedIfSupported);
718 /* Lists 0 to numLists are use for constructing lists in parallel
719 * on the CPU using numLists threads (and then merged into list 0).
721 for (int i = 1; i < numLists; i++)
723 gpuLists_.emplace_back(gmx::PinningPolicy::CannotBePinned);
728 fepLists_.resize(numLists);
730 /* Execute in order to avoid memory interleaving between threads */
731 #pragma omp parallel for num_threads(numLists) schedule(static)
732 for (int i = 0; i < numLists; i++)
736 /* We used to allocate all normal lists locally on each thread
737 * as well. The question is if allocating the object on the
738 * master thread (but all contained list memory thread local)
739 * impacts performance.
741 snew(fepLists_[i], 1);
742 nbnxn_init_pairlist_fep(fepLists_[i]);
744 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
749 /* Print statistics of a pair list, used for debug output */
750 static void print_nblist_statistics(FILE *fp,
751 const NbnxnPairlistCpu &nbl,
752 const Nbnxm::GridSet &gridSet,
755 const Grid &grid = gridSet.grids()[0];
756 const Grid::Dimensions &dims = grid.dimensions();
758 fprintf(fp, "nbl nci %zu ncj %d\n",
759 nbl.ci.size(), nbl.ncjInUse);
760 const int numAtomsJCluster = grid.geometry().numAtomsJCluster;
761 const double numAtomsPerCell = nbl.ncjInUse/static_cast<double>(grid.numCells())*numAtomsJCluster;
762 fprintf(fp, "nbl na_cj %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
763 nbl.na_cj, rl, nbl.ncjInUse, nbl.ncjInUse/static_cast<double>(grid.numCells()),
765 numAtomsPerCell/(0.5*4.0/3.0*M_PI*rl*rl*rl*grid.numCells()*numAtomsJCluster/(dims.gridSize[XX]*dims.gridSize[YY]*dims.gridSize[ZZ])));
767 fprintf(fp, "nbl average j cell list length %.1f\n",
768 0.25*nbl.ncjInUse/std::max(static_cast<double>(nbl.ci.size()), 1.0));
770 int cs[SHIFTS] = { 0 };
772 for (const nbnxn_ci_t &ciEntry : nbl.ci)
774 cs[ciEntry.shift & NBNXN_CI_SHIFT] +=
775 ciEntry.cj_ind_end - ciEntry.cj_ind_start;
777 int j = ciEntry.cj_ind_start;
778 while (j < ciEntry.cj_ind_end &&
779 nbl.cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
785 fprintf(fp, "nbl cell pairs, total: %zu excl: %d %.1f%%\n",
786 nbl.cj.size(), npexcl, 100*npexcl/std::max(static_cast<double>(nbl.cj.size()), 1.0));
787 for (int s = 0; s < SHIFTS; s++)
791 fprintf(fp, "nbl shift %2d ncj %3d\n", s, cs[s]);
796 /* Print statistics of a pair lists, used for debug output */
797 static void print_nblist_statistics(FILE *fp,
798 const NbnxnPairlistGpu &nbl,
799 const Nbnxm::GridSet &gridSet,
802 const Grid &grid = gridSet.grids()[0];
803 const Grid::Dimensions &dims = grid.dimensions();
805 fprintf(fp, "nbl nsci %zu ncj4 %zu nsi %d excl4 %zu\n",
806 nbl.sci.size(), nbl.cj4.size(), nbl.nci_tot, nbl.excl.size());
807 const int numAtomsCluster = grid.geometry().numAtomsICluster;
808 const double numAtomsPerCell = nbl.nci_tot/static_cast<double>(grid.numClusters())*numAtomsCluster;
809 fprintf(fp, "nbl na_c %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
810 nbl.na_ci, rl, nbl.nci_tot, nbl.nci_tot/static_cast<double>(grid.numClusters()),
812 numAtomsPerCell/(0.5*4.0/3.0*M_PI*rl*rl*rl*grid.numClusters()*numAtomsCluster/(dims.gridSize[XX]*dims.gridSize[YY]*dims.gridSize[ZZ])));
817 int c[c_gpuNumClusterPerCell + 1] = { 0 };
818 for (const nbnxn_sci_t &sci : nbl.sci)
821 for (int j4 = sci.cj4_ind_start; j4 < sci.cj4_ind_end; j4++)
823 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
826 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
828 if (nbl.cj4[j4].imei[0].imask & (1U << (j*c_gpuNumClusterPerCell + si)))
839 nsp_max = std::max(nsp_max, nsp);
841 if (!nbl.sci.empty())
843 sum_nsp /= nbl.sci.size();
844 sum_nsp2 /= nbl.sci.size();
846 fprintf(fp, "nbl #cluster-pairs: av %.1f stddev %.1f max %d\n",
847 sum_nsp, std::sqrt(sum_nsp2 - sum_nsp*sum_nsp), nsp_max);
849 if (!nbl.cj4.empty())
851 for (int b = 0; b <= c_gpuNumClusterPerCell; b++)
853 fprintf(fp, "nbl j-list #i-subcell %d %7d %4.1f\n",
854 b, c[b], 100.0*c[b]/size_t {nbl.cj4.size()*c_nbnxnGpuJgroupSize});
859 /* Returns a pointer to the exclusion mask for j-cluster-group \p cj4 and warp \p warp
860 * Generates a new exclusion entry when the j-cluster-group uses
861 * the default all-interaction mask at call time, so the returned mask
862 * can be modified when needed.
864 static nbnxn_excl_t *get_exclusion_mask(NbnxnPairlistGpu *nbl,
868 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
870 /* No exclusions set, make a new list entry */
871 const size_t oldSize = nbl->excl.size();
872 GMX_ASSERT(oldSize >= 1, "We should always have entry [0]");
873 /* Add entry with default values: no exclusions */
874 nbl->excl.resize(oldSize + 1);
875 nbl->cj4[cj4].imei[warp].excl_ind = oldSize;
878 return &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
881 static void set_self_and_newton_excls_supersub(NbnxnPairlistGpu *nbl,
882 int cj4_ind, int sj_offset,
883 int i_cluster_in_cell)
885 nbnxn_excl_t *excl[c_nbnxnGpuClusterpairSplit];
887 /* Here we only set the set self and double pair exclusions */
889 /* Reserve extra elements, so the resize() in get_exclusion_mask()
890 * will not invalidate excl entries in the loop below
892 nbl->excl.reserve(nbl->excl.size() + c_nbnxnGpuClusterpairSplit);
893 for (int w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
895 excl[w] = get_exclusion_mask(nbl, cj4_ind, w);
898 /* Only minor < major bits set */
899 for (int ej = 0; ej < nbl->na_ci; ej++)
902 for (int ei = ej; ei < nbl->na_ci; ei++)
904 excl[w]->pair[(ej & (c_nbnxnGpuJgroupSize-1))*nbl->na_ci + ei] &=
905 ~(1U << (sj_offset*c_gpuNumClusterPerCell + i_cluster_in_cell));
910 /* Returns a diagonal or off-diagonal interaction mask for plain C lists */
911 static unsigned int get_imask(gmx_bool rdiag, int ci, int cj)
913 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
916 /* Returns a diagonal or off-diagonal interaction mask for cj-size=2 */
917 gmx_unused static unsigned int get_imask_simd_j2(gmx_bool rdiag, int ci, int cj)
919 return (rdiag && ci*2 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_0 :
920 (rdiag && ci*2+1 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_1 :
921 NBNXN_INTERACTION_MASK_ALL));
924 /* Returns a diagonal or off-diagonal interaction mask for cj-size=4 */
925 gmx_unused static unsigned int get_imask_simd_j4(gmx_bool rdiag, int ci, int cj)
927 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
930 /* Returns a diagonal or off-diagonal interaction mask for cj-size=8 */
931 gmx_unused static unsigned int get_imask_simd_j8(gmx_bool rdiag, int ci, int cj)
933 return (rdiag && ci == cj*2 ? NBNXN_INTERACTION_MASK_DIAG_J8_0 :
934 (rdiag && ci == cj*2+1 ? NBNXN_INTERACTION_MASK_DIAG_J8_1 :
935 NBNXN_INTERACTION_MASK_ALL));
939 #if GMX_SIMD_REAL_WIDTH == 2
940 #define get_imask_simd_4xn get_imask_simd_j2
942 #if GMX_SIMD_REAL_WIDTH == 4
943 #define get_imask_simd_4xn get_imask_simd_j4
945 #if GMX_SIMD_REAL_WIDTH == 8
946 #define get_imask_simd_4xn get_imask_simd_j8
947 #define get_imask_simd_2xnn get_imask_simd_j4
949 #if GMX_SIMD_REAL_WIDTH == 16
950 #define get_imask_simd_2xnn get_imask_simd_j8
954 /* Plain C code for checking and adding cluster-pairs to the list.
956 * \param[in] gridj The j-grid
957 * \param[in,out] nbl The pair-list to store the cluster pairs in
958 * \param[in] icluster The index of the i-cluster
959 * \param[in] jclusterFirst The first cluster in the j-range
960 * \param[in] jclusterLast The last cluster in the j-range
961 * \param[in] excludeSubDiagonal Exclude atom pairs with i-index > j-index
962 * \param[in] x_j Coordinates for the j-atom, in xyz format
963 * \param[in] rlist2 The squared list cut-off
964 * \param[in] rbb2 The squared cut-off for putting cluster-pairs in the list based on bounding box distance only
965 * \param[in,out] numDistanceChecks The number of distance checks performed
968 makeClusterListSimple(const Grid &jGrid,
969 NbnxnPairlistCpu * nbl,
973 bool excludeSubDiagonal,
974 const real * gmx_restrict x_j,
977 int * gmx_restrict numDistanceChecks)
979 const BoundingBox * gmx_restrict bb_ci = nbl->work->iClusterData.bb.data();
980 const real * gmx_restrict x_ci = nbl->work->iClusterData.x.data();
985 while (!InRange && jclusterFirst <= jclusterLast)
987 real d2 = clusterBoundingBoxDistance2(bb_ci[0], jGrid.jBoundingBoxes()[jclusterFirst]);
988 *numDistanceChecks += 2;
990 /* Check if the distance is within the distance where
991 * we use only the bounding box distance rbb,
992 * or within the cut-off and there is at least one atom pair
993 * within the cut-off.
999 else if (d2 < rlist2)
1001 int cjf_gl = jGrid.cellOffset() + jclusterFirst;
1002 for (int i = 0; i < c_nbnxnCpuIClusterSize && !InRange; i++)
1004 for (int j = 0; j < c_nbnxnCpuIClusterSize; j++)
1006 InRange = InRange ||
1007 (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjf_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+XX]) +
1008 gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjf_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+YY]) +
1009 gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjf_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+ZZ]) < rlist2);
1012 *numDistanceChecks += c_nbnxnCpuIClusterSize*c_nbnxnCpuIClusterSize;
1025 while (!InRange && jclusterLast > jclusterFirst)
1027 real d2 = clusterBoundingBoxDistance2(bb_ci[0], jGrid.jBoundingBoxes()[jclusterLast]);
1028 *numDistanceChecks += 2;
1030 /* Check if the distance is within the distance where
1031 * we use only the bounding box distance rbb,
1032 * or within the cut-off and there is at least one atom pair
1033 * within the cut-off.
1039 else if (d2 < rlist2)
1041 int cjl_gl = jGrid.cellOffset() + jclusterLast;
1042 for (int i = 0; i < c_nbnxnCpuIClusterSize && !InRange; i++)
1044 for (int j = 0; j < c_nbnxnCpuIClusterSize; j++)
1046 InRange = InRange ||
1047 (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjl_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+XX]) +
1048 gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjl_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+YY]) +
1049 gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjl_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+ZZ]) < rlist2);
1052 *numDistanceChecks += c_nbnxnCpuIClusterSize*c_nbnxnCpuIClusterSize;
1060 if (jclusterFirst <= jclusterLast)
1062 for (int jcluster = jclusterFirst; jcluster <= jclusterLast; jcluster++)
1064 /* Store cj and the interaction mask */
1066 cjEntry.cj = jGrid.cellOffset() + jcluster;
1067 cjEntry.excl = get_imask(excludeSubDiagonal, icluster, jcluster);
1068 nbl->cj.push_back(cjEntry);
1070 /* Increase the closing index in the i list */
1071 nbl->ci.back().cj_ind_end = nbl->cj.size();
1075 #ifdef GMX_NBNXN_SIMD_4XN
1076 #include "gromacs/nbnxm/pairlist_simd_4xm.h"
1078 #ifdef GMX_NBNXN_SIMD_2XNN
1079 #include "gromacs/nbnxm/pairlist_simd_2xmm.h"
1082 /* Plain C or SIMD4 code for making a pair list of super-cell sci vs scj.
1083 * Checks bounding box distances and possibly atom pair distances.
1085 static void make_cluster_list_supersub(const Grid &iGrid,
1087 NbnxnPairlistGpu *nbl,
1090 const bool excludeSubDiagonal,
1095 int *numDistanceChecks)
1097 NbnxnPairlistGpuWork &work = *nbl->work;
1100 const float *pbb_ci = work.iSuperClusterData.bbPacked.data();
1102 const BoundingBox *bb_ci = work.iSuperClusterData.bb.data();
1105 assert(c_nbnxnGpuClusterSize == iGrid.geometry().numAtomsICluster);
1106 assert(c_nbnxnGpuClusterSize == jGrid.geometry().numAtomsICluster);
1108 /* We generate the pairlist mainly based on bounding-box distances
1109 * and do atom pair distance based pruning on the GPU.
1110 * Only if a j-group contains a single cluster-pair, we try to prune
1111 * that pair based on atom distances on the CPU to avoid empty j-groups.
1113 #define PRUNE_LIST_CPU_ONE 1
1114 #define PRUNE_LIST_CPU_ALL 0
1116 #if PRUNE_LIST_CPU_ONE
1120 float *d2l = work.distanceBuffer.data();
1122 for (int subc = 0; subc < jGrid.numClustersPerCell()[scj]; subc++)
1124 const int cj4_ind = work.cj_ind/c_nbnxnGpuJgroupSize;
1125 const int cj_offset = work.cj_ind - cj4_ind*c_nbnxnGpuJgroupSize;
1126 const int cj = scj*c_gpuNumClusterPerCell + subc;
1128 const int cj_gl = jGrid.cellOffset()*c_gpuNumClusterPerCell + cj;
1131 if (excludeSubDiagonal && sci == scj)
1137 ci1 = iGrid.numClustersPerCell()[sci];
1141 /* Determine all ci1 bb distances in one call with SIMD4 */
1142 const int offset = packedBoundingBoxesIndex(cj) + (cj & (c_packedBoundingBoxesDimSize - 1));
1143 clusterBoundingBoxDistance2_xxxx_simd4(jGrid.packedBoundingBoxes().data() + offset,
1145 *numDistanceChecks += c_nbnxnGpuClusterSize*2;
1149 unsigned int imask = 0;
1150 /* We use a fixed upper-bound instead of ci1 to help optimization */
1151 for (int ci = 0; ci < c_gpuNumClusterPerCell; ci++)
1159 /* Determine the bb distance between ci and cj */
1160 d2l[ci] = clusterBoundingBoxDistance2(bb_ci[ci], jGrid.jBoundingBoxes()[cj]);
1161 *numDistanceChecks += 2;
1165 #if PRUNE_LIST_CPU_ALL
1166 /* Check if the distance is within the distance where
1167 * we use only the bounding box distance rbb,
1168 * or within the cut-off and there is at least one atom pair
1169 * within the cut-off. This check is very costly.
1171 *numDistanceChecks += c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize;
1174 clusterpair_in_range(work, ci, cj_gl, stride, x, rlist2)))
1176 /* Check if the distance between the two bounding boxes
1177 * in within the pair-list cut-off.
1182 /* Flag this i-subcell to be taken into account */
1183 imask |= (1U << (cj_offset*c_gpuNumClusterPerCell + ci));
1185 #if PRUNE_LIST_CPU_ONE
1193 #if PRUNE_LIST_CPU_ONE
1194 /* If we only found 1 pair, check if any atoms are actually
1195 * within the cut-off, so we could get rid of it.
1197 if (npair == 1 && d2l[ci_last] >= rbb2 &&
1198 !clusterpair_in_range(work, ci_last, cj_gl, stride, x, rlist2))
1200 imask &= ~(1U << (cj_offset*c_gpuNumClusterPerCell + ci_last));
1207 /* We have at least one cluster pair: add a j-entry */
1208 if (static_cast<size_t>(cj4_ind) == nbl->cj4.size())
1210 nbl->cj4.resize(nbl->cj4.size() + 1);
1212 nbnxn_cj4_t *cj4 = &nbl->cj4[cj4_ind];
1214 cj4->cj[cj_offset] = cj_gl;
1216 /* Set the exclusions for the ci==sj entry.
1217 * Here we don't bother to check if this entry is actually flagged,
1218 * as it will nearly always be in the list.
1220 if (excludeSubDiagonal && sci == scj)
1222 set_self_and_newton_excls_supersub(nbl, cj4_ind, cj_offset, subc);
1225 /* Copy the cluster interaction mask to the list */
1226 for (int w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
1228 cj4->imei[w].imask |= imask;
1231 nbl->work->cj_ind++;
1233 /* Keep the count */
1234 nbl->nci_tot += npair;
1236 /* Increase the closing index in i super-cell list */
1237 nbl->sci.back().cj4_ind_end =
1238 (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize;
1243 /* Returns how many contiguous j-clusters we have starting in the i-list */
1244 template <typename CjListType>
1245 static int numContiguousJClusters(const int cjIndexStart,
1246 const int cjIndexEnd,
1247 gmx::ArrayRef<const CjListType> cjList)
1249 const int firstJCluster = nblCj(cjList, cjIndexStart);
1251 int numContiguous = 0;
1253 while (cjIndexStart + numContiguous < cjIndexEnd &&
1254 nblCj(cjList, cjIndexStart + numContiguous) == firstJCluster + numContiguous)
1259 return numContiguous;
1263 * \brief Helper struct for efficient searching for excluded atoms in a j-list
1267 /*! \brief Constructs a j-list range from \p cjList with the given index range */
1268 template <typename CjListType>
1269 JListRanges(int cjIndexStart,
1271 gmx::ArrayRef<const CjListType> cjList);
1273 int cjIndexStart; //!< The start index in the j-list
1274 int cjIndexEnd; //!< The end index in the j-list
1275 int cjFirst; //!< The j-cluster with index cjIndexStart
1276 int cjLast; //!< The j-cluster with index cjIndexEnd-1
1277 int numDirect; //!< Up to cjIndexStart+numDirect the j-clusters are cjFirst + the index offset
1281 template <typename CjListType>
1282 JListRanges::JListRanges(int cjIndexStart,
1284 gmx::ArrayRef<const CjListType> cjList) :
1285 cjIndexStart(cjIndexStart),
1286 cjIndexEnd(cjIndexEnd)
1288 GMX_ASSERT(cjIndexEnd > cjIndexStart, "JListRanges should only be called with non-empty lists");
1290 cjFirst = nblCj(cjList, cjIndexStart);
1291 cjLast = nblCj(cjList, cjIndexEnd - 1);
1293 /* Determine how many contiguous j-cells we have starting
1294 * from the first i-cell. This number can be used to directly
1295 * calculate j-cell indices for excluded atoms.
1297 numDirect = numContiguousJClusters(cjIndexStart, cjIndexEnd, cjList);
1301 /* Return the index of \p jCluster in the given range or -1 when not present
1303 * Note: This code is executed very often and therefore performance is
1304 * important. It should be inlined and fully optimized.
1306 template <typename CjListType>
1308 findJClusterInJList(int jCluster,
1309 const JListRanges &ranges,
1310 gmx::ArrayRef<const CjListType> cjList)
1314 if (jCluster < ranges.cjFirst + ranges.numDirect)
1316 /* We can calculate the index directly using the offset */
1317 index = ranges.cjIndexStart + jCluster - ranges.cjFirst;
1321 /* Search for jCluster using bisection */
1323 int rangeStart = ranges.cjIndexStart + ranges.numDirect;
1324 int rangeEnd = ranges.cjIndexEnd;
1326 while (index == -1 && rangeStart < rangeEnd)
1328 rangeMiddle = (rangeStart + rangeEnd) >> 1;
1330 const int clusterMiddle = nblCj(cjList, rangeMiddle);
1332 if (jCluster == clusterMiddle)
1334 index = rangeMiddle;
1336 else if (jCluster < clusterMiddle)
1338 rangeEnd = rangeMiddle;
1342 rangeStart = rangeMiddle + 1;
1350 // TODO: Get rid of the two functions below by renaming sci to ci (or something better)
1352 /* Return the i-entry in the list we are currently operating on */
1353 static nbnxn_ci_t *getOpenIEntry(NbnxnPairlistCpu *nbl)
1355 return &nbl->ci.back();
1358 /* Return the i-entry in the list we are currently operating on */
1359 static nbnxn_sci_t *getOpenIEntry(NbnxnPairlistGpu *nbl)
1361 return &nbl->sci.back();
1364 /* Set all atom-pair exclusions for a simple type list i-entry
1366 * Set all atom-pair exclusions from the topology stored in exclusions
1367 * as masks in the pair-list for simple list entry iEntry.
1370 setExclusionsForIEntry(const Nbnxm::GridSet &gridSet,
1371 NbnxnPairlistCpu *nbl,
1372 gmx_bool diagRemoved,
1374 const nbnxn_ci_t &iEntry,
1375 const t_blocka &exclusions)
1377 if (iEntry.cj_ind_end == iEntry.cj_ind_start)
1379 /* Empty list: no exclusions */
1383 const JListRanges ranges(iEntry.cj_ind_start, iEntry.cj_ind_end, gmx::makeConstArrayRef(nbl->cj));
1385 const int iCluster = iEntry.ci;
1387 gmx::ArrayRef<const int> cell = gridSet.cells();
1388 gmx::ArrayRef<const int> atomIndices = gridSet.atomIndices();
1390 /* Loop over the atoms in the i-cluster */
1391 for (int i = 0; i < nbl->na_ci; i++)
1393 const int iIndex = iCluster*nbl->na_ci + i;
1394 const int iAtom = atomIndices[iIndex];
1397 /* Loop over the topology-based exclusions for this i-atom */
1398 for (int exclIndex = exclusions.index[iAtom]; exclIndex < exclusions.index[iAtom + 1]; exclIndex++)
1400 const int jAtom = exclusions.a[exclIndex];
1404 /* The self exclusion are already set, save some time */
1408 /* Get the index of the j-atom in the nbnxn atom data */
1409 const int jIndex = cell[jAtom];
1411 /* Without shifts we only calculate interactions j>i
1412 * for one-way pair-lists.
1414 if (diagRemoved && jIndex <= iIndex)
1419 const int jCluster = (jIndex >> na_cj_2log);
1421 /* Could the cluster se be in our list? */
1422 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
1425 findJClusterInJList(jCluster, ranges,
1426 gmx::makeConstArrayRef(nbl->cj));
1430 /* We found an exclusion, clear the corresponding
1433 const int innerJ = jIndex - (jCluster << na_cj_2log);
1435 nbl->cj[index].excl &= ~(1U << ((i << na_cj_2log) + innerJ));
1443 /* Add a new i-entry to the FEP list and copy the i-properties */
1444 static inline void fep_list_new_nri_copy(t_nblist *nlist)
1446 /* Add a new i-entry */
1449 assert(nlist->nri < nlist->maxnri);
1451 /* Duplicate the last i-entry, except for jindex, which continues */
1452 nlist->iinr[nlist->nri] = nlist->iinr[nlist->nri-1];
1453 nlist->shift[nlist->nri] = nlist->shift[nlist->nri-1];
1454 nlist->gid[nlist->nri] = nlist->gid[nlist->nri-1];
1455 nlist->jindex[nlist->nri] = nlist->nrj;
1458 /* Rellocate FEP list for size nl->maxnri, TODO: replace by C++ */
1459 static void reallocate_nblist(t_nblist *nl)
1463 fprintf(debug, "reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, type=%d), maxnri=%d\n",
1464 nl->ielec, nl->ivdw, nl->igeometry, nl->type, nl->maxnri);
1466 srenew(nl->iinr, nl->maxnri);
1467 srenew(nl->gid, nl->maxnri);
1468 srenew(nl->shift, nl->maxnri);
1469 srenew(nl->jindex, nl->maxnri+1);
1472 /* For load balancing of the free-energy lists over threads, we set
1473 * the maximum nrj size of an i-entry to 40. This leads to good
1474 * load balancing in the worst case scenario of a single perturbed
1475 * particle on 16 threads, while not introducing significant overhead.
1476 * Note that half of the perturbed pairs will anyhow end up in very small lists,
1477 * since non perturbed i-particles will see few perturbed j-particles).
1479 const int max_nrj_fep = 40;
1481 /* Exclude the perturbed pairs from the Verlet list. This is only done to avoid
1482 * singularities for overlapping particles (0/0), since the charges and
1483 * LJ parameters have been zeroed in the nbnxn data structure.
1484 * Simultaneously make a group pair list for the perturbed pairs.
1486 static void make_fep_list(gmx::ArrayRef<const int> atomIndices,
1487 const nbnxn_atomdata_t *nbat,
1488 NbnxnPairlistCpu *nbl,
1489 gmx_bool bDiagRemoved,
1491 real gmx_unused shx,
1492 real gmx_unused shy,
1493 real gmx_unused shz,
1494 real gmx_unused rlist_fep2,
1499 int ci, cj_ind_start, cj_ind_end, cja, cjr;
1501 int gid_i = 0, gid_j, gid;
1502 int egp_shift, egp_mask;
1504 int ind_i, ind_j, ai, aj;
1506 gmx_bool bFEP_i, bFEP_i_all;
1508 if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
1516 cj_ind_start = nbl_ci->cj_ind_start;
1517 cj_ind_end = nbl_ci->cj_ind_end;
1519 /* In worst case we have alternating energy groups
1520 * and create #atom-pair lists, which means we need the size
1521 * of a cluster pair (na_ci*na_cj) times the number of cj's.
1523 nri_max = nbl->na_ci*nbl->na_cj*(cj_ind_end - cj_ind_start);
1524 if (nlist->nri + nri_max > nlist->maxnri)
1526 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1527 reallocate_nblist(nlist);
1530 const int numAtomsJCluster = jGrid.geometry().numAtomsJCluster;
1532 const nbnxn_atomdata_t::Params &nbatParams = nbat->params();
1534 const int ngid = nbatParams.nenergrp;
1536 /* TODO: Consider adding a check in grompp and changing this to an assert */
1537 const int numBitsInEnergyGroupIdsForAtomsInJCluster = sizeof(gid_cj)*8;
1538 if (ngid*numAtomsJCluster > numBitsInEnergyGroupIdsForAtomsInJCluster)
1540 gmx_fatal(FARGS, "The Verlet scheme with %dx%d kernels and free-energy only supports up to %zu energy groups",
1541 iGrid.geometry().numAtomsICluster, numAtomsJCluster,
1542 (sizeof(gid_cj)*8)/numAtomsJCluster);
1545 egp_shift = nbatParams.neg_2log;
1546 egp_mask = (1 << egp_shift) - 1;
1548 /* Loop over the atoms in the i sub-cell */
1550 for (int i = 0; i < nbl->na_ci; i++)
1552 ind_i = ci*nbl->na_ci + i;
1553 ai = atomIndices[ind_i];
1557 nlist->jindex[nri+1] = nlist->jindex[nri];
1558 nlist->iinr[nri] = ai;
1559 /* The actual energy group pair index is set later */
1560 nlist->gid[nri] = 0;
1561 nlist->shift[nri] = nbl_ci->shift & NBNXN_CI_SHIFT;
1563 bFEP_i = iGrid.atomIsPerturbed(ci - iGrid.cellOffset(), i);
1565 bFEP_i_all = bFEP_i_all && bFEP_i;
1567 if (nlist->nrj + (cj_ind_end - cj_ind_start)*nbl->na_cj > nlist->maxnrj)
1569 nlist->maxnrj = over_alloc_small(nlist->nrj + (cj_ind_end - cj_ind_start)*nbl->na_cj);
1570 srenew(nlist->jjnr, nlist->maxnrj);
1571 srenew(nlist->excl_fep, nlist->maxnrj);
1576 gid_i = (nbatParams.energrp[ci] >> (egp_shift*i)) & egp_mask;
1579 for (int cj_ind = cj_ind_start; cj_ind < cj_ind_end; cj_ind++)
1581 unsigned int fep_cj;
1583 cja = nbl->cj[cj_ind].cj;
1585 if (numAtomsJCluster == jGrid.geometry().numAtomsICluster)
1587 cjr = cja - jGrid.cellOffset();
1588 fep_cj = jGrid.fepBits(cjr);
1591 gid_cj = nbatParams.energrp[cja];
1594 else if (2*numAtomsJCluster == jGrid.geometry().numAtomsICluster)
1596 cjr = cja - jGrid.cellOffset()*2;
1597 /* Extract half of the ci fep/energrp mask */
1598 fep_cj = (jGrid.fepBits(cjr >> 1) >> ((cjr & 1)*numAtomsJCluster)) & ((1 << numAtomsJCluster) - 1);
1601 gid_cj = nbatParams.energrp[cja >> 1] >> ((cja & 1)*numAtomsJCluster*egp_shift) & ((1 << (numAtomsJCluster*egp_shift)) - 1);
1606 cjr = cja - (jGrid.cellOffset() >> 1);
1607 /* Combine two ci fep masks/energrp */
1608 fep_cj = jGrid.fepBits(cjr*2) + (jGrid.fepBits(cjr*2 + 1) << jGrid.geometry().numAtomsICluster);
1611 gid_cj = nbatParams.energrp[cja*2] + (nbatParams.energrp[cja*2+1] << (jGrid.geometry().numAtomsICluster*egp_shift));
1615 if (bFEP_i || fep_cj != 0)
1617 for (int j = 0; j < nbl->na_cj; j++)
1619 /* Is this interaction perturbed and not excluded? */
1620 ind_j = cja*nbl->na_cj + j;
1621 aj = atomIndices[ind_j];
1623 (bFEP_i || (fep_cj & (1 << j))) &&
1624 (!bDiagRemoved || ind_j >= ind_i))
1628 gid_j = (gid_cj >> (j*egp_shift)) & egp_mask;
1629 gid = GID(gid_i, gid_j, ngid);
1631 if (nlist->nrj > nlist->jindex[nri] &&
1632 nlist->gid[nri] != gid)
1634 /* Energy group pair changed: new list */
1635 fep_list_new_nri_copy(nlist);
1638 nlist->gid[nri] = gid;
1641 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1643 fep_list_new_nri_copy(nlist);
1647 /* Add it to the FEP list */
1648 nlist->jjnr[nlist->nrj] = aj;
1649 nlist->excl_fep[nlist->nrj] = (nbl->cj[cj_ind].excl >> (i*nbl->na_cj + j)) & 1;
1652 /* Exclude it from the normal list.
1653 * Note that the charge has been set to zero,
1654 * but we need to avoid 0/0, as perturbed atoms
1655 * can be on top of each other.
1657 nbl->cj[cj_ind].excl &= ~(1U << (i*nbl->na_cj + j));
1663 if (nlist->nrj > nlist->jindex[nri])
1665 /* Actually add this new, non-empty, list */
1667 nlist->jindex[nlist->nri] = nlist->nrj;
1674 /* All interactions are perturbed, we can skip this entry */
1675 nbl_ci->cj_ind_end = cj_ind_start;
1676 nbl->ncjInUse -= cj_ind_end - cj_ind_start;
1680 /* Return the index of atom a within a cluster */
1681 static inline int cj_mod_cj4(int cj)
1683 return cj & (c_nbnxnGpuJgroupSize - 1);
1686 /* Convert a j-cluster to a cj4 group */
1687 static inline int cj_to_cj4(int cj)
1689 return cj/c_nbnxnGpuJgroupSize;
1692 /* Return the index of an j-atom within a warp */
1693 static inline int a_mod_wj(int a)
1695 return a & (c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit - 1);
1698 /* As make_fep_list above, but for super/sub lists. */
1699 static void make_fep_list(gmx::ArrayRef<const int> atomIndices,
1700 const nbnxn_atomdata_t *nbat,
1701 NbnxnPairlistGpu *nbl,
1702 gmx_bool bDiagRemoved,
1703 const nbnxn_sci_t *nbl_sci,
1714 int ind_i, ind_j, ai, aj;
1718 const nbnxn_cj4_t *cj4;
1720 const int numJClusterGroups = nbl_sci->numJClusterGroups();
1721 if (numJClusterGroups == 0)
1727 const int sci = nbl_sci->sci;
1729 const int cj4_ind_start = nbl_sci->cj4_ind_start;
1730 const int cj4_ind_end = nbl_sci->cj4_ind_end;
1732 /* Here we process one super-cell, max #atoms na_sc, versus a list
1733 * cj4 entries, each with max c_nbnxnGpuJgroupSize cj's, each
1734 * of size na_cj atoms.
1735 * On the GPU we don't support energy groups (yet).
1736 * So for each of the na_sc i-atoms, we need max one FEP list
1737 * for each max_nrj_fep j-atoms.
1739 nri_max = nbl->na_sc*nbl->na_cj*(1 + (numJClusterGroups*c_nbnxnGpuJgroupSize)/max_nrj_fep);
1740 if (nlist->nri + nri_max > nlist->maxnri)
1742 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1743 reallocate_nblist(nlist);
1746 /* Loop over the atoms in the i super-cluster */
1747 for (int c = 0; c < c_gpuNumClusterPerCell; c++)
1749 c_abs = sci*c_gpuNumClusterPerCell + c;
1751 for (int i = 0; i < nbl->na_ci; i++)
1753 ind_i = c_abs*nbl->na_ci + i;
1754 ai = atomIndices[ind_i];
1758 nlist->jindex[nri+1] = nlist->jindex[nri];
1759 nlist->iinr[nri] = ai;
1760 /* With GPUs, energy groups are not supported */
1761 nlist->gid[nri] = 0;
1762 nlist->shift[nri] = nbl_sci->shift & NBNXN_CI_SHIFT;
1764 bFEP_i = iGrid.atomIsPerturbed(c_abs - iGrid.cellOffset()*c_gpuNumClusterPerCell, i);
1766 xi = nbat->x()[ind_i*nbat->xstride+XX] + shx;
1767 yi = nbat->x()[ind_i*nbat->xstride+YY] + shy;
1768 zi = nbat->x()[ind_i*nbat->xstride+ZZ] + shz;
1770 const int nrjMax = nlist->nrj + numJClusterGroups*c_nbnxnGpuJgroupSize*nbl->na_cj;
1771 if (nrjMax > nlist->maxnrj)
1773 nlist->maxnrj = over_alloc_small(nrjMax);
1774 srenew(nlist->jjnr, nlist->maxnrj);
1775 srenew(nlist->excl_fep, nlist->maxnrj);
1778 for (int cj4_ind = cj4_ind_start; cj4_ind < cj4_ind_end; cj4_ind++)
1780 cj4 = &nbl->cj4[cj4_ind];
1782 for (int gcj = 0; gcj < c_nbnxnGpuJgroupSize; gcj++)
1784 if ((cj4->imei[0].imask & (1U << (gcj*c_gpuNumClusterPerCell + c))) == 0)
1786 /* Skip this ci for this cj */
1791 cj4->cj[gcj] - jGrid.cellOffset()*c_gpuNumClusterPerCell;
1793 if (bFEP_i || jGrid.clusterIsPerturbed(cjr))
1795 for (int j = 0; j < nbl->na_cj; j++)
1797 /* Is this interaction perturbed and not excluded? */
1798 ind_j = (jGrid.cellOffset()*c_gpuNumClusterPerCell + cjr)*nbl->na_cj + j;
1799 aj = atomIndices[ind_j];
1801 (bFEP_i || jGrid.atomIsPerturbed(cjr, j)) &&
1802 (!bDiagRemoved || ind_j >= ind_i))
1805 unsigned int excl_bit;
1808 const int jHalf = j/(c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit);
1809 nbnxn_excl_t *excl =
1810 get_exclusion_mask(nbl, cj4_ind, jHalf);
1812 excl_pair = a_mod_wj(j)*nbl->na_ci + i;
1813 excl_bit = (1U << (gcj*c_gpuNumClusterPerCell + c));
1815 dx = nbat->x()[ind_j*nbat->xstride+XX] - xi;
1816 dy = nbat->x()[ind_j*nbat->xstride+YY] - yi;
1817 dz = nbat->x()[ind_j*nbat->xstride+ZZ] - zi;
1819 /* The unpruned GPU list has more than 2/3
1820 * of the atom pairs beyond rlist. Using
1821 * this list will cause a lot of overhead
1822 * in the CPU FEP kernels, especially
1823 * relative to the fast GPU kernels.
1824 * So we prune the FEP list here.
1826 if (dx*dx + dy*dy + dz*dz < rlist_fep2)
1828 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1830 fep_list_new_nri_copy(nlist);
1834 /* Add it to the FEP list */
1835 nlist->jjnr[nlist->nrj] = aj;
1836 nlist->excl_fep[nlist->nrj] = (excl->pair[excl_pair] & excl_bit) ? 1 : 0;
1840 /* Exclude it from the normal list.
1841 * Note that the charge and LJ parameters have
1842 * been set to zero, but we need to avoid 0/0,
1843 * as perturbed atoms can be on top of each other.
1845 excl->pair[excl_pair] &= ~excl_bit;
1849 /* Note that we could mask out this pair in imask
1850 * if all i- and/or all j-particles are perturbed.
1851 * But since the perturbed pairs on the CPU will
1852 * take an order of magnitude more time, the GPU
1853 * will finish before the CPU and there is no gain.
1859 if (nlist->nrj > nlist->jindex[nri])
1861 /* Actually add this new, non-empty, list */
1863 nlist->jindex[nlist->nri] = nlist->nrj;
1870 /* Set all atom-pair exclusions for a GPU type list i-entry
1872 * Sets all atom-pair exclusions from the topology stored in exclusions
1873 * as masks in the pair-list for i-super-cluster list entry iEntry.
1876 setExclusionsForIEntry(const Nbnxm::GridSet &gridSet,
1877 NbnxnPairlistGpu *nbl,
1878 gmx_bool diagRemoved,
1879 int gmx_unused na_cj_2log,
1880 const nbnxn_sci_t &iEntry,
1881 const t_blocka &exclusions)
1883 if (iEntry.numJClusterGroups() == 0)
1889 /* Set the search ranges using start and end j-cluster indices.
1890 * Note that here we can not use cj4_ind_end, since the last cj4
1891 * can be only partially filled, so we use cj_ind.
1893 const JListRanges ranges(iEntry.cj4_ind_start*c_nbnxnGpuJgroupSize,
1895 gmx::makeConstArrayRef(nbl->cj4));
1897 GMX_ASSERT(nbl->na_ci == c_nbnxnGpuClusterSize, "na_ci should match the GPU cluster size");
1898 constexpr int c_clusterSize = c_nbnxnGpuClusterSize;
1899 constexpr int c_superClusterSize = c_nbnxnGpuNumClusterPerSupercluster*c_nbnxnGpuClusterSize;
1901 const int iSuperCluster = iEntry.sci;
1903 gmx::ArrayRef<const int> atomIndices = gridSet.atomIndices();
1904 gmx::ArrayRef<const int> cell = gridSet.cells();
1906 /* Loop over the atoms in the i super-cluster */
1907 for (int i = 0; i < c_superClusterSize; i++)
1909 const int iIndex = iSuperCluster*c_superClusterSize + i;
1910 const int iAtom = atomIndices[iIndex];
1913 const int iCluster = i/c_clusterSize;
1915 /* Loop over the topology-based exclusions for this i-atom */
1916 for (int exclIndex = exclusions.index[iAtom]; exclIndex < exclusions.index[iAtom + 1]; exclIndex++)
1918 const int jAtom = exclusions.a[exclIndex];
1922 /* The self exclusions are already set, save some time */
1926 /* Get the index of the j-atom in the nbnxn atom data */
1927 const int jIndex = cell[jAtom];
1929 /* Without shifts we only calculate interactions j>i
1930 * for one-way pair-lists.
1932 /* NOTE: We would like to use iIndex on the right hand side,
1933 * but that makes this routine 25% slower with gcc6/7.
1934 * Even using c_superClusterSize makes it slower.
1935 * Either of these changes triggers peeling of the exclIndex
1936 * loop, which apparently leads to far less efficient code.
1938 if (diagRemoved && jIndex <= iSuperCluster*nbl->na_sc + i)
1943 const int jCluster = jIndex/c_clusterSize;
1945 /* Check whether the cluster is in our list? */
1946 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
1949 findJClusterInJList(jCluster, ranges,
1950 gmx::makeConstArrayRef(nbl->cj4));
1954 /* We found an exclusion, clear the corresponding
1957 const unsigned int pairMask = (1U << (cj_mod_cj4(index)*c_gpuNumClusterPerCell + iCluster));
1958 /* Check if the i-cluster interacts with the j-cluster */
1959 if (nbl_imask0(nbl, index) & pairMask)
1961 const int innerI = (i & (c_clusterSize - 1));
1962 const int innerJ = (jIndex & (c_clusterSize - 1));
1964 /* Determine which j-half (CUDA warp) we are in */
1965 const int jHalf = innerJ/(c_clusterSize/c_nbnxnGpuClusterpairSplit);
1967 nbnxn_excl_t *interactionMask =
1968 get_exclusion_mask(nbl, cj_to_cj4(index), jHalf);
1970 interactionMask->pair[a_mod_wj(innerJ)*c_clusterSize + innerI] &= ~pairMask;
1979 /* Make a new ci entry at the back of nbl->ci */
1980 static void addNewIEntry(NbnxnPairlistCpu *nbl, int ci, int shift, int flags)
1984 ciEntry.shift = shift;
1985 /* Store the interaction flags along with the shift */
1986 ciEntry.shift |= flags;
1987 ciEntry.cj_ind_start = nbl->cj.size();
1988 ciEntry.cj_ind_end = nbl->cj.size();
1989 nbl->ci.push_back(ciEntry);
1992 /* Make a new sci entry at index nbl->nsci */
1993 static void addNewIEntry(NbnxnPairlistGpu *nbl, int sci, int shift, int gmx_unused flags)
1995 nbnxn_sci_t sciEntry;
1997 sciEntry.shift = shift;
1998 sciEntry.cj4_ind_start = nbl->cj4.size();
1999 sciEntry.cj4_ind_end = nbl->cj4.size();
2001 nbl->sci.push_back(sciEntry);
2004 /* Sort the simple j-list cj on exclusions.
2005 * Entries with exclusions will all be sorted to the beginning of the list.
2007 static void sort_cj_excl(nbnxn_cj_t *cj, int ncj,
2008 NbnxnPairlistCpuWork *work)
2010 work->cj.resize(ncj);
2012 /* Make a list of the j-cells involving exclusions */
2014 for (int j = 0; j < ncj; j++)
2016 if (cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
2018 work->cj[jnew++] = cj[j];
2021 /* Check if there are exclusions at all or not just the first entry */
2022 if (!((jnew == 0) ||
2023 (jnew == 1 && cj[0].excl != NBNXN_INTERACTION_MASK_ALL)))
2025 for (int j = 0; j < ncj; j++)
2027 if (cj[j].excl == NBNXN_INTERACTION_MASK_ALL)
2029 work->cj[jnew++] = cj[j];
2032 for (int j = 0; j < ncj; j++)
2034 cj[j] = work->cj[j];
2039 /* Close this simple list i entry */
2040 static void closeIEntry(NbnxnPairlistCpu *nbl,
2041 int gmx_unused sp_max_av,
2042 gmx_bool gmx_unused progBal,
2043 float gmx_unused nsp_tot_est,
2044 int gmx_unused thread,
2045 int gmx_unused nthread)
2047 nbnxn_ci_t &ciEntry = nbl->ci.back();
2049 /* All content of the new ci entry have already been filled correctly,
2050 * we only need to sort and increase counts or remove the entry when empty.
2052 const int jlen = ciEntry.cj_ind_end - ciEntry.cj_ind_start;
2055 sort_cj_excl(nbl->cj.data() + ciEntry.cj_ind_start, jlen, nbl->work.get());
2057 /* The counts below are used for non-bonded pair/flop counts
2058 * and should therefore match the available kernel setups.
2060 if (!(ciEntry.shift & NBNXN_CI_DO_COUL(0)))
2062 nbl->work->ncj_noq += jlen;
2064 else if ((ciEntry.shift & NBNXN_CI_HALF_LJ(0)) ||
2065 !(ciEntry.shift & NBNXN_CI_DO_LJ(0)))
2067 nbl->work->ncj_hlj += jlen;
2072 /* Entry is empty: remove it */
2077 /* Split sci entry for load balancing on the GPU.
2078 * Splitting ensures we have enough lists to fully utilize the whole GPU.
2079 * With progBal we generate progressively smaller lists, which improves
2080 * load balancing. As we only know the current count on our own thread,
2081 * we will need to estimate the current total amount of i-entries.
2082 * As the lists get concatenated later, this estimate depends
2083 * both on nthread and our own thread index.
2085 static void split_sci_entry(NbnxnPairlistGpu *nbl,
2087 gmx_bool progBal, float nsp_tot_est,
2088 int thread, int nthread)
2096 /* Estimate the total numbers of ci's of the nblist combined
2097 * over all threads using the target number of ci's.
2099 nsp_est = (nsp_tot_est*thread)/nthread + nbl->nci_tot;
2101 /* The first ci blocks should be larger, to avoid overhead.
2102 * The last ci blocks should be smaller, to improve load balancing.
2103 * The factor 3/2 makes the first block 3/2 times the target average
2104 * and ensures that the total number of blocks end up equal to
2105 * that of equally sized blocks of size nsp_target_av.
2107 nsp_max = static_cast<int>(nsp_target_av*(nsp_tot_est*1.5/(nsp_est + nsp_tot_est)));
2111 nsp_max = nsp_target_av;
2114 const int cj4_start = nbl->sci.back().cj4_ind_start;
2115 const int cj4_end = nbl->sci.back().cj4_ind_end;
2116 const int j4len = cj4_end - cj4_start;
2118 if (j4len > 1 && j4len*c_gpuNumClusterPerCell*c_nbnxnGpuJgroupSize > nsp_max)
2120 /* Modify the last ci entry and process the cj4's again */
2126 for (int cj4 = cj4_start; cj4 < cj4_end; cj4++)
2128 int nsp_cj4_p = nsp_cj4;
2129 /* Count the number of cluster pairs in this cj4 group */
2131 for (int p = 0; p < c_gpuNumClusterPerCell*c_nbnxnGpuJgroupSize; p++)
2133 nsp_cj4 += (nbl->cj4[cj4].imei[0].imask >> p) & 1;
2136 /* If adding the current cj4 with nsp_cj4 pairs get us further
2137 * away from our target nsp_max, split the list before this cj4.
2139 if (nsp > 0 && nsp_max - nsp < nsp + nsp_cj4 - nsp_max)
2141 /* Split the list at cj4 */
2142 nbl->sci.back().cj4_ind_end = cj4;
2143 /* Create a new sci entry */
2145 sciNew.sci = nbl->sci.back().sci;
2146 sciNew.shift = nbl->sci.back().shift;
2147 sciNew.cj4_ind_start = cj4;
2148 nbl->sci.push_back(sciNew);
2151 nsp_cj4_e = nsp_cj4_p;
2157 /* Put the remaining cj4's in the last sci entry */
2158 nbl->sci.back().cj4_ind_end = cj4_end;
2160 /* Possibly balance out the last two sci's
2161 * by moving the last cj4 of the second last sci.
2163 if (nsp_sci - nsp_cj4_e >= nsp + nsp_cj4_e)
2165 GMX_ASSERT(nbl->sci.size() >= 2, "We expect at least two elements");
2166 nbl->sci[nbl->sci.size() - 2].cj4_ind_end--;
2167 nbl->sci[nbl->sci.size() - 1].cj4_ind_start--;
2172 /* Clost this super/sub list i entry */
2173 static void closeIEntry(NbnxnPairlistGpu *nbl,
2175 gmx_bool progBal, float nsp_tot_est,
2176 int thread, int nthread)
2178 nbnxn_sci_t &sciEntry = *getOpenIEntry(nbl);
2180 /* All content of the new ci entry have already been filled correctly,
2181 * we only need to, potentially, split or remove the entry when empty.
2183 int j4len = sciEntry.numJClusterGroups();
2186 /* We can only have complete blocks of 4 j-entries in a list,
2187 * so round the count up before closing.
2189 int ncj4 = (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize;
2190 nbl->work->cj_ind = ncj4*c_nbnxnGpuJgroupSize;
2194 /* Measure the size of the new entry and potentially split it */
2195 split_sci_entry(nbl, nsp_max_av, progBal, nsp_tot_est,
2201 /* Entry is empty: remove it */
2202 nbl->sci.pop_back();
2206 /* Syncs the working array before adding another grid pair to the GPU list */
2207 static void sync_work(NbnxnPairlistCpu gmx_unused *nbl)
2211 /* Syncs the working array before adding another grid pair to the GPU list */
2212 static void sync_work(NbnxnPairlistGpu *nbl)
2214 nbl->work->cj_ind = nbl->cj4.size()*c_nbnxnGpuJgroupSize;
2217 /* Clears an NbnxnPairlistCpu data structure */
2218 static void clear_pairlist(NbnxnPairlistCpu *nbl)
2224 nbl->ciOuter.clear();
2225 nbl->cjOuter.clear();
2227 nbl->work->ncj_noq = 0;
2228 nbl->work->ncj_hlj = 0;
2231 /* Clears an NbnxnPairlistGpu data structure */
2232 static void clear_pairlist(NbnxnPairlistGpu *nbl)
2236 nbl->excl.resize(1);
2240 /* Clears a group scheme pair list */
2241 static void clear_pairlist_fep(t_nblist *nl)
2245 if (nl->jindex == nullptr)
2247 snew(nl->jindex, 1);
2252 /* Sets a simple list i-cell bounding box, including PBC shift */
2253 static inline void set_icell_bb_simple(gmx::ArrayRef<const BoundingBox> bb,
2255 real shx, real shy, real shz,
2258 bb_ci->lower.x = bb[ci].lower.x + shx;
2259 bb_ci->lower.y = bb[ci].lower.y + shy;
2260 bb_ci->lower.z = bb[ci].lower.z + shz;
2261 bb_ci->upper.x = bb[ci].upper.x + shx;
2262 bb_ci->upper.y = bb[ci].upper.y + shy;
2263 bb_ci->upper.z = bb[ci].upper.z + shz;
2266 /* Sets a simple list i-cell bounding box, including PBC shift */
2267 static inline void set_icell_bb(const Grid &iGrid,
2269 real shx, real shy, real shz,
2270 NbnxnPairlistCpuWork *work)
2272 set_icell_bb_simple(iGrid.iBoundingBoxes(), ci, shx, shy, shz,
2273 &work->iClusterData.bb[0]);
2277 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2278 static void set_icell_bbxxxx_supersub(gmx::ArrayRef<const float> bb,
2280 real shx, real shy, real shz,
2283 constexpr int cellBBStride = packedBoundingBoxesIndex(c_gpuNumClusterPerCell);
2284 constexpr int pbbStride = c_packedBoundingBoxesDimSize;
2285 const int ia = ci*cellBBStride;
2286 for (int m = 0; m < cellBBStride; m += c_packedBoundingBoxesSize)
2288 for (int i = 0; i < pbbStride; i++)
2290 bb_ci[m + 0*pbbStride + i] = bb[ia + m + 0*pbbStride + i] + shx;
2291 bb_ci[m + 1*pbbStride + i] = bb[ia + m + 1*pbbStride + i] + shy;
2292 bb_ci[m + 2*pbbStride + i] = bb[ia + m + 2*pbbStride + i] + shz;
2293 bb_ci[m + 3*pbbStride + i] = bb[ia + m + 3*pbbStride + i] + shx;
2294 bb_ci[m + 4*pbbStride + i] = bb[ia + m + 4*pbbStride + i] + shy;
2295 bb_ci[m + 5*pbbStride + i] = bb[ia + m + 5*pbbStride + i] + shz;
2301 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2302 gmx_unused static void set_icell_bb_supersub(gmx::ArrayRef<const BoundingBox> bb,
2304 real shx, real shy, real shz,
2307 for (int i = 0; i < c_gpuNumClusterPerCell; i++)
2309 set_icell_bb_simple(bb, ci*c_gpuNumClusterPerCell+i,
2315 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2316 gmx_unused static void set_icell_bb(const Grid &iGrid,
2318 real shx, real shy, real shz,
2319 NbnxnPairlistGpuWork *work)
2322 set_icell_bbxxxx_supersub(iGrid.packedBoundingBoxes(), ci, shx, shy, shz,
2323 work->iSuperClusterData.bbPacked.data());
2325 set_icell_bb_supersub(iGrid.iBoundingBoxes(), ci, shx, shy, shz,
2326 work->iSuperClusterData.bb.data());
2330 /* Copies PBC shifted i-cell atom coordinates x,y,z to working array */
2331 static void icell_set_x_simple(int ci,
2332 real shx, real shy, real shz,
2333 int stride, const real *x,
2334 NbnxnPairlistCpuWork::IClusterData *iClusterData)
2336 const int ia = ci*c_nbnxnCpuIClusterSize;
2338 for (int i = 0; i < c_nbnxnCpuIClusterSize; i++)
2340 iClusterData->x[i*STRIDE_XYZ+XX] = x[(ia+i)*stride+XX] + shx;
2341 iClusterData->x[i*STRIDE_XYZ+YY] = x[(ia+i)*stride+YY] + shy;
2342 iClusterData->x[i*STRIDE_XYZ+ZZ] = x[(ia+i)*stride+ZZ] + shz;
2346 static void icell_set_x(int ci,
2347 real shx, real shy, real shz,
2348 int stride, const real *x,
2349 const ClusterDistanceKernelType kernelType,
2350 NbnxnPairlistCpuWork *work)
2355 #ifdef GMX_NBNXN_SIMD_4XN
2356 case ClusterDistanceKernelType::CpuSimd_4xM:
2357 icell_set_x_simd_4xn(ci, shx, shy, shz, stride, x, work);
2360 #ifdef GMX_NBNXN_SIMD_2XNN
2361 case ClusterDistanceKernelType::CpuSimd_2xMM:
2362 icell_set_x_simd_2xnn(ci, shx, shy, shz, stride, x, work);
2366 case ClusterDistanceKernelType::CpuPlainC:
2367 icell_set_x_simple(ci, shx, shy, shz, stride, x, &work->iClusterData);
2370 GMX_ASSERT(false, "Unhandled case");
2375 /* Copies PBC shifted super-cell atom coordinates x,y,z to working array */
2376 static void icell_set_x(int ci,
2377 real shx, real shy, real shz,
2378 int stride, const real *x,
2379 ClusterDistanceKernelType gmx_unused kernelType,
2380 NbnxnPairlistGpuWork *work)
2382 #if !GMX_SIMD4_HAVE_REAL
2384 real * x_ci = work->iSuperClusterData.x.data();
2386 int ia = ci*c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize;
2387 for (int i = 0; i < c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize; i++)
2389 x_ci[i*DIM + XX] = x[(ia+i)*stride + XX] + shx;
2390 x_ci[i*DIM + YY] = x[(ia+i)*stride + YY] + shy;
2391 x_ci[i*DIM + ZZ] = x[(ia+i)*stride + ZZ] + shz;
2394 #else /* !GMX_SIMD4_HAVE_REAL */
2396 real * x_ci = work->iSuperClusterData.xSimd.data();
2398 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2400 for (int i = 0; i < c_nbnxnGpuClusterSize; i += GMX_SIMD4_WIDTH)
2402 int io = si*c_nbnxnGpuClusterSize + i;
2403 int ia = ci*c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize + io;
2404 for (int j = 0; j < GMX_SIMD4_WIDTH; j++)
2406 x_ci[io*DIM + j + XX*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + XX] + shx;
2407 x_ci[io*DIM + j + YY*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + YY] + shy;
2408 x_ci[io*DIM + j + ZZ*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + ZZ] + shz;
2413 #endif /* !GMX_SIMD4_HAVE_REAL */
2416 static real minimum_subgrid_size_xy(const Grid &grid)
2418 const Grid::Dimensions &dims = grid.dimensions();
2420 if (grid.geometry().isSimple)
2422 return std::min(dims.cellSize[XX], dims.cellSize[YY]);
2426 return std::min(dims.cellSize[XX]/c_gpuNumClusterPerCellX,
2427 dims.cellSize[YY]/c_gpuNumClusterPerCellY);
2431 static real effective_buffer_1x1_vs_MxN(const Grid &iGrid,
2434 const real eff_1x1_buffer_fac_overest = 0.1;
2436 /* Determine an atom-pair list cut-off buffer size for atom pairs,
2437 * to be added to rlist (including buffer) used for MxN.
2438 * This is for converting an MxN list to a 1x1 list. This means we can't
2439 * use the normal buffer estimate, as we have an MxN list in which
2440 * some atom pairs beyond rlist are missing. We want to capture
2441 * the beneficial effect of buffering by extra pairs just outside rlist,
2442 * while removing the useless pairs that are further away from rlist.
2443 * (Also the buffer could have been set manually not using the estimate.)
2444 * This buffer size is an overestimate.
2445 * We add 10% of the smallest grid sub-cell dimensions.
2446 * Note that the z-size differs per cell and we don't use this,
2447 * so we overestimate.
2448 * With PME, the 10% value gives a buffer that is somewhat larger
2449 * than the effective buffer with a tolerance of 0.005 kJ/mol/ps.
2450 * Smaller tolerances or using RF lead to a smaller effective buffer,
2451 * so 10% gives a safe overestimate.
2453 return eff_1x1_buffer_fac_overest*(minimum_subgrid_size_xy(iGrid) +
2454 minimum_subgrid_size_xy(jGrid));
2457 /* Estimates the interaction volume^2 for non-local interactions */
2458 static real nonlocal_vol2(const struct gmx_domdec_zones_t *zones, const rvec ls, real r)
2466 /* Here we simply add up the volumes of 1, 2 or 3 1D decomposition
2467 * not home interaction volume^2. As these volumes are not additive,
2468 * this is an overestimate, but it would only be significant in the limit
2469 * of small cells, where we anyhow need to split the lists into
2470 * as small parts as possible.
2473 for (int z = 0; z < zones->n; z++)
2475 if (zones->shift[z][XX] + zones->shift[z][YY] + zones->shift[z][ZZ] == 1)
2480 for (int d = 0; d < DIM; d++)
2482 if (zones->shift[z][d] == 0)
2486 za *= zones->size[z].x1[d] - zones->size[z].x0[d];
2490 /* 4 octants of a sphere */
2491 vold_est = 0.25*M_PI*r*r*r*r;
2492 /* 4 quarter pie slices on the edges */
2493 vold_est += 4*cl*M_PI/6.0*r*r*r;
2494 /* One rectangular volume on a face */
2495 vold_est += ca*0.5*r*r;
2497 vol2_est_tot += vold_est*za;
2501 return vol2_est_tot;
2504 /* Estimates the average size of a full j-list for super/sub setup */
2505 static void get_nsubpair_target(const Nbnxm::GridSet &gridSet,
2506 const InteractionLocality iloc,
2508 const int min_ci_balanced,
2509 int *nsubpair_target,
2510 float *nsubpair_tot_est)
2512 /* The target value of 36 seems to be the optimum for Kepler.
2513 * Maxwell is less sensitive to the exact value.
2515 const int nsubpair_target_min = 36;
2516 real r_eff_sup, vol_est, nsp_est, nsp_est_nl;
2518 const Grid &grid = gridSet.grids()[0];
2520 /* We don't need to balance list sizes if:
2521 * - We didn't request balancing.
2522 * - The number of grid cells >= the number of lists requested,
2523 * since we will always generate at least #cells lists.
2524 * - We don't have any cells, since then there won't be any lists.
2526 if (min_ci_balanced <= 0 || grid.numCells() >= min_ci_balanced || grid.numCells() == 0)
2528 /* nsubpair_target==0 signals no balancing */
2529 *nsubpair_target = 0;
2530 *nsubpair_tot_est = 0;
2536 const int numAtomsCluster = grid.geometry().numAtomsICluster;
2537 const Grid::Dimensions &dims = grid.dimensions();
2539 ls[XX] = dims.cellSize[XX]/c_gpuNumClusterPerCellX;
2540 ls[YY] = dims.cellSize[YY]/c_gpuNumClusterPerCellY;
2541 ls[ZZ] = numAtomsCluster/(dims.atomDensity*ls[XX]*ls[YY]);
2543 /* The formulas below are a heuristic estimate of the average nsj per si*/
2544 r_eff_sup = rlist + nbnxn_get_rlist_effective_inc(numAtomsCluster, ls);
2546 if (!gridSet.domainSetup().haveMultipleDomains ||
2547 gridSet.domainSetup().zones->n == 1)
2554 gmx::square(dims.atomDensity/numAtomsCluster)*
2555 nonlocal_vol2(gridSet.domainSetup().zones, ls, r_eff_sup);
2558 if (iloc == InteractionLocality::Local)
2560 /* Sub-cell interacts with itself */
2561 vol_est = ls[XX]*ls[YY]*ls[ZZ];
2562 /* 6/2 rectangular volume on the faces */
2563 vol_est += (ls[XX]*ls[YY] + ls[XX]*ls[ZZ] + ls[YY]*ls[ZZ])*r_eff_sup;
2564 /* 12/2 quarter pie slices on the edges */
2565 vol_est += 2*(ls[XX] + ls[YY] + ls[ZZ])*0.25*M_PI*gmx::square(r_eff_sup);
2566 /* 4 octants of a sphere */
2567 vol_est += 0.5*4.0/3.0*M_PI*gmx::power3(r_eff_sup);
2569 /* Estimate the number of cluster pairs as the local number of
2570 * clusters times the volume they interact with times the density.
2572 nsp_est = grid.numClusters()*vol_est*dims.atomDensity/numAtomsCluster;
2574 /* Subtract the non-local pair count */
2575 nsp_est -= nsp_est_nl;
2577 /* For small cut-offs nsp_est will be an underesimate.
2578 * With DD nsp_est_nl is an overestimate so nsp_est can get negative.
2579 * So to avoid too small or negative nsp_est we set a minimum of
2580 * all cells interacting with all 3^3 direct neighbors (3^3-1)/2+1=14.
2581 * This might be a slight overestimate for small non-periodic groups of
2582 * atoms as will occur for a local domain with DD, but for small
2583 * groups of atoms we'll anyhow be limited by nsubpair_target_min,
2584 * so this overestimation will not matter.
2586 nsp_est = std::max(nsp_est, grid.numClusters()*14._real);
2590 fprintf(debug, "nsp_est local %5.1f non-local %5.1f\n",
2591 nsp_est, nsp_est_nl);
2596 nsp_est = nsp_est_nl;
2599 /* Thus the (average) maximum j-list size should be as follows.
2600 * Since there is overhead, we shouldn't make the lists too small
2601 * (and we can't chop up j-groups) so we use a minimum target size of 36.
2603 *nsubpair_target = std::max(nsubpair_target_min,
2604 roundToInt(nsp_est/min_ci_balanced));
2605 *nsubpair_tot_est = static_cast<int>(nsp_est);
2609 fprintf(debug, "nbl nsp estimate %.1f, nsubpair_target %d\n",
2610 nsp_est, *nsubpair_target);
2614 /* Debug list print function */
2615 static void print_nblist_ci_cj(FILE *fp,
2616 const NbnxnPairlistCpu &nbl)
2618 for (const nbnxn_ci_t &ciEntry : nbl.ci)
2620 fprintf(fp, "ci %4d shift %2d ncj %3d\n",
2621 ciEntry.ci, ciEntry.shift,
2622 ciEntry.cj_ind_end - ciEntry.cj_ind_start);
2624 for (int j = ciEntry.cj_ind_start; j < ciEntry.cj_ind_end; j++)
2626 fprintf(fp, " cj %5d imask %x\n",
2633 /* Debug list print function */
2634 static void print_nblist_sci_cj(FILE *fp,
2635 const NbnxnPairlistGpu &nbl)
2637 for (const nbnxn_sci_t &sci : nbl.sci)
2639 fprintf(fp, "ci %4d shift %2d ncj4 %2d\n",
2641 sci.numJClusterGroups());
2644 for (int j4 = sci.cj4_ind_start; j4 < sci.cj4_ind_end; j4++)
2646 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
2648 fprintf(fp, " sj %5d imask %x\n",
2650 nbl.cj4[j4].imei[0].imask);
2651 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2653 if (nbl.cj4[j4].imei[0].imask & (1U << (j*c_gpuNumClusterPerCell + si)))
2660 fprintf(fp, "ci %4d shift %2d ncj4 %2d ncp %3d\n",
2662 sci.numJClusterGroups(),
2667 /* Combine pair lists *nbl generated on multiple threads nblc */
2668 static void combine_nblists(gmx::ArrayRef<const NbnxnPairlistGpu> nbls,
2669 NbnxnPairlistGpu *nblc)
2671 int nsci = nblc->sci.size();
2672 int ncj4 = nblc->cj4.size();
2673 int nexcl = nblc->excl.size();
2674 for (auto &nbl : nbls)
2676 nsci += nbl.sci.size();
2677 ncj4 += nbl.cj4.size();
2678 nexcl += nbl.excl.size();
2681 /* Resize with the final, combined size, so we can fill in parallel */
2682 /* NOTE: For better performance we should use default initialization */
2683 nblc->sci.resize(nsci);
2684 nblc->cj4.resize(ncj4);
2685 nblc->excl.resize(nexcl);
2687 /* Each thread should copy its own data to the combined arrays,
2688 * as otherwise data will go back and forth between different caches.
2690 #if GMX_OPENMP && !(defined __clang_analyzer__)
2691 int nthreads = gmx_omp_nthreads_get(emntPairsearch);
2694 #pragma omp parallel for num_threads(nthreads) schedule(static)
2695 for (int n = 0; n < nbls.ssize(); n++)
2699 /* Determine the offset in the combined data for our thread.
2700 * Note that the original sizes in nblc are lost.
2702 int sci_offset = nsci;
2703 int cj4_offset = ncj4;
2704 int excl_offset = nexcl;
2706 for (int i = n; i < nbls.ssize(); i++)
2708 sci_offset -= nbls[i].sci.size();
2709 cj4_offset -= nbls[i].cj4.size();
2710 excl_offset -= nbls[i].excl.size();
2713 const NbnxnPairlistGpu &nbli = nbls[n];
2715 for (size_t i = 0; i < nbli.sci.size(); i++)
2717 nblc->sci[sci_offset + i] = nbli.sci[i];
2718 nblc->sci[sci_offset + i].cj4_ind_start += cj4_offset;
2719 nblc->sci[sci_offset + i].cj4_ind_end += cj4_offset;
2722 for (size_t j4 = 0; j4 < nbli.cj4.size(); j4++)
2724 nblc->cj4[cj4_offset + j4] = nbli.cj4[j4];
2725 nblc->cj4[cj4_offset + j4].imei[0].excl_ind += excl_offset;
2726 nblc->cj4[cj4_offset + j4].imei[1].excl_ind += excl_offset;
2729 for (size_t j4 = 0; j4 < nbli.excl.size(); j4++)
2731 nblc->excl[excl_offset + j4] = nbli.excl[j4];
2734 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2737 for (auto &nbl : nbls)
2739 nblc->nci_tot += nbl.nci_tot;
2743 static void balance_fep_lists(gmx::ArrayRef<t_nblist *> fepLists,
2744 gmx::ArrayRef<PairsearchWork> work)
2746 const int numLists = fepLists.ssize();
2750 /* Nothing to balance */
2754 /* Count the total i-lists and pairs */
2757 for (auto list : fepLists)
2759 nri_tot += list->nri;
2760 nrj_tot += list->nrj;
2763 const int nrj_target = (nrj_tot + numLists - 1)/numLists;
2765 GMX_ASSERT(gmx_omp_nthreads_get(emntNonbonded) == numLists,
2766 "We should have as many work objects as FEP lists");
2768 #pragma omp parallel for schedule(static) num_threads(numLists)
2769 for (int th = 0; th < numLists; th++)
2773 t_nblist *nbl = work[th].nbl_fep.get();
2775 /* Note that here we allocate for the total size, instead of
2776 * a per-thread esimate (which is hard to obtain).
2778 if (nri_tot > nbl->maxnri)
2780 nbl->maxnri = over_alloc_large(nri_tot);
2781 reallocate_nblist(nbl);
2783 if (nri_tot > nbl->maxnri || nrj_tot > nbl->maxnrj)
2785 nbl->maxnrj = over_alloc_small(nrj_tot);
2786 srenew(nbl->jjnr, nbl->maxnrj);
2787 srenew(nbl->excl_fep, nbl->maxnrj);
2790 clear_pairlist_fep(nbl);
2792 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2795 /* Loop over the source lists and assign and copy i-entries */
2797 t_nblist *nbld = work[th_dest].nbl_fep.get();
2798 for (int th = 0; th < numLists; th++)
2800 t_nblist *nbls = fepLists[th];
2802 for (int i = 0; i < nbls->nri; i++)
2806 /* The number of pairs in this i-entry */
2807 nrj = nbls->jindex[i+1] - nbls->jindex[i];
2809 /* Decide if list th_dest is too large and we should procede
2810 * to the next destination list.
2812 if (th_dest + 1 < numLists && nbld->nrj > 0 &&
2813 nbld->nrj + nrj - nrj_target > nrj_target - nbld->nrj)
2816 nbld = work[th_dest].nbl_fep.get();
2819 nbld->iinr[nbld->nri] = nbls->iinr[i];
2820 nbld->gid[nbld->nri] = nbls->gid[i];
2821 nbld->shift[nbld->nri] = nbls->shift[i];
2823 for (int j = nbls->jindex[i]; j < nbls->jindex[i+1]; j++)
2825 nbld->jjnr[nbld->nrj] = nbls->jjnr[j];
2826 nbld->excl_fep[nbld->nrj] = nbls->excl_fep[j];
2830 nbld->jindex[nbld->nri] = nbld->nrj;
2834 /* Swap the list pointers */
2835 for (int th = 0; th < numLists; th++)
2837 t_nblist *nbl_tmp = work[th].nbl_fep.release();
2838 work[th].nbl_fep.reset(fepLists[th]);
2839 fepLists[th] = nbl_tmp;
2843 fprintf(debug, "nbl_fep[%d] nri %4d nrj %4d\n",
2851 /* Returns the next ci to be processes by our thread */
2852 static gmx_bool next_ci(const Grid &grid,
2853 int nth, int ci_block,
2854 int *ci_x, int *ci_y,
2860 if (*ci_b == ci_block)
2862 /* Jump to the next block assigned to this task */
2863 *ci += (nth - 1)*ci_block;
2867 if (*ci >= grid.numCells())
2872 while (*ci >= grid.firstCellInColumn(*ci_x*grid.dimensions().numCells[YY] + *ci_y + 1))
2875 if (*ci_y == grid.dimensions().numCells[YY])
2885 /* Returns the distance^2 for which we put cell pairs in the list
2886 * without checking atom pair distances. This is usually < rlist^2.
2888 static float boundingbox_only_distance2(const Grid::Dimensions &iGridDims,
2889 const Grid::Dimensions &jGridDims,
2893 /* If the distance between two sub-cell bounding boxes is less
2894 * than this distance, do not check the distance between
2895 * all particle pairs in the sub-cell, since then it is likely
2896 * that the box pair has atom pairs within the cut-off.
2897 * We use the nblist cut-off minus 0.5 times the average x/y diagonal
2898 * spacing of the sub-cells. Around 40% of the checked pairs are pruned.
2899 * Using more than 0.5 gains at most 0.5%.
2900 * If forces are calculated more than twice, the performance gain
2901 * in the force calculation outweighs the cost of checking.
2902 * Note that with subcell lists, the atom-pair distance check
2903 * is only performed when only 1 out of 8 sub-cells in within range,
2904 * this is because the GPU is much faster than the cpu.
2909 bbx = 0.5*(iGridDims.cellSize[XX] + jGridDims.cellSize[XX]);
2910 bby = 0.5*(iGridDims.cellSize[YY] + jGridDims.cellSize[YY]);
2913 bbx /= c_gpuNumClusterPerCellX;
2914 bby /= c_gpuNumClusterPerCellY;
2917 rbb2 = std::max(0.0, rlist - 0.5*std::sqrt(bbx*bbx + bby*bby));
2923 return (float)((1+GMX_FLOAT_EPS)*rbb2);
2927 static int get_ci_block_size(const Grid &iGrid,
2928 const bool haveMultipleDomains,
2931 const int ci_block_enum = 5;
2932 const int ci_block_denom = 11;
2933 const int ci_block_min_atoms = 16;
2936 /* Here we decide how to distribute the blocks over the threads.
2937 * We use prime numbers to try to avoid that the grid size becomes
2938 * a multiple of the number of threads, which would lead to some
2939 * threads getting "inner" pairs and others getting boundary pairs,
2940 * which in turns will lead to load imbalance between threads.
2941 * Set the block size as 5/11/ntask times the average number of cells
2942 * in a y,z slab. This should ensure a quite uniform distribution
2943 * of the grid parts of the different thread along all three grid
2944 * zone boundaries with 3D domain decomposition. At the same time
2945 * the blocks will not become too small.
2947 GMX_ASSERT(iGrid.dimensions().numCells[XX] > 0, "Grid can't be empty");
2948 GMX_ASSERT(numLists > 0, "We need at least one list");
2949 ci_block = (iGrid.numCells()*ci_block_enum)/(ci_block_denom*iGrid.dimensions().numCells[XX]*numLists);
2951 const int numAtomsPerCell = iGrid.geometry().numAtomsPerCell;
2953 /* Ensure the blocks are not too small: avoids cache invalidation */
2954 if (ci_block*numAtomsPerCell < ci_block_min_atoms)
2956 ci_block = (ci_block_min_atoms + numAtomsPerCell - 1)/numAtomsPerCell;
2959 /* Without domain decomposition
2960 * or with less than 3 blocks per task, divide in nth blocks.
2962 if (!haveMultipleDomains || numLists*3*ci_block > iGrid.numCells())
2964 ci_block = (iGrid.numCells() + numLists - 1)/numLists;
2967 if (ci_block > 1 && (numLists - 1)*ci_block >= iGrid.numCells())
2969 /* Some threads have no work. Although reducing the block size
2970 * does not decrease the block count on the first few threads,
2971 * with GPUs better mixing of "upper" cells that have more empty
2972 * clusters results in a somewhat lower max load over all threads.
2973 * Without GPUs the regime of so few atoms per thread is less
2974 * performance relevant, but with 8-wide SIMD the same reasoning
2975 * applies, since the pair list uses 4 i-atom "sub-clusters".
2983 /* Returns the number of bits to right-shift a cluster index to obtain
2984 * the corresponding force buffer flag index.
2986 static int getBufferFlagShift(int numAtomsPerCluster)
2988 int bufferFlagShift = 0;
2989 while ((numAtomsPerCluster << bufferFlagShift) < NBNXN_BUFFERFLAG_SIZE)
2994 return bufferFlagShift;
2997 static bool pairlistIsSimple(const NbnxnPairlistCpu gmx_unused &pairlist)
3002 static bool pairlistIsSimple(const NbnxnPairlistGpu gmx_unused &pairlist)
3008 makeClusterListWrapper(NbnxnPairlistCpu *nbl,
3009 const Grid gmx_unused &iGrid,
3012 const int firstCell,
3014 const bool excludeSubDiagonal,
3015 const nbnxn_atomdata_t *nbat,
3018 const ClusterDistanceKernelType kernelType,
3019 int *numDistanceChecks)
3023 case ClusterDistanceKernelType::CpuPlainC:
3024 makeClusterListSimple(jGrid,
3025 nbl, ci, firstCell, lastCell,
3031 #ifdef GMX_NBNXN_SIMD_4XN
3032 case ClusterDistanceKernelType::CpuSimd_4xM:
3033 makeClusterListSimd4xn(jGrid,
3034 nbl, ci, firstCell, lastCell,
3041 #ifdef GMX_NBNXN_SIMD_2XNN
3042 case ClusterDistanceKernelType::CpuSimd_2xMM:
3043 makeClusterListSimd2xnn(jGrid,
3044 nbl, ci, firstCell, lastCell,
3052 GMX_ASSERT(false, "Unhandled kernel type");
3057 makeClusterListWrapper(NbnxnPairlistGpu *nbl,
3058 const Grid &gmx_unused iGrid,
3061 const int firstCell,
3063 const bool excludeSubDiagonal,
3064 const nbnxn_atomdata_t *nbat,
3067 ClusterDistanceKernelType gmx_unused kernelType,
3068 int *numDistanceChecks)
3070 for (int cj = firstCell; cj <= lastCell; cj++)
3072 make_cluster_list_supersub(iGrid, jGrid,
3075 nbat->xstride, nbat->x().data(),
3081 static int getNumSimpleJClustersInList(const NbnxnPairlistCpu &nbl)
3083 return nbl.cj.size();
3086 static int getNumSimpleJClustersInList(const gmx_unused NbnxnPairlistGpu &nbl)
3091 static void incrementNumSimpleJClustersInList(NbnxnPairlistCpu *nbl,
3094 nbl->ncjInUse += nbl->cj.size() - ncj_old_j;
3097 static void incrementNumSimpleJClustersInList(NbnxnPairlistGpu gmx_unused *nbl,
3098 int gmx_unused ncj_old_j)
3102 static void checkListSizeConsistency(const NbnxnPairlistCpu &nbl,
3103 const bool haveFreeEnergy)
3105 GMX_RELEASE_ASSERT(static_cast<size_t>(nbl.ncjInUse) == nbl.cj.size() || haveFreeEnergy,
3106 "Without free-energy all cj pair-list entries should be in use. "
3107 "Note that subsequent code does not make use of the equality, "
3108 "this check is only here to catch bugs");
3111 static void checkListSizeConsistency(const NbnxnPairlistGpu gmx_unused &nbl,
3112 bool gmx_unused haveFreeEnergy)
3114 /* We currently can not check consistency here */
3117 /* Set the buffer flags for newly added entries in the list */
3118 static void setBufferFlags(const NbnxnPairlistCpu &nbl,
3119 const int ncj_old_j,
3120 const int gridj_flag_shift,
3121 gmx_bitmask_t *gridj_flag,
3124 if (gmx::ssize(nbl.cj) > ncj_old_j)
3126 int cbFirst = nbl.cj[ncj_old_j].cj >> gridj_flag_shift;
3127 int cbLast = nbl.cj.back().cj >> gridj_flag_shift;
3128 for (int cb = cbFirst; cb <= cbLast; cb++)
3130 bitmask_init_bit(&gridj_flag[cb], th);
3135 static void setBufferFlags(const NbnxnPairlistGpu gmx_unused &nbl,
3136 int gmx_unused ncj_old_j,
3137 int gmx_unused gridj_flag_shift,
3138 gmx_bitmask_t gmx_unused *gridj_flag,
3141 GMX_ASSERT(false, "This function should never be called");
3144 /* Generates the part of pair-list nbl assigned to our thread */
3145 template <typename T>
3146 static void nbnxn_make_pairlist_part(const Nbnxm::GridSet &gridSet,
3149 PairsearchWork *work,
3150 const nbnxn_atomdata_t *nbat,
3151 const t_blocka &exclusions,
3153 const PairlistType pairlistType,
3155 gmx_bool bFBufferFlag,
3158 float nsubpair_tot_est,
3167 int ci_b, ci, ci_x, ci_y, ci_xy;
3169 real bx0, bx1, by0, by1, bz0, bz1;
3171 real d2cx, d2z, d2z_cx, d2z_cy, d2zx, d2zxy, d2xy;
3172 int cxf, cxl, cyf, cyf_x, cyl;
3173 int numDistanceChecks;
3174 int gridi_flag_shift = 0, gridj_flag_shift = 0;
3175 gmx_bitmask_t *gridj_flag = nullptr;
3176 int ncj_old_i, ncj_old_j;
3178 if (jGrid.geometry().isSimple != pairlistIsSimple(*nbl) ||
3179 iGrid.geometry().isSimple != pairlistIsSimple(*nbl))
3181 gmx_incons("Grid incompatible with pair-list");
3185 GMX_ASSERT(nbl->na_ci == jGrid.geometry().numAtomsICluster,
3186 "The cluster sizes in the list and grid should match");
3187 nbl->na_cj = JClusterSizePerListType[pairlistType];
3188 na_cj_2log = get_2log(nbl->na_cj);
3194 /* Determine conversion of clusters to flag blocks */
3195 gridi_flag_shift = getBufferFlagShift(nbl->na_ci);
3196 gridj_flag_shift = getBufferFlagShift(nbl->na_cj);
3198 gridj_flag = work->buffer_flags.flag;
3201 gridSet.getBox(box);
3203 const bool haveFep = gridSet.haveFep();
3205 const real rlist2 = nbl->rlist*nbl->rlist;
3207 // Select the cluster pair distance kernel type
3208 const ClusterDistanceKernelType kernelType =
3209 getClusterDistanceKernelType(pairlistType, *nbat);
3211 if (haveFep && !pairlistIsSimple(*nbl))
3213 /* Determine an atom-pair list cut-off distance for FEP atom pairs.
3214 * We should not simply use rlist, since then we would not have
3215 * the small, effective buffering of the NxN lists.
3216 * The buffer is on overestimate, but the resulting cost for pairs
3217 * beyond rlist is neglible compared to the FEP pairs within rlist.
3219 rl_fep2 = nbl->rlist + effective_buffer_1x1_vs_MxN(iGrid, jGrid);
3223 fprintf(debug, "nbl_fep atom-pair rlist %f\n", rl_fep2);
3225 rl_fep2 = rl_fep2*rl_fep2;
3228 const Grid::Dimensions &iGridDims = iGrid.dimensions();
3229 const Grid::Dimensions &jGridDims = jGrid.dimensions();
3231 rbb2 = boundingbox_only_distance2(iGridDims, jGridDims, nbl->rlist, pairlistIsSimple(*nbl));
3235 fprintf(debug, "nbl bounding box only distance %f\n", std::sqrt(rbb2));
3238 const bool isIntraGridList = (&iGrid == &jGrid);
3240 /* Set the shift range */
3241 for (int d = 0; d < DIM; d++)
3243 /* Check if we need periodicity shifts.
3244 * Without PBC or with domain decomposition we don't need them.
3246 if (d >= ePBC2npbcdim(gridSet.domainSetup().ePBC) ||
3247 gridSet.domainSetup().haveMultipleDomainsPerDim[d])
3253 const real listRangeCellToCell =
3254 listRangeForGridCellToGridCell(rlist, iGrid.dimensions(), jGrid.dimensions());
3256 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < listRangeCellToCell)
3266 const bool bSimple = pairlistIsSimple(*nbl);
3267 gmx::ArrayRef<const BoundingBox> bb_i;
3269 gmx::ArrayRef<const float> pbb_i;
3272 bb_i = iGrid.iBoundingBoxes();
3276 pbb_i = iGrid.packedBoundingBoxes();
3279 /* We use the normal bounding box format for both grid types */
3280 bb_i = iGrid.iBoundingBoxes();
3282 gmx::ArrayRef<const BoundingBox1D> bbcz_i = iGrid.zBoundingBoxes();
3283 gmx::ArrayRef<const int> flags_i = iGrid.clusterFlags();
3284 gmx::ArrayRef<const BoundingBox1D> bbcz_j = jGrid.zBoundingBoxes();
3285 int cell0_i = iGrid.cellOffset();
3289 fprintf(debug, "nbl nc_i %d col.av. %.1f ci_block %d\n",
3290 iGrid.numCells(), iGrid.numCells()/static_cast<double>(iGrid.numColumns()), ci_block);
3293 numDistanceChecks = 0;
3295 const real listRangeBBToJCell2 = gmx::square(listRangeForBoundingBoxToGridCell(rlist, jGrid.dimensions()));
3297 /* Initially ci_b and ci to 1 before where we want them to start,
3298 * as they will both be incremented in next_ci.
3301 ci = th*ci_block - 1;
3304 while (next_ci(iGrid, nth, ci_block, &ci_x, &ci_y, &ci_b, &ci))
3306 if (bSimple && flags_i[ci] == 0)
3311 ncj_old_i = getNumSimpleJClustersInList(*nbl);
3314 if (!isIntraGridList && shp[XX] == 0)
3318 bx1 = bb_i[ci].upper.x;
3322 bx1 = iGridDims.lowerCorner[XX] + (ci_x+1)*iGridDims.cellSize[XX];
3324 if (bx1 < jGridDims.lowerCorner[XX])
3326 d2cx = gmx::square(jGridDims.lowerCorner[XX] - bx1);
3328 if (d2cx >= listRangeBBToJCell2)
3335 ci_xy = ci_x*iGridDims.numCells[YY] + ci_y;
3337 /* Loop over shift vectors in three dimensions */
3338 for (int tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
3340 const real shz = tz*box[ZZ][ZZ];
3342 bz0 = bbcz_i[ci].lower + shz;
3343 bz1 = bbcz_i[ci].upper + shz;
3351 d2z = gmx::square(bz1);
3355 d2z = gmx::square(bz0 - box[ZZ][ZZ]);
3358 d2z_cx = d2z + d2cx;
3360 if (d2z_cx >= rlist2)
3365 bz1_frac = bz1/iGrid.numCellsInColumn(ci_xy);
3370 /* The check with bz1_frac close to or larger than 1 comes later */
3372 for (int ty = -shp[YY]; ty <= shp[YY]; ty++)
3374 const real shy = ty*box[YY][YY] + tz*box[ZZ][YY];
3378 by0 = bb_i[ci].lower.y + shy;
3379 by1 = bb_i[ci].upper.y + shy;
3383 by0 = iGridDims.lowerCorner[YY] + (ci_y )*iGridDims.cellSize[YY] + shy;
3384 by1 = iGridDims.lowerCorner[YY] + (ci_y + 1)*iGridDims.cellSize[YY] + shy;
3387 get_cell_range<YY>(by0, by1,
3398 if (by1 < jGridDims.lowerCorner[YY])
3400 d2z_cy += gmx::square(jGridDims.lowerCorner[YY] - by1);
3402 else if (by0 > jGridDims.upperCorner[YY])
3404 d2z_cy += gmx::square(by0 - jGridDims.upperCorner[YY]);
3407 for (int tx = -shp[XX]; tx <= shp[XX]; tx++)
3409 const int shift = XYZ2IS(tx, ty, tz);
3411 const bool excludeSubDiagonal = (isIntraGridList && shift == CENTRAL);
3413 if (c_pbcShiftBackward && isIntraGridList && shift > CENTRAL)
3418 const real shx = tx*box[XX][XX] + ty*box[YY][XX] + tz*box[ZZ][XX];
3422 bx0 = bb_i[ci].lower.x + shx;
3423 bx1 = bb_i[ci].upper.x + shx;
3427 bx0 = iGridDims.lowerCorner[XX] + (ci_x )*iGridDims.cellSize[XX] + shx;
3428 bx1 = iGridDims.lowerCorner[XX] + (ci_x+1)*iGridDims.cellSize[XX] + shx;
3431 get_cell_range<XX>(bx0, bx1,
3441 addNewIEntry(nbl, cell0_i+ci, shift, flags_i[ci]);
3443 if ((!c_pbcShiftBackward || excludeSubDiagonal) &&
3446 /* Leave the pairs with i > j.
3447 * x is the major index, so skip half of it.
3452 set_icell_bb(iGrid, ci, shx, shy, shz,
3455 icell_set_x(cell0_i+ci, shx, shy, shz,
3456 nbat->xstride, nbat->x().data(),
3460 for (int cx = cxf; cx <= cxl; cx++)
3463 if (jGridDims.lowerCorner[XX] + cx*jGridDims.cellSize[XX] > bx1)
3465 d2zx += gmx::square(jGridDims.lowerCorner[XX] + cx*jGridDims.cellSize[XX] - bx1);
3467 else if (jGridDims.lowerCorner[XX] + (cx+1)*jGridDims.cellSize[XX] < bx0)
3469 d2zx += gmx::square(jGridDims.lowerCorner[XX] + (cx+1)*jGridDims.cellSize[XX] - bx0);
3472 if (isIntraGridList &&
3474 (!c_pbcShiftBackward || shift == CENTRAL) &&
3477 /* Leave the pairs with i > j.
3478 * Skip half of y when i and j have the same x.
3487 for (int cy = cyf_x; cy <= cyl; cy++)
3489 const int columnStart = jGrid.firstCellInColumn(cx*jGridDims.numCells[YY] + cy);
3490 const int columnEnd = jGrid.firstCellInColumn(cx*jGridDims.numCells[YY] + cy + 1);
3493 if (jGridDims.lowerCorner[YY] + cy*jGridDims.cellSize[YY] > by1)
3495 d2zxy += gmx::square(jGridDims.lowerCorner[YY] + cy*jGridDims.cellSize[YY] - by1);
3497 else if (jGridDims.lowerCorner[YY] + (cy + 1)*jGridDims.cellSize[YY] < by0)
3499 d2zxy += gmx::square(jGridDims.lowerCorner[YY] + (cy + 1)*jGridDims.cellSize[YY] - by0);
3501 if (columnStart < columnEnd && d2zxy < listRangeBBToJCell2)
3503 /* To improve efficiency in the common case
3504 * of a homogeneous particle distribution,
3505 * we estimate the index of the middle cell
3506 * in range (midCell). We search down and up
3507 * starting from this index.
3509 * Note that the bbcz_j array contains bounds
3510 * for i-clusters, thus for clusters of 4 atoms.
3511 * For the common case where the j-cluster size
3512 * is 8, we could step with a stride of 2,
3513 * but we do not do this because it would
3514 * complicate this code even more.
3516 int midCell = columnStart + static_cast<int>(bz1_frac*(columnEnd - columnStart));
3517 if (midCell >= columnEnd)
3519 midCell = columnEnd - 1;
3524 /* Find the lowest cell that can possibly
3526 * Check if we hit the bottom of the grid,
3527 * if the j-cell is below the i-cell and if so,
3528 * if it is within range.
3530 int downTestCell = midCell;
3531 while (downTestCell >= columnStart &&
3532 (bbcz_j[downTestCell].upper >= bz0 ||
3533 d2xy + gmx::square(bbcz_j[downTestCell].upper - bz0) < rlist2))
3537 int firstCell = downTestCell + 1;
3539 /* Find the highest cell that can possibly
3541 * Check if we hit the top of the grid,
3542 * if the j-cell is above the i-cell and if so,
3543 * if it is within range.
3545 int upTestCell = midCell + 1;
3546 while (upTestCell < columnEnd &&
3547 (bbcz_j[upTestCell].lower <= bz1 ||
3548 d2xy + gmx::square(bbcz_j[upTestCell].lower - bz1) < rlist2))
3552 int lastCell = upTestCell - 1;
3554 #define NBNXN_REFCODE 0
3557 /* Simple reference code, for debugging,
3558 * overrides the more complex code above.
3560 firstCell = columnEnd;
3562 for (int k = columnStart; k < columnEnd; k++)
3564 if (d2xy + gmx::square(bbcz_j[k*NNBSBB_D + 1] - bz0) < rlist2 &&
3569 if (d2xy + gmx::square(bbcz_j[k*NNBSBB_D] - bz1) < rlist2 &&
3578 if (isIntraGridList)
3580 /* We want each atom/cell pair only once,
3581 * only use cj >= ci.
3583 if (!c_pbcShiftBackward || shift == CENTRAL)
3585 firstCell = std::max(firstCell, ci);
3589 if (firstCell <= lastCell)
3591 GMX_ASSERT(firstCell >= columnStart && lastCell < columnEnd, "The range should reside within the current grid column");
3593 /* For f buffer flags with simple lists */
3594 ncj_old_j = getNumSimpleJClustersInList(*nbl);
3596 makeClusterListWrapper(nbl,
3598 jGrid, firstCell, lastCell,
3603 &numDistanceChecks);
3607 setBufferFlags(*nbl, ncj_old_j, gridj_flag_shift,
3611 incrementNumSimpleJClustersInList(nbl, ncj_old_j);
3617 /* Set the exclusions for this ci list */
3618 setExclusionsForIEntry(gridSet,
3622 *getOpenIEntry(nbl),
3627 make_fep_list(gridSet.atomIndices(), nbat, nbl,
3632 iGrid, jGrid, nbl_fep);
3635 /* Close this ci list */
3638 progBal, nsubpair_tot_est,
3644 if (bFBufferFlag && getNumSimpleJClustersInList(*nbl) > ncj_old_i)
3646 bitmask_init_bit(&(work->buffer_flags.flag[(iGrid.cellOffset() + ci) >> gridi_flag_shift]), th);
3650 work->ndistc = numDistanceChecks;
3652 checkListSizeConsistency(*nbl, haveFep);
3656 fprintf(debug, "number of distance checks %d\n", numDistanceChecks);
3658 print_nblist_statistics(debug, *nbl, gridSet, rlist);
3662 fprintf(debug, "nbl FEP list pairs: %d\n", nbl_fep->nrj);
3667 static void reduce_buffer_flags(gmx::ArrayRef<PairsearchWork> searchWork,
3669 const nbnxn_buffer_flags_t *dest)
3671 for (int s = 0; s < nsrc; s++)
3673 gmx_bitmask_t * flag = searchWork[s].buffer_flags.flag;
3675 for (int b = 0; b < dest->nflag; b++)
3677 bitmask_union(&(dest->flag[b]), flag[b]);
3682 static void print_reduction_cost(const nbnxn_buffer_flags_t *flags, int nout)
3684 int nelem, nkeep, ncopy, nred, out;
3685 gmx_bitmask_t mask_0;
3691 bitmask_init_bit(&mask_0, 0);
3692 for (int b = 0; b < flags->nflag; b++)
3694 if (bitmask_is_equal(flags->flag[b], mask_0))
3696 /* Only flag 0 is set, no copy of reduction required */
3700 else if (!bitmask_is_zero(flags->flag[b]))
3703 for (out = 0; out < nout; out++)
3705 if (bitmask_is_set(flags->flag[b], out))
3722 fprintf(debug, "nbnxn reduction: #flag %d #list %d elem %4.2f, keep %4.2f copy %4.2f red %4.2f\n",
3724 nelem/static_cast<double>(flags->nflag),
3725 nkeep/static_cast<double>(flags->nflag),
3726 ncopy/static_cast<double>(flags->nflag),
3727 nred/static_cast<double>(flags->nflag));
3730 /* Copies the list entries from src to dest when cjStart <= *cjGlobal < cjEnd.
3731 * *cjGlobal is updated with the cj count in src.
3732 * When setFlags==true, flag bit t is set in flag for all i and j clusters.
3734 template<bool setFlags>
3735 static void copySelectedListRange(const nbnxn_ci_t * gmx_restrict srcCi,
3736 const NbnxnPairlistCpu * gmx_restrict src,
3737 NbnxnPairlistCpu * gmx_restrict dest,
3738 gmx_bitmask_t *flag,
3739 int iFlagShift, int jFlagShift, int t)
3741 const int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3743 dest->ci.push_back(*srcCi);
3744 dest->ci.back().cj_ind_start = dest->cj.size();
3745 dest->ci.back().cj_ind_end = dest->cj.size() + ncj;
3749 bitmask_init_bit(&flag[srcCi->ci >> iFlagShift], t);
3752 for (int j = srcCi->cj_ind_start; j < srcCi->cj_ind_end; j++)
3754 dest->cj.push_back(src->cj[j]);
3758 /* NOTE: This is relatively expensive, since this
3759 * operation is done for all elements in the list,
3760 * whereas at list generation this is done only
3761 * once for each flag entry.
3763 bitmask_init_bit(&flag[src->cj[j].cj >> jFlagShift], t);
3768 /* This routine re-balances the pairlists such that all are nearly equally
3769 * sized. Only whole i-entries are moved between lists. These are moved
3770 * between the ends of the lists, such that the buffer reduction cost should
3771 * not change significantly.
3772 * Note that all original reduction flags are currently kept. This can lead
3773 * to reduction of parts of the force buffer that could be avoided. But since
3774 * the original lists are quite balanced, this will only give minor overhead.
3776 static void rebalanceSimpleLists(gmx::ArrayRef<const NbnxnPairlistCpu> srcSet,
3777 gmx::ArrayRef<NbnxnPairlistCpu> destSet,
3778 gmx::ArrayRef<PairsearchWork> searchWork)
3781 for (auto &src : srcSet)
3783 ncjTotal += src.ncjInUse;
3785 const int numLists = srcSet.ssize();
3786 const int ncjTarget = (ncjTotal + numLists - 1)/numLists;
3788 #pragma omp parallel num_threads(numLists)
3790 int t = gmx_omp_get_thread_num();
3792 int cjStart = ncjTarget* t;
3793 int cjEnd = ncjTarget*(t + 1);
3795 /* The destination pair-list for task/thread t */
3796 NbnxnPairlistCpu &dest = destSet[t];
3798 clear_pairlist(&dest);
3799 dest.na_cj = srcSet[0].na_cj;
3801 /* Note that the flags in the work struct (still) contain flags
3802 * for all entries that are present in srcSet->nbl[t].
3804 gmx_bitmask_t *flag = searchWork[t].buffer_flags.flag;
3806 int iFlagShift = getBufferFlagShift(dest.na_ci);
3807 int jFlagShift = getBufferFlagShift(dest.na_cj);
3810 for (int s = 0; s < numLists && cjGlobal < cjEnd; s++)
3812 const NbnxnPairlistCpu *src = &srcSet[s];
3814 if (cjGlobal + src->ncjInUse > cjStart)
3816 for (gmx::index i = 0; i < gmx::ssize(src->ci) && cjGlobal < cjEnd; i++)
3818 const nbnxn_ci_t *srcCi = &src->ci[i];
3819 int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3820 if (cjGlobal >= cjStart)
3822 /* If the source list is not our own, we need to set
3823 * extra flags (the template bool parameter).
3827 copySelectedListRange
3830 flag, iFlagShift, jFlagShift, t);
3834 copySelectedListRange
3837 &dest, flag, iFlagShift, jFlagShift, t);
3845 cjGlobal += src->ncjInUse;
3849 dest.ncjInUse = dest.cj.size();
3853 int ncjTotalNew = 0;
3854 for (auto &dest : destSet)
3856 ncjTotalNew += dest.ncjInUse;
3858 GMX_RELEASE_ASSERT(ncjTotalNew == ncjTotal, "The total size of the lists before and after rebalancing should match");
3862 /* Returns if the pairlists are so imbalanced that it is worth rebalancing. */
3863 static bool checkRebalanceSimpleLists(gmx::ArrayRef<const NbnxnPairlistCpu> lists)
3865 int numLists = lists.ssize();
3868 for (int s = 0; s < numLists; s++)
3870 ncjMax = std::max(ncjMax, lists[s].ncjInUse);
3871 ncjTotal += lists[s].ncjInUse;
3875 fprintf(debug, "Pair-list ncjMax %d ncjTotal %d\n", ncjMax, ncjTotal);
3877 /* The rebalancing adds 3% extra time to the search. Heuristically we
3878 * determined that under common conditions the non-bonded kernel balance
3879 * improvement will outweigh this when the imbalance is more than 3%.
3880 * But this will, obviously, depend on search vs kernel time and nstlist.
3882 const real rebalanceTolerance = 1.03;
3884 return numLists*ncjMax > ncjTotal*rebalanceTolerance;
3887 /* Perform a count (linear) sort to sort the smaller lists to the end.
3888 * This avoids load imbalance on the GPU, as large lists will be
3889 * scheduled and executed first and the smaller lists later.
3890 * Load balancing between multi-processors only happens at the end
3891 * and there smaller lists lead to more effective load balancing.
3892 * The sorting is done on the cj4 count, not on the actual pair counts.
3893 * Not only does this make the sort faster, but it also results in
3894 * better load balancing than using a list sorted on exact load.
3895 * This function swaps the pointer in the pair list to avoid a copy operation.
3897 static void sort_sci(NbnxnPairlistGpu *nbl)
3899 if (nbl->cj4.size() <= nbl->sci.size())
3901 /* nsci = 0 or all sci have size 1, sorting won't change the order */
3905 NbnxnPairlistGpuWork &work = *nbl->work;
3907 /* We will distinguish differences up to double the average */
3908 const int m = (2*nbl->cj4.size())/nbl->sci.size();
3910 /* Resize work.sci_sort so we can sort into it */
3911 work.sci_sort.resize(nbl->sci.size());
3913 std::vector<int> &sort = work.sortBuffer;
3914 /* Set up m + 1 entries in sort, initialized at 0 */
3916 sort.resize(m + 1, 0);
3917 /* Count the entries of each size */
3918 for (const nbnxn_sci_t &sci : nbl->sci)
3920 int i = std::min(m, sci.numJClusterGroups());
3923 /* Calculate the offset for each count */
3926 for (int i = m - 1; i >= 0; i--)
3929 sort[i] = sort[i + 1] + s0;
3933 /* Sort entries directly into place */
3934 gmx::ArrayRef<nbnxn_sci_t> sci_sort = work.sci_sort;
3935 for (const nbnxn_sci_t &sci : nbl->sci)
3937 int i = std::min(m, sci.numJClusterGroups());
3938 sci_sort[sort[i]++] = sci;
3941 /* Swap the sci pointers so we use the new, sorted list */
3942 std::swap(nbl->sci, work.sci_sort);
3945 //! Prepares CPU lists produced by the search for dynamic pruning
3946 static void prepareListsForDynamicPruning(gmx::ArrayRef<NbnxnPairlistCpu> lists);
3949 PairlistSet::constructPairlists(const Nbnxm::GridSet &gridSet,
3950 gmx::ArrayRef<PairsearchWork> searchWork,
3951 nbnxn_atomdata_t *nbat,
3952 const t_blocka *excl,
3953 const int minimumIlistCountForGpuBalancing,
3955 SearchCycleCounting *searchCycleCounting)
3957 const real rlist = params_.rlistOuter;
3959 int nsubpair_target;
3960 float nsubpair_tot_est;
3963 int np_tot, np_noq, np_hlj, nap;
3965 const int numLists = (isCpuType_ ? cpuLists_.size() : gpuLists_.size());
3969 fprintf(debug, "ns making %d nblists\n", numLists);
3972 nbat->bUseBufferFlags = (nbat->out.size() > 1);
3973 /* We should re-init the flags before making the first list */
3974 if (nbat->bUseBufferFlags && locality_ == InteractionLocality::Local)
3976 init_buffer_flags(&nbat->buffer_flags, nbat->numAtoms());
3980 if (locality_ == InteractionLocality::Local)
3982 /* Only zone (grid) 0 vs 0 */
3987 nzi = gridSet.domainSetup().zones->nizone;
3990 if (!isCpuType_ && minimumIlistCountForGpuBalancing > 0)
3992 get_nsubpair_target(gridSet, locality_, rlist, minimumIlistCountForGpuBalancing,
3993 &nsubpair_target, &nsubpair_tot_est);
3997 nsubpair_target = 0;
3998 nsubpair_tot_est = 0;
4001 /* Clear all pair-lists */
4002 for (int th = 0; th < numLists; th++)
4006 clear_pairlist(&cpuLists_[th]);
4010 clear_pairlist(&gpuLists_[th]);
4013 if (params_.haveFep)
4015 clear_pairlist_fep(fepLists_[th]);
4019 const gmx_domdec_zones_t *ddZones = gridSet.domainSetup().zones;
4021 for (int zi = 0; zi < nzi; zi++)
4023 const Grid &iGrid = gridSet.grids()[zi];
4027 if (locality_ == InteractionLocality::Local)
4034 zj0 = ddZones->izone[zi].j0;
4035 zj1 = ddZones->izone[zi].j1;
4041 for (int zj = zj0; zj < zj1; zj++)
4043 const Grid &jGrid = gridSet.grids()[zj];
4047 fprintf(debug, "ns search grid %d vs %d\n", zi, zj);
4050 searchCycleCounting->start(enbsCCsearch);
4052 ci_block = get_ci_block_size(iGrid, gridSet.domainSetup().haveMultipleDomains, numLists);
4054 /* With GPU: generate progressively smaller lists for
4055 * load balancing for local only or non-local with 2 zones.
4057 progBal = (locality_ == InteractionLocality::Local || ddZones->n <= 2);
4059 #pragma omp parallel for num_threads(numLists) schedule(static)
4060 for (int th = 0; th < numLists; th++)
4064 /* Re-init the thread-local work flag data before making
4065 * the first list (not an elegant conditional).
4067 if (nbat->bUseBufferFlags && ((zi == 0 && zj == 0)))
4069 init_buffer_flags(&searchWork[th].buffer_flags, nbat->numAtoms());
4072 if (combineLists_ && th > 0)
4074 GMX_ASSERT(!isCpuType_, "Can only combine GPU lists");
4076 clear_pairlist(&gpuLists_[th]);
4079 PairsearchWork &work = searchWork[th];
4081 work.cycleCounter.start();
4083 t_nblist *fepListPtr = (fepLists_.empty() ? nullptr : fepLists_[th]);
4085 /* Divide the i cells equally over the pairlists */
4088 nbnxn_make_pairlist_part(gridSet, iGrid, jGrid,
4091 params_.pairlistType,
4093 nbat->bUseBufferFlags,
4095 progBal, nsubpair_tot_est,
4102 nbnxn_make_pairlist_part(gridSet, iGrid, jGrid,
4105 params_.pairlistType,
4107 nbat->bUseBufferFlags,
4109 progBal, nsubpair_tot_est,
4115 work.cycleCounter.stop();
4117 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4119 searchCycleCounting->stop(enbsCCsearch);
4124 for (int th = 0; th < numLists; th++)
4126 inc_nrnb(nrnb, eNR_NBNXN_DIST2, searchWork[th].ndistc);
4130 const NbnxnPairlistCpu &nbl = cpuLists_[th];
4131 np_tot += nbl.cj.size();
4132 np_noq += nbl.work->ncj_noq;
4133 np_hlj += nbl.work->ncj_hlj;
4137 const NbnxnPairlistGpu &nbl = gpuLists_[th];
4138 /* This count ignores potential subsequent pair pruning */
4139 np_tot += nbl.nci_tot;
4144 nap = cpuLists_[0].na_ci*cpuLists_[0].na_cj;
4148 nap = gmx::square(gpuLists_[0].na_ci);
4150 natpair_ljq_ = (np_tot - np_noq)*nap - np_hlj*nap/2;
4151 natpair_lj_ = np_noq*nap;
4152 natpair_q_ = np_hlj*nap/2;
4154 if (combineLists_ && numLists > 1)
4156 GMX_ASSERT(!isCpuType_, "Can only combine GPU lists");
4158 searchCycleCounting->start(enbsCCcombine);
4160 combine_nblists(gmx::constArrayRefFromArray(&gpuLists_[1], numLists - 1),
4163 searchCycleCounting->stop(enbsCCcombine);
4170 if (numLists > 1 && checkRebalanceSimpleLists(cpuLists_))
4172 rebalanceSimpleLists(cpuLists_, cpuListsWork_, searchWork);
4174 /* Swap the sets of pair lists */
4175 cpuLists_.swap(cpuListsWork_);
4180 /* Sort the entries on size, large ones first */
4181 if (combineLists_ || gpuLists_.size() == 1)
4183 sort_sci(&gpuLists_[0]);
4187 #pragma omp parallel for num_threads(numLists) schedule(static)
4188 for (int th = 0; th < numLists; th++)
4192 sort_sci(&gpuLists_[th]);
4194 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4199 if (nbat->bUseBufferFlags)
4201 reduce_buffer_flags(searchWork, numLists, &nbat->buffer_flags);
4204 if (gridSet.haveFep())
4206 /* Balance the free-energy lists over all the threads */
4207 balance_fep_lists(fepLists_, searchWork);
4212 /* This is a fresh list, so not pruned, stored using ci.
4213 * ciOuter is invalid at this point.
4215 GMX_ASSERT(cpuLists_[0].ciOuter.empty(), "ciOuter is invalid so it should be empty");
4218 /* If we have more than one list, they either got rebalancing (CPU)
4219 * or combined (GPU), so we should dump the final result to debug.
4223 if (isCpuType_ && cpuLists_.size() > 1)
4225 for (auto &cpuList : cpuLists_)
4227 print_nblist_statistics(debug, cpuList, gridSet, rlist);
4230 else if (!isCpuType_ && gpuLists_.size() > 1)
4232 print_nblist_statistics(debug, gpuLists_[0], gridSet, rlist);
4242 for (auto &cpuList : cpuLists_)
4244 print_nblist_ci_cj(debug, cpuList);
4249 print_nblist_sci_cj(debug, gpuLists_[0]);
4253 if (nbat->bUseBufferFlags)
4255 print_reduction_cost(&nbat->buffer_flags, numLists);
4259 if (params_.useDynamicPruning && isCpuType_)
4261 prepareListsForDynamicPruning(cpuLists_);
4266 PairlistSets::construct(const InteractionLocality iLocality,
4267 PairSearch *pairSearch,
4268 nbnxn_atomdata_t *nbat,
4269 const t_blocka *excl,
4273 pairlistSet(iLocality).constructPairlists(pairSearch->gridSet(), pairSearch->work(),
4274 nbat, excl, minimumIlistCountForGpuBalancing_,
4275 nrnb, &pairSearch->cycleCounting_);
4277 if (iLocality == Nbnxm::InteractionLocality::Local)
4279 outerListCreationStep_ = step;
4283 GMX_RELEASE_ASSERT(outerListCreationStep_ == step,
4284 "Outer list should be created at the same step as the inner list");
4287 /* Special performance logging stuff (env.var. GMX_NBNXN_CYCLE) */
4288 if (iLocality == InteractionLocality::Local)
4290 pairSearch->cycleCounting_.searchCount_++;
4292 if (pairSearch->cycleCounting_.recordCycles_ &&
4293 (!pairSearch->gridSet().domainSetup().haveMultipleDomains || iLocality == InteractionLocality::NonLocal) &&
4294 pairSearch->cycleCounting_.searchCount_ % 100 == 0)
4296 pairSearch->cycleCounting_.printCycles(stderr, pairSearch->work());
4301 nonbonded_verlet_t::constructPairlist(const Nbnxm::InteractionLocality iLocality,
4302 const t_blocka *excl,
4306 pairlistSets_->construct(iLocality, pairSearch_.get(), nbat.get(), excl,
4311 /* Launch the transfer of the pairlist to the GPU.
4313 * NOTE: The launch overhead is currently not timed separately
4315 Nbnxm::gpu_init_pairlist(gpu_nbv,
4316 pairlistSets().pairlistSet(iLocality).gpuList(),
4321 static void prepareListsForDynamicPruning(gmx::ArrayRef<NbnxnPairlistCpu> lists)
4323 /* TODO: Restructure the lists so we have actual outer and inner
4324 * list objects so we can set a single pointer instead of
4325 * swapping several pointers.
4328 for (auto &list : lists)
4330 /* The search produced a list in ci/cj.
4331 * Swap the list pointers so we get the outer list is ciOuter,cjOuter
4332 * and we can prune that to get an inner list in ci/cj.
4334 GMX_RELEASE_ASSERT(list.ciOuter.empty() && list.cjOuter.empty(),
4335 "The outer lists should be empty before preparation");
4337 std::swap(list.ci, list.ciOuter);
4338 std::swap(list.cj, list.cjOuter);