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/mdlib/ns.h"
55 #include "gromacs/mdtypes/group.h"
56 #include "gromacs/mdtypes/md_enums.h"
57 #include "gromacs/nbnxm/atomdata.h"
58 #include "gromacs/nbnxm/gpu_data_mgmt.h"
59 #include "gromacs/nbnxm/nbnxm.h"
60 #include "gromacs/nbnxm/nbnxm_geometry.h"
61 #include "gromacs/nbnxm/nbnxm_simd.h"
62 #include "gromacs/nbnxm/pairlistset.h"
63 #include "gromacs/pbcutil/ishift.h"
64 #include "gromacs/pbcutil/pbc.h"
65 #include "gromacs/simd/simd.h"
66 #include "gromacs/simd/vector_operations.h"
67 #include "gromacs/topology/block.h"
68 #include "gromacs/utility/exceptions.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/gmxomp.h"
71 #include "gromacs/utility/smalloc.h"
75 #include "pairlistwork.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;
94 static void nbs_cycle_clear(nbnxn_cycle_t *cc)
96 for (int i = 0; i < enbsCCnr; i++)
103 static double Mcyc_av(const nbnxn_cycle_t *cc)
105 return static_cast<double>(cc->c)*1e-6/cc->count;
108 static void nbs_cycle_print(FILE *fp, const nbnxn_search *nbs)
111 fprintf(fp, "ns %4d grid %4.1f search %4.1f red.f %5.3f",
112 nbs->cc[enbsCCgrid].count,
113 Mcyc_av(&nbs->cc[enbsCCgrid]),
114 Mcyc_av(&nbs->cc[enbsCCsearch]),
115 Mcyc_av(&nbs->cc[enbsCCreducef]));
117 if (nbs->work.size() > 1)
119 if (nbs->cc[enbsCCcombine].count > 0)
121 fprintf(fp, " comb %5.2f",
122 Mcyc_av(&nbs->cc[enbsCCcombine]));
124 fprintf(fp, " s. th");
125 for (const nbnxn_search_work_t &work : nbs->work)
127 fprintf(fp, " %4.1f",
128 Mcyc_av(&work.cc[enbsCCsearch]));
134 /* Layout for the nonbonded NxN pair lists */
135 enum class NbnxnLayout
137 NoSimd4x4, // i-cluster size 4, j-cluster size 4
138 Simd4xN, // i-cluster size 4, j-cluster size SIMD width
139 Simd2xNN, // i-cluster size 4, j-cluster size half SIMD width
140 Gpu8x8x8 // i-cluster size 8, j-cluster size 8 + super-clustering
144 /* Returns the j-cluster size */
145 template <NbnxnLayout layout>
146 static constexpr int jClusterSize()
148 static_assert(layout == NbnxnLayout::NoSimd4x4 || layout == NbnxnLayout::Simd4xN || layout == NbnxnLayout::Simd2xNN, "Currently jClusterSize only supports CPU layouts");
150 return layout == NbnxnLayout::Simd4xN ? GMX_SIMD_REAL_WIDTH : (layout == NbnxnLayout::Simd2xNN ? GMX_SIMD_REAL_WIDTH/2 : c_nbnxnCpuIClusterSize);
153 /*! \brief Returns the j-cluster index given the i-cluster index.
155 * \tparam jClusterSize The number of atoms in a j-cluster
156 * \tparam jSubClusterIndex The j-sub-cluster index (0/1), used when size(j-cluster) < size(i-cluster)
157 * \param[in] ci The i-cluster index
159 template <int jClusterSize, int jSubClusterIndex>
160 static inline int cjFromCi(int ci)
162 static_assert(jClusterSize == c_nbnxnCpuIClusterSize/2 || jClusterSize == c_nbnxnCpuIClusterSize || jClusterSize == c_nbnxnCpuIClusterSize*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
164 static_assert(jSubClusterIndex == 0 || jSubClusterIndex == 1,
165 "Only sub-cluster indices 0 and 1 are supported");
167 if (jClusterSize == c_nbnxnCpuIClusterSize/2)
169 if (jSubClusterIndex == 0)
175 return ((ci + 1) << 1) - 1;
178 else if (jClusterSize == c_nbnxnCpuIClusterSize)
188 /*! \brief Returns the j-cluster index given the i-cluster index.
190 * \tparam layout The pair-list layout
191 * \tparam jSubClusterIndex The j-sub-cluster index (0/1), used when size(j-cluster) < size(i-cluster)
192 * \param[in] ci The i-cluster index
194 template <NbnxnLayout layout, int jSubClusterIndex>
195 static inline int cjFromCi(int ci)
197 constexpr int clusterSize = jClusterSize<layout>();
199 return cjFromCi<clusterSize, jSubClusterIndex>(ci);
202 /* Returns the nbnxn coordinate data index given the i-cluster index */
203 template <NbnxnLayout layout>
204 static inline int xIndexFromCi(int ci)
206 constexpr int clusterSize = jClusterSize<layout>();
208 static_assert(clusterSize == c_nbnxnCpuIClusterSize/2 || clusterSize == c_nbnxnCpuIClusterSize || clusterSize == c_nbnxnCpuIClusterSize*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
210 if (clusterSize <= c_nbnxnCpuIClusterSize)
212 /* Coordinates are stored packed in groups of 4 */
217 /* Coordinates packed in 8, i-cluster size is half the packing width */
218 return (ci >> 1)*STRIDE_P8 + (ci & 1)*(c_packX8 >> 1);
222 /* Returns the nbnxn coordinate data index given the j-cluster index */
223 template <NbnxnLayout layout>
224 static inline int xIndexFromCj(int cj)
226 constexpr int clusterSize = jClusterSize<layout>();
228 static_assert(clusterSize == c_nbnxnCpuIClusterSize/2 || clusterSize == c_nbnxnCpuIClusterSize || clusterSize == c_nbnxnCpuIClusterSize*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
230 if (clusterSize == c_nbnxnCpuIClusterSize/2)
232 /* Coordinates are stored packed in groups of 4 */
233 return (cj >> 1)*STRIDE_P4 + (cj & 1)*(c_packX4 >> 1);
235 else if (clusterSize == c_nbnxnCpuIClusterSize)
237 /* Coordinates are stored packed in groups of 4 */
242 /* Coordinates are stored packed in groups of 8 */
248 /* Initializes a single nbnxn_pairlist_t data structure */
249 static void nbnxn_init_pairlist_fep(t_nblist *nl)
251 nl->type = GMX_NBLIST_INTERACTION_FREE_ENERGY;
252 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
253 /* The interaction functions are set in the free energy kernel fuction */
266 nl->jindex = nullptr;
268 nl->excl_fep = nullptr;
272 static void free_nblist(t_nblist *nl)
282 nbnxn_search_work_t::nbnxn_search_work_t() :
285 buffer_flags({0, nullptr, 0}),
287 nbl_fep(new t_nblist),
290 nbnxn_init_pairlist_fep(nbl_fep.get());
295 nbnxn_search_work_t::~nbnxn_search_work_t()
297 sfree(buffer_flags.flag);
299 free_nblist(nbl_fep.get());
302 nbnxn_search::nbnxn_search(const int ePBC,
303 const ivec *n_dd_cells,
304 const gmx_domdec_zones_t *zones,
305 const PairlistType pairlistType,
307 const int maxNumThreads) :
310 gridSet_(n_dd_cells, pairlistType, haveFep, maxNumThreads),
315 DomDec = n_dd_cells != nullptr;
318 for (int d = 0; d < DIM; d++)
320 if ((*n_dd_cells)[d] > 1)
327 /* Initialize detailed nbsearch cycle counting */
328 print_cycles = (getenv("GMX_NBNXN_CYCLE") != nullptr);
332 static void init_buffer_flags(nbnxn_buffer_flags_t *flags,
335 flags->nflag = (natoms + NBNXN_BUFFERFLAG_SIZE - 1)/NBNXN_BUFFERFLAG_SIZE;
336 if (flags->nflag > flags->flag_nalloc)
338 flags->flag_nalloc = over_alloc_large(flags->nflag);
339 srenew(flags->flag, flags->flag_nalloc);
341 for (int b = 0; b < flags->nflag; b++)
343 bitmask_clear(&(flags->flag[b]));
347 /* Returns the pair-list cutoff between a bounding box and a grid cell given an atom-to-atom pair-list cutoff
349 * Given a cutoff distance between atoms, this functions returns the cutoff
350 * distance2 between a bounding box of a group of atoms and a grid cell.
351 * Since atoms can be geometrically outside of the cell they have been
352 * assigned to (when atom groups instead of individual atoms are assigned
353 * to cells), this distance returned can be larger than the input.
356 listRangeForBoundingBoxToGridCell(real rlist,
357 const Grid::Dimensions &gridDims)
359 return rlist + gridDims.maxAtomGroupRadius;
362 /* Returns the pair-list cutoff between a grid cells given an atom-to-atom pair-list cutoff
364 * Given a cutoff distance between atoms, this functions returns the cutoff
365 * distance2 between two grid cells.
366 * Since atoms can be geometrically outside of the cell they have been
367 * assigned to (when atom groups instead of individual atoms are assigned
368 * to cells), this distance returned can be larger than the input.
371 listRangeForGridCellToGridCell(real rlist,
372 const Grid::Dimensions &iGridDims,
373 const Grid::Dimensions &jGridDims)
375 return rlist + iGridDims.maxAtomGroupRadius + jGridDims.maxAtomGroupRadius;
378 /* Determines the cell range along one dimension that
379 * the bounding box b0 - b1 sees.
382 static void get_cell_range(real b0, real b1,
383 const Grid::Dimensions &jGridDims,
384 real d2, real rlist, int *cf, int *cl)
386 real listRangeBBToCell2 = gmx::square(listRangeForBoundingBoxToGridCell(rlist, jGridDims));
387 real distanceInCells = (b0 - jGridDims.lowerCorner[dim])*jGridDims.invCellSize[dim];
388 *cf = std::max(static_cast<int>(distanceInCells), 0);
391 d2 + gmx::square((b0 - jGridDims.lowerCorner[dim]) - (*cf - 1 + 1)*jGridDims.cellSize[dim]) < listRangeBBToCell2)
396 *cl = std::min(static_cast<int>((b1 - jGridDims.lowerCorner[dim])*jGridDims.invCellSize[dim]), jGridDims.numCells[dim] - 1);
397 while (*cl < jGridDims.numCells[dim] - 1 &&
398 d2 + gmx::square((*cl + 1)*jGridDims.cellSize[dim] - (b1 - jGridDims.lowerCorner[dim])) < listRangeBBToCell2)
404 /* Reference code calculating the distance^2 between two bounding boxes */
406 static float box_dist2(float bx0, float bx1, float by0,
407 float by1, float bz0, float bz1,
408 const BoundingBox *bb)
411 float dl, dh, dm, dm0;
415 dl = bx0 - bb->upper.x;
416 dh = bb->lower.x - bx1;
417 dm = std::max(dl, dh);
418 dm0 = std::max(dm, 0.0f);
421 dl = by0 - bb->upper.y;
422 dh = bb->lower.y - by1;
423 dm = std::max(dl, dh);
424 dm0 = std::max(dm, 0.0f);
427 dl = bz0 - bb->upper.z;
428 dh = bb->lower.z - bz1;
429 dm = std::max(dl, dh);
430 dm0 = std::max(dm, 0.0f);
437 #if !NBNXN_SEARCH_BB_SIMD4
439 /*! \brief Plain C code calculating the distance^2 between two bounding boxes in xyz0 format
441 * \param[in] bb_i First bounding box
442 * \param[in] bb_j Second bounding box
444 static float clusterBoundingBoxDistance2(const BoundingBox &bb_i,
445 const BoundingBox &bb_j)
447 float dl = bb_i.lower.x - bb_j.upper.x;
448 float dh = bb_j.lower.x - bb_i.upper.x;
449 float dm = std::max(dl, dh);
450 float dm0 = std::max(dm, 0.0f);
453 dl = bb_i.lower.y - bb_j.upper.y;
454 dh = bb_j.lower.y - bb_i.upper.y;
455 dm = std::max(dl, dh);
456 dm0 = std::max(dm, 0.0f);
459 dl = bb_i.lower.z - bb_j.upper.z;
460 dh = bb_j.lower.z - bb_i.upper.z;
461 dm = std::max(dl, dh);
462 dm0 = std::max(dm, 0.0f);
468 #else /* NBNXN_SEARCH_BB_SIMD4 */
470 /*! \brief 4-wide SIMD code calculating the distance^2 between two bounding boxes in xyz0 format
472 * \param[in] bb_i First bounding box, should be aligned for 4-wide SIMD
473 * \param[in] bb_j Second bounding box, should be aligned for 4-wide SIMD
475 static float clusterBoundingBoxDistance2(const BoundingBox &bb_i,
476 const BoundingBox &bb_j)
478 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
481 const Simd4Float bb_i_S0 = load4(bb_i.lower.ptr());
482 const Simd4Float bb_i_S1 = load4(bb_i.upper.ptr());
483 const Simd4Float bb_j_S0 = load4(bb_j.lower.ptr());
484 const Simd4Float bb_j_S1 = load4(bb_j.upper.ptr());
486 const Simd4Float dl_S = bb_i_S0 - bb_j_S1;
487 const Simd4Float dh_S = bb_j_S0 - bb_i_S1;
489 const Simd4Float dm_S = max(dl_S, dh_S);
490 const Simd4Float dm0_S = max(dm_S, simd4SetZeroF());
492 return dotProduct(dm0_S, dm0_S);
495 /* Calculate bb bounding distances of bb_i[si,...,si+3] and store them in d2 */
496 template <int boundingBoxStart>
497 static inline void gmx_simdcall
498 clusterBoundingBoxDistance2_xxxx_simd4_inner(const float *bb_i,
500 const Simd4Float xj_l,
501 const Simd4Float yj_l,
502 const Simd4Float zj_l,
503 const Simd4Float xj_h,
504 const Simd4Float yj_h,
505 const Simd4Float zj_h)
507 constexpr int stride = c_packedBoundingBoxesDimSize;
509 const int shi = boundingBoxStart*Nbnxm::c_numBoundingBoxBounds1D*DIM;
511 const Simd4Float zero = setZero();
513 const Simd4Float xi_l = load4(bb_i + shi + 0*stride);
514 const Simd4Float yi_l = load4(bb_i + shi + 1*stride);
515 const Simd4Float zi_l = load4(bb_i + shi + 2*stride);
516 const Simd4Float xi_h = load4(bb_i + shi + 3*stride);
517 const Simd4Float yi_h = load4(bb_i + shi + 4*stride);
518 const Simd4Float zi_h = load4(bb_i + shi + 5*stride);
520 const Simd4Float dx_0 = xi_l - xj_h;
521 const Simd4Float dy_0 = yi_l - yj_h;
522 const Simd4Float dz_0 = zi_l - zj_h;
524 const Simd4Float dx_1 = xj_l - xi_h;
525 const Simd4Float dy_1 = yj_l - yi_h;
526 const Simd4Float dz_1 = zj_l - zi_h;
528 const Simd4Float mx = max(dx_0, dx_1);
529 const Simd4Float my = max(dy_0, dy_1);
530 const Simd4Float mz = max(dz_0, dz_1);
532 const Simd4Float m0x = max(mx, zero);
533 const Simd4Float m0y = max(my, zero);
534 const Simd4Float m0z = max(mz, zero);
536 const Simd4Float d2x = m0x * m0x;
537 const Simd4Float d2y = m0y * m0y;
538 const Simd4Float d2z = m0z * m0z;
540 const Simd4Float d2s = d2x + d2y;
541 const Simd4Float d2t = d2s + d2z;
543 store4(d2 + boundingBoxStart, d2t);
546 /* 4-wide SIMD code for nsi bb distances for bb format xxxxyyyyzzzz */
548 clusterBoundingBoxDistance2_xxxx_simd4(const float *bb_j,
553 constexpr int stride = c_packedBoundingBoxesDimSize;
555 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
558 const Simd4Float xj_l = Simd4Float(bb_j[0*stride]);
559 const Simd4Float yj_l = Simd4Float(bb_j[1*stride]);
560 const Simd4Float zj_l = Simd4Float(bb_j[2*stride]);
561 const Simd4Float xj_h = Simd4Float(bb_j[3*stride]);
562 const Simd4Float yj_h = Simd4Float(bb_j[4*stride]);
563 const Simd4Float zj_h = Simd4Float(bb_j[5*stride]);
565 /* Here we "loop" over si (0,stride) from 0 to nsi with step stride.
566 * But as we know the number of iterations is 1 or 2, we unroll manually.
568 clusterBoundingBoxDistance2_xxxx_simd4_inner<0>(bb_i, d2,
573 clusterBoundingBoxDistance2_xxxx_simd4_inner<stride>(bb_i, d2,
579 #endif /* NBNXN_SEARCH_BB_SIMD4 */
582 /* Returns if any atom pair from two clusters is within distance sqrt(rlist2) */
583 static inline gmx_bool
584 clusterpair_in_range(const NbnxnPairlistGpuWork &work,
586 int csj, int stride, const real *x_j,
589 #if !GMX_SIMD4_HAVE_REAL
592 * All coordinates are stored as xyzxyz...
595 const real *x_i = work.iSuperClusterData.x.data();
597 for (int i = 0; i < c_nbnxnGpuClusterSize; i++)
599 int i0 = (si*c_nbnxnGpuClusterSize + i)*DIM;
600 for (int j = 0; j < c_nbnxnGpuClusterSize; j++)
602 int j0 = (csj*c_nbnxnGpuClusterSize + j)*stride;
604 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]);
615 #else /* !GMX_SIMD4_HAVE_REAL */
617 /* 4-wide SIMD version.
618 * The coordinates x_i are stored as xxxxyyyy..., x_j is stored xyzxyz...
619 * Using 8-wide AVX(2) is not faster on Intel Sandy Bridge and Haswell.
621 static_assert(c_nbnxnGpuClusterSize == 8 || c_nbnxnGpuClusterSize == 4,
622 "A cluster is hard-coded to 4/8 atoms.");
624 Simd4Real rc2_S = Simd4Real(rlist2);
626 const real *x_i = work.iSuperClusterData.xSimd.data();
628 int dim_stride = c_nbnxnGpuClusterSize*DIM;
629 Simd4Real ix_S0 = load4(x_i + si*dim_stride + 0*GMX_SIMD4_WIDTH);
630 Simd4Real iy_S0 = load4(x_i + si*dim_stride + 1*GMX_SIMD4_WIDTH);
631 Simd4Real iz_S0 = load4(x_i + si*dim_stride + 2*GMX_SIMD4_WIDTH);
633 Simd4Real ix_S1, iy_S1, iz_S1;
634 if (c_nbnxnGpuClusterSize == 8)
636 ix_S1 = load4(x_i + si*dim_stride + 3*GMX_SIMD4_WIDTH);
637 iy_S1 = load4(x_i + si*dim_stride + 4*GMX_SIMD4_WIDTH);
638 iz_S1 = load4(x_i + si*dim_stride + 5*GMX_SIMD4_WIDTH);
640 /* We loop from the outer to the inner particles to maximize
641 * the chance that we find a pair in range quickly and return.
643 int j0 = csj*c_nbnxnGpuClusterSize;
644 int j1 = j0 + c_nbnxnGpuClusterSize - 1;
647 Simd4Real jx0_S, jy0_S, jz0_S;
648 Simd4Real jx1_S, jy1_S, jz1_S;
650 Simd4Real dx_S0, dy_S0, dz_S0;
651 Simd4Real dx_S1, dy_S1, dz_S1;
652 Simd4Real dx_S2, dy_S2, dz_S2;
653 Simd4Real dx_S3, dy_S3, dz_S3;
664 Simd4Bool wco_any_S01, wco_any_S23, wco_any_S;
666 jx0_S = Simd4Real(x_j[j0*stride+0]);
667 jy0_S = Simd4Real(x_j[j0*stride+1]);
668 jz0_S = Simd4Real(x_j[j0*stride+2]);
670 jx1_S = Simd4Real(x_j[j1*stride+0]);
671 jy1_S = Simd4Real(x_j[j1*stride+1]);
672 jz1_S = Simd4Real(x_j[j1*stride+2]);
674 /* Calculate distance */
675 dx_S0 = ix_S0 - jx0_S;
676 dy_S0 = iy_S0 - jy0_S;
677 dz_S0 = iz_S0 - jz0_S;
678 dx_S2 = ix_S0 - jx1_S;
679 dy_S2 = iy_S0 - jy1_S;
680 dz_S2 = iz_S0 - jz1_S;
681 if (c_nbnxnGpuClusterSize == 8)
683 dx_S1 = ix_S1 - jx0_S;
684 dy_S1 = iy_S1 - jy0_S;
685 dz_S1 = iz_S1 - jz0_S;
686 dx_S3 = ix_S1 - jx1_S;
687 dy_S3 = iy_S1 - jy1_S;
688 dz_S3 = iz_S1 - jz1_S;
691 /* rsq = dx*dx+dy*dy+dz*dz */
692 rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
693 rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
694 if (c_nbnxnGpuClusterSize == 8)
696 rsq_S1 = norm2(dx_S1, dy_S1, dz_S1);
697 rsq_S3 = norm2(dx_S3, dy_S3, dz_S3);
700 wco_S0 = (rsq_S0 < rc2_S);
701 wco_S2 = (rsq_S2 < rc2_S);
702 if (c_nbnxnGpuClusterSize == 8)
704 wco_S1 = (rsq_S1 < rc2_S);
705 wco_S3 = (rsq_S3 < rc2_S);
707 if (c_nbnxnGpuClusterSize == 8)
709 wco_any_S01 = wco_S0 || wco_S1;
710 wco_any_S23 = wco_S2 || wco_S3;
711 wco_any_S = wco_any_S01 || wco_any_S23;
715 wco_any_S = wco_S0 || wco_S2;
718 if (anyTrue(wco_any_S))
729 #endif /* !GMX_SIMD4_HAVE_REAL */
732 /* Returns the j-cluster index for index cjIndex in a cj list */
733 static inline int nblCj(gmx::ArrayRef<const nbnxn_cj_t> cjList,
736 return cjList[cjIndex].cj;
739 /* Returns the j-cluster index for index cjIndex in a cj4 list */
740 static inline int nblCj(gmx::ArrayRef<const nbnxn_cj4_t> cj4List,
743 return cj4List[cjIndex/c_nbnxnGpuJgroupSize].cj[cjIndex & (c_nbnxnGpuJgroupSize - 1)];
746 /* Returns the i-interaction mask of the j sub-cell for index cj_ind */
747 static unsigned int nbl_imask0(const NbnxnPairlistGpu *nbl, int cj_ind)
749 return nbl->cj4[cj_ind/c_nbnxnGpuJgroupSize].imei[0].imask;
752 /* Initializes a single NbnxnPairlistCpu data structure */
753 static void nbnxn_init_pairlist(NbnxnPairlistCpu *nbl)
755 nbl->na_ci = c_nbnxnCpuIClusterSize;
758 nbl->ciOuter.clear();
761 nbl->cjOuter.clear();
764 nbl->work = new NbnxnPairlistCpuWork();
767 NbnxnPairlistGpu::NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy) :
768 na_ci(c_nbnxnGpuClusterSize),
769 na_cj(c_nbnxnGpuClusterSize),
770 na_sc(c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize),
772 sci({}, {pinningPolicy}),
773 cj4({}, {pinningPolicy}),
774 excl({}, {pinningPolicy}),
777 static_assert(c_nbnxnGpuNumClusterPerSupercluster == c_gpuNumClusterPerCell,
778 "The search code assumes that the a super-cluster matches a search grid cell");
780 static_assert(sizeof(cj4[0].imei[0].imask)*8 >= c_nbnxnGpuJgroupSize*c_gpuNumClusterPerCell,
781 "The i super-cluster cluster interaction mask does not contain a sufficient number of bits");
783 static_assert(sizeof(excl[0])*8 >= c_nbnxnGpuJgroupSize*c_gpuNumClusterPerCell, "The GPU exclusion mask does not contain a sufficient number of bits");
785 // We always want a first entry without any exclusions
788 work = new NbnxnPairlistGpuWork();
791 void nbnxn_init_pairlist_set(nbnxn_pairlist_set_t *nbl_list)
794 (nbl_list->params.pairlistType == PairlistType::Simple4x2 ||
795 nbl_list->params.pairlistType == PairlistType::Simple4x4 ||
796 nbl_list->params.pairlistType == PairlistType::Simple4x8);
797 // Currently GPU lists are always combined
798 nbl_list->bCombined = !nbl_list->bSimple;
800 nbl_list->nnbl = gmx_omp_nthreads_get(emntNonbonded);
802 if (!nbl_list->bCombined &&
803 nbl_list->nnbl > NBNXN_BUFFERFLAG_MAX_THREADS)
805 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.",
806 nbl_list->nnbl, NBNXN_BUFFERFLAG_MAX_THREADS, NBNXN_BUFFERFLAG_MAX_THREADS);
809 if (nbl_list->bSimple)
811 snew(nbl_list->nbl, nbl_list->nnbl);
812 if (nbl_list->nnbl > 1)
814 snew(nbl_list->nbl_work, nbl_list->nnbl);
819 snew(nbl_list->nblGpu, nbl_list->nnbl);
821 nbl_list->nbl_fep.resize(nbl_list->nnbl);
822 /* Execute in order to avoid memory interleaving between threads */
823 #pragma omp parallel for num_threads(nbl_list->nnbl) schedule(static)
824 for (int i = 0; i < nbl_list->nnbl; i++)
828 /* Allocate the nblist data structure locally on each thread
829 * to optimize memory access for NUMA architectures.
831 if (nbl_list->bSimple)
833 nbl_list->nbl[i] = new NbnxnPairlistCpu();
835 nbnxn_init_pairlist(nbl_list->nbl[i]);
836 if (nbl_list->nnbl > 1)
838 nbl_list->nbl_work[i] = new NbnxnPairlistCpu();
839 nbnxn_init_pairlist(nbl_list->nbl_work[i]);
844 /* Only list 0 is used on the GPU, use normal allocation for i>0 */
845 auto pinningPolicy = (i == 0 ? gmx::PinningPolicy::PinnedIfSupported : gmx::PinningPolicy::CannotBePinned);
847 nbl_list->nblGpu[i] = new NbnxnPairlistGpu(pinningPolicy);
850 snew(nbl_list->nbl_fep[i], 1);
851 nbnxn_init_pairlist_fep(nbl_list->nbl_fep[i]);
853 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
857 /* Print statistics of a pair list, used for debug output */
858 static void print_nblist_statistics(FILE *fp, const NbnxnPairlistCpu *nbl,
859 const nbnxn_search *nbs, real rl)
861 const Grid &grid = nbs->gridSet().grids()[0];
862 const Grid::Dimensions &dims = grid.dimensions();
864 fprintf(fp, "nbl nci %zu ncj %d\n",
865 nbl->ci.size(), nbl->ncjInUse);
866 const int numAtomsJCluster = grid.geometry().numAtomsJCluster;
867 const double numAtomsPerCell = nbl->ncjInUse/static_cast<double>(grid.numCells())*numAtomsJCluster;
868 fprintf(fp, "nbl na_cj %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
869 nbl->na_cj, rl, nbl->ncjInUse, nbl->ncjInUse/static_cast<double>(grid.numCells()),
871 numAtomsPerCell/(0.5*4.0/3.0*M_PI*rl*rl*rl*grid.numCells()*numAtomsJCluster/(dims.gridSize[XX]*dims.gridSize[YY]*dims.gridSize[ZZ])));
873 fprintf(fp, "nbl average j cell list length %.1f\n",
874 0.25*nbl->ncjInUse/std::max(static_cast<double>(nbl->ci.size()), 1.0));
876 int cs[SHIFTS] = { 0 };
878 for (const nbnxn_ci_t &ciEntry : nbl->ci)
880 cs[ciEntry.shift & NBNXN_CI_SHIFT] +=
881 ciEntry.cj_ind_end - ciEntry.cj_ind_start;
883 int j = ciEntry.cj_ind_start;
884 while (j < ciEntry.cj_ind_end &&
885 nbl->cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
891 fprintf(fp, "nbl cell pairs, total: %zu excl: %d %.1f%%\n",
892 nbl->cj.size(), npexcl, 100*npexcl/std::max(static_cast<double>(nbl->cj.size()), 1.0));
893 for (int s = 0; s < SHIFTS; s++)
897 fprintf(fp, "nbl shift %2d ncj %3d\n", s, cs[s]);
902 /* Print statistics of a pair lists, used for debug output */
903 static void print_nblist_statistics(FILE *fp, const NbnxnPairlistGpu *nbl,
904 const nbnxn_search *nbs, real rl)
906 const Grid &grid = nbs->gridSet().grids()[0];
907 const Grid::Dimensions &dims = grid.dimensions();
909 fprintf(fp, "nbl nsci %zu ncj4 %zu nsi %d excl4 %zu\n",
910 nbl->sci.size(), nbl->cj4.size(), nbl->nci_tot, nbl->excl.size());
911 const int numAtomsCluster = grid.geometry().numAtomsICluster;
912 const double numAtomsPerCell = nbl->nci_tot/static_cast<double>(grid.numClusters())*numAtomsCluster;
913 fprintf(fp, "nbl na_c %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
914 nbl->na_ci, rl, nbl->nci_tot, nbl->nci_tot/static_cast<double>(grid.numClusters()),
916 numAtomsPerCell/(0.5*4.0/3.0*M_PI*rl*rl*rl*grid.numClusters()*numAtomsCluster/(dims.gridSize[XX]*dims.gridSize[YY]*dims.gridSize[ZZ])));
921 int c[c_gpuNumClusterPerCell + 1] = { 0 };
922 for (const nbnxn_sci_t &sci : nbl->sci)
925 for (int j4 = sci.cj4_ind_start; j4 < sci.cj4_ind_end; j4++)
927 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
930 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
932 if (nbl->cj4[j4].imei[0].imask & (1U << (j*c_gpuNumClusterPerCell + si)))
943 nsp_max = std::max(nsp_max, nsp);
945 if (!nbl->sci.empty())
947 sum_nsp /= nbl->sci.size();
948 sum_nsp2 /= nbl->sci.size();
950 fprintf(fp, "nbl #cluster-pairs: av %.1f stddev %.1f max %d\n",
951 sum_nsp, std::sqrt(sum_nsp2 - sum_nsp*sum_nsp), nsp_max);
953 if (!nbl->cj4.empty())
955 for (int b = 0; b <= c_gpuNumClusterPerCell; b++)
957 fprintf(fp, "nbl j-list #i-subcell %d %7d %4.1f\n",
958 b, c[b], 100.0*c[b]/size_t {nbl->cj4.size()*c_nbnxnGpuJgroupSize});
963 /* Returns a pointer to the exclusion mask for j-cluster-group \p cj4 and warp \p warp
964 * Generates a new exclusion entry when the j-cluster-group uses
965 * the default all-interaction mask at call time, so the returned mask
966 * can be modified when needed.
968 static nbnxn_excl_t *get_exclusion_mask(NbnxnPairlistGpu *nbl,
972 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
974 /* No exclusions set, make a new list entry */
975 const size_t oldSize = nbl->excl.size();
976 GMX_ASSERT(oldSize >= 1, "We should always have entry [0]");
977 /* Add entry with default values: no exclusions */
978 nbl->excl.resize(oldSize + 1);
979 nbl->cj4[cj4].imei[warp].excl_ind = oldSize;
982 return &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
985 static void set_self_and_newton_excls_supersub(NbnxnPairlistGpu *nbl,
986 int cj4_ind, int sj_offset,
987 int i_cluster_in_cell)
989 nbnxn_excl_t *excl[c_nbnxnGpuClusterpairSplit];
991 /* Here we only set the set self and double pair exclusions */
993 /* Reserve extra elements, so the resize() in get_exclusion_mask()
994 * will not invalidate excl entries in the loop below
996 nbl->excl.reserve(nbl->excl.size() + c_nbnxnGpuClusterpairSplit);
997 for (int w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
999 excl[w] = get_exclusion_mask(nbl, cj4_ind, w);
1002 /* Only minor < major bits set */
1003 for (int ej = 0; ej < nbl->na_ci; ej++)
1006 for (int ei = ej; ei < nbl->na_ci; ei++)
1008 excl[w]->pair[(ej & (c_nbnxnGpuJgroupSize-1))*nbl->na_ci + ei] &=
1009 ~(1U << (sj_offset*c_gpuNumClusterPerCell + i_cluster_in_cell));
1014 /* Returns a diagonal or off-diagonal interaction mask for plain C lists */
1015 static unsigned int get_imask(gmx_bool rdiag, int ci, int cj)
1017 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
1020 /* Returns a diagonal or off-diagonal interaction mask for cj-size=2 */
1021 gmx_unused static unsigned int get_imask_simd_j2(gmx_bool rdiag, int ci, int cj)
1023 return (rdiag && ci*2 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_0 :
1024 (rdiag && ci*2+1 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_1 :
1025 NBNXN_INTERACTION_MASK_ALL));
1028 /* Returns a diagonal or off-diagonal interaction mask for cj-size=4 */
1029 gmx_unused static unsigned int get_imask_simd_j4(gmx_bool rdiag, int ci, int cj)
1031 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
1034 /* Returns a diagonal or off-diagonal interaction mask for cj-size=8 */
1035 gmx_unused static unsigned int get_imask_simd_j8(gmx_bool rdiag, int ci, int cj)
1037 return (rdiag && ci == cj*2 ? NBNXN_INTERACTION_MASK_DIAG_J8_0 :
1038 (rdiag && ci == cj*2+1 ? NBNXN_INTERACTION_MASK_DIAG_J8_1 :
1039 NBNXN_INTERACTION_MASK_ALL));
1043 #if GMX_SIMD_REAL_WIDTH == 2
1044 #define get_imask_simd_4xn get_imask_simd_j2
1046 #if GMX_SIMD_REAL_WIDTH == 4
1047 #define get_imask_simd_4xn get_imask_simd_j4
1049 #if GMX_SIMD_REAL_WIDTH == 8
1050 #define get_imask_simd_4xn get_imask_simd_j8
1051 #define get_imask_simd_2xnn get_imask_simd_j4
1053 #if GMX_SIMD_REAL_WIDTH == 16
1054 #define get_imask_simd_2xnn get_imask_simd_j8
1058 /* Plain C code for checking and adding cluster-pairs to the list.
1060 * \param[in] gridj The j-grid
1061 * \param[in,out] nbl The pair-list to store the cluster pairs in
1062 * \param[in] icluster The index of the i-cluster
1063 * \param[in] jclusterFirst The first cluster in the j-range
1064 * \param[in] jclusterLast The last cluster in the j-range
1065 * \param[in] excludeSubDiagonal Exclude atom pairs with i-index > j-index
1066 * \param[in] x_j Coordinates for the j-atom, in xyz format
1067 * \param[in] rlist2 The squared list cut-off
1068 * \param[in] rbb2 The squared cut-off for putting cluster-pairs in the list based on bounding box distance only
1069 * \param[in,out] numDistanceChecks The number of distance checks performed
1072 makeClusterListSimple(const Grid &jGrid,
1073 NbnxnPairlistCpu * nbl,
1077 bool excludeSubDiagonal,
1078 const real * gmx_restrict x_j,
1081 int * gmx_restrict numDistanceChecks)
1083 const BoundingBox * gmx_restrict bb_ci = nbl->work->iClusterData.bb.data();
1084 const real * gmx_restrict x_ci = nbl->work->iClusterData.x.data();
1089 while (!InRange && jclusterFirst <= jclusterLast)
1091 real d2 = clusterBoundingBoxDistance2(bb_ci[0], jGrid.jBoundingBoxes()[jclusterFirst]);
1092 *numDistanceChecks += 2;
1094 /* Check if the distance is within the distance where
1095 * we use only the bounding box distance rbb,
1096 * or within the cut-off and there is at least one atom pair
1097 * within the cut-off.
1103 else if (d2 < rlist2)
1105 int cjf_gl = jGrid.cellOffset() + jclusterFirst;
1106 for (int i = 0; i < c_nbnxnCpuIClusterSize && !InRange; i++)
1108 for (int j = 0; j < c_nbnxnCpuIClusterSize; j++)
1110 InRange = InRange ||
1111 (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjf_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+XX]) +
1112 gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjf_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+YY]) +
1113 gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjf_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+ZZ]) < rlist2);
1116 *numDistanceChecks += c_nbnxnCpuIClusterSize*c_nbnxnCpuIClusterSize;
1129 while (!InRange && jclusterLast > jclusterFirst)
1131 real d2 = clusterBoundingBoxDistance2(bb_ci[0], jGrid.jBoundingBoxes()[jclusterLast]);
1132 *numDistanceChecks += 2;
1134 /* Check if the distance is within the distance where
1135 * we use only the bounding box distance rbb,
1136 * or within the cut-off and there is at least one atom pair
1137 * within the cut-off.
1143 else if (d2 < rlist2)
1145 int cjl_gl = jGrid.cellOffset() + jclusterLast;
1146 for (int i = 0; i < c_nbnxnCpuIClusterSize && !InRange; i++)
1148 for (int j = 0; j < c_nbnxnCpuIClusterSize; j++)
1150 InRange = InRange ||
1151 (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjl_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+XX]) +
1152 gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjl_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+YY]) +
1153 gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjl_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+ZZ]) < rlist2);
1156 *numDistanceChecks += c_nbnxnCpuIClusterSize*c_nbnxnCpuIClusterSize;
1164 if (jclusterFirst <= jclusterLast)
1166 for (int jcluster = jclusterFirst; jcluster <= jclusterLast; jcluster++)
1168 /* Store cj and the interaction mask */
1170 cjEntry.cj = jGrid.cellOffset() + jcluster;
1171 cjEntry.excl = get_imask(excludeSubDiagonal, icluster, jcluster);
1172 nbl->cj.push_back(cjEntry);
1174 /* Increase the closing index in the i list */
1175 nbl->ci.back().cj_ind_end = nbl->cj.size();
1179 #ifdef GMX_NBNXN_SIMD_4XN
1180 #include "gromacs/nbnxm/pairlist_simd_4xm.h"
1182 #ifdef GMX_NBNXN_SIMD_2XNN
1183 #include "gromacs/nbnxm/pairlist_simd_2xmm.h"
1186 /* Plain C or SIMD4 code for making a pair list of super-cell sci vs scj.
1187 * Checks bounding box distances and possibly atom pair distances.
1189 static void make_cluster_list_supersub(const Grid &iGrid,
1191 NbnxnPairlistGpu *nbl,
1194 const bool excludeSubDiagonal,
1199 int *numDistanceChecks)
1201 NbnxnPairlistGpuWork &work = *nbl->work;
1204 const float *pbb_ci = work.iSuperClusterData.bbPacked.data();
1206 const BoundingBox *bb_ci = work.iSuperClusterData.bb.data();
1209 assert(c_nbnxnGpuClusterSize == iGrid.geometry().numAtomsICluster);
1210 assert(c_nbnxnGpuClusterSize == jGrid.geometry().numAtomsICluster);
1212 /* We generate the pairlist mainly based on bounding-box distances
1213 * and do atom pair distance based pruning on the GPU.
1214 * Only if a j-group contains a single cluster-pair, we try to prune
1215 * that pair based on atom distances on the CPU to avoid empty j-groups.
1217 #define PRUNE_LIST_CPU_ONE 1
1218 #define PRUNE_LIST_CPU_ALL 0
1220 #if PRUNE_LIST_CPU_ONE
1224 float *d2l = work.distanceBuffer.data();
1226 for (int subc = 0; subc < jGrid.numClustersPerCell()[scj]; subc++)
1228 const int cj4_ind = work.cj_ind/c_nbnxnGpuJgroupSize;
1229 const int cj_offset = work.cj_ind - cj4_ind*c_nbnxnGpuJgroupSize;
1230 const int cj = scj*c_gpuNumClusterPerCell + subc;
1232 const int cj_gl = jGrid.cellOffset()*c_gpuNumClusterPerCell + cj;
1235 if (excludeSubDiagonal && sci == scj)
1241 ci1 = iGrid.numClustersPerCell()[sci];
1245 /* Determine all ci1 bb distances in one call with SIMD4 */
1246 const int offset = packedBoundingBoxesIndex(cj) + (cj & (c_packedBoundingBoxesDimSize - 1));
1247 clusterBoundingBoxDistance2_xxxx_simd4(jGrid.packedBoundingBoxes().data() + offset,
1249 *numDistanceChecks += c_nbnxnGpuClusterSize*2;
1253 unsigned int imask = 0;
1254 /* We use a fixed upper-bound instead of ci1 to help optimization */
1255 for (int ci = 0; ci < c_gpuNumClusterPerCell; ci++)
1263 /* Determine the bb distance between ci and cj */
1264 d2l[ci] = clusterBoundingBoxDistance2(bb_ci[ci], jGrid.jBoundingBoxes()[cj]);
1265 *numDistanceChecks += 2;
1269 #if PRUNE_LIST_CPU_ALL
1270 /* Check if the distance is within the distance where
1271 * we use only the bounding box distance rbb,
1272 * or within the cut-off and there is at least one atom pair
1273 * within the cut-off. This check is very costly.
1275 *numDistanceChecks += c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize;
1278 clusterpair_in_range(work, ci, cj_gl, stride, x, rlist2)))
1280 /* Check if the distance between the two bounding boxes
1281 * in within the pair-list cut-off.
1286 /* Flag this i-subcell to be taken into account */
1287 imask |= (1U << (cj_offset*c_gpuNumClusterPerCell + ci));
1289 #if PRUNE_LIST_CPU_ONE
1297 #if PRUNE_LIST_CPU_ONE
1298 /* If we only found 1 pair, check if any atoms are actually
1299 * within the cut-off, so we could get rid of it.
1301 if (npair == 1 && d2l[ci_last] >= rbb2 &&
1302 !clusterpair_in_range(work, ci_last, cj_gl, stride, x, rlist2))
1304 imask &= ~(1U << (cj_offset*c_gpuNumClusterPerCell + ci_last));
1311 /* We have at least one cluster pair: add a j-entry */
1312 if (static_cast<size_t>(cj4_ind) == nbl->cj4.size())
1314 nbl->cj4.resize(nbl->cj4.size() + 1);
1316 nbnxn_cj4_t *cj4 = &nbl->cj4[cj4_ind];
1318 cj4->cj[cj_offset] = cj_gl;
1320 /* Set the exclusions for the ci==sj entry.
1321 * Here we don't bother to check if this entry is actually flagged,
1322 * as it will nearly always be in the list.
1324 if (excludeSubDiagonal && sci == scj)
1326 set_self_and_newton_excls_supersub(nbl, cj4_ind, cj_offset, subc);
1329 /* Copy the cluster interaction mask to the list */
1330 for (int w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
1332 cj4->imei[w].imask |= imask;
1335 nbl->work->cj_ind++;
1337 /* Keep the count */
1338 nbl->nci_tot += npair;
1340 /* Increase the closing index in i super-cell list */
1341 nbl->sci.back().cj4_ind_end =
1342 (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize;
1347 /* Returns how many contiguous j-clusters we have starting in the i-list */
1348 template <typename CjListType>
1349 static int numContiguousJClusters(const int cjIndexStart,
1350 const int cjIndexEnd,
1351 gmx::ArrayRef<const CjListType> cjList)
1353 const int firstJCluster = nblCj(cjList, cjIndexStart);
1355 int numContiguous = 0;
1357 while (cjIndexStart + numContiguous < cjIndexEnd &&
1358 nblCj(cjList, cjIndexStart + numContiguous) == firstJCluster + numContiguous)
1363 return numContiguous;
1367 * \brief Helper struct for efficient searching for excluded atoms in a j-list
1371 /*! \brief Constructs a j-list range from \p cjList with the given index range */
1372 template <typename CjListType>
1373 JListRanges(int cjIndexStart,
1375 gmx::ArrayRef<const CjListType> cjList);
1377 int cjIndexStart; //!< The start index in the j-list
1378 int cjIndexEnd; //!< The end index in the j-list
1379 int cjFirst; //!< The j-cluster with index cjIndexStart
1380 int cjLast; //!< The j-cluster with index cjIndexEnd-1
1381 int numDirect; //!< Up to cjIndexStart+numDirect the j-clusters are cjFirst + the index offset
1385 template <typename CjListType>
1386 JListRanges::JListRanges(int cjIndexStart,
1388 gmx::ArrayRef<const CjListType> cjList) :
1389 cjIndexStart(cjIndexStart),
1390 cjIndexEnd(cjIndexEnd)
1392 GMX_ASSERT(cjIndexEnd > cjIndexStart, "JListRanges should only be called with non-empty lists");
1394 cjFirst = nblCj(cjList, cjIndexStart);
1395 cjLast = nblCj(cjList, cjIndexEnd - 1);
1397 /* Determine how many contiguous j-cells we have starting
1398 * from the first i-cell. This number can be used to directly
1399 * calculate j-cell indices for excluded atoms.
1401 numDirect = numContiguousJClusters(cjIndexStart, cjIndexEnd, cjList);
1405 /* Return the index of \p jCluster in the given range or -1 when not present
1407 * Note: This code is executed very often and therefore performance is
1408 * important. It should be inlined and fully optimized.
1410 template <typename CjListType>
1412 findJClusterInJList(int jCluster,
1413 const JListRanges &ranges,
1414 gmx::ArrayRef<const CjListType> cjList)
1418 if (jCluster < ranges.cjFirst + ranges.numDirect)
1420 /* We can calculate the index directly using the offset */
1421 index = ranges.cjIndexStart + jCluster - ranges.cjFirst;
1425 /* Search for jCluster using bisection */
1427 int rangeStart = ranges.cjIndexStart + ranges.numDirect;
1428 int rangeEnd = ranges.cjIndexEnd;
1430 while (index == -1 && rangeStart < rangeEnd)
1432 rangeMiddle = (rangeStart + rangeEnd) >> 1;
1434 const int clusterMiddle = nblCj(cjList, rangeMiddle);
1436 if (jCluster == clusterMiddle)
1438 index = rangeMiddle;
1440 else if (jCluster < clusterMiddle)
1442 rangeEnd = rangeMiddle;
1446 rangeStart = rangeMiddle + 1;
1454 // TODO: Get rid of the two functions below by renaming sci to ci (or something better)
1456 /* Return the i-entry in the list we are currently operating on */
1457 static nbnxn_ci_t *getOpenIEntry(NbnxnPairlistCpu *nbl)
1459 return &nbl->ci.back();
1462 /* Return the i-entry in the list we are currently operating on */
1463 static nbnxn_sci_t *getOpenIEntry(NbnxnPairlistGpu *nbl)
1465 return &nbl->sci.back();
1468 /* Set all atom-pair exclusions for a simple type list i-entry
1470 * Set all atom-pair exclusions from the topology stored in exclusions
1471 * as masks in the pair-list for simple list entry iEntry.
1474 setExclusionsForIEntry(const Nbnxm::GridSet &gridSet,
1475 NbnxnPairlistCpu *nbl,
1476 gmx_bool diagRemoved,
1478 const nbnxn_ci_t &iEntry,
1479 const t_blocka &exclusions)
1481 if (iEntry.cj_ind_end == iEntry.cj_ind_start)
1483 /* Empty list: no exclusions */
1487 const JListRanges ranges(iEntry.cj_ind_start, iEntry.cj_ind_end, gmx::makeConstArrayRef(nbl->cj));
1489 const int iCluster = iEntry.ci;
1491 gmx::ArrayRef<const int> cell = gridSet.cells();
1492 gmx::ArrayRef<const int> atomIndices = gridSet.atomIndices();
1494 /* Loop over the atoms in the i-cluster */
1495 for (int i = 0; i < nbl->na_ci; i++)
1497 const int iIndex = iCluster*nbl->na_ci + i;
1498 const int iAtom = atomIndices[iIndex];
1501 /* Loop over the topology-based exclusions for this i-atom */
1502 for (int exclIndex = exclusions.index[iAtom]; exclIndex < exclusions.index[iAtom + 1]; exclIndex++)
1504 const int jAtom = exclusions.a[exclIndex];
1508 /* The self exclusion are already set, save some time */
1512 /* Get the index of the j-atom in the nbnxn atom data */
1513 const int jIndex = cell[jAtom];
1515 /* Without shifts we only calculate interactions j>i
1516 * for one-way pair-lists.
1518 if (diagRemoved && jIndex <= iIndex)
1523 const int jCluster = (jIndex >> na_cj_2log);
1525 /* Could the cluster se be in our list? */
1526 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
1529 findJClusterInJList(jCluster, ranges,
1530 gmx::makeConstArrayRef(nbl->cj));
1534 /* We found an exclusion, clear the corresponding
1537 const int innerJ = jIndex - (jCluster << na_cj_2log);
1539 nbl->cj[index].excl &= ~(1U << ((i << na_cj_2log) + innerJ));
1547 /* Add a new i-entry to the FEP list and copy the i-properties */
1548 static inline void fep_list_new_nri_copy(t_nblist *nlist)
1550 /* Add a new i-entry */
1553 assert(nlist->nri < nlist->maxnri);
1555 /* Duplicate the last i-entry, except for jindex, which continues */
1556 nlist->iinr[nlist->nri] = nlist->iinr[nlist->nri-1];
1557 nlist->shift[nlist->nri] = nlist->shift[nlist->nri-1];
1558 nlist->gid[nlist->nri] = nlist->gid[nlist->nri-1];
1559 nlist->jindex[nlist->nri] = nlist->nrj;
1562 /* For load balancing of the free-energy lists over threads, we set
1563 * the maximum nrj size of an i-entry to 40. This leads to good
1564 * load balancing in the worst case scenario of a single perturbed
1565 * particle on 16 threads, while not introducing significant overhead.
1566 * Note that half of the perturbed pairs will anyhow end up in very small lists,
1567 * since non perturbed i-particles will see few perturbed j-particles).
1569 const int max_nrj_fep = 40;
1571 /* Exclude the perturbed pairs from the Verlet list. This is only done to avoid
1572 * singularities for overlapping particles (0/0), since the charges and
1573 * LJ parameters have been zeroed in the nbnxn data structure.
1574 * Simultaneously make a group pair list for the perturbed pairs.
1576 static void make_fep_list(gmx::ArrayRef<const int> atomIndices,
1577 const nbnxn_atomdata_t *nbat,
1578 NbnxnPairlistCpu *nbl,
1579 gmx_bool bDiagRemoved,
1581 real gmx_unused shx,
1582 real gmx_unused shy,
1583 real gmx_unused shz,
1584 real gmx_unused rlist_fep2,
1589 int ci, cj_ind_start, cj_ind_end, cja, cjr;
1591 int gid_i = 0, gid_j, gid;
1592 int egp_shift, egp_mask;
1594 int ind_i, ind_j, ai, aj;
1596 gmx_bool bFEP_i, bFEP_i_all;
1598 if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
1606 cj_ind_start = nbl_ci->cj_ind_start;
1607 cj_ind_end = nbl_ci->cj_ind_end;
1609 /* In worst case we have alternating energy groups
1610 * and create #atom-pair lists, which means we need the size
1611 * of a cluster pair (na_ci*na_cj) times the number of cj's.
1613 nri_max = nbl->na_ci*nbl->na_cj*(cj_ind_end - cj_ind_start);
1614 if (nlist->nri + nri_max > nlist->maxnri)
1616 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1617 reallocate_nblist(nlist);
1620 const int numAtomsJCluster = jGrid.geometry().numAtomsJCluster;
1622 const nbnxn_atomdata_t::Params &nbatParams = nbat->params();
1624 const int ngid = nbatParams.nenergrp;
1626 /* TODO: Consider adding a check in grompp and changing this to an assert */
1627 const int numBitsInEnergyGroupIdsForAtomsInJCluster = sizeof(gid_cj)*8;
1628 if (ngid*numAtomsJCluster > numBitsInEnergyGroupIdsForAtomsInJCluster)
1630 gmx_fatal(FARGS, "The Verlet scheme with %dx%d kernels and free-energy only supports up to %zu energy groups",
1631 iGrid.geometry().numAtomsICluster, numAtomsJCluster,
1632 (sizeof(gid_cj)*8)/numAtomsJCluster);
1635 egp_shift = nbatParams.neg_2log;
1636 egp_mask = (1 << egp_shift) - 1;
1638 /* Loop over the atoms in the i sub-cell */
1640 for (int i = 0; i < nbl->na_ci; i++)
1642 ind_i = ci*nbl->na_ci + i;
1643 ai = atomIndices[ind_i];
1647 nlist->jindex[nri+1] = nlist->jindex[nri];
1648 nlist->iinr[nri] = ai;
1649 /* The actual energy group pair index is set later */
1650 nlist->gid[nri] = 0;
1651 nlist->shift[nri] = nbl_ci->shift & NBNXN_CI_SHIFT;
1653 bFEP_i = iGrid.atomIsPerturbed(ci - iGrid.cellOffset(), i);
1655 bFEP_i_all = bFEP_i_all && bFEP_i;
1657 if (nlist->nrj + (cj_ind_end - cj_ind_start)*nbl->na_cj > nlist->maxnrj)
1659 nlist->maxnrj = over_alloc_small(nlist->nrj + (cj_ind_end - cj_ind_start)*nbl->na_cj);
1660 srenew(nlist->jjnr, nlist->maxnrj);
1661 srenew(nlist->excl_fep, nlist->maxnrj);
1666 gid_i = (nbatParams.energrp[ci] >> (egp_shift*i)) & egp_mask;
1669 for (int cj_ind = cj_ind_start; cj_ind < cj_ind_end; cj_ind++)
1671 unsigned int fep_cj;
1673 cja = nbl->cj[cj_ind].cj;
1675 if (numAtomsJCluster == jGrid.geometry().numAtomsICluster)
1677 cjr = cja - jGrid.cellOffset();
1678 fep_cj = jGrid.fepBits(cjr);
1681 gid_cj = nbatParams.energrp[cja];
1684 else if (2*numAtomsJCluster == jGrid.geometry().numAtomsICluster)
1686 cjr = cja - jGrid.cellOffset()*2;
1687 /* Extract half of the ci fep/energrp mask */
1688 fep_cj = (jGrid.fepBits(cjr >> 1) >> ((cjr & 1)*numAtomsJCluster)) & ((1 << numAtomsJCluster) - 1);
1691 gid_cj = nbatParams.energrp[cja >> 1] >> ((cja & 1)*numAtomsJCluster*egp_shift) & ((1 << (numAtomsJCluster*egp_shift)) - 1);
1696 cjr = cja - (jGrid.cellOffset() >> 1);
1697 /* Combine two ci fep masks/energrp */
1698 fep_cj = jGrid.fepBits(cjr*2) + (jGrid.fepBits(cjr*2 + 1) << jGrid.geometry().numAtomsICluster);
1701 gid_cj = nbatParams.energrp[cja*2] + (nbatParams.energrp[cja*2+1] << (jGrid.geometry().numAtomsICluster*egp_shift));
1705 if (bFEP_i || fep_cj != 0)
1707 for (int j = 0; j < nbl->na_cj; j++)
1709 /* Is this interaction perturbed and not excluded? */
1710 ind_j = cja*nbl->na_cj + j;
1711 aj = atomIndices[ind_j];
1713 (bFEP_i || (fep_cj & (1 << j))) &&
1714 (!bDiagRemoved || ind_j >= ind_i))
1718 gid_j = (gid_cj >> (j*egp_shift)) & egp_mask;
1719 gid = GID(gid_i, gid_j, ngid);
1721 if (nlist->nrj > nlist->jindex[nri] &&
1722 nlist->gid[nri] != gid)
1724 /* Energy group pair changed: new list */
1725 fep_list_new_nri_copy(nlist);
1728 nlist->gid[nri] = gid;
1731 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1733 fep_list_new_nri_copy(nlist);
1737 /* Add it to the FEP list */
1738 nlist->jjnr[nlist->nrj] = aj;
1739 nlist->excl_fep[nlist->nrj] = (nbl->cj[cj_ind].excl >> (i*nbl->na_cj + j)) & 1;
1742 /* Exclude it from the normal list.
1743 * Note that the charge has been set to zero,
1744 * but we need to avoid 0/0, as perturbed atoms
1745 * can be on top of each other.
1747 nbl->cj[cj_ind].excl &= ~(1U << (i*nbl->na_cj + j));
1753 if (nlist->nrj > nlist->jindex[nri])
1755 /* Actually add this new, non-empty, list */
1757 nlist->jindex[nlist->nri] = nlist->nrj;
1764 /* All interactions are perturbed, we can skip this entry */
1765 nbl_ci->cj_ind_end = cj_ind_start;
1766 nbl->ncjInUse -= cj_ind_end - cj_ind_start;
1770 /* Return the index of atom a within a cluster */
1771 static inline int cj_mod_cj4(int cj)
1773 return cj & (c_nbnxnGpuJgroupSize - 1);
1776 /* Convert a j-cluster to a cj4 group */
1777 static inline int cj_to_cj4(int cj)
1779 return cj/c_nbnxnGpuJgroupSize;
1782 /* Return the index of an j-atom within a warp */
1783 static inline int a_mod_wj(int a)
1785 return a & (c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit - 1);
1788 /* As make_fep_list above, but for super/sub lists. */
1789 static void make_fep_list(gmx::ArrayRef<const int> atomIndices,
1790 const nbnxn_atomdata_t *nbat,
1791 NbnxnPairlistGpu *nbl,
1792 gmx_bool bDiagRemoved,
1793 const nbnxn_sci_t *nbl_sci,
1804 int ind_i, ind_j, ai, aj;
1808 const nbnxn_cj4_t *cj4;
1810 const int numJClusterGroups = nbl_sci->numJClusterGroups();
1811 if (numJClusterGroups == 0)
1817 const int sci = nbl_sci->sci;
1819 const int cj4_ind_start = nbl_sci->cj4_ind_start;
1820 const int cj4_ind_end = nbl_sci->cj4_ind_end;
1822 /* Here we process one super-cell, max #atoms na_sc, versus a list
1823 * cj4 entries, each with max c_nbnxnGpuJgroupSize cj's, each
1824 * of size na_cj atoms.
1825 * On the GPU we don't support energy groups (yet).
1826 * So for each of the na_sc i-atoms, we need max one FEP list
1827 * for each max_nrj_fep j-atoms.
1829 nri_max = nbl->na_sc*nbl->na_cj*(1 + (numJClusterGroups*c_nbnxnGpuJgroupSize)/max_nrj_fep);
1830 if (nlist->nri + nri_max > nlist->maxnri)
1832 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1833 reallocate_nblist(nlist);
1836 /* Loop over the atoms in the i super-cluster */
1837 for (int c = 0; c < c_gpuNumClusterPerCell; c++)
1839 c_abs = sci*c_gpuNumClusterPerCell + c;
1841 for (int i = 0; i < nbl->na_ci; i++)
1843 ind_i = c_abs*nbl->na_ci + i;
1844 ai = atomIndices[ind_i];
1848 nlist->jindex[nri+1] = nlist->jindex[nri];
1849 nlist->iinr[nri] = ai;
1850 /* With GPUs, energy groups are not supported */
1851 nlist->gid[nri] = 0;
1852 nlist->shift[nri] = nbl_sci->shift & NBNXN_CI_SHIFT;
1854 bFEP_i = iGrid.atomIsPerturbed(c_abs - iGrid.cellOffset()*c_gpuNumClusterPerCell, i);
1856 xi = nbat->x()[ind_i*nbat->xstride+XX] + shx;
1857 yi = nbat->x()[ind_i*nbat->xstride+YY] + shy;
1858 zi = nbat->x()[ind_i*nbat->xstride+ZZ] + shz;
1860 const int nrjMax = nlist->nrj + numJClusterGroups*c_nbnxnGpuJgroupSize*nbl->na_cj;
1861 if (nrjMax > nlist->maxnrj)
1863 nlist->maxnrj = over_alloc_small(nrjMax);
1864 srenew(nlist->jjnr, nlist->maxnrj);
1865 srenew(nlist->excl_fep, nlist->maxnrj);
1868 for (int cj4_ind = cj4_ind_start; cj4_ind < cj4_ind_end; cj4_ind++)
1870 cj4 = &nbl->cj4[cj4_ind];
1872 for (int gcj = 0; gcj < c_nbnxnGpuJgroupSize; gcj++)
1874 if ((cj4->imei[0].imask & (1U << (gcj*c_gpuNumClusterPerCell + c))) == 0)
1876 /* Skip this ci for this cj */
1881 cj4->cj[gcj] - jGrid.cellOffset()*c_gpuNumClusterPerCell;
1883 if (bFEP_i || jGrid.clusterIsPerturbed(cjr))
1885 for (int j = 0; j < nbl->na_cj; j++)
1887 /* Is this interaction perturbed and not excluded? */
1888 ind_j = (jGrid.cellOffset()*c_gpuNumClusterPerCell + cjr)*nbl->na_cj + j;
1889 aj = atomIndices[ind_j];
1891 (bFEP_i || jGrid.atomIsPerturbed(cjr, j)) &&
1892 (!bDiagRemoved || ind_j >= ind_i))
1895 unsigned int excl_bit;
1898 const int jHalf = j/(c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit);
1899 nbnxn_excl_t *excl =
1900 get_exclusion_mask(nbl, cj4_ind, jHalf);
1902 excl_pair = a_mod_wj(j)*nbl->na_ci + i;
1903 excl_bit = (1U << (gcj*c_gpuNumClusterPerCell + c));
1905 dx = nbat->x()[ind_j*nbat->xstride+XX] - xi;
1906 dy = nbat->x()[ind_j*nbat->xstride+YY] - yi;
1907 dz = nbat->x()[ind_j*nbat->xstride+ZZ] - zi;
1909 /* The unpruned GPU list has more than 2/3
1910 * of the atom pairs beyond rlist. Using
1911 * this list will cause a lot of overhead
1912 * in the CPU FEP kernels, especially
1913 * relative to the fast GPU kernels.
1914 * So we prune the FEP list here.
1916 if (dx*dx + dy*dy + dz*dz < rlist_fep2)
1918 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1920 fep_list_new_nri_copy(nlist);
1924 /* Add it to the FEP list */
1925 nlist->jjnr[nlist->nrj] = aj;
1926 nlist->excl_fep[nlist->nrj] = (excl->pair[excl_pair] & excl_bit) ? 1 : 0;
1930 /* Exclude it from the normal list.
1931 * Note that the charge and LJ parameters have
1932 * been set to zero, but we need to avoid 0/0,
1933 * as perturbed atoms can be on top of each other.
1935 excl->pair[excl_pair] &= ~excl_bit;
1939 /* Note that we could mask out this pair in imask
1940 * if all i- and/or all j-particles are perturbed.
1941 * But since the perturbed pairs on the CPU will
1942 * take an order of magnitude more time, the GPU
1943 * will finish before the CPU and there is no gain.
1949 if (nlist->nrj > nlist->jindex[nri])
1951 /* Actually add this new, non-empty, list */
1953 nlist->jindex[nlist->nri] = nlist->nrj;
1960 /* Set all atom-pair exclusions for a GPU type list i-entry
1962 * Sets all atom-pair exclusions from the topology stored in exclusions
1963 * as masks in the pair-list for i-super-cluster list entry iEntry.
1966 setExclusionsForIEntry(const Nbnxm::GridSet &gridSet,
1967 NbnxnPairlistGpu *nbl,
1968 gmx_bool diagRemoved,
1969 int gmx_unused na_cj_2log,
1970 const nbnxn_sci_t &iEntry,
1971 const t_blocka &exclusions)
1973 if (iEntry.numJClusterGroups() == 0)
1979 /* Set the search ranges using start and end j-cluster indices.
1980 * Note that here we can not use cj4_ind_end, since the last cj4
1981 * can be only partially filled, so we use cj_ind.
1983 const JListRanges ranges(iEntry.cj4_ind_start*c_nbnxnGpuJgroupSize,
1985 gmx::makeConstArrayRef(nbl->cj4));
1987 GMX_ASSERT(nbl->na_ci == c_nbnxnGpuClusterSize, "na_ci should match the GPU cluster size");
1988 constexpr int c_clusterSize = c_nbnxnGpuClusterSize;
1989 constexpr int c_superClusterSize = c_nbnxnGpuNumClusterPerSupercluster*c_nbnxnGpuClusterSize;
1991 const int iSuperCluster = iEntry.sci;
1993 gmx::ArrayRef<const int> atomIndices = gridSet.atomIndices();
1994 gmx::ArrayRef<const int> cell = gridSet.cells();
1996 /* Loop over the atoms in the i super-cluster */
1997 for (int i = 0; i < c_superClusterSize; i++)
1999 const int iIndex = iSuperCluster*c_superClusterSize + i;
2000 const int iAtom = atomIndices[iIndex];
2003 const int iCluster = i/c_clusterSize;
2005 /* Loop over the topology-based exclusions for this i-atom */
2006 for (int exclIndex = exclusions.index[iAtom]; exclIndex < exclusions.index[iAtom + 1]; exclIndex++)
2008 const int jAtom = exclusions.a[exclIndex];
2012 /* The self exclusions are already set, save some time */
2016 /* Get the index of the j-atom in the nbnxn atom data */
2017 const int jIndex = cell[jAtom];
2019 /* Without shifts we only calculate interactions j>i
2020 * for one-way pair-lists.
2022 /* NOTE: We would like to use iIndex on the right hand side,
2023 * but that makes this routine 25% slower with gcc6/7.
2024 * Even using c_superClusterSize makes it slower.
2025 * Either of these changes triggers peeling of the exclIndex
2026 * loop, which apparently leads to far less efficient code.
2028 if (diagRemoved && jIndex <= iSuperCluster*nbl->na_sc + i)
2033 const int jCluster = jIndex/c_clusterSize;
2035 /* Check whether the cluster is in our list? */
2036 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
2039 findJClusterInJList(jCluster, ranges,
2040 gmx::makeConstArrayRef(nbl->cj4));
2044 /* We found an exclusion, clear the corresponding
2047 const unsigned int pairMask = (1U << (cj_mod_cj4(index)*c_gpuNumClusterPerCell + iCluster));
2048 /* Check if the i-cluster interacts with the j-cluster */
2049 if (nbl_imask0(nbl, index) & pairMask)
2051 const int innerI = (i & (c_clusterSize - 1));
2052 const int innerJ = (jIndex & (c_clusterSize - 1));
2054 /* Determine which j-half (CUDA warp) we are in */
2055 const int jHalf = innerJ/(c_clusterSize/c_nbnxnGpuClusterpairSplit);
2057 nbnxn_excl_t *interactionMask =
2058 get_exclusion_mask(nbl, cj_to_cj4(index), jHalf);
2060 interactionMask->pair[a_mod_wj(innerJ)*c_clusterSize + innerI] &= ~pairMask;
2069 /* Make a new ci entry at the back of nbl->ci */
2070 static void addNewIEntry(NbnxnPairlistCpu *nbl, int ci, int shift, int flags)
2074 ciEntry.shift = shift;
2075 /* Store the interaction flags along with the shift */
2076 ciEntry.shift |= flags;
2077 ciEntry.cj_ind_start = nbl->cj.size();
2078 ciEntry.cj_ind_end = nbl->cj.size();
2079 nbl->ci.push_back(ciEntry);
2082 /* Make a new sci entry at index nbl->nsci */
2083 static void addNewIEntry(NbnxnPairlistGpu *nbl, int sci, int shift, int gmx_unused flags)
2085 nbnxn_sci_t sciEntry;
2087 sciEntry.shift = shift;
2088 sciEntry.cj4_ind_start = nbl->cj4.size();
2089 sciEntry.cj4_ind_end = nbl->cj4.size();
2091 nbl->sci.push_back(sciEntry);
2094 /* Sort the simple j-list cj on exclusions.
2095 * Entries with exclusions will all be sorted to the beginning of the list.
2097 static void sort_cj_excl(nbnxn_cj_t *cj, int ncj,
2098 NbnxnPairlistCpuWork *work)
2100 work->cj.resize(ncj);
2102 /* Make a list of the j-cells involving exclusions */
2104 for (int j = 0; j < ncj; j++)
2106 if (cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
2108 work->cj[jnew++] = cj[j];
2111 /* Check if there are exclusions at all or not just the first entry */
2112 if (!((jnew == 0) ||
2113 (jnew == 1 && cj[0].excl != NBNXN_INTERACTION_MASK_ALL)))
2115 for (int j = 0; j < ncj; j++)
2117 if (cj[j].excl == NBNXN_INTERACTION_MASK_ALL)
2119 work->cj[jnew++] = cj[j];
2122 for (int j = 0; j < ncj; j++)
2124 cj[j] = work->cj[j];
2129 /* Close this simple list i entry */
2130 static void closeIEntry(NbnxnPairlistCpu *nbl,
2131 int gmx_unused sp_max_av,
2132 gmx_bool gmx_unused progBal,
2133 float gmx_unused nsp_tot_est,
2134 int gmx_unused thread,
2135 int gmx_unused nthread)
2137 nbnxn_ci_t &ciEntry = nbl->ci.back();
2139 /* All content of the new ci entry have already been filled correctly,
2140 * we only need to sort and increase counts or remove the entry when empty.
2142 const int jlen = ciEntry.cj_ind_end - ciEntry.cj_ind_start;
2145 sort_cj_excl(nbl->cj.data() + ciEntry.cj_ind_start, jlen, nbl->work);
2147 /* The counts below are used for non-bonded pair/flop counts
2148 * and should therefore match the available kernel setups.
2150 if (!(ciEntry.shift & NBNXN_CI_DO_COUL(0)))
2152 nbl->work->ncj_noq += jlen;
2154 else if ((ciEntry.shift & NBNXN_CI_HALF_LJ(0)) ||
2155 !(ciEntry.shift & NBNXN_CI_DO_LJ(0)))
2157 nbl->work->ncj_hlj += jlen;
2162 /* Entry is empty: remove it */
2167 /* Split sci entry for load balancing on the GPU.
2168 * Splitting ensures we have enough lists to fully utilize the whole GPU.
2169 * With progBal we generate progressively smaller lists, which improves
2170 * load balancing. As we only know the current count on our own thread,
2171 * we will need to estimate the current total amount of i-entries.
2172 * As the lists get concatenated later, this estimate depends
2173 * both on nthread and our own thread index.
2175 static void split_sci_entry(NbnxnPairlistGpu *nbl,
2177 gmx_bool progBal, float nsp_tot_est,
2178 int thread, int nthread)
2186 /* Estimate the total numbers of ci's of the nblist combined
2187 * over all threads using the target number of ci's.
2189 nsp_est = (nsp_tot_est*thread)/nthread + nbl->nci_tot;
2191 /* The first ci blocks should be larger, to avoid overhead.
2192 * The last ci blocks should be smaller, to improve load balancing.
2193 * The factor 3/2 makes the first block 3/2 times the target average
2194 * and ensures that the total number of blocks end up equal to
2195 * that of equally sized blocks of size nsp_target_av.
2197 nsp_max = static_cast<int>(nsp_target_av*(nsp_tot_est*1.5/(nsp_est + nsp_tot_est)));
2201 nsp_max = nsp_target_av;
2204 const int cj4_start = nbl->sci.back().cj4_ind_start;
2205 const int cj4_end = nbl->sci.back().cj4_ind_end;
2206 const int j4len = cj4_end - cj4_start;
2208 if (j4len > 1 && j4len*c_gpuNumClusterPerCell*c_nbnxnGpuJgroupSize > nsp_max)
2210 /* Modify the last ci entry and process the cj4's again */
2216 for (int cj4 = cj4_start; cj4 < cj4_end; cj4++)
2218 int nsp_cj4_p = nsp_cj4;
2219 /* Count the number of cluster pairs in this cj4 group */
2221 for (int p = 0; p < c_gpuNumClusterPerCell*c_nbnxnGpuJgroupSize; p++)
2223 nsp_cj4 += (nbl->cj4[cj4].imei[0].imask >> p) & 1;
2226 /* If adding the current cj4 with nsp_cj4 pairs get us further
2227 * away from our target nsp_max, split the list before this cj4.
2229 if (nsp > 0 && nsp_max - nsp < nsp + nsp_cj4 - nsp_max)
2231 /* Split the list at cj4 */
2232 nbl->sci.back().cj4_ind_end = cj4;
2233 /* Create a new sci entry */
2235 sciNew.sci = nbl->sci.back().sci;
2236 sciNew.shift = nbl->sci.back().shift;
2237 sciNew.cj4_ind_start = cj4;
2238 nbl->sci.push_back(sciNew);
2241 nsp_cj4_e = nsp_cj4_p;
2247 /* Put the remaining cj4's in the last sci entry */
2248 nbl->sci.back().cj4_ind_end = cj4_end;
2250 /* Possibly balance out the last two sci's
2251 * by moving the last cj4 of the second last sci.
2253 if (nsp_sci - nsp_cj4_e >= nsp + nsp_cj4_e)
2255 GMX_ASSERT(nbl->sci.size() >= 2, "We expect at least two elements");
2256 nbl->sci[nbl->sci.size() - 2].cj4_ind_end--;
2257 nbl->sci[nbl->sci.size() - 1].cj4_ind_start--;
2262 /* Clost this super/sub list i entry */
2263 static void closeIEntry(NbnxnPairlistGpu *nbl,
2265 gmx_bool progBal, float nsp_tot_est,
2266 int thread, int nthread)
2268 nbnxn_sci_t &sciEntry = *getOpenIEntry(nbl);
2270 /* All content of the new ci entry have already been filled correctly,
2271 * we only need to, potentially, split or remove the entry when empty.
2273 int j4len = sciEntry.numJClusterGroups();
2276 /* We can only have complete blocks of 4 j-entries in a list,
2277 * so round the count up before closing.
2279 int ncj4 = (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize;
2280 nbl->work->cj_ind = ncj4*c_nbnxnGpuJgroupSize;
2284 /* Measure the size of the new entry and potentially split it */
2285 split_sci_entry(nbl, nsp_max_av, progBal, nsp_tot_est,
2291 /* Entry is empty: remove it */
2292 nbl->sci.pop_back();
2296 /* Syncs the working array before adding another grid pair to the GPU list */
2297 static void sync_work(NbnxnPairlistCpu gmx_unused *nbl)
2301 /* Syncs the working array before adding another grid pair to the GPU list */
2302 static void sync_work(NbnxnPairlistGpu *nbl)
2304 nbl->work->cj_ind = nbl->cj4.size()*c_nbnxnGpuJgroupSize;
2307 /* Clears an NbnxnPairlistCpu data structure */
2308 static void clear_pairlist(NbnxnPairlistCpu *nbl)
2314 nbl->ciOuter.clear();
2315 nbl->cjOuter.clear();
2317 nbl->work->ncj_noq = 0;
2318 nbl->work->ncj_hlj = 0;
2321 /* Clears an NbnxnPairlistGpu data structure */
2322 static void clear_pairlist(NbnxnPairlistGpu *nbl)
2326 nbl->excl.resize(1);
2330 /* Clears a group scheme pair list */
2331 static void clear_pairlist_fep(t_nblist *nl)
2335 if (nl->jindex == nullptr)
2337 snew(nl->jindex, 1);
2342 /* Sets a simple list i-cell bounding box, including PBC shift */
2343 static inline void set_icell_bb_simple(gmx::ArrayRef<const BoundingBox> bb,
2345 real shx, real shy, real shz,
2348 bb_ci->lower.x = bb[ci].lower.x + shx;
2349 bb_ci->lower.y = bb[ci].lower.y + shy;
2350 bb_ci->lower.z = bb[ci].lower.z + shz;
2351 bb_ci->upper.x = bb[ci].upper.x + shx;
2352 bb_ci->upper.y = bb[ci].upper.y + shy;
2353 bb_ci->upper.z = bb[ci].upper.z + shz;
2356 /* Sets a simple list i-cell bounding box, including PBC shift */
2357 static inline void set_icell_bb(const Grid &iGrid,
2359 real shx, real shy, real shz,
2360 NbnxnPairlistCpuWork *work)
2362 set_icell_bb_simple(iGrid.iBoundingBoxes(), ci, shx, shy, shz,
2363 &work->iClusterData.bb[0]);
2367 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2368 static void set_icell_bbxxxx_supersub(gmx::ArrayRef<const float> bb,
2370 real shx, real shy, real shz,
2373 constexpr int cellBBStride = packedBoundingBoxesIndex(c_gpuNumClusterPerCell);
2374 constexpr int pbbStride = c_packedBoundingBoxesDimSize;
2375 const int ia = ci*cellBBStride;
2376 for (int m = 0; m < cellBBStride; m += c_packedBoundingBoxesSize)
2378 for (int i = 0; i < pbbStride; i++)
2380 bb_ci[m + 0*pbbStride + i] = bb[ia + m + 0*pbbStride + i] + shx;
2381 bb_ci[m + 1*pbbStride + i] = bb[ia + m + 1*pbbStride + i] + shy;
2382 bb_ci[m + 2*pbbStride + i] = bb[ia + m + 2*pbbStride + i] + shz;
2383 bb_ci[m + 3*pbbStride + i] = bb[ia + m + 3*pbbStride + i] + shx;
2384 bb_ci[m + 4*pbbStride + i] = bb[ia + m + 4*pbbStride + i] + shy;
2385 bb_ci[m + 5*pbbStride + i] = bb[ia + m + 5*pbbStride + i] + shz;
2391 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2392 gmx_unused static void set_icell_bb_supersub(gmx::ArrayRef<const BoundingBox> bb,
2394 real shx, real shy, real shz,
2397 for (int i = 0; i < c_gpuNumClusterPerCell; i++)
2399 set_icell_bb_simple(bb, ci*c_gpuNumClusterPerCell+i,
2405 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2406 gmx_unused static void set_icell_bb(const Grid &iGrid,
2408 real shx, real shy, real shz,
2409 NbnxnPairlistGpuWork *work)
2412 set_icell_bbxxxx_supersub(iGrid.packedBoundingBoxes(), ci, shx, shy, shz,
2413 work->iSuperClusterData.bbPacked.data());
2415 set_icell_bb_supersub(iGrid.iBoundingBoxes(), ci, shx, shy, shz,
2416 work->iSuperClusterData.bb.data());
2420 /* Copies PBC shifted i-cell atom coordinates x,y,z to working array */
2421 static void icell_set_x_simple(int ci,
2422 real shx, real shy, real shz,
2423 int stride, const real *x,
2424 NbnxnPairlistCpuWork::IClusterData *iClusterData)
2426 const int ia = ci*c_nbnxnCpuIClusterSize;
2428 for (int i = 0; i < c_nbnxnCpuIClusterSize; i++)
2430 iClusterData->x[i*STRIDE_XYZ+XX] = x[(ia+i)*stride+XX] + shx;
2431 iClusterData->x[i*STRIDE_XYZ+YY] = x[(ia+i)*stride+YY] + shy;
2432 iClusterData->x[i*STRIDE_XYZ+ZZ] = x[(ia+i)*stride+ZZ] + shz;
2436 static void icell_set_x(int ci,
2437 real shx, real shy, real shz,
2438 int stride, const real *x,
2439 const Nbnxm::KernelType kernelType,
2440 NbnxnPairlistCpuWork *work)
2445 #ifdef GMX_NBNXN_SIMD_4XN
2446 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
2447 icell_set_x_simd_4xn(ci, shx, shy, shz, stride, x, work);
2450 #ifdef GMX_NBNXN_SIMD_2XNN
2451 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
2452 icell_set_x_simd_2xnn(ci, shx, shy, shz, stride, x, work);
2456 case Nbnxm::KernelType::Cpu4x4_PlainC:
2457 icell_set_x_simple(ci, shx, shy, shz, stride, x, &work->iClusterData);
2460 GMX_ASSERT(false, "Unhandled case");
2465 /* Copies PBC shifted super-cell atom coordinates x,y,z to working array */
2466 static void icell_set_x(int ci,
2467 real shx, real shy, real shz,
2468 int stride, const real *x,
2469 Nbnxm::KernelType gmx_unused kernelType,
2470 NbnxnPairlistGpuWork *work)
2472 #if !GMX_SIMD4_HAVE_REAL
2474 real * x_ci = work->iSuperClusterData.x.data();
2476 int ia = ci*c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize;
2477 for (int i = 0; i < c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize; i++)
2479 x_ci[i*DIM + XX] = x[(ia+i)*stride + XX] + shx;
2480 x_ci[i*DIM + YY] = x[(ia+i)*stride + YY] + shy;
2481 x_ci[i*DIM + ZZ] = x[(ia+i)*stride + ZZ] + shz;
2484 #else /* !GMX_SIMD4_HAVE_REAL */
2486 real * x_ci = work->iSuperClusterData.xSimd.data();
2488 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2490 for (int i = 0; i < c_nbnxnGpuClusterSize; i += GMX_SIMD4_WIDTH)
2492 int io = si*c_nbnxnGpuClusterSize + i;
2493 int ia = ci*c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize + io;
2494 for (int j = 0; j < GMX_SIMD4_WIDTH; j++)
2496 x_ci[io*DIM + j + XX*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + XX] + shx;
2497 x_ci[io*DIM + j + YY*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + YY] + shy;
2498 x_ci[io*DIM + j + ZZ*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + ZZ] + shz;
2503 #endif /* !GMX_SIMD4_HAVE_REAL */
2506 static real minimum_subgrid_size_xy(const Grid &grid)
2508 const Grid::Dimensions &dims = grid.dimensions();
2510 if (grid.geometry().isSimple)
2512 return std::min(dims.cellSize[XX], dims.cellSize[YY]);
2516 return std::min(dims.cellSize[XX]/c_gpuNumClusterPerCellX,
2517 dims.cellSize[YY]/c_gpuNumClusterPerCellY);
2521 static real effective_buffer_1x1_vs_MxN(const Grid &iGrid,
2524 const real eff_1x1_buffer_fac_overest = 0.1;
2526 /* Determine an atom-pair list cut-off buffer size for atom pairs,
2527 * to be added to rlist (including buffer) used for MxN.
2528 * This is for converting an MxN list to a 1x1 list. This means we can't
2529 * use the normal buffer estimate, as we have an MxN list in which
2530 * some atom pairs beyond rlist are missing. We want to capture
2531 * the beneficial effect of buffering by extra pairs just outside rlist,
2532 * while removing the useless pairs that are further away from rlist.
2533 * (Also the buffer could have been set manually not using the estimate.)
2534 * This buffer size is an overestimate.
2535 * We add 10% of the smallest grid sub-cell dimensions.
2536 * Note that the z-size differs per cell and we don't use this,
2537 * so we overestimate.
2538 * With PME, the 10% value gives a buffer that is somewhat larger
2539 * than the effective buffer with a tolerance of 0.005 kJ/mol/ps.
2540 * Smaller tolerances or using RF lead to a smaller effective buffer,
2541 * so 10% gives a safe overestimate.
2543 return eff_1x1_buffer_fac_overest*(minimum_subgrid_size_xy(iGrid) +
2544 minimum_subgrid_size_xy(jGrid));
2547 /* Estimates the interaction volume^2 for non-local interactions */
2548 static real nonlocal_vol2(const struct gmx_domdec_zones_t *zones, const rvec ls, real r)
2556 /* Here we simply add up the volumes of 1, 2 or 3 1D decomposition
2557 * not home interaction volume^2. As these volumes are not additive,
2558 * this is an overestimate, but it would only be significant in the limit
2559 * of small cells, where we anyhow need to split the lists into
2560 * as small parts as possible.
2563 for (int z = 0; z < zones->n; z++)
2565 if (zones->shift[z][XX] + zones->shift[z][YY] + zones->shift[z][ZZ] == 1)
2570 for (int d = 0; d < DIM; d++)
2572 if (zones->shift[z][d] == 0)
2576 za *= zones->size[z].x1[d] - zones->size[z].x0[d];
2580 /* 4 octants of a sphere */
2581 vold_est = 0.25*M_PI*r*r*r*r;
2582 /* 4 quarter pie slices on the edges */
2583 vold_est += 4*cl*M_PI/6.0*r*r*r;
2584 /* One rectangular volume on a face */
2585 vold_est += ca*0.5*r*r;
2587 vol2_est_tot += vold_est*za;
2591 return vol2_est_tot;
2594 /* Estimates the average size of a full j-list for super/sub setup */
2595 static void get_nsubpair_target(const nbnxn_search *nbs,
2596 const InteractionLocality iloc,
2598 const int min_ci_balanced,
2599 int *nsubpair_target,
2600 float *nsubpair_tot_est)
2602 /* The target value of 36 seems to be the optimum for Kepler.
2603 * Maxwell is less sensitive to the exact value.
2605 const int nsubpair_target_min = 36;
2606 real r_eff_sup, vol_est, nsp_est, nsp_est_nl;
2608 const Grid &grid = nbs->gridSet().grids()[0];
2610 /* We don't need to balance list sizes if:
2611 * - We didn't request balancing.
2612 * - The number of grid cells >= the number of lists requested,
2613 * since we will always generate at least #cells lists.
2614 * - We don't have any cells, since then there won't be any lists.
2616 if (min_ci_balanced <= 0 || grid.numCells() >= min_ci_balanced || grid.numCells() == 0)
2618 /* nsubpair_target==0 signals no balancing */
2619 *nsubpair_target = 0;
2620 *nsubpair_tot_est = 0;
2626 const int numAtomsCluster = grid.geometry().numAtomsICluster;
2627 const Grid::Dimensions &dims = grid.dimensions();
2629 ls[XX] = dims.cellSize[XX]/c_gpuNumClusterPerCellX;
2630 ls[YY] = dims.cellSize[YY]/c_gpuNumClusterPerCellY;
2631 ls[ZZ] = numAtomsCluster/(dims.atomDensity*ls[XX]*ls[YY]);
2633 /* The formulas below are a heuristic estimate of the average nsj per si*/
2634 r_eff_sup = rlist + nbnxn_get_rlist_effective_inc(numAtomsCluster, ls);
2636 if (!nbs->DomDec || nbs->zones->n == 1)
2643 gmx::square(dims.atomDensity/numAtomsCluster)*
2644 nonlocal_vol2(nbs->zones, ls, r_eff_sup);
2647 if (iloc == InteractionLocality::Local)
2649 /* Sub-cell interacts with itself */
2650 vol_est = ls[XX]*ls[YY]*ls[ZZ];
2651 /* 6/2 rectangular volume on the faces */
2652 vol_est += (ls[XX]*ls[YY] + ls[XX]*ls[ZZ] + ls[YY]*ls[ZZ])*r_eff_sup;
2653 /* 12/2 quarter pie slices on the edges */
2654 vol_est += 2*(ls[XX] + ls[YY] + ls[ZZ])*0.25*M_PI*gmx::square(r_eff_sup);
2655 /* 4 octants of a sphere */
2656 vol_est += 0.5*4.0/3.0*M_PI*gmx::power3(r_eff_sup);
2658 /* Estimate the number of cluster pairs as the local number of
2659 * clusters times the volume they interact with times the density.
2661 nsp_est = grid.numClusters()*vol_est*dims.atomDensity/numAtomsCluster;
2663 /* Subtract the non-local pair count */
2664 nsp_est -= nsp_est_nl;
2666 /* For small cut-offs nsp_est will be an underesimate.
2667 * With DD nsp_est_nl is an overestimate so nsp_est can get negative.
2668 * So to avoid too small or negative nsp_est we set a minimum of
2669 * all cells interacting with all 3^3 direct neighbors (3^3-1)/2+1=14.
2670 * This might be a slight overestimate for small non-periodic groups of
2671 * atoms as will occur for a local domain with DD, but for small
2672 * groups of atoms we'll anyhow be limited by nsubpair_target_min,
2673 * so this overestimation will not matter.
2675 nsp_est = std::max(nsp_est, grid.numClusters()*14._real);
2679 fprintf(debug, "nsp_est local %5.1f non-local %5.1f\n",
2680 nsp_est, nsp_est_nl);
2685 nsp_est = nsp_est_nl;
2688 /* Thus the (average) maximum j-list size should be as follows.
2689 * Since there is overhead, we shouldn't make the lists too small
2690 * (and we can't chop up j-groups) so we use a minimum target size of 36.
2692 *nsubpair_target = std::max(nsubpair_target_min,
2693 roundToInt(nsp_est/min_ci_balanced));
2694 *nsubpair_tot_est = static_cast<int>(nsp_est);
2698 fprintf(debug, "nbl nsp estimate %.1f, nsubpair_target %d\n",
2699 nsp_est, *nsubpair_target);
2703 /* Debug list print function */
2704 static void print_nblist_ci_cj(FILE *fp, const NbnxnPairlistCpu *nbl)
2706 for (const nbnxn_ci_t &ciEntry : nbl->ci)
2708 fprintf(fp, "ci %4d shift %2d ncj %3d\n",
2709 ciEntry.ci, ciEntry.shift,
2710 ciEntry.cj_ind_end - ciEntry.cj_ind_start);
2712 for (int j = ciEntry.cj_ind_start; j < ciEntry.cj_ind_end; j++)
2714 fprintf(fp, " cj %5d imask %x\n",
2721 /* Debug list print function */
2722 static void print_nblist_sci_cj(FILE *fp, const NbnxnPairlistGpu *nbl)
2724 for (const nbnxn_sci_t &sci : nbl->sci)
2726 fprintf(fp, "ci %4d shift %2d ncj4 %2d\n",
2728 sci.numJClusterGroups());
2731 for (int j4 = sci.cj4_ind_start; j4 < sci.cj4_ind_end; j4++)
2733 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
2735 fprintf(fp, " sj %5d imask %x\n",
2737 nbl->cj4[j4].imei[0].imask);
2738 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2740 if (nbl->cj4[j4].imei[0].imask & (1U << (j*c_gpuNumClusterPerCell + si)))
2747 fprintf(fp, "ci %4d shift %2d ncj4 %2d ncp %3d\n",
2749 sci.numJClusterGroups(),
2754 /* Combine pair lists *nbl generated on multiple threads nblc */
2755 static void combine_nblists(int nnbl, NbnxnPairlistGpu **nbl,
2756 NbnxnPairlistGpu *nblc)
2758 int nsci = nblc->sci.size();
2759 int ncj4 = nblc->cj4.size();
2760 int nexcl = nblc->excl.size();
2761 for (int i = 0; i < nnbl; i++)
2763 nsci += nbl[i]->sci.size();
2764 ncj4 += nbl[i]->cj4.size();
2765 nexcl += nbl[i]->excl.size();
2768 /* Resize with the final, combined size, so we can fill in parallel */
2769 /* NOTE: For better performance we should use default initialization */
2770 nblc->sci.resize(nsci);
2771 nblc->cj4.resize(ncj4);
2772 nblc->excl.resize(nexcl);
2774 /* Each thread should copy its own data to the combined arrays,
2775 * as otherwise data will go back and forth between different caches.
2777 #if GMX_OPENMP && !(defined __clang_analyzer__)
2778 int nthreads = gmx_omp_nthreads_get(emntPairsearch);
2781 #pragma omp parallel for num_threads(nthreads) schedule(static)
2782 for (int n = 0; n < nnbl; n++)
2786 /* Determine the offset in the combined data for our thread.
2787 * Note that the original sizes in nblc are lost.
2789 int sci_offset = nsci;
2790 int cj4_offset = ncj4;
2791 int excl_offset = nexcl;
2793 for (int i = n; i < nnbl; i++)
2795 sci_offset -= nbl[i]->sci.size();
2796 cj4_offset -= nbl[i]->cj4.size();
2797 excl_offset -= nbl[i]->excl.size();
2800 const NbnxnPairlistGpu &nbli = *nbl[n];
2802 for (size_t i = 0; i < nbli.sci.size(); i++)
2804 nblc->sci[sci_offset + i] = nbli.sci[i];
2805 nblc->sci[sci_offset + i].cj4_ind_start += cj4_offset;
2806 nblc->sci[sci_offset + i].cj4_ind_end += cj4_offset;
2809 for (size_t j4 = 0; j4 < nbli.cj4.size(); j4++)
2811 nblc->cj4[cj4_offset + j4] = nbli.cj4[j4];
2812 nblc->cj4[cj4_offset + j4].imei[0].excl_ind += excl_offset;
2813 nblc->cj4[cj4_offset + j4].imei[1].excl_ind += excl_offset;
2816 for (size_t j4 = 0; j4 < nbli.excl.size(); j4++)
2818 nblc->excl[excl_offset + j4] = nbli.excl[j4];
2821 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2824 for (int n = 0; n < nnbl; n++)
2826 nblc->nci_tot += nbl[n]->nci_tot;
2830 static void balance_fep_lists(const nbnxn_search *nbs,
2831 nbnxn_pairlist_set_t *nbl_lists)
2834 int nri_tot, nrj_tot, nrj_target;
2838 nnbl = nbl_lists->nnbl;
2842 /* Nothing to balance */
2846 /* Count the total i-lists and pairs */
2849 for (int th = 0; th < nnbl; th++)
2851 nri_tot += nbl_lists->nbl_fep[th]->nri;
2852 nrj_tot += nbl_lists->nbl_fep[th]->nrj;
2855 nrj_target = (nrj_tot + nnbl - 1)/nnbl;
2857 assert(gmx_omp_nthreads_get(emntNonbonded) == nnbl);
2859 #pragma omp parallel for schedule(static) num_threads(nnbl)
2860 for (int th = 0; th < nnbl; th++)
2864 t_nblist *nbl = nbs->work[th].nbl_fep.get();
2866 /* Note that here we allocate for the total size, instead of
2867 * a per-thread esimate (which is hard to obtain).
2869 if (nri_tot > nbl->maxnri)
2871 nbl->maxnri = over_alloc_large(nri_tot);
2872 reallocate_nblist(nbl);
2874 if (nri_tot > nbl->maxnri || nrj_tot > nbl->maxnrj)
2876 nbl->maxnrj = over_alloc_small(nrj_tot);
2877 srenew(nbl->jjnr, nbl->maxnrj);
2878 srenew(nbl->excl_fep, nbl->maxnrj);
2881 clear_pairlist_fep(nbl);
2883 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2886 /* Loop over the source lists and assign and copy i-entries */
2888 nbld = nbs->work[th_dest].nbl_fep.get();
2889 for (int th = 0; th < nnbl; th++)
2893 nbls = nbl_lists->nbl_fep[th];
2895 for (int i = 0; i < nbls->nri; i++)
2899 /* The number of pairs in this i-entry */
2900 nrj = nbls->jindex[i+1] - nbls->jindex[i];
2902 /* Decide if list th_dest is too large and we should procede
2903 * to the next destination list.
2905 if (th_dest+1 < nnbl && nbld->nrj > 0 &&
2906 nbld->nrj + nrj - nrj_target > nrj_target - nbld->nrj)
2909 nbld = nbs->work[th_dest].nbl_fep.get();
2912 nbld->iinr[nbld->nri] = nbls->iinr[i];
2913 nbld->gid[nbld->nri] = nbls->gid[i];
2914 nbld->shift[nbld->nri] = nbls->shift[i];
2916 for (int j = nbls->jindex[i]; j < nbls->jindex[i+1]; j++)
2918 nbld->jjnr[nbld->nrj] = nbls->jjnr[j];
2919 nbld->excl_fep[nbld->nrj] = nbls->excl_fep[j];
2923 nbld->jindex[nbld->nri] = nbld->nrj;
2927 /* Swap the list pointers */
2928 for (int th = 0; th < nnbl; th++)
2930 t_nblist *nbl_tmp = nbs->work[th].nbl_fep.release();
2931 nbs->work[th].nbl_fep.reset(nbl_lists->nbl_fep[th]);
2932 nbl_lists->nbl_fep[th] = nbl_tmp;
2936 fprintf(debug, "nbl_fep[%d] nri %4d nrj %4d\n",
2938 nbl_lists->nbl_fep[th]->nri,
2939 nbl_lists->nbl_fep[th]->nrj);
2944 /* Returns the next ci to be processes by our thread */
2945 static gmx_bool next_ci(const Grid &grid,
2946 int nth, int ci_block,
2947 int *ci_x, int *ci_y,
2953 if (*ci_b == ci_block)
2955 /* Jump to the next block assigned to this task */
2956 *ci += (nth - 1)*ci_block;
2960 if (*ci >= grid.numCells())
2965 while (*ci >= grid.firstCellInColumn(*ci_x*grid.dimensions().numCells[YY] + *ci_y + 1))
2968 if (*ci_y == grid.dimensions().numCells[YY])
2978 /* Returns the distance^2 for which we put cell pairs in the list
2979 * without checking atom pair distances. This is usually < rlist^2.
2981 static float boundingbox_only_distance2(const Grid::Dimensions &iGridDims,
2982 const Grid::Dimensions &jGridDims,
2986 /* If the distance between two sub-cell bounding boxes is less
2987 * than this distance, do not check the distance between
2988 * all particle pairs in the sub-cell, since then it is likely
2989 * that the box pair has atom pairs within the cut-off.
2990 * We use the nblist cut-off minus 0.5 times the average x/y diagonal
2991 * spacing of the sub-cells. Around 40% of the checked pairs are pruned.
2992 * Using more than 0.5 gains at most 0.5%.
2993 * If forces are calculated more than twice, the performance gain
2994 * in the force calculation outweighs the cost of checking.
2995 * Note that with subcell lists, the atom-pair distance check
2996 * is only performed when only 1 out of 8 sub-cells in within range,
2997 * this is because the GPU is much faster than the cpu.
3002 bbx = 0.5*(iGridDims.cellSize[XX] + jGridDims.cellSize[XX]);
3003 bby = 0.5*(iGridDims.cellSize[YY] + jGridDims.cellSize[YY]);
3006 bbx /= c_gpuNumClusterPerCellX;
3007 bby /= c_gpuNumClusterPerCellY;
3010 rbb2 = std::max(0.0, rlist - 0.5*std::sqrt(bbx*bbx + bby*bby));
3016 return (float)((1+GMX_FLOAT_EPS)*rbb2);
3020 static int get_ci_block_size(const Grid &iGrid,
3021 gmx_bool bDomDec, int nth)
3023 const int ci_block_enum = 5;
3024 const int ci_block_denom = 11;
3025 const int ci_block_min_atoms = 16;
3028 /* Here we decide how to distribute the blocks over the threads.
3029 * We use prime numbers to try to avoid that the grid size becomes
3030 * a multiple of the number of threads, which would lead to some
3031 * threads getting "inner" pairs and others getting boundary pairs,
3032 * which in turns will lead to load imbalance between threads.
3033 * Set the block size as 5/11/ntask times the average number of cells
3034 * in a y,z slab. This should ensure a quite uniform distribution
3035 * of the grid parts of the different thread along all three grid
3036 * zone boundaries with 3D domain decomposition. At the same time
3037 * the blocks will not become too small.
3039 ci_block = (iGrid.numCells()*ci_block_enum)/(ci_block_denom*iGrid.dimensions().numCells[XX]*nth);
3041 const int numAtomsPerCell = iGrid.geometry().numAtomsPerCell;
3043 /* Ensure the blocks are not too small: avoids cache invalidation */
3044 if (ci_block*numAtomsPerCell < ci_block_min_atoms)
3046 ci_block = (ci_block_min_atoms + numAtomsPerCell - 1)/numAtomsPerCell;
3049 /* Without domain decomposition
3050 * or with less than 3 blocks per task, divide in nth blocks.
3052 if (!bDomDec || nth*3*ci_block > iGrid.numCells())
3054 ci_block = (iGrid.numCells() + nth - 1)/nth;
3057 if (ci_block > 1 && (nth - 1)*ci_block >= iGrid.numCells())
3059 /* Some threads have no work. Although reducing the block size
3060 * does not decrease the block count on the first few threads,
3061 * with GPUs better mixing of "upper" cells that have more empty
3062 * clusters results in a somewhat lower max load over all threads.
3063 * Without GPUs the regime of so few atoms per thread is less
3064 * performance relevant, but with 8-wide SIMD the same reasoning
3065 * applies, since the pair list uses 4 i-atom "sub-clusters".
3073 /* Returns the number of bits to right-shift a cluster index to obtain
3074 * the corresponding force buffer flag index.
3076 static int getBufferFlagShift(int numAtomsPerCluster)
3078 int bufferFlagShift = 0;
3079 while ((numAtomsPerCluster << bufferFlagShift) < NBNXN_BUFFERFLAG_SIZE)
3084 return bufferFlagShift;
3087 static bool pairlistIsSimple(const NbnxnPairlistCpu gmx_unused &pairlist)
3092 static bool pairlistIsSimple(const NbnxnPairlistGpu gmx_unused &pairlist)
3097 static void makeClusterListWrapper(NbnxnPairlistCpu *nbl,
3098 const Grid gmx_unused &iGrid,
3101 const int firstCell,
3103 const bool excludeSubDiagonal,
3104 const nbnxn_atomdata_t *nbat,
3107 const Nbnxm::KernelType kernelType,
3108 int *numDistanceChecks)
3112 case Nbnxm::KernelType::Cpu4x4_PlainC:
3113 makeClusterListSimple(jGrid,
3114 nbl, ci, firstCell, lastCell,
3120 #ifdef GMX_NBNXN_SIMD_4XN
3121 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
3122 makeClusterListSimd4xn(jGrid,
3123 nbl, ci, firstCell, lastCell,
3130 #ifdef GMX_NBNXN_SIMD_2XNN
3131 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
3132 makeClusterListSimd2xnn(jGrid,
3133 nbl, ci, firstCell, lastCell,
3141 GMX_ASSERT(false, "Unhandled kernel type");
3145 static void makeClusterListWrapper(NbnxnPairlistGpu *nbl,
3146 const Grid &gmx_unused iGrid,
3149 const int firstCell,
3151 const bool excludeSubDiagonal,
3152 const nbnxn_atomdata_t *nbat,
3155 Nbnxm::KernelType gmx_unused kernelType,
3156 int *numDistanceChecks)
3158 for (int cj = firstCell; cj <= lastCell; cj++)
3160 make_cluster_list_supersub(iGrid, jGrid,
3163 nbat->xstride, nbat->x().data(),
3169 static int getNumSimpleJClustersInList(const NbnxnPairlistCpu &nbl)
3171 return nbl.cj.size();
3174 static int getNumSimpleJClustersInList(const gmx_unused NbnxnPairlistGpu &nbl)
3179 static void incrementNumSimpleJClustersInList(NbnxnPairlistCpu *nbl,
3182 nbl->ncjInUse += nbl->cj.size() - ncj_old_j;
3185 static void incrementNumSimpleJClustersInList(NbnxnPairlistGpu gmx_unused *nbl,
3186 int gmx_unused ncj_old_j)
3190 static void checkListSizeConsistency(const NbnxnPairlistCpu &nbl,
3191 const bool haveFreeEnergy)
3193 GMX_RELEASE_ASSERT(static_cast<size_t>(nbl.ncjInUse) == nbl.cj.size() || haveFreeEnergy,
3194 "Without free-energy all cj pair-list entries should be in use. "
3195 "Note that subsequent code does not make use of the equality, "
3196 "this check is only here to catch bugs");
3199 static void checkListSizeConsistency(const NbnxnPairlistGpu gmx_unused &nbl,
3200 bool gmx_unused haveFreeEnergy)
3202 /* We currently can not check consistency here */
3205 /* Set the buffer flags for newly added entries in the list */
3206 static void setBufferFlags(const NbnxnPairlistCpu &nbl,
3207 const int ncj_old_j,
3208 const int gridj_flag_shift,
3209 gmx_bitmask_t *gridj_flag,
3212 if (gmx::ssize(nbl.cj) > ncj_old_j)
3214 int cbFirst = nbl.cj[ncj_old_j].cj >> gridj_flag_shift;
3215 int cbLast = nbl.cj.back().cj >> gridj_flag_shift;
3216 for (int cb = cbFirst; cb <= cbLast; cb++)
3218 bitmask_init_bit(&gridj_flag[cb], th);
3223 static void setBufferFlags(const NbnxnPairlistGpu gmx_unused &nbl,
3224 int gmx_unused ncj_old_j,
3225 int gmx_unused gridj_flag_shift,
3226 gmx_bitmask_t gmx_unused *gridj_flag,
3229 GMX_ASSERT(false, "This function should never be called");
3232 /* Generates the part of pair-list nbl assigned to our thread */
3233 template <typename T>
3234 static void nbnxn_make_pairlist_part(const nbnxn_search *nbs,
3237 nbnxn_search_work_t *work,
3238 const nbnxn_atomdata_t *nbat,
3239 const t_blocka &exclusions,
3241 const Nbnxm::KernelType kernelType,
3243 gmx_bool bFBufferFlag,
3246 float nsubpair_tot_est,
3255 int ci_b, ci, ci_x, ci_y, ci_xy;
3257 real bx0, bx1, by0, by1, bz0, bz1;
3259 real d2cx, d2z, d2z_cx, d2z_cy, d2zx, d2zxy, d2xy;
3260 int cxf, cxl, cyf, cyf_x, cyl;
3261 int numDistanceChecks;
3262 int gridi_flag_shift = 0, gridj_flag_shift = 0;
3263 gmx_bitmask_t *gridj_flag = nullptr;
3264 int ncj_old_i, ncj_old_j;
3266 nbs_cycle_start(&work->cc[enbsCCsearch]);
3268 if (jGrid.geometry().isSimple != pairlistIsSimple(*nbl) ||
3269 iGrid.geometry().isSimple != pairlistIsSimple(*nbl))
3271 gmx_incons("Grid incompatible with pair-list");
3275 GMX_ASSERT(nbl->na_ci == jGrid.geometry().numAtomsICluster,
3276 "The cluster sizes in the list and grid should match");
3277 nbl->na_cj = Nbnxm::JClusterSizePerKernelType[kernelType];
3278 na_cj_2log = get_2log(nbl->na_cj);
3284 /* Determine conversion of clusters to flag blocks */
3285 gridi_flag_shift = getBufferFlagShift(nbl->na_ci);
3286 gridj_flag_shift = getBufferFlagShift(nbl->na_cj);
3288 gridj_flag = work->buffer_flags.flag;
3291 const Nbnxm::GridSet &gridSet = nbs->gridSet();
3293 gridSet.getBox(box);
3295 const bool haveFep = gridSet.haveFep();
3297 const real rlist2 = nbl->rlist*nbl->rlist;
3299 if (haveFep && !pairlistIsSimple(*nbl))
3301 /* Determine an atom-pair list cut-off distance for FEP atom pairs.
3302 * We should not simply use rlist, since then we would not have
3303 * the small, effective buffering of the NxN lists.
3304 * The buffer is on overestimate, but the resulting cost for pairs
3305 * beyond rlist is neglible compared to the FEP pairs within rlist.
3307 rl_fep2 = nbl->rlist + effective_buffer_1x1_vs_MxN(iGrid, jGrid);
3311 fprintf(debug, "nbl_fep atom-pair rlist %f\n", rl_fep2);
3313 rl_fep2 = rl_fep2*rl_fep2;
3316 const Grid::Dimensions &iGridDims = iGrid.dimensions();
3317 const Grid::Dimensions &jGridDims = jGrid.dimensions();
3319 rbb2 = boundingbox_only_distance2(iGridDims, jGridDims, nbl->rlist, pairlistIsSimple(*nbl));
3323 fprintf(debug, "nbl bounding box only distance %f\n", std::sqrt(rbb2));
3326 const bool isIntraGridList = (&iGrid == &jGrid);
3328 /* Set the shift range */
3329 for (int d = 0; d < DIM; d++)
3331 /* Check if we need periodicity shifts.
3332 * Without PBC or with domain decomposition we don't need them.
3334 if (d >= ePBC2npbcdim(nbs->ePBC) || nbs->dd_dim[d])
3340 const real listRangeCellToCell =
3341 listRangeForGridCellToGridCell(rlist, iGrid.dimensions(), jGrid.dimensions());
3343 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < listRangeCellToCell)
3353 const bool bSimple = pairlistIsSimple(*nbl);
3354 gmx::ArrayRef<const BoundingBox> bb_i;
3356 gmx::ArrayRef<const float> pbb_i;
3359 bb_i = iGrid.iBoundingBoxes();
3363 pbb_i = iGrid.packedBoundingBoxes();
3366 /* We use the normal bounding box format for both grid types */
3367 bb_i = iGrid.iBoundingBoxes();
3369 gmx::ArrayRef<const BoundingBox1D> bbcz_i = iGrid.zBoundingBoxes();
3370 gmx::ArrayRef<const int> flags_i = iGrid.clusterFlags();
3371 gmx::ArrayRef<const BoundingBox1D> bbcz_j = jGrid.zBoundingBoxes();
3372 int cell0_i = iGrid.cellOffset();
3376 fprintf(debug, "nbl nc_i %d col.av. %.1f ci_block %d\n",
3377 iGrid.numCells(), iGrid.numCells()/static_cast<double>(iGrid.numColumns()), ci_block);
3380 numDistanceChecks = 0;
3382 const real listRangeBBToJCell2 = gmx::square(listRangeForBoundingBoxToGridCell(rlist, jGrid.dimensions()));
3384 /* Initially ci_b and ci to 1 before where we want them to start,
3385 * as they will both be incremented in next_ci.
3388 ci = th*ci_block - 1;
3391 while (next_ci(iGrid, nth, ci_block, &ci_x, &ci_y, &ci_b, &ci))
3393 if (bSimple && flags_i[ci] == 0)
3398 ncj_old_i = getNumSimpleJClustersInList(*nbl);
3401 if (!isIntraGridList && shp[XX] == 0)
3405 bx1 = bb_i[ci].upper.x;
3409 bx1 = iGridDims.lowerCorner[XX] + (ci_x+1)*iGridDims.cellSize[XX];
3411 if (bx1 < jGridDims.lowerCorner[XX])
3413 d2cx = gmx::square(jGridDims.lowerCorner[XX] - bx1);
3415 if (d2cx >= listRangeBBToJCell2)
3422 ci_xy = ci_x*iGridDims.numCells[YY] + ci_y;
3424 /* Loop over shift vectors in three dimensions */
3425 for (int tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
3427 const real shz = tz*box[ZZ][ZZ];
3429 bz0 = bbcz_i[ci].lower + shz;
3430 bz1 = bbcz_i[ci].upper + shz;
3438 d2z = gmx::square(bz1);
3442 d2z = gmx::square(bz0 - box[ZZ][ZZ]);
3445 d2z_cx = d2z + d2cx;
3447 if (d2z_cx >= rlist2)
3452 bz1_frac = bz1/iGrid.numCellsInColumn(ci_xy);
3457 /* The check with bz1_frac close to or larger than 1 comes later */
3459 for (int ty = -shp[YY]; ty <= shp[YY]; ty++)
3461 const real shy = ty*box[YY][YY] + tz*box[ZZ][YY];
3465 by0 = bb_i[ci].lower.y + shy;
3466 by1 = bb_i[ci].upper.y + shy;
3470 by0 = iGridDims.lowerCorner[YY] + (ci_y )*iGridDims.cellSize[YY] + shy;
3471 by1 = iGridDims.lowerCorner[YY] + (ci_y + 1)*iGridDims.cellSize[YY] + shy;
3474 get_cell_range<YY>(by0, by1,
3485 if (by1 < jGridDims.lowerCorner[YY])
3487 d2z_cy += gmx::square(jGridDims.lowerCorner[YY] - by1);
3489 else if (by0 > jGridDims.upperCorner[YY])
3491 d2z_cy += gmx::square(by0 - jGridDims.upperCorner[YY]);
3494 for (int tx = -shp[XX]; tx <= shp[XX]; tx++)
3496 const int shift = XYZ2IS(tx, ty, tz);
3498 const bool excludeSubDiagonal = (isIntraGridList && shift == CENTRAL);
3500 if (c_pbcShiftBackward && isIntraGridList && shift > CENTRAL)
3505 const real shx = tx*box[XX][XX] + ty*box[YY][XX] + tz*box[ZZ][XX];
3509 bx0 = bb_i[ci].lower.x + shx;
3510 bx1 = bb_i[ci].upper.x + shx;
3514 bx0 = iGridDims.lowerCorner[XX] + (ci_x )*iGridDims.cellSize[XX] + shx;
3515 bx1 = iGridDims.lowerCorner[XX] + (ci_x+1)*iGridDims.cellSize[XX] + shx;
3518 get_cell_range<XX>(bx0, bx1,
3528 addNewIEntry(nbl, cell0_i+ci, shift, flags_i[ci]);
3530 if ((!c_pbcShiftBackward || excludeSubDiagonal) &&
3533 /* Leave the pairs with i > j.
3534 * x is the major index, so skip half of it.
3539 set_icell_bb(iGrid, ci, shx, shy, shz,
3542 icell_set_x(cell0_i+ci, shx, shy, shz,
3543 nbat->xstride, nbat->x().data(),
3547 for (int cx = cxf; cx <= cxl; cx++)
3550 if (jGridDims.lowerCorner[XX] + cx*jGridDims.cellSize[XX] > bx1)
3552 d2zx += gmx::square(jGridDims.lowerCorner[XX] + cx*jGridDims.cellSize[XX] - bx1);
3554 else if (jGridDims.lowerCorner[XX] + (cx+1)*jGridDims.cellSize[XX] < bx0)
3556 d2zx += gmx::square(jGridDims.lowerCorner[XX] + (cx+1)*jGridDims.cellSize[XX] - bx0);
3559 if (isIntraGridList &&
3561 (!c_pbcShiftBackward || shift == CENTRAL) &&
3564 /* Leave the pairs with i > j.
3565 * Skip half of y when i and j have the same x.
3574 for (int cy = cyf_x; cy <= cyl; cy++)
3576 const int columnStart = jGrid.firstCellInColumn(cx*jGridDims.numCells[YY] + cy);
3577 const int columnEnd = jGrid.firstCellInColumn(cx*jGridDims.numCells[YY] + cy + 1);
3580 if (jGridDims.lowerCorner[YY] + cy*jGridDims.cellSize[YY] > by1)
3582 d2zxy += gmx::square(jGridDims.lowerCorner[YY] + cy*jGridDims.cellSize[YY] - by1);
3584 else if (jGridDims.lowerCorner[YY] + (cy + 1)*jGridDims.cellSize[YY] < by0)
3586 d2zxy += gmx::square(jGridDims.lowerCorner[YY] + (cy + 1)*jGridDims.cellSize[YY] - by0);
3588 if (columnStart < columnEnd && d2zxy < listRangeBBToJCell2)
3590 /* To improve efficiency in the common case
3591 * of a homogeneous particle distribution,
3592 * we estimate the index of the middle cell
3593 * in range (midCell). We search down and up
3594 * starting from this index.
3596 * Note that the bbcz_j array contains bounds
3597 * for i-clusters, thus for clusters of 4 atoms.
3598 * For the common case where the j-cluster size
3599 * is 8, we could step with a stride of 2,
3600 * but we do not do this because it would
3601 * complicate this code even more.
3603 int midCell = columnStart + static_cast<int>(bz1_frac*(columnEnd - columnStart));
3604 if (midCell >= columnEnd)
3606 midCell = columnEnd - 1;
3611 /* Find the lowest cell that can possibly
3613 * Check if we hit the bottom of the grid,
3614 * if the j-cell is below the i-cell and if so,
3615 * if it is within range.
3617 int downTestCell = midCell;
3618 while (downTestCell >= columnStart &&
3619 (bbcz_j[downTestCell].upper >= bz0 ||
3620 d2xy + gmx::square(bbcz_j[downTestCell].upper - bz0) < rlist2))
3624 int firstCell = downTestCell + 1;
3626 /* Find the highest cell that can possibly
3628 * Check if we hit the top of the grid,
3629 * if the j-cell is above the i-cell and if so,
3630 * if it is within range.
3632 int upTestCell = midCell + 1;
3633 while (upTestCell < columnEnd &&
3634 (bbcz_j[upTestCell].lower <= bz1 ||
3635 d2xy + gmx::square(bbcz_j[upTestCell].lower - bz1) < rlist2))
3639 int lastCell = upTestCell - 1;
3641 #define NBNXN_REFCODE 0
3644 /* Simple reference code, for debugging,
3645 * overrides the more complex code above.
3647 firstCell = columnEnd;
3649 for (int k = columnStart; k < columnEnd; k++)
3651 if (d2xy + gmx::square(bbcz_j[k*NNBSBB_D + 1] - bz0) < rlist2 &&
3656 if (d2xy + gmx::square(bbcz_j[k*NNBSBB_D] - bz1) < rlist2 &&
3665 if (isIntraGridList)
3667 /* We want each atom/cell pair only once,
3668 * only use cj >= ci.
3670 if (!c_pbcShiftBackward || shift == CENTRAL)
3672 firstCell = std::max(firstCell, ci);
3676 if (firstCell <= lastCell)
3678 GMX_ASSERT(firstCell >= columnStart && lastCell < columnEnd, "The range should reside within the current grid column");
3680 /* For f buffer flags with simple lists */
3681 ncj_old_j = getNumSimpleJClustersInList(*nbl);
3683 makeClusterListWrapper(nbl,
3685 jGrid, firstCell, lastCell,
3690 &numDistanceChecks);
3694 setBufferFlags(*nbl, ncj_old_j, gridj_flag_shift,
3698 incrementNumSimpleJClustersInList(nbl, ncj_old_j);
3704 /* Set the exclusions for this ci list */
3705 setExclusionsForIEntry(gridSet,
3709 *getOpenIEntry(nbl),
3714 make_fep_list(gridSet.atomIndices(), nbat, nbl,
3719 iGrid, jGrid, nbl_fep);
3722 /* Close this ci list */
3725 progBal, nsubpair_tot_est,
3731 if (bFBufferFlag && getNumSimpleJClustersInList(*nbl) > ncj_old_i)
3733 bitmask_init_bit(&(work->buffer_flags.flag[(iGrid.cellOffset() + ci) >> gridi_flag_shift]), th);
3737 work->ndistc = numDistanceChecks;
3739 nbs_cycle_stop(&work->cc[enbsCCsearch]);
3741 checkListSizeConsistency(*nbl, haveFep);
3745 fprintf(debug, "number of distance checks %d\n", numDistanceChecks);
3747 print_nblist_statistics(debug, nbl, nbs, rlist);
3751 fprintf(debug, "nbl FEP list pairs: %d\n", nbl_fep->nrj);
3756 static void reduce_buffer_flags(const nbnxn_search *nbs,
3758 const nbnxn_buffer_flags_t *dest)
3760 for (int s = 0; s < nsrc; s++)
3762 gmx_bitmask_t * flag = nbs->work[s].buffer_flags.flag;
3764 for (int b = 0; b < dest->nflag; b++)
3766 bitmask_union(&(dest->flag[b]), flag[b]);
3771 static void print_reduction_cost(const nbnxn_buffer_flags_t *flags, int nout)
3773 int nelem, nkeep, ncopy, nred, out;
3774 gmx_bitmask_t mask_0;
3780 bitmask_init_bit(&mask_0, 0);
3781 for (int b = 0; b < flags->nflag; b++)
3783 if (bitmask_is_equal(flags->flag[b], mask_0))
3785 /* Only flag 0 is set, no copy of reduction required */
3789 else if (!bitmask_is_zero(flags->flag[b]))
3792 for (out = 0; out < nout; out++)
3794 if (bitmask_is_set(flags->flag[b], out))
3811 fprintf(debug, "nbnxn reduction: #flag %d #list %d elem %4.2f, keep %4.2f copy %4.2f red %4.2f\n",
3813 nelem/static_cast<double>(flags->nflag),
3814 nkeep/static_cast<double>(flags->nflag),
3815 ncopy/static_cast<double>(flags->nflag),
3816 nred/static_cast<double>(flags->nflag));
3819 /* Copies the list entries from src to dest when cjStart <= *cjGlobal < cjEnd.
3820 * *cjGlobal is updated with the cj count in src.
3821 * When setFlags==true, flag bit t is set in flag for all i and j clusters.
3823 template<bool setFlags>
3824 static void copySelectedListRange(const nbnxn_ci_t * gmx_restrict srcCi,
3825 const NbnxnPairlistCpu * gmx_restrict src,
3826 NbnxnPairlistCpu * gmx_restrict dest,
3827 gmx_bitmask_t *flag,
3828 int iFlagShift, int jFlagShift, int t)
3830 const int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3832 dest->ci.push_back(*srcCi);
3833 dest->ci.back().cj_ind_start = dest->cj.size();
3834 dest->ci.back().cj_ind_end = dest->cj.size() + ncj;
3838 bitmask_init_bit(&flag[srcCi->ci >> iFlagShift], t);
3841 for (int j = srcCi->cj_ind_start; j < srcCi->cj_ind_end; j++)
3843 dest->cj.push_back(src->cj[j]);
3847 /* NOTE: This is relatively expensive, since this
3848 * operation is done for all elements in the list,
3849 * whereas at list generation this is done only
3850 * once for each flag entry.
3852 bitmask_init_bit(&flag[src->cj[j].cj >> jFlagShift], t);
3857 /* This routine re-balances the pairlists such that all are nearly equally
3858 * sized. Only whole i-entries are moved between lists. These are moved
3859 * between the ends of the lists, such that the buffer reduction cost should
3860 * not change significantly.
3861 * Note that all original reduction flags are currently kept. This can lead
3862 * to reduction of parts of the force buffer that could be avoided. But since
3863 * the original lists are quite balanced, this will only give minor overhead.
3865 static void rebalanceSimpleLists(int numLists,
3866 NbnxnPairlistCpu * const * const srcSet,
3867 NbnxnPairlistCpu **destSet,
3868 gmx::ArrayRef<nbnxn_search_work_t> searchWork)
3871 for (int s = 0; s < numLists; s++)
3873 ncjTotal += srcSet[s]->ncjInUse;
3875 int ncjTarget = (ncjTotal + numLists - 1)/numLists;
3877 #pragma omp parallel num_threads(numLists)
3879 int t = gmx_omp_get_thread_num();
3881 int cjStart = ncjTarget* t;
3882 int cjEnd = ncjTarget*(t + 1);
3884 /* The destination pair-list for task/thread t */
3885 NbnxnPairlistCpu *dest = destSet[t];
3887 clear_pairlist(dest);
3888 dest->na_cj = srcSet[0]->na_cj;
3890 /* Note that the flags in the work struct (still) contain flags
3891 * for all entries that are present in srcSet->nbl[t].
3893 gmx_bitmask_t *flag = searchWork[t].buffer_flags.flag;
3895 int iFlagShift = getBufferFlagShift(dest->na_ci);
3896 int jFlagShift = getBufferFlagShift(dest->na_cj);
3899 for (int s = 0; s < numLists && cjGlobal < cjEnd; s++)
3901 const NbnxnPairlistCpu *src = srcSet[s];
3903 if (cjGlobal + src->ncjInUse > cjStart)
3905 for (gmx::index i = 0; i < gmx::ssize(src->ci) && cjGlobal < cjEnd; i++)
3907 const nbnxn_ci_t *srcCi = &src->ci[i];
3908 int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3909 if (cjGlobal >= cjStart)
3911 /* If the source list is not our own, we need to set
3912 * extra flags (the template bool parameter).
3916 copySelectedListRange
3919 flag, iFlagShift, jFlagShift, t);
3923 copySelectedListRange
3926 dest, flag, iFlagShift, jFlagShift, t);
3934 cjGlobal += src->ncjInUse;
3938 dest->ncjInUse = dest->cj.size();
3942 int ncjTotalNew = 0;
3943 for (int s = 0; s < numLists; s++)
3945 ncjTotalNew += destSet[s]->ncjInUse;
3947 GMX_RELEASE_ASSERT(ncjTotalNew == ncjTotal, "The total size of the lists before and after rebalancing should match");
3951 /* Returns if the pairlists are so imbalanced that it is worth rebalancing. */
3952 static bool checkRebalanceSimpleLists(const nbnxn_pairlist_set_t *listSet)
3954 int numLists = listSet->nnbl;
3957 for (int s = 0; s < numLists; s++)
3959 ncjMax = std::max(ncjMax, listSet->nbl[s]->ncjInUse);
3960 ncjTotal += listSet->nbl[s]->ncjInUse;
3964 fprintf(debug, "Pair-list ncjMax %d ncjTotal %d\n", ncjMax, ncjTotal);
3966 /* The rebalancing adds 3% extra time to the search. Heuristically we
3967 * determined that under common conditions the non-bonded kernel balance
3968 * improvement will outweigh this when the imbalance is more than 3%.
3969 * But this will, obviously, depend on search vs kernel time and nstlist.
3971 const real rebalanceTolerance = 1.03;
3973 return numLists*ncjMax > ncjTotal*rebalanceTolerance;
3976 /* Perform a count (linear) sort to sort the smaller lists to the end.
3977 * This avoids load imbalance on the GPU, as large lists will be
3978 * scheduled and executed first and the smaller lists later.
3979 * Load balancing between multi-processors only happens at the end
3980 * and there smaller lists lead to more effective load balancing.
3981 * The sorting is done on the cj4 count, not on the actual pair counts.
3982 * Not only does this make the sort faster, but it also results in
3983 * better load balancing than using a list sorted on exact load.
3984 * This function swaps the pointer in the pair list to avoid a copy operation.
3986 static void sort_sci(NbnxnPairlistGpu *nbl)
3988 if (nbl->cj4.size() <= nbl->sci.size())
3990 /* nsci = 0 or all sci have size 1, sorting won't change the order */
3994 NbnxnPairlistGpuWork &work = *nbl->work;
3996 /* We will distinguish differences up to double the average */
3997 const int m = (2*nbl->cj4.size())/nbl->sci.size();
3999 /* Resize work.sci_sort so we can sort into it */
4000 work.sci_sort.resize(nbl->sci.size());
4002 std::vector<int> &sort = work.sortBuffer;
4003 /* Set up m + 1 entries in sort, initialized at 0 */
4005 sort.resize(m + 1, 0);
4006 /* Count the entries of each size */
4007 for (const nbnxn_sci_t &sci : nbl->sci)
4009 int i = std::min(m, sci.numJClusterGroups());
4012 /* Calculate the offset for each count */
4015 for (int i = m - 1; i >= 0; i--)
4018 sort[i] = sort[i + 1] + s0;
4022 /* Sort entries directly into place */
4023 gmx::ArrayRef<nbnxn_sci_t> sci_sort = work.sci_sort;
4024 for (const nbnxn_sci_t &sci : nbl->sci)
4026 int i = std::min(m, sci.numJClusterGroups());
4027 sci_sort[sort[i]++] = sci;
4030 /* Swap the sci pointers so we use the new, sorted list */
4031 std::swap(nbl->sci, work.sci_sort);
4035 nonbonded_verlet_t::PairlistSets::construct(const InteractionLocality iLocality,
4037 nbnxn_atomdata_t *nbat,
4038 const t_blocka *excl,
4039 const Nbnxm::KernelType kernelType,
4043 nbnxn_pairlist_set_t *nbl_list = &pairlistSet(iLocality);
4045 const real rlist = nbl_list->params.rlistOuter;
4047 int nsubpair_target;
4048 float nsubpair_tot_est;
4051 gmx_bool CombineNBLists;
4053 int np_tot, np_noq, np_hlj, nap;
4055 nnbl = nbl_list->nnbl;
4056 CombineNBLists = nbl_list->bCombined;
4060 fprintf(debug, "ns making %d nblists\n", nnbl);
4063 nbat->bUseBufferFlags = (nbat->out.size() > 1);
4064 /* We should re-init the flags before making the first list */
4065 if (nbat->bUseBufferFlags && iLocality == InteractionLocality::Local)
4067 init_buffer_flags(&nbat->buffer_flags, nbat->numAtoms());
4071 if (iLocality == InteractionLocality::Local)
4073 /* Only zone (grid) 0 vs 0 */
4078 nzi = nbs->zones->nizone;
4081 if (!nbl_list->bSimple && minimumIlistCountForGpuBalancing_ > 0)
4083 get_nsubpair_target(nbs, iLocality, rlist, minimumIlistCountForGpuBalancing_,
4084 &nsubpair_target, &nsubpair_tot_est);
4088 nsubpair_target = 0;
4089 nsubpair_tot_est = 0;
4092 /* Clear all pair-lists */
4093 for (int th = 0; th < nnbl; th++)
4095 if (nbl_list->bSimple)
4097 clear_pairlist(nbl_list->nbl[th]);
4101 clear_pairlist(nbl_list->nblGpu[th]);
4104 if (nbs->gridSet().haveFep())
4106 clear_pairlist_fep(nbl_list->nbl_fep[th]);
4110 for (int zi = 0; zi < nzi; zi++)
4112 const Grid &iGrid = nbs->gridSet().grids()[zi];
4116 if (iLocality == InteractionLocality::Local)
4123 zj0 = nbs->zones->izone[zi].j0;
4124 zj1 = nbs->zones->izone[zi].j1;
4130 for (int zj = zj0; zj < zj1; zj++)
4132 const Grid &jGrid = nbs->gridSet().grids()[zj];
4136 fprintf(debug, "ns search grid %d vs %d\n", zi, zj);
4139 nbs_cycle_start(&nbs->cc[enbsCCsearch]);
4141 ci_block = get_ci_block_size(iGrid, nbs->DomDec, nnbl);
4143 /* With GPU: generate progressively smaller lists for
4144 * load balancing for local only or non-local with 2 zones.
4146 progBal = (iLocality == InteractionLocality::Local || nbs->zones->n <= 2);
4148 #pragma omp parallel for num_threads(nnbl) schedule(static)
4149 for (int th = 0; th < nnbl; th++)
4153 /* Re-init the thread-local work flag data before making
4154 * the first list (not an elegant conditional).
4156 if (nbat->bUseBufferFlags && ((zi == 0 && zj == 0)))
4158 init_buffer_flags(&nbs->work[th].buffer_flags, nbat->numAtoms());
4161 if (CombineNBLists && th > 0)
4163 GMX_ASSERT(!nbl_list->bSimple, "Can only combine GPU lists");
4165 clear_pairlist(nbl_list->nblGpu[th]);
4168 /* Divide the i super cell equally over the nblists */
4169 if (nbl_list->bSimple)
4171 nbnxn_make_pairlist_part(nbs, iGrid, jGrid,
4172 &nbs->work[th], nbat, *excl,
4176 nbat->bUseBufferFlags,
4178 progBal, nsubpair_tot_est,
4181 nbl_list->nbl_fep[th]);
4185 nbnxn_make_pairlist_part(nbs, iGrid, jGrid,
4186 &nbs->work[th], nbat, *excl,
4190 nbat->bUseBufferFlags,
4192 progBal, nsubpair_tot_est,
4194 nbl_list->nblGpu[th],
4195 nbl_list->nbl_fep[th]);
4198 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4200 nbs_cycle_stop(&nbs->cc[enbsCCsearch]);
4205 for (int th = 0; th < nnbl; th++)
4207 inc_nrnb(nrnb, eNR_NBNXN_DIST2, nbs->work[th].ndistc);
4209 if (nbl_list->bSimple)
4211 NbnxnPairlistCpu *nbl = nbl_list->nbl[th];
4212 np_tot += nbl->cj.size();
4213 np_noq += nbl->work->ncj_noq;
4214 np_hlj += nbl->work->ncj_hlj;
4218 NbnxnPairlistGpu *nbl = nbl_list->nblGpu[th];
4219 /* This count ignores potential subsequent pair pruning */
4220 np_tot += nbl->nci_tot;
4223 if (nbl_list->bSimple)
4225 nap = nbl_list->nbl[0]->na_ci*nbl_list->nbl[0]->na_cj;
4229 nap = gmx::square(nbl_list->nblGpu[0]->na_ci);
4231 nbl_list->natpair_ljq = (np_tot - np_noq)*nap - np_hlj*nap/2;
4232 nbl_list->natpair_lj = np_noq*nap;
4233 nbl_list->natpair_q = np_hlj*nap/2;
4235 if (CombineNBLists && nnbl > 1)
4237 GMX_ASSERT(!nbl_list->bSimple, "Can only combine GPU lists");
4238 NbnxnPairlistGpu **nbl = nbl_list->nblGpu;
4240 nbs_cycle_start(&nbs->cc[enbsCCcombine]);
4242 combine_nblists(nnbl-1, nbl+1, nbl[0]);
4244 nbs_cycle_stop(&nbs->cc[enbsCCcombine]);
4249 if (nbl_list->bSimple)
4251 if (nnbl > 1 && checkRebalanceSimpleLists(nbl_list))
4253 rebalanceSimpleLists(nbl_list->nnbl, nbl_list->nbl, nbl_list->nbl_work, nbs->work);
4255 /* Swap the pointer of the sets of pair lists */
4256 NbnxnPairlistCpu **tmp = nbl_list->nbl;
4257 nbl_list->nbl = nbl_list->nbl_work;
4258 nbl_list->nbl_work = tmp;
4263 /* Sort the entries on size, large ones first */
4264 if (CombineNBLists || nnbl == 1)
4266 sort_sci(nbl_list->nblGpu[0]);
4270 #pragma omp parallel for num_threads(nnbl) schedule(static)
4271 for (int th = 0; th < nnbl; th++)
4275 sort_sci(nbl_list->nblGpu[th]);
4277 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4282 if (nbat->bUseBufferFlags)
4284 reduce_buffer_flags(nbs, nbl_list->nnbl, &nbat->buffer_flags);
4287 if (nbs->gridSet().haveFep())
4289 /* Balance the free-energy lists over all the threads */
4290 balance_fep_lists(nbs, nbl_list);
4293 if (nbl_list->bSimple)
4295 /* This is a fresh list, so not pruned, stored using ci.
4296 * ciOuter is invalid at this point.
4298 GMX_ASSERT(nbl_list->nbl[0]->ciOuter.empty(), "ciOuter is invalid so it should be empty");
4301 if (iLocality == Nbnxm::InteractionLocality::Local)
4303 outerListCreationStep_ = step;
4307 GMX_RELEASE_ASSERT(outerListCreationStep_ == step,
4308 "Outer list should be created at the same step as the inner list");
4311 /* Special performance logging stuff (env.var. GMX_NBNXN_CYCLE) */
4312 if (iLocality == InteractionLocality::Local)
4314 nbs->search_count++;
4316 if (nbs->print_cycles &&
4317 (!nbs->DomDec || iLocality == InteractionLocality::NonLocal) &&
4318 nbs->search_count % 100 == 0)
4320 nbs_cycle_print(stderr, nbs);
4323 /* If we have more than one list, they either got rebalancing (CPU)
4324 * or combined (GPU), so we should dump the final result to debug.
4326 if (debug && nbl_list->nnbl > 1)
4328 if (nbl_list->bSimple)
4330 for (int t = 0; t < nbl_list->nnbl; t++)
4332 print_nblist_statistics(debug, nbl_list->nbl[t], nbs, rlist);
4337 print_nblist_statistics(debug, nbl_list->nblGpu[0], nbs, rlist);
4345 if (nbl_list->bSimple)
4347 for (int t = 0; t < nbl_list->nnbl; t++)
4349 print_nblist_ci_cj(debug, nbl_list->nbl[t]);
4354 print_nblist_sci_cj(debug, nbl_list->nblGpu[0]);
4358 if (nbat->bUseBufferFlags)
4360 print_reduction_cost(&nbat->buffer_flags, nbl_list->nnbl);
4364 if (params_.useDynamicPruning && nbl_list->bSimple)
4366 nbnxnPrepareListForDynamicPruning(nbl_list);
4371 nonbonded_verlet_t::constructPairlist(const Nbnxm::InteractionLocality iLocality,
4372 const t_blocka *excl,
4376 pairlistSets_->construct(iLocality, nbs.get(), nbat.get(), excl,
4377 kernelSetup_.kernelType,
4382 /* Launch the transfer of the pairlist to the GPU.
4384 * NOTE: The launch overhead is currently not timed separately
4386 Nbnxm::gpu_init_pairlist(gpu_nbv,
4387 pairlistSets().pairlistSet(iLocality).nblGpu[0],
4392 void nbnxnPrepareListForDynamicPruning(nbnxn_pairlist_set_t *listSet)
4394 GMX_RELEASE_ASSERT(listSet->bSimple, "Should only be called for simple lists");
4396 /* TODO: Restructure the lists so we have actual outer and inner
4397 * list objects so we can set a single pointer instead of
4398 * swapping several pointers.
4401 for (int i = 0; i < listSet->nnbl; i++)
4403 NbnxnPairlistCpu &list = *listSet->nbl[i];
4405 /* The search produced a list in ci/cj.
4406 * Swap the list pointers so we get the outer list is ciOuter,cjOuter
4407 * and we can prune that to get an inner list in ci/cj.
4409 GMX_RELEASE_ASSERT(list.ciOuter.empty() && list.cjOuter.empty(),
4410 "The outer lists should be empty before preparation");
4412 std::swap(list.ci, list.ciOuter);
4413 std::swap(list.cj, list.cjOuter);