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 nbnxn_bb_t = Nbnxm::BoundingBox; // TODO: Remove when refactoring this file
81 using Grid = Nbnxm::Grid; // TODO: Remove when refactoring this file
83 // Convience alias for partial Nbnxn namespace usage
84 using InteractionLocality = Nbnxm::InteractionLocality;
86 /* We shift the i-particles backward for PBC.
87 * This leads to more conditionals than shifting forward.
88 * We do this to get more balanced pair lists.
90 constexpr bool c_pbcShiftBackward = true;
93 static void nbs_cycle_clear(nbnxn_cycle_t *cc)
95 for (int i = 0; i < enbsCCnr; i++)
102 static double Mcyc_av(const nbnxn_cycle_t *cc)
104 return static_cast<double>(cc->c)*1e-6/cc->count;
107 static void nbs_cycle_print(FILE *fp, const nbnxn_search *nbs)
110 fprintf(fp, "ns %4d grid %4.1f search %4.1f red.f %5.3f",
111 nbs->cc[enbsCCgrid].count,
112 Mcyc_av(&nbs->cc[enbsCCgrid]),
113 Mcyc_av(&nbs->cc[enbsCCsearch]),
114 Mcyc_av(&nbs->cc[enbsCCreducef]));
116 if (nbs->work.size() > 1)
118 if (nbs->cc[enbsCCcombine].count > 0)
120 fprintf(fp, " comb %5.2f",
121 Mcyc_av(&nbs->cc[enbsCCcombine]));
123 fprintf(fp, " s. th");
124 for (const nbnxn_search_work_t &work : nbs->work)
126 fprintf(fp, " %4.1f",
127 Mcyc_av(&work.cc[enbsCCsearch]));
133 /* Layout for the nonbonded NxN pair lists */
134 enum class NbnxnLayout
136 NoSimd4x4, // i-cluster size 4, j-cluster size 4
137 Simd4xN, // i-cluster size 4, j-cluster size SIMD width
138 Simd2xNN, // i-cluster size 4, j-cluster size half SIMD width
139 Gpu8x8x8 // i-cluster size 8, j-cluster size 8 + super-clustering
143 /* Returns the j-cluster size */
144 template <NbnxnLayout layout>
145 static constexpr int jClusterSize()
147 static_assert(layout == NbnxnLayout::NoSimd4x4 || layout == NbnxnLayout::Simd4xN || layout == NbnxnLayout::Simd2xNN, "Currently jClusterSize only supports CPU layouts");
149 return layout == NbnxnLayout::Simd4xN ? GMX_SIMD_REAL_WIDTH : (layout == NbnxnLayout::Simd2xNN ? GMX_SIMD_REAL_WIDTH/2 : c_nbnxnCpuIClusterSize);
152 /*! \brief Returns the j-cluster index given the i-cluster index.
154 * \tparam jClusterSize The number of atoms in a j-cluster
155 * \tparam jSubClusterIndex The j-sub-cluster index (0/1), used when size(j-cluster) < size(i-cluster)
156 * \param[in] ci The i-cluster index
158 template <int jClusterSize, int jSubClusterIndex>
159 static inline int cjFromCi(int ci)
161 static_assert(jClusterSize == c_nbnxnCpuIClusterSize/2 || jClusterSize == c_nbnxnCpuIClusterSize || jClusterSize == c_nbnxnCpuIClusterSize*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
163 static_assert(jSubClusterIndex == 0 || jSubClusterIndex == 1,
164 "Only sub-cluster indices 0 and 1 are supported");
166 if (jClusterSize == c_nbnxnCpuIClusterSize/2)
168 if (jSubClusterIndex == 0)
174 return ((ci + 1) << 1) - 1;
177 else if (jClusterSize == c_nbnxnCpuIClusterSize)
187 /*! \brief Returns the j-cluster index given the i-cluster index.
189 * \tparam layout The pair-list layout
190 * \tparam jSubClusterIndex The j-sub-cluster index (0/1), used when size(j-cluster) < size(i-cluster)
191 * \param[in] ci The i-cluster index
193 template <NbnxnLayout layout, int jSubClusterIndex>
194 static inline int cjFromCi(int ci)
196 constexpr int clusterSize = jClusterSize<layout>();
198 return cjFromCi<clusterSize, jSubClusterIndex>(ci);
201 /* Returns the nbnxn coordinate data index given the i-cluster index */
202 template <NbnxnLayout layout>
203 static inline int xIndexFromCi(int ci)
205 constexpr int clusterSize = jClusterSize<layout>();
207 static_assert(clusterSize == c_nbnxnCpuIClusterSize/2 || clusterSize == c_nbnxnCpuIClusterSize || clusterSize == c_nbnxnCpuIClusterSize*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
209 if (clusterSize <= c_nbnxnCpuIClusterSize)
211 /* Coordinates are stored packed in groups of 4 */
216 /* Coordinates packed in 8, i-cluster size is half the packing width */
217 return (ci >> 1)*STRIDE_P8 + (ci & 1)*(c_packX8 >> 1);
221 /* Returns the nbnxn coordinate data index given the j-cluster index */
222 template <NbnxnLayout layout>
223 static inline int xIndexFromCj(int cj)
225 constexpr int clusterSize = jClusterSize<layout>();
227 static_assert(clusterSize == c_nbnxnCpuIClusterSize/2 || clusterSize == c_nbnxnCpuIClusterSize || clusterSize == c_nbnxnCpuIClusterSize*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
229 if (clusterSize == c_nbnxnCpuIClusterSize/2)
231 /* Coordinates are stored packed in groups of 4 */
232 return (cj >> 1)*STRIDE_P4 + (cj & 1)*(c_packX4 >> 1);
234 else if (clusterSize == c_nbnxnCpuIClusterSize)
236 /* Coordinates are stored packed in groups of 4 */
241 /* Coordinates are stored packed in groups of 8 */
247 /* Initializes a single nbnxn_pairlist_t data structure */
248 static void nbnxn_init_pairlist_fep(t_nblist *nl)
250 nl->type = GMX_NBLIST_INTERACTION_FREE_ENERGY;
251 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
252 /* The interaction functions are set in the free energy kernel fuction */
265 nl->jindex = nullptr;
267 nl->excl_fep = nullptr;
271 static void free_nblist(t_nblist *nl)
281 nbnxn_search_work_t::nbnxn_search_work_t() :
284 buffer_flags({0, nullptr, 0}),
286 nbl_fep(new t_nblist),
289 nbnxn_init_pairlist_fep(nbl_fep.get());
294 nbnxn_search_work_t::~nbnxn_search_work_t()
296 sfree(buffer_flags.flag);
298 free_nblist(nbl_fep.get());
301 nbnxn_search::nbnxn_search(int ePBC,
302 const ivec *n_dd_cells,
303 const gmx_domdec_zones_t *zones,
304 const PairlistType pairlistType,
315 // The correct value will be set during the gridding
319 DomDec = n_dd_cells != nullptr;
322 for (int d = 0; d < DIM; d++)
324 if ((*n_dd_cells)[d] > 1)
327 /* Each grid matches a DD zone */
333 grid.resize(numGrids, pairlistType);
335 /* Initialize detailed nbsearch cycle counting */
336 print_cycles = (getenv("GMX_NBNXN_CYCLE") != nullptr);
340 static void init_buffer_flags(nbnxn_buffer_flags_t *flags,
343 flags->nflag = (natoms + NBNXN_BUFFERFLAG_SIZE - 1)/NBNXN_BUFFERFLAG_SIZE;
344 if (flags->nflag > flags->flag_nalloc)
346 flags->flag_nalloc = over_alloc_large(flags->nflag);
347 srenew(flags->flag, flags->flag_nalloc);
349 for (int b = 0; b < flags->nflag; b++)
351 bitmask_clear(&(flags->flag[b]));
355 /* Returns the pair-list cutoff between a bounding box and a grid cell given an atom-to-atom pair-list cutoff
357 * Given a cutoff distance between atoms, this functions returns the cutoff
358 * distance2 between a bounding box of a group of atoms and a grid cell.
359 * Since atoms can be geometrically outside of the cell they have been
360 * assigned to (when atom groups instead of individual atoms are assigned
361 * to cells), this distance returned can be larger than the input.
364 listRangeForBoundingBoxToGridCell(real rlist,
365 const Grid::Dimensions &gridDims)
367 return rlist + gridDims.maxAtomGroupRadius;
370 /* Returns the pair-list cutoff between a grid cells given an atom-to-atom pair-list cutoff
372 * Given a cutoff distance between atoms, this functions returns the cutoff
373 * distance2 between two grid cells.
374 * Since atoms can be geometrically outside of the cell they have been
375 * assigned to (when atom groups instead of individual atoms are assigned
376 * to cells), this distance returned can be larger than the input.
379 listRangeForGridCellToGridCell(real rlist,
380 const Grid::Dimensions &iGridDims,
381 const Grid::Dimensions &jGridDims)
383 return rlist + iGridDims.maxAtomGroupRadius + jGridDims.maxAtomGroupRadius;
386 /* Determines the cell range along one dimension that
387 * the bounding box b0 - b1 sees.
390 static void get_cell_range(real b0, real b1,
391 const Grid::Dimensions &jGridDims,
392 real d2, real rlist, int *cf, int *cl)
394 real listRangeBBToCell2 = gmx::square(listRangeForBoundingBoxToGridCell(rlist, jGridDims));
395 real distanceInCells = (b0 - jGridDims.lowerCorner[dim])*jGridDims.invCellSize[dim];
396 *cf = std::max(static_cast<int>(distanceInCells), 0);
399 d2 + gmx::square((b0 - jGridDims.lowerCorner[dim]) - (*cf - 1 + 1)*jGridDims.cellSize[dim]) < listRangeBBToCell2)
404 *cl = std::min(static_cast<int>((b1 - jGridDims.lowerCorner[dim])*jGridDims.invCellSize[dim]), jGridDims.numCells[dim] - 1);
405 while (*cl < jGridDims.numCells[dim] - 1 &&
406 d2 + gmx::square((*cl + 1)*jGridDims.cellSize[dim] - (b1 - jGridDims.lowerCorner[dim])) < listRangeBBToCell2)
412 /* Reference code calculating the distance^2 between two bounding boxes */
414 static float box_dist2(float bx0, float bx1, float by0,
415 float by1, float bz0, float bz1,
416 const nbnxn_bb_t *bb)
419 float dl, dh, dm, dm0;
423 dl = bx0 - bb->upper.x;
424 dh = bb->lower.x - bx1;
425 dm = std::max(dl, dh);
426 dm0 = std::max(dm, 0.0f);
429 dl = by0 - bb->upper.y;
430 dh = bb->lower.y - by1;
431 dm = std::max(dl, dh);
432 dm0 = std::max(dm, 0.0f);
435 dl = bz0 - bb->upper.z;
436 dh = bb->lower.z - bz1;
437 dm = std::max(dl, dh);
438 dm0 = std::max(dm, 0.0f);
445 /* Plain C code calculating the distance^2 between two bounding boxes */
446 static float subc_bb_dist2(int si,
447 const nbnxn_bb_t *bb_i_ci,
449 gmx::ArrayRef<const nbnxn_bb_t> bb_j_all)
451 const nbnxn_bb_t &bb_i = bb_i_ci[si];
452 const nbnxn_bb_t &bb_j = bb_j_all[csj];
455 float dl, dh, dm, dm0;
457 dl = bb_i.lower.x - bb_j.upper.x;
458 dh = bb_j.lower.x - bb_i.upper.x;
459 dm = std::max(dl, dh);
460 dm0 = std::max(dm, 0.0f);
463 dl = bb_i.lower.y - bb_j.upper.y;
464 dh = bb_j.lower.y - bb_i.upper.y;
465 dm = std::max(dl, dh);
466 dm0 = std::max(dm, 0.0f);
469 dl = bb_i.lower.z - bb_j.upper.z;
470 dh = bb_j.lower.z - bb_i.upper.z;
471 dm = std::max(dl, dh);
472 dm0 = std::max(dm, 0.0f);
478 #if NBNXN_SEARCH_BB_SIMD4
480 /* 4-wide SIMD code for bb distance for bb format xyz0 */
481 static float subc_bb_dist2_simd4(int si,
482 const nbnxn_bb_t *bb_i_ci,
484 gmx::ArrayRef<const nbnxn_bb_t> bb_j_all)
486 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
489 Simd4Float bb_i_S0, bb_i_S1;
490 Simd4Float bb_j_S0, bb_j_S1;
496 bb_i_S0 = load4(bb_i_ci[si].lower.ptr());
497 bb_i_S1 = load4(bb_i_ci[si].upper.ptr());
498 bb_j_S0 = load4(bb_j_all[csj].lower.ptr());
499 bb_j_S1 = load4(bb_j_all[csj].upper.ptr());
501 dl_S = bb_i_S0 - bb_j_S1;
502 dh_S = bb_j_S0 - bb_i_S1;
504 dm_S = max(dl_S, dh_S);
505 dm0_S = max(dm_S, simd4SetZeroF());
507 return dotProduct(dm0_S, dm0_S);
510 /* Calculate bb bounding distances of bb_i[si,...,si+3] and store them in d2 */
511 #define SUBC_BB_DIST2_SIMD4_XXXX_INNER(si, bb_i, d2) \
515 Simd4Float dx_0, dy_0, dz_0; \
516 Simd4Float dx_1, dy_1, dz_1; \
518 Simd4Float mx, my, mz; \
519 Simd4Float m0x, m0y, m0z; \
521 Simd4Float d2x, d2y, d2z; \
522 Simd4Float d2s, d2t; \
524 shi = (si)*NNBSBB_D*DIM; \
526 xi_l = load4((bb_i)+shi+0*STRIDE_PBB); \
527 yi_l = load4((bb_i)+shi+1*STRIDE_PBB); \
528 zi_l = load4((bb_i)+shi+2*STRIDE_PBB); \
529 xi_h = load4((bb_i)+shi+3*STRIDE_PBB); \
530 yi_h = load4((bb_i)+shi+4*STRIDE_PBB); \
531 zi_h = load4((bb_i)+shi+5*STRIDE_PBB); \
533 dx_0 = xi_l - xj_h; \
534 dy_0 = yi_l - yj_h; \
535 dz_0 = zi_l - zj_h; \
537 dx_1 = xj_l - xi_h; \
538 dy_1 = yj_l - yi_h; \
539 dz_1 = zj_l - zi_h; \
541 mx = max(dx_0, dx_1); \
542 my = max(dy_0, dy_1); \
543 mz = max(dz_0, dz_1); \
545 m0x = max(mx, zero); \
546 m0y = max(my, zero); \
547 m0z = max(mz, zero); \
556 store4((d2)+(si), d2t); \
559 /* 4-wide SIMD code for nsi bb distances for bb format xxxxyyyyzzzz */
560 static void subc_bb_dist2_simd4_xxxx(const float *bb_j,
561 int nsi, const float *bb_i,
564 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
567 Simd4Float xj_l, yj_l, zj_l;
568 Simd4Float xj_h, yj_h, zj_h;
569 Simd4Float xi_l, yi_l, zi_l;
570 Simd4Float xi_h, yi_h, zi_h;
576 xj_l = Simd4Float(bb_j[0*STRIDE_PBB]);
577 yj_l = Simd4Float(bb_j[1*STRIDE_PBB]);
578 zj_l = Simd4Float(bb_j[2*STRIDE_PBB]);
579 xj_h = Simd4Float(bb_j[3*STRIDE_PBB]);
580 yj_h = Simd4Float(bb_j[4*STRIDE_PBB]);
581 zj_h = Simd4Float(bb_j[5*STRIDE_PBB]);
583 /* Here we "loop" over si (0,STRIDE_PBB) from 0 to nsi with step STRIDE_PBB.
584 * But as we know the number of iterations is 1 or 2, we unroll manually.
586 SUBC_BB_DIST2_SIMD4_XXXX_INNER(0, bb_i, d2);
587 if (STRIDE_PBB < nsi)
589 SUBC_BB_DIST2_SIMD4_XXXX_INNER(STRIDE_PBB, bb_i, d2);
593 #endif /* NBNXN_SEARCH_BB_SIMD4 */
596 /* Returns if any atom pair from two clusters is within distance sqrt(rlist2) */
597 static inline gmx_bool
598 clusterpair_in_range(const NbnxnPairlistGpuWork &work,
600 int csj, int stride, const real *x_j,
603 #if !GMX_SIMD4_HAVE_REAL
606 * All coordinates are stored as xyzxyz...
609 const real *x_i = work.iSuperClusterData.x.data();
611 for (int i = 0; i < c_nbnxnGpuClusterSize; i++)
613 int i0 = (si*c_nbnxnGpuClusterSize + i)*DIM;
614 for (int j = 0; j < c_nbnxnGpuClusterSize; j++)
616 int j0 = (csj*c_nbnxnGpuClusterSize + j)*stride;
618 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]);
629 #else /* !GMX_SIMD4_HAVE_REAL */
631 /* 4-wide SIMD version.
632 * The coordinates x_i are stored as xxxxyyyy..., x_j is stored xyzxyz...
633 * Using 8-wide AVX(2) is not faster on Intel Sandy Bridge and Haswell.
635 static_assert(c_nbnxnGpuClusterSize == 8 || c_nbnxnGpuClusterSize == 4,
636 "A cluster is hard-coded to 4/8 atoms.");
638 Simd4Real rc2_S = Simd4Real(rlist2);
640 const real *x_i = work.iSuperClusterData.xSimd.data();
642 int dim_stride = c_nbnxnGpuClusterSize*DIM;
643 Simd4Real ix_S0 = load4(x_i + si*dim_stride + 0*GMX_SIMD4_WIDTH);
644 Simd4Real iy_S0 = load4(x_i + si*dim_stride + 1*GMX_SIMD4_WIDTH);
645 Simd4Real iz_S0 = load4(x_i + si*dim_stride + 2*GMX_SIMD4_WIDTH);
647 Simd4Real ix_S1, iy_S1, iz_S1;
648 if (c_nbnxnGpuClusterSize == 8)
650 ix_S1 = load4(x_i + si*dim_stride + 3*GMX_SIMD4_WIDTH);
651 iy_S1 = load4(x_i + si*dim_stride + 4*GMX_SIMD4_WIDTH);
652 iz_S1 = load4(x_i + si*dim_stride + 5*GMX_SIMD4_WIDTH);
654 /* We loop from the outer to the inner particles to maximize
655 * the chance that we find a pair in range quickly and return.
657 int j0 = csj*c_nbnxnGpuClusterSize;
658 int j1 = j0 + c_nbnxnGpuClusterSize - 1;
661 Simd4Real jx0_S, jy0_S, jz0_S;
662 Simd4Real jx1_S, jy1_S, jz1_S;
664 Simd4Real dx_S0, dy_S0, dz_S0;
665 Simd4Real dx_S1, dy_S1, dz_S1;
666 Simd4Real dx_S2, dy_S2, dz_S2;
667 Simd4Real dx_S3, dy_S3, dz_S3;
678 Simd4Bool wco_any_S01, wco_any_S23, wco_any_S;
680 jx0_S = Simd4Real(x_j[j0*stride+0]);
681 jy0_S = Simd4Real(x_j[j0*stride+1]);
682 jz0_S = Simd4Real(x_j[j0*stride+2]);
684 jx1_S = Simd4Real(x_j[j1*stride+0]);
685 jy1_S = Simd4Real(x_j[j1*stride+1]);
686 jz1_S = Simd4Real(x_j[j1*stride+2]);
688 /* Calculate distance */
689 dx_S0 = ix_S0 - jx0_S;
690 dy_S0 = iy_S0 - jy0_S;
691 dz_S0 = iz_S0 - jz0_S;
692 dx_S2 = ix_S0 - jx1_S;
693 dy_S2 = iy_S0 - jy1_S;
694 dz_S2 = iz_S0 - jz1_S;
695 if (c_nbnxnGpuClusterSize == 8)
697 dx_S1 = ix_S1 - jx0_S;
698 dy_S1 = iy_S1 - jy0_S;
699 dz_S1 = iz_S1 - jz0_S;
700 dx_S3 = ix_S1 - jx1_S;
701 dy_S3 = iy_S1 - jy1_S;
702 dz_S3 = iz_S1 - jz1_S;
705 /* rsq = dx*dx+dy*dy+dz*dz */
706 rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
707 rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
708 if (c_nbnxnGpuClusterSize == 8)
710 rsq_S1 = norm2(dx_S1, dy_S1, dz_S1);
711 rsq_S3 = norm2(dx_S3, dy_S3, dz_S3);
714 wco_S0 = (rsq_S0 < rc2_S);
715 wco_S2 = (rsq_S2 < rc2_S);
716 if (c_nbnxnGpuClusterSize == 8)
718 wco_S1 = (rsq_S1 < rc2_S);
719 wco_S3 = (rsq_S3 < rc2_S);
721 if (c_nbnxnGpuClusterSize == 8)
723 wco_any_S01 = wco_S0 || wco_S1;
724 wco_any_S23 = wco_S2 || wco_S3;
725 wco_any_S = wco_any_S01 || wco_any_S23;
729 wco_any_S = wco_S0 || wco_S2;
732 if (anyTrue(wco_any_S))
743 #endif /* !GMX_SIMD4_HAVE_REAL */
746 /* Returns the j-cluster index for index cjIndex in a cj list */
747 static inline int nblCj(gmx::ArrayRef<const nbnxn_cj_t> cjList,
750 return cjList[cjIndex].cj;
753 /* Returns the j-cluster index for index cjIndex in a cj4 list */
754 static inline int nblCj(gmx::ArrayRef<const nbnxn_cj4_t> cj4List,
757 return cj4List[cjIndex/c_nbnxnGpuJgroupSize].cj[cjIndex & (c_nbnxnGpuJgroupSize - 1)];
760 /* Returns the i-interaction mask of the j sub-cell for index cj_ind */
761 static unsigned int nbl_imask0(const NbnxnPairlistGpu *nbl, int cj_ind)
763 return nbl->cj4[cj_ind/c_nbnxnGpuJgroupSize].imei[0].imask;
766 /* Initializes a single NbnxnPairlistCpu data structure */
767 static void nbnxn_init_pairlist(NbnxnPairlistCpu *nbl)
769 nbl->na_ci = c_nbnxnCpuIClusterSize;
772 nbl->ciOuter.clear();
775 nbl->cjOuter.clear();
778 nbl->work = new NbnxnPairlistCpuWork();
781 NbnxnPairlistGpu::NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy) :
782 na_ci(c_nbnxnGpuClusterSize),
783 na_cj(c_nbnxnGpuClusterSize),
784 na_sc(c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize),
786 sci({}, {pinningPolicy}),
787 cj4({}, {pinningPolicy}),
788 excl({}, {pinningPolicy}),
791 static_assert(c_nbnxnGpuNumClusterPerSupercluster == c_gpuNumClusterPerCell,
792 "The search code assumes that the a super-cluster matches a search grid cell");
794 static_assert(sizeof(cj4[0].imei[0].imask)*8 >= c_nbnxnGpuJgroupSize*c_gpuNumClusterPerCell,
795 "The i super-cluster cluster interaction mask does not contain a sufficient number of bits");
797 static_assert(sizeof(excl[0])*8 >= c_nbnxnGpuJgroupSize*c_gpuNumClusterPerCell, "The GPU exclusion mask does not contain a sufficient number of bits");
799 // We always want a first entry without any exclusions
802 work = new NbnxnPairlistGpuWork();
805 void nbnxn_init_pairlist_set(nbnxn_pairlist_set_t *nbl_list)
808 (nbl_list->params.pairlistType == PairlistType::Simple4x2 ||
809 nbl_list->params.pairlistType == PairlistType::Simple4x4 ||
810 nbl_list->params.pairlistType == PairlistType::Simple4x8);
811 // Currently GPU lists are always combined
812 nbl_list->bCombined = !nbl_list->bSimple;
814 nbl_list->nnbl = gmx_omp_nthreads_get(emntNonbonded);
816 if (!nbl_list->bCombined &&
817 nbl_list->nnbl > NBNXN_BUFFERFLAG_MAX_THREADS)
819 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.",
820 nbl_list->nnbl, NBNXN_BUFFERFLAG_MAX_THREADS, NBNXN_BUFFERFLAG_MAX_THREADS);
823 if (nbl_list->bSimple)
825 snew(nbl_list->nbl, nbl_list->nnbl);
826 if (nbl_list->nnbl > 1)
828 snew(nbl_list->nbl_work, nbl_list->nnbl);
833 snew(nbl_list->nblGpu, nbl_list->nnbl);
835 nbl_list->nbl_fep.resize(nbl_list->nnbl);
836 /* Execute in order to avoid memory interleaving between threads */
837 #pragma omp parallel for num_threads(nbl_list->nnbl) schedule(static)
838 for (int i = 0; i < nbl_list->nnbl; i++)
842 /* Allocate the nblist data structure locally on each thread
843 * to optimize memory access for NUMA architectures.
845 if (nbl_list->bSimple)
847 nbl_list->nbl[i] = new NbnxnPairlistCpu();
849 nbnxn_init_pairlist(nbl_list->nbl[i]);
850 if (nbl_list->nnbl > 1)
852 nbl_list->nbl_work[i] = new NbnxnPairlistCpu();
853 nbnxn_init_pairlist(nbl_list->nbl_work[i]);
858 /* Only list 0 is used on the GPU, use normal allocation for i>0 */
859 auto pinningPolicy = (i == 0 ? gmx::PinningPolicy::PinnedIfSupported : gmx::PinningPolicy::CannotBePinned);
861 nbl_list->nblGpu[i] = new NbnxnPairlistGpu(pinningPolicy);
864 snew(nbl_list->nbl_fep[i], 1);
865 nbnxn_init_pairlist_fep(nbl_list->nbl_fep[i]);
867 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
871 /* Print statistics of a pair list, used for debug output */
872 static void print_nblist_statistics(FILE *fp, const NbnxnPairlistCpu *nbl,
873 const nbnxn_search *nbs, real rl)
875 const Grid &grid = nbs->grid[0];
876 const Grid::Dimensions &dims = grid.dimensions();
878 fprintf(fp, "nbl nci %zu ncj %d\n",
879 nbl->ci.size(), nbl->ncjInUse);
880 const int numAtomsJCluster = grid.geometry().numAtomsJCluster;
881 const double numAtomsPerCell = nbl->ncjInUse/static_cast<double>(grid.numCells())*numAtomsJCluster;
882 fprintf(fp, "nbl na_cj %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
883 nbl->na_cj, rl, nbl->ncjInUse, nbl->ncjInUse/static_cast<double>(grid.numCells()),
885 numAtomsPerCell/(0.5*4.0/3.0*M_PI*rl*rl*rl*grid.numCells()*numAtomsJCluster/(dims.gridSize[XX]*dims.gridSize[YY]*dims.gridSize[ZZ])));
887 fprintf(fp, "nbl average j cell list length %.1f\n",
888 0.25*nbl->ncjInUse/std::max(static_cast<double>(nbl->ci.size()), 1.0));
890 int cs[SHIFTS] = { 0 };
892 for (const nbnxn_ci_t &ciEntry : nbl->ci)
894 cs[ciEntry.shift & NBNXN_CI_SHIFT] +=
895 ciEntry.cj_ind_end - ciEntry.cj_ind_start;
897 int j = ciEntry.cj_ind_start;
898 while (j < ciEntry.cj_ind_end &&
899 nbl->cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
905 fprintf(fp, "nbl cell pairs, total: %zu excl: %d %.1f%%\n",
906 nbl->cj.size(), npexcl, 100*npexcl/std::max(static_cast<double>(nbl->cj.size()), 1.0));
907 for (int s = 0; s < SHIFTS; s++)
911 fprintf(fp, "nbl shift %2d ncj %3d\n", s, cs[s]);
916 /* Print statistics of a pair lists, used for debug output */
917 static void print_nblist_statistics(FILE *fp, const NbnxnPairlistGpu *nbl,
918 const nbnxn_search *nbs, real rl)
920 const Grid &grid = nbs->grid[0];
921 const Grid::Dimensions &dims = grid.dimensions();
923 fprintf(fp, "nbl nsci %zu ncj4 %zu nsi %d excl4 %zu\n",
924 nbl->sci.size(), nbl->cj4.size(), nbl->nci_tot, nbl->excl.size());
925 const int numAtomsCluster = grid.geometry().numAtomsICluster;
926 const double numAtomsPerCell = nbl->nci_tot/static_cast<double>(grid.numClusters())*numAtomsCluster;
927 fprintf(fp, "nbl na_c %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
928 nbl->na_ci, rl, nbl->nci_tot, nbl->nci_tot/static_cast<double>(grid.numClusters()),
930 numAtomsPerCell/(0.5*4.0/3.0*M_PI*rl*rl*rl*grid.numClusters()*numAtomsCluster/(dims.gridSize[XX]*dims.gridSize[YY]*dims.gridSize[ZZ])));
935 int c[c_gpuNumClusterPerCell + 1] = { 0 };
936 for (const nbnxn_sci_t &sci : nbl->sci)
939 for (int j4 = sci.cj4_ind_start; j4 < sci.cj4_ind_end; j4++)
941 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
944 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
946 if (nbl->cj4[j4].imei[0].imask & (1U << (j*c_gpuNumClusterPerCell + si)))
957 nsp_max = std::max(nsp_max, nsp);
959 if (!nbl->sci.empty())
961 sum_nsp /= nbl->sci.size();
962 sum_nsp2 /= nbl->sci.size();
964 fprintf(fp, "nbl #cluster-pairs: av %.1f stddev %.1f max %d\n",
965 sum_nsp, std::sqrt(sum_nsp2 - sum_nsp*sum_nsp), nsp_max);
967 if (!nbl->cj4.empty())
969 for (int b = 0; b <= c_gpuNumClusterPerCell; b++)
971 fprintf(fp, "nbl j-list #i-subcell %d %7d %4.1f\n",
972 b, c[b], 100.0*c[b]/size_t {nbl->cj4.size()*c_nbnxnGpuJgroupSize});
977 /* Returns a pointer to the exclusion mask for j-cluster-group \p cj4 and warp \p warp
978 * Generates a new exclusion entry when the j-cluster-group uses
979 * the default all-interaction mask at call time, so the returned mask
980 * can be modified when needed.
982 static nbnxn_excl_t *get_exclusion_mask(NbnxnPairlistGpu *nbl,
986 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
988 /* No exclusions set, make a new list entry */
989 const size_t oldSize = nbl->excl.size();
990 GMX_ASSERT(oldSize >= 1, "We should always have entry [0]");
991 /* Add entry with default values: no exclusions */
992 nbl->excl.resize(oldSize + 1);
993 nbl->cj4[cj4].imei[warp].excl_ind = oldSize;
996 return &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
999 static void set_self_and_newton_excls_supersub(NbnxnPairlistGpu *nbl,
1000 int cj4_ind, int sj_offset,
1001 int i_cluster_in_cell)
1003 nbnxn_excl_t *excl[c_nbnxnGpuClusterpairSplit];
1005 /* Here we only set the set self and double pair exclusions */
1007 /* Reserve extra elements, so the resize() in get_exclusion_mask()
1008 * will not invalidate excl entries in the loop below
1010 nbl->excl.reserve(nbl->excl.size() + c_nbnxnGpuClusterpairSplit);
1011 for (int w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
1013 excl[w] = get_exclusion_mask(nbl, cj4_ind, w);
1016 /* Only minor < major bits set */
1017 for (int ej = 0; ej < nbl->na_ci; ej++)
1020 for (int ei = ej; ei < nbl->na_ci; ei++)
1022 excl[w]->pair[(ej & (c_nbnxnGpuJgroupSize-1))*nbl->na_ci + ei] &=
1023 ~(1U << (sj_offset*c_gpuNumClusterPerCell + i_cluster_in_cell));
1028 /* Returns a diagonal or off-diagonal interaction mask for plain C lists */
1029 static unsigned int get_imask(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=2 */
1035 gmx_unused static unsigned int get_imask_simd_j2(gmx_bool rdiag, int ci, int cj)
1037 return (rdiag && ci*2 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_0 :
1038 (rdiag && ci*2+1 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_1 :
1039 NBNXN_INTERACTION_MASK_ALL));
1042 /* Returns a diagonal or off-diagonal interaction mask for cj-size=4 */
1043 gmx_unused static unsigned int get_imask_simd_j4(gmx_bool rdiag, int ci, int cj)
1045 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
1048 /* Returns a diagonal or off-diagonal interaction mask for cj-size=8 */
1049 gmx_unused static unsigned int get_imask_simd_j8(gmx_bool rdiag, int ci, int cj)
1051 return (rdiag && ci == cj*2 ? NBNXN_INTERACTION_MASK_DIAG_J8_0 :
1052 (rdiag && ci == cj*2+1 ? NBNXN_INTERACTION_MASK_DIAG_J8_1 :
1053 NBNXN_INTERACTION_MASK_ALL));
1057 #if GMX_SIMD_REAL_WIDTH == 2
1058 #define get_imask_simd_4xn get_imask_simd_j2
1060 #if GMX_SIMD_REAL_WIDTH == 4
1061 #define get_imask_simd_4xn get_imask_simd_j4
1063 #if GMX_SIMD_REAL_WIDTH == 8
1064 #define get_imask_simd_4xn get_imask_simd_j8
1065 #define get_imask_simd_2xnn get_imask_simd_j4
1067 #if GMX_SIMD_REAL_WIDTH == 16
1068 #define get_imask_simd_2xnn get_imask_simd_j8
1072 /* Plain C code for checking and adding cluster-pairs to the list.
1074 * \param[in] gridj The j-grid
1075 * \param[in,out] nbl The pair-list to store the cluster pairs in
1076 * \param[in] icluster The index of the i-cluster
1077 * \param[in] jclusterFirst The first cluster in the j-range
1078 * \param[in] jclusterLast The last cluster in the j-range
1079 * \param[in] excludeSubDiagonal Exclude atom pairs with i-index > j-index
1080 * \param[in] x_j Coordinates for the j-atom, in xyz format
1081 * \param[in] rlist2 The squared list cut-off
1082 * \param[in] rbb2 The squared cut-off for putting cluster-pairs in the list based on bounding box distance only
1083 * \param[in,out] numDistanceChecks The number of distance checks performed
1086 makeClusterListSimple(const Grid &jGrid,
1087 NbnxnPairlistCpu * nbl,
1091 bool excludeSubDiagonal,
1092 const real * gmx_restrict x_j,
1095 int * gmx_restrict numDistanceChecks)
1097 const nbnxn_bb_t * gmx_restrict bb_ci = nbl->work->iClusterData.bb.data();
1098 const real * gmx_restrict x_ci = nbl->work->iClusterData.x.data();
1103 while (!InRange && jclusterFirst <= jclusterLast)
1105 real d2 = subc_bb_dist2(0, bb_ci, jclusterFirst, jGrid.jBoundingBoxes());
1106 *numDistanceChecks += 2;
1108 /* Check if the distance is within the distance where
1109 * we use only the bounding box distance rbb,
1110 * or within the cut-off and there is at least one atom pair
1111 * within the cut-off.
1117 else if (d2 < rlist2)
1119 int cjf_gl = jGrid.cellOffset() + jclusterFirst;
1120 for (int i = 0; i < c_nbnxnCpuIClusterSize && !InRange; i++)
1122 for (int j = 0; j < c_nbnxnCpuIClusterSize; j++)
1124 InRange = InRange ||
1125 (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjf_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+XX]) +
1126 gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjf_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+YY]) +
1127 gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjf_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+ZZ]) < rlist2);
1130 *numDistanceChecks += c_nbnxnCpuIClusterSize*c_nbnxnCpuIClusterSize;
1143 while (!InRange && jclusterLast > jclusterFirst)
1145 real d2 = subc_bb_dist2(0, bb_ci, jclusterLast, jGrid.jBoundingBoxes());
1146 *numDistanceChecks += 2;
1148 /* Check if the distance is within the distance where
1149 * we use only the bounding box distance rbb,
1150 * or within the cut-off and there is at least one atom pair
1151 * within the cut-off.
1157 else if (d2 < rlist2)
1159 int cjl_gl = jGrid.cellOffset() + jclusterLast;
1160 for (int i = 0; i < c_nbnxnCpuIClusterSize && !InRange; i++)
1162 for (int j = 0; j < c_nbnxnCpuIClusterSize; j++)
1164 InRange = InRange ||
1165 (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjl_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+XX]) +
1166 gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjl_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+YY]) +
1167 gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjl_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+ZZ]) < rlist2);
1170 *numDistanceChecks += c_nbnxnCpuIClusterSize*c_nbnxnCpuIClusterSize;
1178 if (jclusterFirst <= jclusterLast)
1180 for (int jcluster = jclusterFirst; jcluster <= jclusterLast; jcluster++)
1182 /* Store cj and the interaction mask */
1184 cjEntry.cj = jGrid.cellOffset() + jcluster;
1185 cjEntry.excl = get_imask(excludeSubDiagonal, icluster, jcluster);
1186 nbl->cj.push_back(cjEntry);
1188 /* Increase the closing index in the i list */
1189 nbl->ci.back().cj_ind_end = nbl->cj.size();
1193 #ifdef GMX_NBNXN_SIMD_4XN
1194 #include "gromacs/nbnxm/pairlist_simd_4xm.h"
1196 #ifdef GMX_NBNXN_SIMD_2XNN
1197 #include "gromacs/nbnxm/pairlist_simd_2xmm.h"
1200 /* Plain C or SIMD4 code for making a pair list of super-cell sci vs scj.
1201 * Checks bounding box distances and possibly atom pair distances.
1203 static void make_cluster_list_supersub(const Grid &iGrid,
1205 NbnxnPairlistGpu *nbl,
1208 const bool excludeSubDiagonal,
1213 int *numDistanceChecks)
1215 NbnxnPairlistGpuWork &work = *nbl->work;
1218 const float *pbb_ci = work.iSuperClusterData.bbPacked.data();
1220 const nbnxn_bb_t *bb_ci = work.iSuperClusterData.bb.data();
1223 assert(c_nbnxnGpuClusterSize == iGrid.geometry().numAtomsICluster);
1224 assert(c_nbnxnGpuClusterSize == jGrid.geometry().numAtomsICluster);
1226 /* We generate the pairlist mainly based on bounding-box distances
1227 * and do atom pair distance based pruning on the GPU.
1228 * Only if a j-group contains a single cluster-pair, we try to prune
1229 * that pair based on atom distances on the CPU to avoid empty j-groups.
1231 #define PRUNE_LIST_CPU_ONE 1
1232 #define PRUNE_LIST_CPU_ALL 0
1234 #if PRUNE_LIST_CPU_ONE
1238 float *d2l = work.distanceBuffer.data();
1240 for (int subc = 0; subc < jGrid.numClustersPerCell()[scj]; subc++)
1242 const int cj4_ind = work.cj_ind/c_nbnxnGpuJgroupSize;
1243 const int cj_offset = work.cj_ind - cj4_ind*c_nbnxnGpuJgroupSize;
1244 const int cj = scj*c_gpuNumClusterPerCell + subc;
1246 const int cj_gl = jGrid.cellOffset()*c_gpuNumClusterPerCell + cj;
1249 if (excludeSubDiagonal && sci == scj)
1255 ci1 = iGrid.numClustersPerCell()[sci];
1259 /* Determine all ci1 bb distances in one call with SIMD4 */
1260 subc_bb_dist2_simd4_xxxx(jGrid.packedBoundingBoxes().data() + (cj >> STRIDE_PBB_2LOG)*NNBSBB_XXXX + (cj & (STRIDE_PBB-1)),
1262 *numDistanceChecks += c_nbnxnGpuClusterSize*2;
1266 unsigned int imask = 0;
1267 /* We use a fixed upper-bound instead of ci1 to help optimization */
1268 for (int ci = 0; ci < c_gpuNumClusterPerCell; ci++)
1276 /* Determine the bb distance between ci and cj */
1277 d2l[ci] = subc_bb_dist2(ci, bb_ci, cj, jGrid.jBoundingBoxes());
1278 *numDistanceChecks += 2;
1282 #if PRUNE_LIST_CPU_ALL
1283 /* Check if the distance is within the distance where
1284 * we use only the bounding box distance rbb,
1285 * or within the cut-off and there is at least one atom pair
1286 * within the cut-off. This check is very costly.
1288 *numDistanceChecks += c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize;
1291 clusterpair_in_range(work, ci, cj_gl, stride, x, rlist2)))
1293 /* Check if the distance between the two bounding boxes
1294 * in within the pair-list cut-off.
1299 /* Flag this i-subcell to be taken into account */
1300 imask |= (1U << (cj_offset*c_gpuNumClusterPerCell + ci));
1302 #if PRUNE_LIST_CPU_ONE
1310 #if PRUNE_LIST_CPU_ONE
1311 /* If we only found 1 pair, check if any atoms are actually
1312 * within the cut-off, so we could get rid of it.
1314 if (npair == 1 && d2l[ci_last] >= rbb2 &&
1315 !clusterpair_in_range(work, ci_last, cj_gl, stride, x, rlist2))
1317 imask &= ~(1U << (cj_offset*c_gpuNumClusterPerCell + ci_last));
1324 /* We have at least one cluster pair: add a j-entry */
1325 if (static_cast<size_t>(cj4_ind) == nbl->cj4.size())
1327 nbl->cj4.resize(nbl->cj4.size() + 1);
1329 nbnxn_cj4_t *cj4 = &nbl->cj4[cj4_ind];
1331 cj4->cj[cj_offset] = cj_gl;
1333 /* Set the exclusions for the ci==sj entry.
1334 * Here we don't bother to check if this entry is actually flagged,
1335 * as it will nearly always be in the list.
1337 if (excludeSubDiagonal && sci == scj)
1339 set_self_and_newton_excls_supersub(nbl, cj4_ind, cj_offset, subc);
1342 /* Copy the cluster interaction mask to the list */
1343 for (int w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
1345 cj4->imei[w].imask |= imask;
1348 nbl->work->cj_ind++;
1350 /* Keep the count */
1351 nbl->nci_tot += npair;
1353 /* Increase the closing index in i super-cell list */
1354 nbl->sci.back().cj4_ind_end =
1355 (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize;
1360 /* Returns how many contiguous j-clusters we have starting in the i-list */
1361 template <typename CjListType>
1362 static int numContiguousJClusters(const int cjIndexStart,
1363 const int cjIndexEnd,
1364 gmx::ArrayRef<const CjListType> cjList)
1366 const int firstJCluster = nblCj(cjList, cjIndexStart);
1368 int numContiguous = 0;
1370 while (cjIndexStart + numContiguous < cjIndexEnd &&
1371 nblCj(cjList, cjIndexStart + numContiguous) == firstJCluster + numContiguous)
1376 return numContiguous;
1380 * \brief Helper struct for efficient searching for excluded atoms in a j-list
1384 /*! \brief Constructs a j-list range from \p cjList with the given index range */
1385 template <typename CjListType>
1386 JListRanges(int cjIndexStart,
1388 gmx::ArrayRef<const CjListType> cjList);
1390 int cjIndexStart; //!< The start index in the j-list
1391 int cjIndexEnd; //!< The end index in the j-list
1392 int cjFirst; //!< The j-cluster with index cjIndexStart
1393 int cjLast; //!< The j-cluster with index cjIndexEnd-1
1394 int numDirect; //!< Up to cjIndexStart+numDirect the j-clusters are cjFirst + the index offset
1398 template <typename CjListType>
1399 JListRanges::JListRanges(int cjIndexStart,
1401 gmx::ArrayRef<const CjListType> cjList) :
1402 cjIndexStart(cjIndexStart),
1403 cjIndexEnd(cjIndexEnd)
1405 GMX_ASSERT(cjIndexEnd > cjIndexStart, "JListRanges should only be called with non-empty lists");
1407 cjFirst = nblCj(cjList, cjIndexStart);
1408 cjLast = nblCj(cjList, cjIndexEnd - 1);
1410 /* Determine how many contiguous j-cells we have starting
1411 * from the first i-cell. This number can be used to directly
1412 * calculate j-cell indices for excluded atoms.
1414 numDirect = numContiguousJClusters(cjIndexStart, cjIndexEnd, cjList);
1418 /* Return the index of \p jCluster in the given range or -1 when not present
1420 * Note: This code is executed very often and therefore performance is
1421 * important. It should be inlined and fully optimized.
1423 template <typename CjListType>
1425 findJClusterInJList(int jCluster,
1426 const JListRanges &ranges,
1427 gmx::ArrayRef<const CjListType> cjList)
1431 if (jCluster < ranges.cjFirst + ranges.numDirect)
1433 /* We can calculate the index directly using the offset */
1434 index = ranges.cjIndexStart + jCluster - ranges.cjFirst;
1438 /* Search for jCluster using bisection */
1440 int rangeStart = ranges.cjIndexStart + ranges.numDirect;
1441 int rangeEnd = ranges.cjIndexEnd;
1443 while (index == -1 && rangeStart < rangeEnd)
1445 rangeMiddle = (rangeStart + rangeEnd) >> 1;
1447 const int clusterMiddle = nblCj(cjList, rangeMiddle);
1449 if (jCluster == clusterMiddle)
1451 index = rangeMiddle;
1453 else if (jCluster < clusterMiddle)
1455 rangeEnd = rangeMiddle;
1459 rangeStart = rangeMiddle + 1;
1467 // TODO: Get rid of the two functions below by renaming sci to ci (or something better)
1469 /* Return the i-entry in the list we are currently operating on */
1470 static nbnxn_ci_t *getOpenIEntry(NbnxnPairlistCpu *nbl)
1472 return &nbl->ci.back();
1475 /* Return the i-entry in the list we are currently operating on */
1476 static nbnxn_sci_t *getOpenIEntry(NbnxnPairlistGpu *nbl)
1478 return &nbl->sci.back();
1481 /* Set all atom-pair exclusions for a simple type list i-entry
1483 * Set all atom-pair exclusions from the topology stored in exclusions
1484 * as masks in the pair-list for simple list entry iEntry.
1487 setExclusionsForIEntry(const nbnxn_search *nbs,
1488 NbnxnPairlistCpu *nbl,
1489 gmx_bool diagRemoved,
1491 const nbnxn_ci_t &iEntry,
1492 const t_blocka &exclusions)
1494 if (iEntry.cj_ind_end == iEntry.cj_ind_start)
1496 /* Empty list: no exclusions */
1500 const JListRanges ranges(iEntry.cj_ind_start, iEntry.cj_ind_end, gmx::makeConstArrayRef(nbl->cj));
1502 const int iCluster = iEntry.ci;
1504 gmx::ArrayRef<const int> cell = nbs->cell;
1506 /* Loop over the atoms in the i-cluster */
1507 for (int i = 0; i < nbl->na_ci; i++)
1509 const int iIndex = iCluster*nbl->na_ci + i;
1510 const int iAtom = nbs->a[iIndex];
1513 /* Loop over the topology-based exclusions for this i-atom */
1514 for (int exclIndex = exclusions.index[iAtom]; exclIndex < exclusions.index[iAtom + 1]; exclIndex++)
1516 const int jAtom = exclusions.a[exclIndex];
1520 /* The self exclusion are already set, save some time */
1524 /* Get the index of the j-atom in the nbnxn atom data */
1525 const int jIndex = cell[jAtom];
1527 /* Without shifts we only calculate interactions j>i
1528 * for one-way pair-lists.
1530 if (diagRemoved && jIndex <= iIndex)
1535 const int jCluster = (jIndex >> na_cj_2log);
1537 /* Could the cluster se be in our list? */
1538 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
1541 findJClusterInJList(jCluster, ranges,
1542 gmx::makeConstArrayRef(nbl->cj));
1546 /* We found an exclusion, clear the corresponding
1549 const int innerJ = jIndex - (jCluster << na_cj_2log);
1551 nbl->cj[index].excl &= ~(1U << ((i << na_cj_2log) + innerJ));
1559 /* Add a new i-entry to the FEP list and copy the i-properties */
1560 static inline void fep_list_new_nri_copy(t_nblist *nlist)
1562 /* Add a new i-entry */
1565 assert(nlist->nri < nlist->maxnri);
1567 /* Duplicate the last i-entry, except for jindex, which continues */
1568 nlist->iinr[nlist->nri] = nlist->iinr[nlist->nri-1];
1569 nlist->shift[nlist->nri] = nlist->shift[nlist->nri-1];
1570 nlist->gid[nlist->nri] = nlist->gid[nlist->nri-1];
1571 nlist->jindex[nlist->nri] = nlist->nrj;
1574 /* For load balancing of the free-energy lists over threads, we set
1575 * the maximum nrj size of an i-entry to 40. This leads to good
1576 * load balancing in the worst case scenario of a single perturbed
1577 * particle on 16 threads, while not introducing significant overhead.
1578 * Note that half of the perturbed pairs will anyhow end up in very small lists,
1579 * since non perturbed i-particles will see few perturbed j-particles).
1581 const int max_nrj_fep = 40;
1583 /* Exclude the perturbed pairs from the Verlet list. This is only done to avoid
1584 * singularities for overlapping particles (0/0), since the charges and
1585 * LJ parameters have been zeroed in the nbnxn data structure.
1586 * Simultaneously make a group pair list for the perturbed pairs.
1588 static void make_fep_list(const nbnxn_search *nbs,
1589 const nbnxn_atomdata_t *nbat,
1590 NbnxnPairlistCpu *nbl,
1591 gmx_bool bDiagRemoved,
1593 real gmx_unused shx,
1594 real gmx_unused shy,
1595 real gmx_unused shz,
1596 real gmx_unused rlist_fep2,
1601 int ci, cj_ind_start, cj_ind_end, cja, cjr;
1603 int gid_i = 0, gid_j, gid;
1604 int egp_shift, egp_mask;
1606 int ind_i, ind_j, ai, aj;
1608 gmx_bool bFEP_i, bFEP_i_all;
1610 if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
1618 cj_ind_start = nbl_ci->cj_ind_start;
1619 cj_ind_end = nbl_ci->cj_ind_end;
1621 /* In worst case we have alternating energy groups
1622 * and create #atom-pair lists, which means we need the size
1623 * of a cluster pair (na_ci*na_cj) times the number of cj's.
1625 nri_max = nbl->na_ci*nbl->na_cj*(cj_ind_end - cj_ind_start);
1626 if (nlist->nri + nri_max > nlist->maxnri)
1628 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1629 reallocate_nblist(nlist);
1632 const int numAtomsJCluster = jGrid.geometry().numAtomsJCluster;
1634 const nbnxn_atomdata_t::Params &nbatParams = nbat->params();
1636 const int ngid = nbatParams.nenergrp;
1638 /* TODO: Consider adding a check in grompp and changing this to an assert */
1639 const int numBitsInEnergyGroupIdsForAtomsInJCluster = sizeof(gid_cj)*8;
1640 if (ngid*numAtomsJCluster > numBitsInEnergyGroupIdsForAtomsInJCluster)
1642 gmx_fatal(FARGS, "The Verlet scheme with %dx%d kernels and free-energy only supports up to %zu energy groups",
1643 iGrid.geometry().numAtomsICluster, numAtomsJCluster,
1644 (sizeof(gid_cj)*8)/numAtomsJCluster);
1647 egp_shift = nbatParams.neg_2log;
1648 egp_mask = (1 << egp_shift) - 1;
1650 /* Loop over the atoms in the i sub-cell */
1652 for (int i = 0; i < nbl->na_ci; i++)
1654 ind_i = ci*nbl->na_ci + i;
1659 nlist->jindex[nri+1] = nlist->jindex[nri];
1660 nlist->iinr[nri] = ai;
1661 /* The actual energy group pair index is set later */
1662 nlist->gid[nri] = 0;
1663 nlist->shift[nri] = nbl_ci->shift & NBNXN_CI_SHIFT;
1665 bFEP_i = iGrid.atomIsPerturbed(ci - iGrid.cellOffset(), i);
1667 bFEP_i_all = bFEP_i_all && bFEP_i;
1669 if (nlist->nrj + (cj_ind_end - cj_ind_start)*nbl->na_cj > nlist->maxnrj)
1671 nlist->maxnrj = over_alloc_small(nlist->nrj + (cj_ind_end - cj_ind_start)*nbl->na_cj);
1672 srenew(nlist->jjnr, nlist->maxnrj);
1673 srenew(nlist->excl_fep, nlist->maxnrj);
1678 gid_i = (nbatParams.energrp[ci] >> (egp_shift*i)) & egp_mask;
1681 for (int cj_ind = cj_ind_start; cj_ind < cj_ind_end; cj_ind++)
1683 unsigned int fep_cj;
1685 cja = nbl->cj[cj_ind].cj;
1687 if (numAtomsJCluster == jGrid.geometry().numAtomsICluster)
1689 cjr = cja - jGrid.cellOffset();
1690 fep_cj = jGrid.fepBits(cjr);
1693 gid_cj = nbatParams.energrp[cja];
1696 else if (2*numAtomsJCluster == jGrid.geometry().numAtomsICluster)
1698 cjr = cja - jGrid.cellOffset()*2;
1699 /* Extract half of the ci fep/energrp mask */
1700 fep_cj = (jGrid.fepBits(cjr >> 1) >> ((cjr & 1)*numAtomsJCluster)) & ((1 << numAtomsJCluster) - 1);
1703 gid_cj = nbatParams.energrp[cja >> 1] >> ((cja & 1)*numAtomsJCluster*egp_shift) & ((1 << (numAtomsJCluster*egp_shift)) - 1);
1708 cjr = cja - (jGrid.cellOffset() >> 1);
1709 /* Combine two ci fep masks/energrp */
1710 fep_cj = jGrid.fepBits(cjr*2) + (jGrid.fepBits(cjr*2 + 1) << jGrid.geometry().numAtomsICluster);
1713 gid_cj = nbatParams.energrp[cja*2] + (nbatParams.energrp[cja*2+1] << (jGrid.geometry().numAtomsICluster*egp_shift));
1717 if (bFEP_i || fep_cj != 0)
1719 for (int j = 0; j < nbl->na_cj; j++)
1721 /* Is this interaction perturbed and not excluded? */
1722 ind_j = cja*nbl->na_cj + j;
1725 (bFEP_i || (fep_cj & (1 << j))) &&
1726 (!bDiagRemoved || ind_j >= ind_i))
1730 gid_j = (gid_cj >> (j*egp_shift)) & egp_mask;
1731 gid = GID(gid_i, gid_j, ngid);
1733 if (nlist->nrj > nlist->jindex[nri] &&
1734 nlist->gid[nri] != gid)
1736 /* Energy group pair changed: new list */
1737 fep_list_new_nri_copy(nlist);
1740 nlist->gid[nri] = gid;
1743 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1745 fep_list_new_nri_copy(nlist);
1749 /* Add it to the FEP list */
1750 nlist->jjnr[nlist->nrj] = aj;
1751 nlist->excl_fep[nlist->nrj] = (nbl->cj[cj_ind].excl >> (i*nbl->na_cj + j)) & 1;
1754 /* Exclude it from the normal list.
1755 * Note that the charge has been set to zero,
1756 * but we need to avoid 0/0, as perturbed atoms
1757 * can be on top of each other.
1759 nbl->cj[cj_ind].excl &= ~(1U << (i*nbl->na_cj + j));
1765 if (nlist->nrj > nlist->jindex[nri])
1767 /* Actually add this new, non-empty, list */
1769 nlist->jindex[nlist->nri] = nlist->nrj;
1776 /* All interactions are perturbed, we can skip this entry */
1777 nbl_ci->cj_ind_end = cj_ind_start;
1778 nbl->ncjInUse -= cj_ind_end - cj_ind_start;
1782 /* Return the index of atom a within a cluster */
1783 static inline int cj_mod_cj4(int cj)
1785 return cj & (c_nbnxnGpuJgroupSize - 1);
1788 /* Convert a j-cluster to a cj4 group */
1789 static inline int cj_to_cj4(int cj)
1791 return cj/c_nbnxnGpuJgroupSize;
1794 /* Return the index of an j-atom within a warp */
1795 static inline int a_mod_wj(int a)
1797 return a & (c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit - 1);
1800 /* As make_fep_list above, but for super/sub lists. */
1801 static void make_fep_list(const nbnxn_search *nbs,
1802 const nbnxn_atomdata_t *nbat,
1803 NbnxnPairlistGpu *nbl,
1804 gmx_bool bDiagRemoved,
1805 const nbnxn_sci_t *nbl_sci,
1816 int ind_i, ind_j, ai, aj;
1820 const nbnxn_cj4_t *cj4;
1822 const int numJClusterGroups = nbl_sci->numJClusterGroups();
1823 if (numJClusterGroups == 0)
1829 const int sci = nbl_sci->sci;
1831 const int cj4_ind_start = nbl_sci->cj4_ind_start;
1832 const int cj4_ind_end = nbl_sci->cj4_ind_end;
1834 /* Here we process one super-cell, max #atoms na_sc, versus a list
1835 * cj4 entries, each with max c_nbnxnGpuJgroupSize cj's, each
1836 * of size na_cj atoms.
1837 * On the GPU we don't support energy groups (yet).
1838 * So for each of the na_sc i-atoms, we need max one FEP list
1839 * for each max_nrj_fep j-atoms.
1841 nri_max = nbl->na_sc*nbl->na_cj*(1 + (numJClusterGroups*c_nbnxnGpuJgroupSize)/max_nrj_fep);
1842 if (nlist->nri + nri_max > nlist->maxnri)
1844 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1845 reallocate_nblist(nlist);
1848 /* Loop over the atoms in the i super-cluster */
1849 for (int c = 0; c < c_gpuNumClusterPerCell; c++)
1851 c_abs = sci*c_gpuNumClusterPerCell + c;
1853 for (int i = 0; i < nbl->na_ci; i++)
1855 ind_i = c_abs*nbl->na_ci + i;
1860 nlist->jindex[nri+1] = nlist->jindex[nri];
1861 nlist->iinr[nri] = ai;
1862 /* With GPUs, energy groups are not supported */
1863 nlist->gid[nri] = 0;
1864 nlist->shift[nri] = nbl_sci->shift & NBNXN_CI_SHIFT;
1866 bFEP_i = iGrid.atomIsPerturbed(c_abs - iGrid.cellOffset()*c_gpuNumClusterPerCell, i);
1868 xi = nbat->x()[ind_i*nbat->xstride+XX] + shx;
1869 yi = nbat->x()[ind_i*nbat->xstride+YY] + shy;
1870 zi = nbat->x()[ind_i*nbat->xstride+ZZ] + shz;
1872 const int nrjMax = nlist->nrj + numJClusterGroups*c_nbnxnGpuJgroupSize*nbl->na_cj;
1873 if (nrjMax > nlist->maxnrj)
1875 nlist->maxnrj = over_alloc_small(nrjMax);
1876 srenew(nlist->jjnr, nlist->maxnrj);
1877 srenew(nlist->excl_fep, nlist->maxnrj);
1880 for (int cj4_ind = cj4_ind_start; cj4_ind < cj4_ind_end; cj4_ind++)
1882 cj4 = &nbl->cj4[cj4_ind];
1884 for (int gcj = 0; gcj < c_nbnxnGpuJgroupSize; gcj++)
1886 if ((cj4->imei[0].imask & (1U << (gcj*c_gpuNumClusterPerCell + c))) == 0)
1888 /* Skip this ci for this cj */
1893 cj4->cj[gcj] - jGrid.cellOffset()*c_gpuNumClusterPerCell;
1895 if (bFEP_i || jGrid.clusterIsPerturbed(cjr))
1897 for (int j = 0; j < nbl->na_cj; j++)
1899 /* Is this interaction perturbed and not excluded? */
1900 ind_j = (jGrid.cellOffset()*c_gpuNumClusterPerCell + cjr)*nbl->na_cj + j;
1903 (bFEP_i || jGrid.atomIsPerturbed(cjr, j)) &&
1904 (!bDiagRemoved || ind_j >= ind_i))
1907 unsigned int excl_bit;
1910 const int jHalf = j/(c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit);
1911 nbnxn_excl_t *excl =
1912 get_exclusion_mask(nbl, cj4_ind, jHalf);
1914 excl_pair = a_mod_wj(j)*nbl->na_ci + i;
1915 excl_bit = (1U << (gcj*c_gpuNumClusterPerCell + c));
1917 dx = nbat->x()[ind_j*nbat->xstride+XX] - xi;
1918 dy = nbat->x()[ind_j*nbat->xstride+YY] - yi;
1919 dz = nbat->x()[ind_j*nbat->xstride+ZZ] - zi;
1921 /* The unpruned GPU list has more than 2/3
1922 * of the atom pairs beyond rlist. Using
1923 * this list will cause a lot of overhead
1924 * in the CPU FEP kernels, especially
1925 * relative to the fast GPU kernels.
1926 * So we prune the FEP list here.
1928 if (dx*dx + dy*dy + dz*dz < rlist_fep2)
1930 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1932 fep_list_new_nri_copy(nlist);
1936 /* Add it to the FEP list */
1937 nlist->jjnr[nlist->nrj] = aj;
1938 nlist->excl_fep[nlist->nrj] = (excl->pair[excl_pair] & excl_bit) ? 1 : 0;
1942 /* Exclude it from the normal list.
1943 * Note that the charge and LJ parameters have
1944 * been set to zero, but we need to avoid 0/0,
1945 * as perturbed atoms can be on top of each other.
1947 excl->pair[excl_pair] &= ~excl_bit;
1951 /* Note that we could mask out this pair in imask
1952 * if all i- and/or all j-particles are perturbed.
1953 * But since the perturbed pairs on the CPU will
1954 * take an order of magnitude more time, the GPU
1955 * will finish before the CPU and there is no gain.
1961 if (nlist->nrj > nlist->jindex[nri])
1963 /* Actually add this new, non-empty, list */
1965 nlist->jindex[nlist->nri] = nlist->nrj;
1972 /* Set all atom-pair exclusions for a GPU type list i-entry
1974 * Sets all atom-pair exclusions from the topology stored in exclusions
1975 * as masks in the pair-list for i-super-cluster list entry iEntry.
1978 setExclusionsForIEntry(const nbnxn_search *nbs,
1979 NbnxnPairlistGpu *nbl,
1980 gmx_bool diagRemoved,
1981 int gmx_unused na_cj_2log,
1982 const nbnxn_sci_t &iEntry,
1983 const t_blocka &exclusions)
1985 if (iEntry.numJClusterGroups() == 0)
1991 /* Set the search ranges using start and end j-cluster indices.
1992 * Note that here we can not use cj4_ind_end, since the last cj4
1993 * can be only partially filled, so we use cj_ind.
1995 const JListRanges ranges(iEntry.cj4_ind_start*c_nbnxnGpuJgroupSize,
1997 gmx::makeConstArrayRef(nbl->cj4));
1999 GMX_ASSERT(nbl->na_ci == c_nbnxnGpuClusterSize, "na_ci should match the GPU cluster size");
2000 constexpr int c_clusterSize = c_nbnxnGpuClusterSize;
2001 constexpr int c_superClusterSize = c_nbnxnGpuNumClusterPerSupercluster*c_nbnxnGpuClusterSize;
2003 const int iSuperCluster = iEntry.sci;
2005 gmx::ArrayRef<const int> cell = nbs->cell;
2007 /* Loop over the atoms in the i super-cluster */
2008 for (int i = 0; i < c_superClusterSize; i++)
2010 const int iIndex = iSuperCluster*c_superClusterSize + i;
2011 const int iAtom = nbs->a[iIndex];
2014 const int iCluster = i/c_clusterSize;
2016 /* Loop over the topology-based exclusions for this i-atom */
2017 for (int exclIndex = exclusions.index[iAtom]; exclIndex < exclusions.index[iAtom + 1]; exclIndex++)
2019 const int jAtom = exclusions.a[exclIndex];
2023 /* The self exclusions are already set, save some time */
2027 /* Get the index of the j-atom in the nbnxn atom data */
2028 const int jIndex = cell[jAtom];
2030 /* Without shifts we only calculate interactions j>i
2031 * for one-way pair-lists.
2033 /* NOTE: We would like to use iIndex on the right hand side,
2034 * but that makes this routine 25% slower with gcc6/7.
2035 * Even using c_superClusterSize makes it slower.
2036 * Either of these changes triggers peeling of the exclIndex
2037 * loop, which apparently leads to far less efficient code.
2039 if (diagRemoved && jIndex <= iSuperCluster*nbl->na_sc + i)
2044 const int jCluster = jIndex/c_clusterSize;
2046 /* Check whether the cluster is in our list? */
2047 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
2050 findJClusterInJList(jCluster, ranges,
2051 gmx::makeConstArrayRef(nbl->cj4));
2055 /* We found an exclusion, clear the corresponding
2058 const unsigned int pairMask = (1U << (cj_mod_cj4(index)*c_gpuNumClusterPerCell + iCluster));
2059 /* Check if the i-cluster interacts with the j-cluster */
2060 if (nbl_imask0(nbl, index) & pairMask)
2062 const int innerI = (i & (c_clusterSize - 1));
2063 const int innerJ = (jIndex & (c_clusterSize - 1));
2065 /* Determine which j-half (CUDA warp) we are in */
2066 const int jHalf = innerJ/(c_clusterSize/c_nbnxnGpuClusterpairSplit);
2068 nbnxn_excl_t *interactionMask =
2069 get_exclusion_mask(nbl, cj_to_cj4(index), jHalf);
2071 interactionMask->pair[a_mod_wj(innerJ)*c_clusterSize + innerI] &= ~pairMask;
2080 /* Make a new ci entry at the back of nbl->ci */
2081 static void addNewIEntry(NbnxnPairlistCpu *nbl, int ci, int shift, int flags)
2085 ciEntry.shift = shift;
2086 /* Store the interaction flags along with the shift */
2087 ciEntry.shift |= flags;
2088 ciEntry.cj_ind_start = nbl->cj.size();
2089 ciEntry.cj_ind_end = nbl->cj.size();
2090 nbl->ci.push_back(ciEntry);
2093 /* Make a new sci entry at index nbl->nsci */
2094 static void addNewIEntry(NbnxnPairlistGpu *nbl, int sci, int shift, int gmx_unused flags)
2096 nbnxn_sci_t sciEntry;
2098 sciEntry.shift = shift;
2099 sciEntry.cj4_ind_start = nbl->cj4.size();
2100 sciEntry.cj4_ind_end = nbl->cj4.size();
2102 nbl->sci.push_back(sciEntry);
2105 /* Sort the simple j-list cj on exclusions.
2106 * Entries with exclusions will all be sorted to the beginning of the list.
2108 static void sort_cj_excl(nbnxn_cj_t *cj, int ncj,
2109 NbnxnPairlistCpuWork *work)
2111 work->cj.resize(ncj);
2113 /* Make a list of the j-cells involving exclusions */
2115 for (int j = 0; j < ncj; j++)
2117 if (cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
2119 work->cj[jnew++] = cj[j];
2122 /* Check if there are exclusions at all or not just the first entry */
2123 if (!((jnew == 0) ||
2124 (jnew == 1 && cj[0].excl != NBNXN_INTERACTION_MASK_ALL)))
2126 for (int j = 0; j < ncj; j++)
2128 if (cj[j].excl == NBNXN_INTERACTION_MASK_ALL)
2130 work->cj[jnew++] = cj[j];
2133 for (int j = 0; j < ncj; j++)
2135 cj[j] = work->cj[j];
2140 /* Close this simple list i entry */
2141 static void closeIEntry(NbnxnPairlistCpu *nbl,
2142 int gmx_unused sp_max_av,
2143 gmx_bool gmx_unused progBal,
2144 float gmx_unused nsp_tot_est,
2145 int gmx_unused thread,
2146 int gmx_unused nthread)
2148 nbnxn_ci_t &ciEntry = nbl->ci.back();
2150 /* All content of the new ci entry have already been filled correctly,
2151 * we only need to sort and increase counts or remove the entry when empty.
2153 const int jlen = ciEntry.cj_ind_end - ciEntry.cj_ind_start;
2156 sort_cj_excl(nbl->cj.data() + ciEntry.cj_ind_start, jlen, nbl->work);
2158 /* The counts below are used for non-bonded pair/flop counts
2159 * and should therefore match the available kernel setups.
2161 if (!(ciEntry.shift & NBNXN_CI_DO_COUL(0)))
2163 nbl->work->ncj_noq += jlen;
2165 else if ((ciEntry.shift & NBNXN_CI_HALF_LJ(0)) ||
2166 !(ciEntry.shift & NBNXN_CI_DO_LJ(0)))
2168 nbl->work->ncj_hlj += jlen;
2173 /* Entry is empty: remove it */
2178 /* Split sci entry for load balancing on the GPU.
2179 * Splitting ensures we have enough lists to fully utilize the whole GPU.
2180 * With progBal we generate progressively smaller lists, which improves
2181 * load balancing. As we only know the current count on our own thread,
2182 * we will need to estimate the current total amount of i-entries.
2183 * As the lists get concatenated later, this estimate depends
2184 * both on nthread and our own thread index.
2186 static void split_sci_entry(NbnxnPairlistGpu *nbl,
2188 gmx_bool progBal, float nsp_tot_est,
2189 int thread, int nthread)
2197 /* Estimate the total numbers of ci's of the nblist combined
2198 * over all threads using the target number of ci's.
2200 nsp_est = (nsp_tot_est*thread)/nthread + nbl->nci_tot;
2202 /* The first ci blocks should be larger, to avoid overhead.
2203 * The last ci blocks should be smaller, to improve load balancing.
2204 * The factor 3/2 makes the first block 3/2 times the target average
2205 * and ensures that the total number of blocks end up equal to
2206 * that of equally sized blocks of size nsp_target_av.
2208 nsp_max = static_cast<int>(nsp_target_av*(nsp_tot_est*1.5/(nsp_est + nsp_tot_est)));
2212 nsp_max = nsp_target_av;
2215 const int cj4_start = nbl->sci.back().cj4_ind_start;
2216 const int cj4_end = nbl->sci.back().cj4_ind_end;
2217 const int j4len = cj4_end - cj4_start;
2219 if (j4len > 1 && j4len*c_gpuNumClusterPerCell*c_nbnxnGpuJgroupSize > nsp_max)
2221 /* Modify the last ci entry and process the cj4's again */
2227 for (int cj4 = cj4_start; cj4 < cj4_end; cj4++)
2229 int nsp_cj4_p = nsp_cj4;
2230 /* Count the number of cluster pairs in this cj4 group */
2232 for (int p = 0; p < c_gpuNumClusterPerCell*c_nbnxnGpuJgroupSize; p++)
2234 nsp_cj4 += (nbl->cj4[cj4].imei[0].imask >> p) & 1;
2237 /* If adding the current cj4 with nsp_cj4 pairs get us further
2238 * away from our target nsp_max, split the list before this cj4.
2240 if (nsp > 0 && nsp_max - nsp < nsp + nsp_cj4 - nsp_max)
2242 /* Split the list at cj4 */
2243 nbl->sci.back().cj4_ind_end = cj4;
2244 /* Create a new sci entry */
2246 sciNew.sci = nbl->sci.back().sci;
2247 sciNew.shift = nbl->sci.back().shift;
2248 sciNew.cj4_ind_start = cj4;
2249 nbl->sci.push_back(sciNew);
2252 nsp_cj4_e = nsp_cj4_p;
2258 /* Put the remaining cj4's in the last sci entry */
2259 nbl->sci.back().cj4_ind_end = cj4_end;
2261 /* Possibly balance out the last two sci's
2262 * by moving the last cj4 of the second last sci.
2264 if (nsp_sci - nsp_cj4_e >= nsp + nsp_cj4_e)
2266 GMX_ASSERT(nbl->sci.size() >= 2, "We expect at least two elements");
2267 nbl->sci[nbl->sci.size() - 2].cj4_ind_end--;
2268 nbl->sci[nbl->sci.size() - 1].cj4_ind_start--;
2273 /* Clost this super/sub list i entry */
2274 static void closeIEntry(NbnxnPairlistGpu *nbl,
2276 gmx_bool progBal, float nsp_tot_est,
2277 int thread, int nthread)
2279 nbnxn_sci_t &sciEntry = *getOpenIEntry(nbl);
2281 /* All content of the new ci entry have already been filled correctly,
2282 * we only need to, potentially, split or remove the entry when empty.
2284 int j4len = sciEntry.numJClusterGroups();
2287 /* We can only have complete blocks of 4 j-entries in a list,
2288 * so round the count up before closing.
2290 int ncj4 = (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize;
2291 nbl->work->cj_ind = ncj4*c_nbnxnGpuJgroupSize;
2295 /* Measure the size of the new entry and potentially split it */
2296 split_sci_entry(nbl, nsp_max_av, progBal, nsp_tot_est,
2302 /* Entry is empty: remove it */
2303 nbl->sci.pop_back();
2307 /* Syncs the working array before adding another grid pair to the GPU list */
2308 static void sync_work(NbnxnPairlistCpu gmx_unused *nbl)
2312 /* Syncs the working array before adding another grid pair to the GPU list */
2313 static void sync_work(NbnxnPairlistGpu *nbl)
2315 nbl->work->cj_ind = nbl->cj4.size()*c_nbnxnGpuJgroupSize;
2318 /* Clears an NbnxnPairlistCpu data structure */
2319 static void clear_pairlist(NbnxnPairlistCpu *nbl)
2325 nbl->ciOuter.clear();
2326 nbl->cjOuter.clear();
2328 nbl->work->ncj_noq = 0;
2329 nbl->work->ncj_hlj = 0;
2332 /* Clears an NbnxnPairlistGpu data structure */
2333 static void clear_pairlist(NbnxnPairlistGpu *nbl)
2337 nbl->excl.resize(1);
2341 /* Clears a group scheme pair list */
2342 static void clear_pairlist_fep(t_nblist *nl)
2346 if (nl->jindex == nullptr)
2348 snew(nl->jindex, 1);
2353 /* Sets a simple list i-cell bounding box, including PBC shift */
2354 static inline void set_icell_bb_simple(gmx::ArrayRef<const nbnxn_bb_t> bb,
2356 real shx, real shy, real shz,
2359 bb_ci->lower.x = bb[ci].lower.x + shx;
2360 bb_ci->lower.y = bb[ci].lower.y + shy;
2361 bb_ci->lower.z = bb[ci].lower.z + shz;
2362 bb_ci->upper.x = bb[ci].upper.x + shx;
2363 bb_ci->upper.y = bb[ci].upper.y + shy;
2364 bb_ci->upper.z = bb[ci].upper.z + shz;
2367 /* Sets a simple list i-cell bounding box, including PBC shift */
2368 static inline void set_icell_bb(const Grid &iGrid,
2370 real shx, real shy, real shz,
2371 NbnxnPairlistCpuWork *work)
2373 set_icell_bb_simple(iGrid.iBoundingBoxes(), ci, shx, shy, shz,
2374 &work->iClusterData.bb[0]);
2378 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2379 static void set_icell_bbxxxx_supersub(gmx::ArrayRef<const float> bb,
2381 real shx, real shy, real shz,
2384 int ia = ci*(c_gpuNumClusterPerCell >> STRIDE_PBB_2LOG)*NNBSBB_XXXX;
2385 for (int m = 0; m < (c_gpuNumClusterPerCell >> STRIDE_PBB_2LOG)*NNBSBB_XXXX; m += NNBSBB_XXXX)
2387 for (int i = 0; i < STRIDE_PBB; i++)
2389 bb_ci[m+0*STRIDE_PBB+i] = bb[ia+m+0*STRIDE_PBB+i] + shx;
2390 bb_ci[m+1*STRIDE_PBB+i] = bb[ia+m+1*STRIDE_PBB+i] + shy;
2391 bb_ci[m+2*STRIDE_PBB+i] = bb[ia+m+2*STRIDE_PBB+i] + shz;
2392 bb_ci[m+3*STRIDE_PBB+i] = bb[ia+m+3*STRIDE_PBB+i] + shx;
2393 bb_ci[m+4*STRIDE_PBB+i] = bb[ia+m+4*STRIDE_PBB+i] + shy;
2394 bb_ci[m+5*STRIDE_PBB+i] = bb[ia+m+5*STRIDE_PBB+i] + shz;
2400 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2401 gmx_unused static void set_icell_bb_supersub(gmx::ArrayRef<const nbnxn_bb_t> bb,
2403 real shx, real shy, real shz,
2406 for (int i = 0; i < c_gpuNumClusterPerCell; i++)
2408 set_icell_bb_simple(bb, ci*c_gpuNumClusterPerCell+i,
2414 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2415 gmx_unused static void set_icell_bb(const Grid &iGrid,
2417 real shx, real shy, real shz,
2418 NbnxnPairlistGpuWork *work)
2421 set_icell_bbxxxx_supersub(iGrid.packedBoundingBoxes(), ci, shx, shy, shz,
2422 work->iSuperClusterData.bbPacked.data());
2424 set_icell_bb_supersub(iGrid.iBoundingBoxes(), ci, shx, shy, shz,
2425 work->iSuperClusterData.bb.data());
2429 /* Copies PBC shifted i-cell atom coordinates x,y,z to working array */
2430 static void icell_set_x_simple(int ci,
2431 real shx, real shy, real shz,
2432 int stride, const real *x,
2433 NbnxnPairlistCpuWork::IClusterData *iClusterData)
2435 const int ia = ci*c_nbnxnCpuIClusterSize;
2437 for (int i = 0; i < c_nbnxnCpuIClusterSize; i++)
2439 iClusterData->x[i*STRIDE_XYZ+XX] = x[(ia+i)*stride+XX] + shx;
2440 iClusterData->x[i*STRIDE_XYZ+YY] = x[(ia+i)*stride+YY] + shy;
2441 iClusterData->x[i*STRIDE_XYZ+ZZ] = x[(ia+i)*stride+ZZ] + shz;
2445 static void icell_set_x(int ci,
2446 real shx, real shy, real shz,
2447 int stride, const real *x,
2448 const Nbnxm::KernelType kernelType,
2449 NbnxnPairlistCpuWork *work)
2454 #ifdef GMX_NBNXN_SIMD_4XN
2455 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
2456 icell_set_x_simd_4xn(ci, shx, shy, shz, stride, x, work);
2459 #ifdef GMX_NBNXN_SIMD_2XNN
2460 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
2461 icell_set_x_simd_2xnn(ci, shx, shy, shz, stride, x, work);
2465 case Nbnxm::KernelType::Cpu4x4_PlainC:
2466 icell_set_x_simple(ci, shx, shy, shz, stride, x, &work->iClusterData);
2469 GMX_ASSERT(false, "Unhandled case");
2474 /* Copies PBC shifted super-cell atom coordinates x,y,z to working array */
2475 static void icell_set_x(int ci,
2476 real shx, real shy, real shz,
2477 int stride, const real *x,
2478 Nbnxm::KernelType gmx_unused kernelType,
2479 NbnxnPairlistGpuWork *work)
2481 #if !GMX_SIMD4_HAVE_REAL
2483 real * x_ci = work->iSuperClusterData.x.data();
2485 int ia = ci*c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize;
2486 for (int i = 0; i < c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize; i++)
2488 x_ci[i*DIM + XX] = x[(ia+i)*stride + XX] + shx;
2489 x_ci[i*DIM + YY] = x[(ia+i)*stride + YY] + shy;
2490 x_ci[i*DIM + ZZ] = x[(ia+i)*stride + ZZ] + shz;
2493 #else /* !GMX_SIMD4_HAVE_REAL */
2495 real * x_ci = work->iSuperClusterData.xSimd.data();
2497 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2499 for (int i = 0; i < c_nbnxnGpuClusterSize; i += GMX_SIMD4_WIDTH)
2501 int io = si*c_nbnxnGpuClusterSize + i;
2502 int ia = ci*c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize + io;
2503 for (int j = 0; j < GMX_SIMD4_WIDTH; j++)
2505 x_ci[io*DIM + j + XX*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + XX] + shx;
2506 x_ci[io*DIM + j + YY*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + YY] + shy;
2507 x_ci[io*DIM + j + ZZ*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + ZZ] + shz;
2512 #endif /* !GMX_SIMD4_HAVE_REAL */
2515 static real minimum_subgrid_size_xy(const Grid &grid)
2517 const Grid::Dimensions &dims = grid.dimensions();
2519 if (grid.geometry().isSimple)
2521 return std::min(dims.cellSize[XX], dims.cellSize[YY]);
2525 return std::min(dims.cellSize[XX]/c_gpuNumClusterPerCellX,
2526 dims.cellSize[YY]/c_gpuNumClusterPerCellY);
2530 static real effective_buffer_1x1_vs_MxN(const Grid &iGrid,
2533 const real eff_1x1_buffer_fac_overest = 0.1;
2535 /* Determine an atom-pair list cut-off buffer size for atom pairs,
2536 * to be added to rlist (including buffer) used for MxN.
2537 * This is for converting an MxN list to a 1x1 list. This means we can't
2538 * use the normal buffer estimate, as we have an MxN list in which
2539 * some atom pairs beyond rlist are missing. We want to capture
2540 * the beneficial effect of buffering by extra pairs just outside rlist,
2541 * while removing the useless pairs that are further away from rlist.
2542 * (Also the buffer could have been set manually not using the estimate.)
2543 * This buffer size is an overestimate.
2544 * We add 10% of the smallest grid sub-cell dimensions.
2545 * Note that the z-size differs per cell and we don't use this,
2546 * so we overestimate.
2547 * With PME, the 10% value gives a buffer that is somewhat larger
2548 * than the effective buffer with a tolerance of 0.005 kJ/mol/ps.
2549 * Smaller tolerances or using RF lead to a smaller effective buffer,
2550 * so 10% gives a safe overestimate.
2552 return eff_1x1_buffer_fac_overest*(minimum_subgrid_size_xy(iGrid) +
2553 minimum_subgrid_size_xy(jGrid));
2556 /* Estimates the interaction volume^2 for non-local interactions */
2557 static real nonlocal_vol2(const struct gmx_domdec_zones_t *zones, const rvec ls, real r)
2565 /* Here we simply add up the volumes of 1, 2 or 3 1D decomposition
2566 * not home interaction volume^2. As these volumes are not additive,
2567 * this is an overestimate, but it would only be significant in the limit
2568 * of small cells, where we anyhow need to split the lists into
2569 * as small parts as possible.
2572 for (int z = 0; z < zones->n; z++)
2574 if (zones->shift[z][XX] + zones->shift[z][YY] + zones->shift[z][ZZ] == 1)
2579 for (int d = 0; d < DIM; d++)
2581 if (zones->shift[z][d] == 0)
2585 za *= zones->size[z].x1[d] - zones->size[z].x0[d];
2589 /* 4 octants of a sphere */
2590 vold_est = 0.25*M_PI*r*r*r*r;
2591 /* 4 quarter pie slices on the edges */
2592 vold_est += 4*cl*M_PI/6.0*r*r*r;
2593 /* One rectangular volume on a face */
2594 vold_est += ca*0.5*r*r;
2596 vol2_est_tot += vold_est*za;
2600 return vol2_est_tot;
2603 /* Estimates the average size of a full j-list for super/sub setup */
2604 static void get_nsubpair_target(const nbnxn_search *nbs,
2605 const InteractionLocality iloc,
2607 const int min_ci_balanced,
2608 int *nsubpair_target,
2609 float *nsubpair_tot_est)
2611 /* The target value of 36 seems to be the optimum for Kepler.
2612 * Maxwell is less sensitive to the exact value.
2614 const int nsubpair_target_min = 36;
2615 real r_eff_sup, vol_est, nsp_est, nsp_est_nl;
2617 const Grid &grid = nbs->grid[0];
2619 /* We don't need to balance list sizes if:
2620 * - We didn't request balancing.
2621 * - The number of grid cells >= the number of lists requested,
2622 * since we will always generate at least #cells lists.
2623 * - We don't have any cells, since then there won't be any lists.
2625 if (min_ci_balanced <= 0 || grid.numCells() >= min_ci_balanced || grid.numCells() == 0)
2627 /* nsubpair_target==0 signals no balancing */
2628 *nsubpair_target = 0;
2629 *nsubpair_tot_est = 0;
2635 const int numAtomsCluster = grid.geometry().numAtomsICluster;
2636 const Grid::Dimensions &dims = grid.dimensions();
2638 ls[XX] = dims.cellSize[XX]/c_gpuNumClusterPerCellX;
2639 ls[YY] = dims.cellSize[YY]/c_gpuNumClusterPerCellY;
2640 ls[ZZ] = numAtomsCluster/(dims.atomDensity*ls[XX]*ls[YY]);
2642 /* The formulas below are a heuristic estimate of the average nsj per si*/
2643 r_eff_sup = rlist + nbnxn_get_rlist_effective_inc(numAtomsCluster, ls);
2645 if (!nbs->DomDec || nbs->zones->n == 1)
2652 gmx::square(dims.atomDensity/numAtomsCluster)*
2653 nonlocal_vol2(nbs->zones, ls, r_eff_sup);
2656 if (iloc == InteractionLocality::Local)
2658 /* Sub-cell interacts with itself */
2659 vol_est = ls[XX]*ls[YY]*ls[ZZ];
2660 /* 6/2 rectangular volume on the faces */
2661 vol_est += (ls[XX]*ls[YY] + ls[XX]*ls[ZZ] + ls[YY]*ls[ZZ])*r_eff_sup;
2662 /* 12/2 quarter pie slices on the edges */
2663 vol_est += 2*(ls[XX] + ls[YY] + ls[ZZ])*0.25*M_PI*gmx::square(r_eff_sup);
2664 /* 4 octants of a sphere */
2665 vol_est += 0.5*4.0/3.0*M_PI*gmx::power3(r_eff_sup);
2667 /* Estimate the number of cluster pairs as the local number of
2668 * clusters times the volume they interact with times the density.
2670 nsp_est = grid.numClusters()*vol_est*dims.atomDensity/numAtomsCluster;
2672 /* Subtract the non-local pair count */
2673 nsp_est -= nsp_est_nl;
2675 /* For small cut-offs nsp_est will be an underesimate.
2676 * With DD nsp_est_nl is an overestimate so nsp_est can get negative.
2677 * So to avoid too small or negative nsp_est we set a minimum of
2678 * all cells interacting with all 3^3 direct neighbors (3^3-1)/2+1=14.
2679 * This might be a slight overestimate for small non-periodic groups of
2680 * atoms as will occur for a local domain with DD, but for small
2681 * groups of atoms we'll anyhow be limited by nsubpair_target_min,
2682 * so this overestimation will not matter.
2684 nsp_est = std::max(nsp_est, grid.numClusters()*14._real);
2688 fprintf(debug, "nsp_est local %5.1f non-local %5.1f\n",
2689 nsp_est, nsp_est_nl);
2694 nsp_est = nsp_est_nl;
2697 /* Thus the (average) maximum j-list size should be as follows.
2698 * Since there is overhead, we shouldn't make the lists too small
2699 * (and we can't chop up j-groups) so we use a minimum target size of 36.
2701 *nsubpair_target = std::max(nsubpair_target_min,
2702 roundToInt(nsp_est/min_ci_balanced));
2703 *nsubpair_tot_est = static_cast<int>(nsp_est);
2707 fprintf(debug, "nbl nsp estimate %.1f, nsubpair_target %d\n",
2708 nsp_est, *nsubpair_target);
2712 /* Debug list print function */
2713 static void print_nblist_ci_cj(FILE *fp, const NbnxnPairlistCpu *nbl)
2715 for (const nbnxn_ci_t &ciEntry : nbl->ci)
2717 fprintf(fp, "ci %4d shift %2d ncj %3d\n",
2718 ciEntry.ci, ciEntry.shift,
2719 ciEntry.cj_ind_end - ciEntry.cj_ind_start);
2721 for (int j = ciEntry.cj_ind_start; j < ciEntry.cj_ind_end; j++)
2723 fprintf(fp, " cj %5d imask %x\n",
2730 /* Debug list print function */
2731 static void print_nblist_sci_cj(FILE *fp, const NbnxnPairlistGpu *nbl)
2733 for (const nbnxn_sci_t &sci : nbl->sci)
2735 fprintf(fp, "ci %4d shift %2d ncj4 %2d\n",
2737 sci.numJClusterGroups());
2740 for (int j4 = sci.cj4_ind_start; j4 < sci.cj4_ind_end; j4++)
2742 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
2744 fprintf(fp, " sj %5d imask %x\n",
2746 nbl->cj4[j4].imei[0].imask);
2747 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2749 if (nbl->cj4[j4].imei[0].imask & (1U << (j*c_gpuNumClusterPerCell + si)))
2756 fprintf(fp, "ci %4d shift %2d ncj4 %2d ncp %3d\n",
2758 sci.numJClusterGroups(),
2763 /* Combine pair lists *nbl generated on multiple threads nblc */
2764 static void combine_nblists(int nnbl, NbnxnPairlistGpu **nbl,
2765 NbnxnPairlistGpu *nblc)
2767 int nsci = nblc->sci.size();
2768 int ncj4 = nblc->cj4.size();
2769 int nexcl = nblc->excl.size();
2770 for (int i = 0; i < nnbl; i++)
2772 nsci += nbl[i]->sci.size();
2773 ncj4 += nbl[i]->cj4.size();
2774 nexcl += nbl[i]->excl.size();
2777 /* Resize with the final, combined size, so we can fill in parallel */
2778 /* NOTE: For better performance we should use default initialization */
2779 nblc->sci.resize(nsci);
2780 nblc->cj4.resize(ncj4);
2781 nblc->excl.resize(nexcl);
2783 /* Each thread should copy its own data to the combined arrays,
2784 * as otherwise data will go back and forth between different caches.
2786 #if GMX_OPENMP && !(defined __clang_analyzer__)
2787 int nthreads = gmx_omp_nthreads_get(emntPairsearch);
2790 #pragma omp parallel for num_threads(nthreads) schedule(static)
2791 for (int n = 0; n < nnbl; n++)
2795 /* Determine the offset in the combined data for our thread.
2796 * Note that the original sizes in nblc are lost.
2798 int sci_offset = nsci;
2799 int cj4_offset = ncj4;
2800 int excl_offset = nexcl;
2802 for (int i = n; i < nnbl; i++)
2804 sci_offset -= nbl[i]->sci.size();
2805 cj4_offset -= nbl[i]->cj4.size();
2806 excl_offset -= nbl[i]->excl.size();
2809 const NbnxnPairlistGpu &nbli = *nbl[n];
2811 for (size_t i = 0; i < nbli.sci.size(); i++)
2813 nblc->sci[sci_offset + i] = nbli.sci[i];
2814 nblc->sci[sci_offset + i].cj4_ind_start += cj4_offset;
2815 nblc->sci[sci_offset + i].cj4_ind_end += cj4_offset;
2818 for (size_t j4 = 0; j4 < nbli.cj4.size(); j4++)
2820 nblc->cj4[cj4_offset + j4] = nbli.cj4[j4];
2821 nblc->cj4[cj4_offset + j4].imei[0].excl_ind += excl_offset;
2822 nblc->cj4[cj4_offset + j4].imei[1].excl_ind += excl_offset;
2825 for (size_t j4 = 0; j4 < nbli.excl.size(); j4++)
2827 nblc->excl[excl_offset + j4] = nbli.excl[j4];
2830 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2833 for (int n = 0; n < nnbl; n++)
2835 nblc->nci_tot += nbl[n]->nci_tot;
2839 static void balance_fep_lists(const nbnxn_search *nbs,
2840 nbnxn_pairlist_set_t *nbl_lists)
2843 int nri_tot, nrj_tot, nrj_target;
2847 nnbl = nbl_lists->nnbl;
2851 /* Nothing to balance */
2855 /* Count the total i-lists and pairs */
2858 for (int th = 0; th < nnbl; th++)
2860 nri_tot += nbl_lists->nbl_fep[th]->nri;
2861 nrj_tot += nbl_lists->nbl_fep[th]->nrj;
2864 nrj_target = (nrj_tot + nnbl - 1)/nnbl;
2866 assert(gmx_omp_nthreads_get(emntNonbonded) == nnbl);
2868 #pragma omp parallel for schedule(static) num_threads(nnbl)
2869 for (int th = 0; th < nnbl; th++)
2873 t_nblist *nbl = nbs->work[th].nbl_fep.get();
2875 /* Note that here we allocate for the total size, instead of
2876 * a per-thread esimate (which is hard to obtain).
2878 if (nri_tot > nbl->maxnri)
2880 nbl->maxnri = over_alloc_large(nri_tot);
2881 reallocate_nblist(nbl);
2883 if (nri_tot > nbl->maxnri || nrj_tot > nbl->maxnrj)
2885 nbl->maxnrj = over_alloc_small(nrj_tot);
2886 srenew(nbl->jjnr, nbl->maxnrj);
2887 srenew(nbl->excl_fep, nbl->maxnrj);
2890 clear_pairlist_fep(nbl);
2892 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2895 /* Loop over the source lists and assign and copy i-entries */
2897 nbld = nbs->work[th_dest].nbl_fep.get();
2898 for (int th = 0; th < nnbl; th++)
2902 nbls = nbl_lists->nbl_fep[th];
2904 for (int i = 0; i < nbls->nri; i++)
2908 /* The number of pairs in this i-entry */
2909 nrj = nbls->jindex[i+1] - nbls->jindex[i];
2911 /* Decide if list th_dest is too large and we should procede
2912 * to the next destination list.
2914 if (th_dest+1 < nnbl && nbld->nrj > 0 &&
2915 nbld->nrj + nrj - nrj_target > nrj_target - nbld->nrj)
2918 nbld = nbs->work[th_dest].nbl_fep.get();
2921 nbld->iinr[nbld->nri] = nbls->iinr[i];
2922 nbld->gid[nbld->nri] = nbls->gid[i];
2923 nbld->shift[nbld->nri] = nbls->shift[i];
2925 for (int j = nbls->jindex[i]; j < nbls->jindex[i+1]; j++)
2927 nbld->jjnr[nbld->nrj] = nbls->jjnr[j];
2928 nbld->excl_fep[nbld->nrj] = nbls->excl_fep[j];
2932 nbld->jindex[nbld->nri] = nbld->nrj;
2936 /* Swap the list pointers */
2937 for (int th = 0; th < nnbl; th++)
2939 t_nblist *nbl_tmp = nbs->work[th].nbl_fep.release();
2940 nbs->work[th].nbl_fep.reset(nbl_lists->nbl_fep[th]);
2941 nbl_lists->nbl_fep[th] = nbl_tmp;
2945 fprintf(debug, "nbl_fep[%d] nri %4d nrj %4d\n",
2947 nbl_lists->nbl_fep[th]->nri,
2948 nbl_lists->nbl_fep[th]->nrj);
2953 /* Returns the next ci to be processes by our thread */
2954 static gmx_bool next_ci(const Grid &grid,
2955 int nth, int ci_block,
2956 int *ci_x, int *ci_y,
2962 if (*ci_b == ci_block)
2964 /* Jump to the next block assigned to this task */
2965 *ci += (nth - 1)*ci_block;
2969 if (*ci >= grid.numCells())
2974 while (*ci >= grid.firstCellInColumn(*ci_x*grid.dimensions().numCells[YY] + *ci_y + 1))
2977 if (*ci_y == grid.dimensions().numCells[YY])
2987 /* Returns the distance^2 for which we put cell pairs in the list
2988 * without checking atom pair distances. This is usually < rlist^2.
2990 static float boundingbox_only_distance2(const Grid::Dimensions &iGridDims,
2991 const Grid::Dimensions &jGridDims,
2995 /* If the distance between two sub-cell bounding boxes is less
2996 * than this distance, do not check the distance between
2997 * all particle pairs in the sub-cell, since then it is likely
2998 * that the box pair has atom pairs within the cut-off.
2999 * We use the nblist cut-off minus 0.5 times the average x/y diagonal
3000 * spacing of the sub-cells. Around 40% of the checked pairs are pruned.
3001 * Using more than 0.5 gains at most 0.5%.
3002 * If forces are calculated more than twice, the performance gain
3003 * in the force calculation outweighs the cost of checking.
3004 * Note that with subcell lists, the atom-pair distance check
3005 * is only performed when only 1 out of 8 sub-cells in within range,
3006 * this is because the GPU is much faster than the cpu.
3011 bbx = 0.5*(iGridDims.cellSize[XX] + jGridDims.cellSize[XX]);
3012 bby = 0.5*(iGridDims.cellSize[YY] + jGridDims.cellSize[YY]);
3015 bbx /= c_gpuNumClusterPerCellX;
3016 bby /= c_gpuNumClusterPerCellY;
3019 rbb2 = std::max(0.0, rlist - 0.5*std::sqrt(bbx*bbx + bby*bby));
3025 return (float)((1+GMX_FLOAT_EPS)*rbb2);
3029 static int get_ci_block_size(const Grid &iGrid,
3030 gmx_bool bDomDec, int nth)
3032 const int ci_block_enum = 5;
3033 const int ci_block_denom = 11;
3034 const int ci_block_min_atoms = 16;
3037 /* Here we decide how to distribute the blocks over the threads.
3038 * We use prime numbers to try to avoid that the grid size becomes
3039 * a multiple of the number of threads, which would lead to some
3040 * threads getting "inner" pairs and others getting boundary pairs,
3041 * which in turns will lead to load imbalance between threads.
3042 * Set the block size as 5/11/ntask times the average number of cells
3043 * in a y,z slab. This should ensure a quite uniform distribution
3044 * of the grid parts of the different thread along all three grid
3045 * zone boundaries with 3D domain decomposition. At the same time
3046 * the blocks will not become too small.
3048 ci_block = (iGrid.numCells()*ci_block_enum)/(ci_block_denom*iGrid.dimensions().numCells[XX]*nth);
3050 const int numAtomsPerCell = iGrid.geometry().numAtomsPerCell;
3052 /* Ensure the blocks are not too small: avoids cache invalidation */
3053 if (ci_block*numAtomsPerCell < ci_block_min_atoms)
3055 ci_block = (ci_block_min_atoms + numAtomsPerCell - 1)/numAtomsPerCell;
3058 /* Without domain decomposition
3059 * or with less than 3 blocks per task, divide in nth blocks.
3061 if (!bDomDec || nth*3*ci_block > iGrid.numCells())
3063 ci_block = (iGrid.numCells() + nth - 1)/nth;
3066 if (ci_block > 1 && (nth - 1)*ci_block >= iGrid.numCells())
3068 /* Some threads have no work. Although reducing the block size
3069 * does not decrease the block count on the first few threads,
3070 * with GPUs better mixing of "upper" cells that have more empty
3071 * clusters results in a somewhat lower max load over all threads.
3072 * Without GPUs the regime of so few atoms per thread is less
3073 * performance relevant, but with 8-wide SIMD the same reasoning
3074 * applies, since the pair list uses 4 i-atom "sub-clusters".
3082 /* Returns the number of bits to right-shift a cluster index to obtain
3083 * the corresponding force buffer flag index.
3085 static int getBufferFlagShift(int numAtomsPerCluster)
3087 int bufferFlagShift = 0;
3088 while ((numAtomsPerCluster << bufferFlagShift) < NBNXN_BUFFERFLAG_SIZE)
3093 return bufferFlagShift;
3096 static bool pairlistIsSimple(const NbnxnPairlistCpu gmx_unused &pairlist)
3101 static bool pairlistIsSimple(const NbnxnPairlistGpu gmx_unused &pairlist)
3106 static void makeClusterListWrapper(NbnxnPairlistCpu *nbl,
3107 const Grid gmx_unused &iGrid,
3110 const int firstCell,
3112 const bool excludeSubDiagonal,
3113 const nbnxn_atomdata_t *nbat,
3116 const Nbnxm::KernelType kernelType,
3117 int *numDistanceChecks)
3121 case Nbnxm::KernelType::Cpu4x4_PlainC:
3122 makeClusterListSimple(jGrid,
3123 nbl, ci, firstCell, lastCell,
3129 #ifdef GMX_NBNXN_SIMD_4XN
3130 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
3131 makeClusterListSimd4xn(jGrid,
3132 nbl, ci, firstCell, lastCell,
3139 #ifdef GMX_NBNXN_SIMD_2XNN
3140 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
3141 makeClusterListSimd2xnn(jGrid,
3142 nbl, ci, firstCell, lastCell,
3150 GMX_ASSERT(false, "Unhandled kernel type");
3154 static void makeClusterListWrapper(NbnxnPairlistGpu *nbl,
3155 const Grid &gmx_unused iGrid,
3158 const int firstCell,
3160 const bool excludeSubDiagonal,
3161 const nbnxn_atomdata_t *nbat,
3164 Nbnxm::KernelType gmx_unused kernelType,
3165 int *numDistanceChecks)
3167 for (int cj = firstCell; cj <= lastCell; cj++)
3169 make_cluster_list_supersub(iGrid, jGrid,
3172 nbat->xstride, nbat->x().data(),
3178 static int getNumSimpleJClustersInList(const NbnxnPairlistCpu &nbl)
3180 return nbl.cj.size();
3183 static int getNumSimpleJClustersInList(const gmx_unused NbnxnPairlistGpu &nbl)
3188 static void incrementNumSimpleJClustersInList(NbnxnPairlistCpu *nbl,
3191 nbl->ncjInUse += nbl->cj.size() - ncj_old_j;
3194 static void incrementNumSimpleJClustersInList(NbnxnPairlistGpu gmx_unused *nbl,
3195 int gmx_unused ncj_old_j)
3199 static void checkListSizeConsistency(const NbnxnPairlistCpu &nbl,
3200 const bool haveFreeEnergy)
3202 GMX_RELEASE_ASSERT(static_cast<size_t>(nbl.ncjInUse) == nbl.cj.size() || haveFreeEnergy,
3203 "Without free-energy all cj pair-list entries should be in use. "
3204 "Note that subsequent code does not make use of the equality, "
3205 "this check is only here to catch bugs");
3208 static void checkListSizeConsistency(const NbnxnPairlistGpu gmx_unused &nbl,
3209 bool gmx_unused haveFreeEnergy)
3211 /* We currently can not check consistency here */
3214 /* Set the buffer flags for newly added entries in the list */
3215 static void setBufferFlags(const NbnxnPairlistCpu &nbl,
3216 const int ncj_old_j,
3217 const int gridj_flag_shift,
3218 gmx_bitmask_t *gridj_flag,
3221 if (gmx::ssize(nbl.cj) > ncj_old_j)
3223 int cbFirst = nbl.cj[ncj_old_j].cj >> gridj_flag_shift;
3224 int cbLast = nbl.cj.back().cj >> gridj_flag_shift;
3225 for (int cb = cbFirst; cb <= cbLast; cb++)
3227 bitmask_init_bit(&gridj_flag[cb], th);
3232 static void setBufferFlags(const NbnxnPairlistGpu gmx_unused &nbl,
3233 int gmx_unused ncj_old_j,
3234 int gmx_unused gridj_flag_shift,
3235 gmx_bitmask_t gmx_unused *gridj_flag,
3238 GMX_ASSERT(false, "This function should never be called");
3241 /* Generates the part of pair-list nbl assigned to our thread */
3242 template <typename T>
3243 static void nbnxn_make_pairlist_part(const nbnxn_search *nbs,
3246 nbnxn_search_work_t *work,
3247 const nbnxn_atomdata_t *nbat,
3248 const t_blocka &exclusions,
3250 const Nbnxm::KernelType kernelType,
3252 gmx_bool bFBufferFlag,
3255 float nsubpair_tot_est,
3262 real rlist2, rl_fep2 = 0;
3264 int ci_b, ci, ci_x, ci_y, ci_xy;
3266 real bx0, bx1, by0, by1, bz0, bz1;
3268 real d2cx, d2z, d2z_cx, d2z_cy, d2zx, d2zxy, d2xy;
3269 int cxf, cxl, cyf, cyf_x, cyl;
3270 int numDistanceChecks;
3271 int gridi_flag_shift = 0, gridj_flag_shift = 0;
3272 gmx_bitmask_t *gridj_flag = nullptr;
3273 int ncj_old_i, ncj_old_j;
3275 nbs_cycle_start(&work->cc[enbsCCsearch]);
3277 if (jGrid.geometry().isSimple != pairlistIsSimple(*nbl) ||
3278 iGrid.geometry().isSimple != pairlistIsSimple(*nbl))
3280 gmx_incons("Grid incompatible with pair-list");
3284 GMX_ASSERT(nbl->na_ci == jGrid.geometry().numAtomsICluster,
3285 "The cluster sizes in the list and grid should match");
3286 nbl->na_cj = Nbnxm::JClusterSizePerKernelType[kernelType];
3287 na_cj_2log = get_2log(nbl->na_cj);
3293 /* Determine conversion of clusters to flag blocks */
3294 gridi_flag_shift = getBufferFlagShift(nbl->na_ci);
3295 gridj_flag_shift = getBufferFlagShift(nbl->na_cj);
3297 gridj_flag = work->buffer_flags.flag;
3300 copy_mat(nbs->box, box);
3302 rlist2 = nbl->rlist*nbl->rlist;
3304 if (nbs->bFEP && !pairlistIsSimple(*nbl))
3306 /* Determine an atom-pair list cut-off distance for FEP atom pairs.
3307 * We should not simply use rlist, since then we would not have
3308 * the small, effective buffering of the NxN lists.
3309 * The buffer is on overestimate, but the resulting cost for pairs
3310 * beyond rlist is neglible compared to the FEP pairs within rlist.
3312 rl_fep2 = nbl->rlist + effective_buffer_1x1_vs_MxN(iGrid, jGrid);
3316 fprintf(debug, "nbl_fep atom-pair rlist %f\n", rl_fep2);
3318 rl_fep2 = rl_fep2*rl_fep2;
3321 const Grid::Dimensions &iGridDims = iGrid.dimensions();
3322 const Grid::Dimensions &jGridDims = jGrid.dimensions();
3324 rbb2 = boundingbox_only_distance2(iGridDims, jGridDims, nbl->rlist, pairlistIsSimple(*nbl));
3328 fprintf(debug, "nbl bounding box only distance %f\n", std::sqrt(rbb2));
3331 const bool isIntraGridList = (&iGrid == &jGrid);
3333 /* Set the shift range */
3334 for (int d = 0; d < DIM; d++)
3336 /* Check if we need periodicity shifts.
3337 * Without PBC or with domain decomposition we don't need them.
3339 if (d >= ePBC2npbcdim(nbs->ePBC) || nbs->dd_dim[d])
3345 const real listRangeCellToCell =
3346 listRangeForGridCellToGridCell(rlist, iGrid.dimensions(), jGrid.dimensions());
3348 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < listRangeCellToCell)
3358 const bool bSimple = pairlistIsSimple(*nbl);
3359 gmx::ArrayRef<const nbnxn_bb_t> bb_i;
3361 gmx::ArrayRef<const float> pbb_i;
3364 bb_i = iGrid.iBoundingBoxes();
3368 pbb_i = iGrid.packedBoundingBoxes();
3371 /* We use the normal bounding box format for both grid types */
3372 bb_i = iGrid.iBoundingBoxes();
3374 gmx::ArrayRef<const float> bbcz_i = iGrid.zBoundingBoxes();
3375 gmx::ArrayRef<const int> flags_i = iGrid.clusterFlags();
3376 gmx::ArrayRef<const float> bbcz_j = jGrid.zBoundingBoxes();
3377 int cell0_i = iGrid.cellOffset();
3381 fprintf(debug, "nbl nc_i %d col.av. %.1f ci_block %d\n",
3382 iGrid.numCells(), iGrid.numCells()/static_cast<double>(iGrid.numColumns()), ci_block);
3385 numDistanceChecks = 0;
3387 const real listRangeBBToJCell2 = gmx::square(listRangeForBoundingBoxToGridCell(rlist, jGrid.dimensions()));
3389 /* Initially ci_b and ci to 1 before where we want them to start,
3390 * as they will both be incremented in next_ci.
3393 ci = th*ci_block - 1;
3396 while (next_ci(iGrid, nth, ci_block, &ci_x, &ci_y, &ci_b, &ci))
3398 if (bSimple && flags_i[ci] == 0)
3403 ncj_old_i = getNumSimpleJClustersInList(*nbl);
3406 if (!isIntraGridList && shp[XX] == 0)
3410 bx1 = bb_i[ci].upper.x;
3414 bx1 = iGridDims.lowerCorner[XX] + (ci_x+1)*iGridDims.cellSize[XX];
3416 if (bx1 < jGridDims.lowerCorner[XX])
3418 d2cx = gmx::square(jGridDims.lowerCorner[XX] - bx1);
3420 if (d2cx >= listRangeBBToJCell2)
3427 ci_xy = ci_x*iGridDims.numCells[YY] + ci_y;
3429 /* Loop over shift vectors in three dimensions */
3430 for (int tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
3432 const real shz = tz*box[ZZ][ZZ];
3434 bz0 = bbcz_i[ci*NNBSBB_D ] + shz;
3435 bz1 = bbcz_i[ci*NNBSBB_D+1] + shz;
3443 d2z = gmx::square(bz1);
3447 d2z = gmx::square(bz0 - box[ZZ][ZZ]);
3450 d2z_cx = d2z + d2cx;
3452 if (d2z_cx >= rlist2)
3457 bz1_frac = bz1/iGrid.numCellsInColumn(ci_xy);
3462 /* The check with bz1_frac close to or larger than 1 comes later */
3464 for (int ty = -shp[YY]; ty <= shp[YY]; ty++)
3466 const real shy = ty*box[YY][YY] + tz*box[ZZ][YY];
3470 by0 = bb_i[ci].lower.y + shy;
3471 by1 = bb_i[ci].upper.y + shy;
3475 by0 = iGridDims.lowerCorner[YY] + (ci_y )*iGridDims.cellSize[YY] + shy;
3476 by1 = iGridDims.lowerCorner[YY] + (ci_y + 1)*iGridDims.cellSize[YY] + shy;
3479 get_cell_range<YY>(by0, by1,
3490 if (by1 < jGridDims.lowerCorner[YY])
3492 d2z_cy += gmx::square(jGridDims.lowerCorner[YY] - by1);
3494 else if (by0 > jGridDims.upperCorner[YY])
3496 d2z_cy += gmx::square(by0 - jGridDims.upperCorner[YY]);
3499 for (int tx = -shp[XX]; tx <= shp[XX]; tx++)
3501 const int shift = XYZ2IS(tx, ty, tz);
3503 const bool excludeSubDiagonal = (isIntraGridList && shift == CENTRAL);
3505 if (c_pbcShiftBackward && isIntraGridList && shift > CENTRAL)
3510 const real shx = tx*box[XX][XX] + ty*box[YY][XX] + tz*box[ZZ][XX];
3514 bx0 = bb_i[ci].lower.x + shx;
3515 bx1 = bb_i[ci].upper.x + shx;
3519 bx0 = iGridDims.lowerCorner[XX] + (ci_x )*iGridDims.cellSize[XX] + shx;
3520 bx1 = iGridDims.lowerCorner[XX] + (ci_x+1)*iGridDims.cellSize[XX] + shx;
3523 get_cell_range<XX>(bx0, bx1,
3533 addNewIEntry(nbl, cell0_i+ci, shift, flags_i[ci]);
3535 if ((!c_pbcShiftBackward || excludeSubDiagonal) &&
3538 /* Leave the pairs with i > j.
3539 * x is the major index, so skip half of it.
3544 set_icell_bb(iGrid, ci, shx, shy, shz,
3547 icell_set_x(cell0_i+ci, shx, shy, shz,
3548 nbat->xstride, nbat->x().data(),
3552 for (int cx = cxf; cx <= cxl; cx++)
3555 if (jGridDims.lowerCorner[XX] + cx*jGridDims.cellSize[XX] > bx1)
3557 d2zx += gmx::square(jGridDims.lowerCorner[XX] + cx*jGridDims.cellSize[XX] - bx1);
3559 else if (jGridDims.lowerCorner[XX] + (cx+1)*jGridDims.cellSize[XX] < bx0)
3561 d2zx += gmx::square(jGridDims.lowerCorner[XX] + (cx+1)*jGridDims.cellSize[XX] - bx0);
3564 if (isIntraGridList &&
3566 (!c_pbcShiftBackward || shift == CENTRAL) &&
3569 /* Leave the pairs with i > j.
3570 * Skip half of y when i and j have the same x.
3579 for (int cy = cyf_x; cy <= cyl; cy++)
3581 const int columnStart = jGrid.firstCellInColumn(cx*jGridDims.numCells[YY] + cy);
3582 const int columnEnd = jGrid.firstCellInColumn(cx*jGridDims.numCells[YY] + cy + 1);
3585 if (jGridDims.lowerCorner[YY] + cy*jGridDims.cellSize[YY] > by1)
3587 d2zxy += gmx::square(jGridDims.lowerCorner[YY] + cy*jGridDims.cellSize[YY] - by1);
3589 else if (jGridDims.lowerCorner[YY] + (cy + 1)*jGridDims.cellSize[YY] < by0)
3591 d2zxy += gmx::square(jGridDims.lowerCorner[YY] + (cy + 1)*jGridDims.cellSize[YY] - by0);
3593 if (columnStart < columnEnd && d2zxy < listRangeBBToJCell2)
3595 /* To improve efficiency in the common case
3596 * of a homogeneous particle distribution,
3597 * we estimate the index of the middle cell
3598 * in range (midCell). We search down and up
3599 * starting from this index.
3601 * Note that the bbcz_j array contains bounds
3602 * for i-clusters, thus for clusters of 4 atoms.
3603 * For the common case where the j-cluster size
3604 * is 8, we could step with a stride of 2,
3605 * but we do not do this because it would
3606 * complicate this code even more.
3608 int midCell = columnStart + static_cast<int>(bz1_frac*(columnEnd - columnStart));
3609 if (midCell >= columnEnd)
3611 midCell = columnEnd - 1;
3616 /* Find the lowest cell that can possibly
3618 * Check if we hit the bottom of the grid,
3619 * if the j-cell is below the i-cell and if so,
3620 * if it is within range.
3622 int downTestCell = midCell;
3623 while (downTestCell >= columnStart &&
3624 (bbcz_j[downTestCell*NNBSBB_D + 1] >= bz0 ||
3625 d2xy + gmx::square(bbcz_j[downTestCell*NNBSBB_D + 1] - bz0) < rlist2))
3629 int firstCell = downTestCell + 1;
3631 /* Find the highest cell that can possibly
3633 * Check if we hit the top of the grid,
3634 * if the j-cell is above the i-cell and if so,
3635 * if it is within range.
3637 int upTestCell = midCell + 1;
3638 while (upTestCell < columnEnd &&
3639 (bbcz_j[upTestCell*NNBSBB_D] <= bz1 ||
3640 d2xy + gmx::square(bbcz_j[upTestCell*NNBSBB_D] - bz1) < rlist2))
3644 int lastCell = upTestCell - 1;
3646 #define NBNXN_REFCODE 0
3649 /* Simple reference code, for debugging,
3650 * overrides the more complex code above.
3652 firstCell = columnEnd;
3654 for (int k = columnStart; k < columnEnd; k++)
3656 if (d2xy + gmx::square(bbcz_j[k*NNBSBB_D + 1] - bz0) < rlist2 &&
3661 if (d2xy + gmx::square(bbcz_j[k*NNBSBB_D] - bz1) < rlist2 &&
3670 if (isIntraGridList)
3672 /* We want each atom/cell pair only once,
3673 * only use cj >= ci.
3675 if (!c_pbcShiftBackward || shift == CENTRAL)
3677 firstCell = std::max(firstCell, ci);
3681 if (firstCell <= lastCell)
3683 GMX_ASSERT(firstCell >= columnStart && lastCell < columnEnd, "The range should reside within the current grid column");
3685 /* For f buffer flags with simple lists */
3686 ncj_old_j = getNumSimpleJClustersInList(*nbl);
3688 makeClusterListWrapper(nbl,
3690 jGrid, firstCell, lastCell,
3695 &numDistanceChecks);
3699 setBufferFlags(*nbl, ncj_old_j, gridj_flag_shift,
3703 incrementNumSimpleJClustersInList(nbl, ncj_old_j);
3709 /* Set the exclusions for this ci list */
3710 setExclusionsForIEntry(nbs,
3714 *getOpenIEntry(nbl),
3719 make_fep_list(nbs, nbat, nbl,
3724 iGrid, jGrid, nbl_fep);
3727 /* Close this ci list */
3730 progBal, nsubpair_tot_est,
3736 if (bFBufferFlag && getNumSimpleJClustersInList(*nbl) > ncj_old_i)
3738 bitmask_init_bit(&(work->buffer_flags.flag[(iGrid.cellOffset() + ci) >> gridi_flag_shift]), th);
3742 work->ndistc = numDistanceChecks;
3744 nbs_cycle_stop(&work->cc[enbsCCsearch]);
3746 checkListSizeConsistency(*nbl, nbs->bFEP);
3750 fprintf(debug, "number of distance checks %d\n", numDistanceChecks);
3752 print_nblist_statistics(debug, nbl, nbs, rlist);
3756 fprintf(debug, "nbl FEP list pairs: %d\n", nbl_fep->nrj);
3761 static void reduce_buffer_flags(const nbnxn_search *nbs,
3763 const nbnxn_buffer_flags_t *dest)
3765 for (int s = 0; s < nsrc; s++)
3767 gmx_bitmask_t * flag = nbs->work[s].buffer_flags.flag;
3769 for (int b = 0; b < dest->nflag; b++)
3771 bitmask_union(&(dest->flag[b]), flag[b]);
3776 static void print_reduction_cost(const nbnxn_buffer_flags_t *flags, int nout)
3778 int nelem, nkeep, ncopy, nred, out;
3779 gmx_bitmask_t mask_0;
3785 bitmask_init_bit(&mask_0, 0);
3786 for (int b = 0; b < flags->nflag; b++)
3788 if (bitmask_is_equal(flags->flag[b], mask_0))
3790 /* Only flag 0 is set, no copy of reduction required */
3794 else if (!bitmask_is_zero(flags->flag[b]))
3797 for (out = 0; out < nout; out++)
3799 if (bitmask_is_set(flags->flag[b], out))
3816 fprintf(debug, "nbnxn reduction: #flag %d #list %d elem %4.2f, keep %4.2f copy %4.2f red %4.2f\n",
3818 nelem/static_cast<double>(flags->nflag),
3819 nkeep/static_cast<double>(flags->nflag),
3820 ncopy/static_cast<double>(flags->nflag),
3821 nred/static_cast<double>(flags->nflag));
3824 /* Copies the list entries from src to dest when cjStart <= *cjGlobal < cjEnd.
3825 * *cjGlobal is updated with the cj count in src.
3826 * When setFlags==true, flag bit t is set in flag for all i and j clusters.
3828 template<bool setFlags>
3829 static void copySelectedListRange(const nbnxn_ci_t * gmx_restrict srcCi,
3830 const NbnxnPairlistCpu * gmx_restrict src,
3831 NbnxnPairlistCpu * gmx_restrict dest,
3832 gmx_bitmask_t *flag,
3833 int iFlagShift, int jFlagShift, int t)
3835 const int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3837 dest->ci.push_back(*srcCi);
3838 dest->ci.back().cj_ind_start = dest->cj.size();
3839 dest->ci.back().cj_ind_end = dest->cj.size() + ncj;
3843 bitmask_init_bit(&flag[srcCi->ci >> iFlagShift], t);
3846 for (int j = srcCi->cj_ind_start; j < srcCi->cj_ind_end; j++)
3848 dest->cj.push_back(src->cj[j]);
3852 /* NOTE: This is relatively expensive, since this
3853 * operation is done for all elements in the list,
3854 * whereas at list generation this is done only
3855 * once for each flag entry.
3857 bitmask_init_bit(&flag[src->cj[j].cj >> jFlagShift], t);
3862 /* This routine re-balances the pairlists such that all are nearly equally
3863 * sized. Only whole i-entries are moved between lists. These are moved
3864 * between the ends of the lists, such that the buffer reduction cost should
3865 * not change significantly.
3866 * Note that all original reduction flags are currently kept. This can lead
3867 * to reduction of parts of the force buffer that could be avoided. But since
3868 * the original lists are quite balanced, this will only give minor overhead.
3870 static void rebalanceSimpleLists(int numLists,
3871 NbnxnPairlistCpu * const * const srcSet,
3872 NbnxnPairlistCpu **destSet,
3873 gmx::ArrayRef<nbnxn_search_work_t> searchWork)
3876 for (int s = 0; s < numLists; s++)
3878 ncjTotal += srcSet[s]->ncjInUse;
3880 int ncjTarget = (ncjTotal + numLists - 1)/numLists;
3882 #pragma omp parallel num_threads(numLists)
3884 int t = gmx_omp_get_thread_num();
3886 int cjStart = ncjTarget* t;
3887 int cjEnd = ncjTarget*(t + 1);
3889 /* The destination pair-list for task/thread t */
3890 NbnxnPairlistCpu *dest = destSet[t];
3892 clear_pairlist(dest);
3893 dest->na_cj = srcSet[0]->na_cj;
3895 /* Note that the flags in the work struct (still) contain flags
3896 * for all entries that are present in srcSet->nbl[t].
3898 gmx_bitmask_t *flag = searchWork[t].buffer_flags.flag;
3900 int iFlagShift = getBufferFlagShift(dest->na_ci);
3901 int jFlagShift = getBufferFlagShift(dest->na_cj);
3904 for (int s = 0; s < numLists && cjGlobal < cjEnd; s++)
3906 const NbnxnPairlistCpu *src = srcSet[s];
3908 if (cjGlobal + src->ncjInUse > cjStart)
3910 for (gmx::index i = 0; i < gmx::ssize(src->ci) && cjGlobal < cjEnd; i++)
3912 const nbnxn_ci_t *srcCi = &src->ci[i];
3913 int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3914 if (cjGlobal >= cjStart)
3916 /* If the source list is not our own, we need to set
3917 * extra flags (the template bool parameter).
3921 copySelectedListRange
3924 flag, iFlagShift, jFlagShift, t);
3928 copySelectedListRange
3931 dest, flag, iFlagShift, jFlagShift, t);
3939 cjGlobal += src->ncjInUse;
3943 dest->ncjInUse = dest->cj.size();
3947 int ncjTotalNew = 0;
3948 for (int s = 0; s < numLists; s++)
3950 ncjTotalNew += destSet[s]->ncjInUse;
3952 GMX_RELEASE_ASSERT(ncjTotalNew == ncjTotal, "The total size of the lists before and after rebalancing should match");
3956 /* Returns if the pairlists are so imbalanced that it is worth rebalancing. */
3957 static bool checkRebalanceSimpleLists(const nbnxn_pairlist_set_t *listSet)
3959 int numLists = listSet->nnbl;
3962 for (int s = 0; s < numLists; s++)
3964 ncjMax = std::max(ncjMax, listSet->nbl[s]->ncjInUse);
3965 ncjTotal += listSet->nbl[s]->ncjInUse;
3969 fprintf(debug, "Pair-list ncjMax %d ncjTotal %d\n", ncjMax, ncjTotal);
3971 /* The rebalancing adds 3% extra time to the search. Heuristically we
3972 * determined that under common conditions the non-bonded kernel balance
3973 * improvement will outweigh this when the imbalance is more than 3%.
3974 * But this will, obviously, depend on search vs kernel time and nstlist.
3976 const real rebalanceTolerance = 1.03;
3978 return numLists*ncjMax > ncjTotal*rebalanceTolerance;
3981 /* Perform a count (linear) sort to sort the smaller lists to the end.
3982 * This avoids load imbalance on the GPU, as large lists will be
3983 * scheduled and executed first and the smaller lists later.
3984 * Load balancing between multi-processors only happens at the end
3985 * and there smaller lists lead to more effective load balancing.
3986 * The sorting is done on the cj4 count, not on the actual pair counts.
3987 * Not only does this make the sort faster, but it also results in
3988 * better load balancing than using a list sorted on exact load.
3989 * This function swaps the pointer in the pair list to avoid a copy operation.
3991 static void sort_sci(NbnxnPairlistGpu *nbl)
3993 if (nbl->cj4.size() <= nbl->sci.size())
3995 /* nsci = 0 or all sci have size 1, sorting won't change the order */
3999 NbnxnPairlistGpuWork &work = *nbl->work;
4001 /* We will distinguish differences up to double the average */
4002 const int m = (2*nbl->cj4.size())/nbl->sci.size();
4004 /* Resize work.sci_sort so we can sort into it */
4005 work.sci_sort.resize(nbl->sci.size());
4007 std::vector<int> &sort = work.sortBuffer;
4008 /* Set up m + 1 entries in sort, initialized at 0 */
4010 sort.resize(m + 1, 0);
4011 /* Count the entries of each size */
4012 for (const nbnxn_sci_t &sci : nbl->sci)
4014 int i = std::min(m, sci.numJClusterGroups());
4017 /* Calculate the offset for each count */
4020 for (int i = m - 1; i >= 0; i--)
4023 sort[i] = sort[i + 1] + s0;
4027 /* Sort entries directly into place */
4028 gmx::ArrayRef<nbnxn_sci_t> sci_sort = work.sci_sort;
4029 for (const nbnxn_sci_t &sci : nbl->sci)
4031 int i = std::min(m, sci.numJClusterGroups());
4032 sci_sort[sort[i]++] = sci;
4035 /* Swap the sci pointers so we use the new, sorted list */
4036 std::swap(nbl->sci, work.sci_sort);
4040 nonbonded_verlet_t::PairlistSets::construct(const InteractionLocality iLocality,
4042 nbnxn_atomdata_t *nbat,
4043 const t_blocka *excl,
4044 const Nbnxm::KernelType kernelType,
4048 nbnxn_pairlist_set_t *nbl_list = &pairlistSet(iLocality);
4050 const real rlist = nbl_list->params.rlistOuter;
4052 int nsubpair_target;
4053 float nsubpair_tot_est;
4056 gmx_bool CombineNBLists;
4058 int np_tot, np_noq, np_hlj, nap;
4060 nnbl = nbl_list->nnbl;
4061 CombineNBLists = nbl_list->bCombined;
4065 fprintf(debug, "ns making %d nblists\n", nnbl);
4068 nbat->bUseBufferFlags = (nbat->out.size() > 1);
4069 /* We should re-init the flags before making the first list */
4070 if (nbat->bUseBufferFlags && iLocality == InteractionLocality::Local)
4072 init_buffer_flags(&nbat->buffer_flags, nbat->numAtoms());
4076 if (iLocality == InteractionLocality::Local)
4078 /* Only zone (grid) 0 vs 0 */
4083 nzi = nbs->zones->nizone;
4086 if (!nbl_list->bSimple && minimumIlistCountForGpuBalancing_ > 0)
4088 get_nsubpair_target(nbs, iLocality, rlist, minimumIlistCountForGpuBalancing_,
4089 &nsubpair_target, &nsubpair_tot_est);
4093 nsubpair_target = 0;
4094 nsubpair_tot_est = 0;
4097 /* Clear all pair-lists */
4098 for (int th = 0; th < nnbl; th++)
4100 if (nbl_list->bSimple)
4102 clear_pairlist(nbl_list->nbl[th]);
4106 clear_pairlist(nbl_list->nblGpu[th]);
4111 clear_pairlist_fep(nbl_list->nbl_fep[th]);
4115 for (int zi = 0; zi < nzi; zi++)
4117 const Grid &iGrid = nbs->grid[zi];
4121 if (iLocality == InteractionLocality::Local)
4128 zj0 = nbs->zones->izone[zi].j0;
4129 zj1 = nbs->zones->izone[zi].j1;
4135 for (int zj = zj0; zj < zj1; zj++)
4137 const Grid &jGrid = nbs->grid[zj];
4141 fprintf(debug, "ns search grid %d vs %d\n", zi, zj);
4144 nbs_cycle_start(&nbs->cc[enbsCCsearch]);
4146 ci_block = get_ci_block_size(iGrid, nbs->DomDec, nnbl);
4148 /* With GPU: generate progressively smaller lists for
4149 * load balancing for local only or non-local with 2 zones.
4151 progBal = (iLocality == InteractionLocality::Local || nbs->zones->n <= 2);
4153 #pragma omp parallel for num_threads(nnbl) schedule(static)
4154 for (int th = 0; th < nnbl; th++)
4158 /* Re-init the thread-local work flag data before making
4159 * the first list (not an elegant conditional).
4161 if (nbat->bUseBufferFlags && ((zi == 0 && zj == 0)))
4163 init_buffer_flags(&nbs->work[th].buffer_flags, nbat->numAtoms());
4166 if (CombineNBLists && th > 0)
4168 GMX_ASSERT(!nbl_list->bSimple, "Can only combine GPU lists");
4170 clear_pairlist(nbl_list->nblGpu[th]);
4173 /* Divide the i super cell equally over the nblists */
4174 if (nbl_list->bSimple)
4176 nbnxn_make_pairlist_part(nbs, iGrid, jGrid,
4177 &nbs->work[th], nbat, *excl,
4181 nbat->bUseBufferFlags,
4183 progBal, nsubpair_tot_est,
4186 nbl_list->nbl_fep[th]);
4190 nbnxn_make_pairlist_part(nbs, iGrid, jGrid,
4191 &nbs->work[th], nbat, *excl,
4195 nbat->bUseBufferFlags,
4197 progBal, nsubpair_tot_est,
4199 nbl_list->nblGpu[th],
4200 nbl_list->nbl_fep[th]);
4203 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4205 nbs_cycle_stop(&nbs->cc[enbsCCsearch]);
4210 for (int th = 0; th < nnbl; th++)
4212 inc_nrnb(nrnb, eNR_NBNXN_DIST2, nbs->work[th].ndistc);
4214 if (nbl_list->bSimple)
4216 NbnxnPairlistCpu *nbl = nbl_list->nbl[th];
4217 np_tot += nbl->cj.size();
4218 np_noq += nbl->work->ncj_noq;
4219 np_hlj += nbl->work->ncj_hlj;
4223 NbnxnPairlistGpu *nbl = nbl_list->nblGpu[th];
4224 /* This count ignores potential subsequent pair pruning */
4225 np_tot += nbl->nci_tot;
4228 if (nbl_list->bSimple)
4230 nap = nbl_list->nbl[0]->na_ci*nbl_list->nbl[0]->na_cj;
4234 nap = gmx::square(nbl_list->nblGpu[0]->na_ci);
4236 nbl_list->natpair_ljq = (np_tot - np_noq)*nap - np_hlj*nap/2;
4237 nbl_list->natpair_lj = np_noq*nap;
4238 nbl_list->natpair_q = np_hlj*nap/2;
4240 if (CombineNBLists && nnbl > 1)
4242 GMX_ASSERT(!nbl_list->bSimple, "Can only combine GPU lists");
4243 NbnxnPairlistGpu **nbl = nbl_list->nblGpu;
4245 nbs_cycle_start(&nbs->cc[enbsCCcombine]);
4247 combine_nblists(nnbl-1, nbl+1, nbl[0]);
4249 nbs_cycle_stop(&nbs->cc[enbsCCcombine]);
4254 if (nbl_list->bSimple)
4256 if (nnbl > 1 && checkRebalanceSimpleLists(nbl_list))
4258 rebalanceSimpleLists(nbl_list->nnbl, nbl_list->nbl, nbl_list->nbl_work, nbs->work);
4260 /* Swap the pointer of the sets of pair lists */
4261 NbnxnPairlistCpu **tmp = nbl_list->nbl;
4262 nbl_list->nbl = nbl_list->nbl_work;
4263 nbl_list->nbl_work = tmp;
4268 /* Sort the entries on size, large ones first */
4269 if (CombineNBLists || nnbl == 1)
4271 sort_sci(nbl_list->nblGpu[0]);
4275 #pragma omp parallel for num_threads(nnbl) schedule(static)
4276 for (int th = 0; th < nnbl; th++)
4280 sort_sci(nbl_list->nblGpu[th]);
4282 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4287 if (nbat->bUseBufferFlags)
4289 reduce_buffer_flags(nbs, nbl_list->nnbl, &nbat->buffer_flags);
4294 /* Balance the free-energy lists over all the threads */
4295 balance_fep_lists(nbs, nbl_list);
4298 if (nbl_list->bSimple)
4300 /* This is a fresh list, so not pruned, stored using ci.
4301 * ciOuter is invalid at this point.
4303 GMX_ASSERT(nbl_list->nbl[0]->ciOuter.empty(), "ciOuter is invalid so it should be empty");
4306 if (iLocality == Nbnxm::InteractionLocality::Local)
4308 outerListCreationStep_ = step;
4312 GMX_RELEASE_ASSERT(outerListCreationStep_ == step,
4313 "Outer list should be created at the same step as the inner list");
4316 /* Special performance logging stuff (env.var. GMX_NBNXN_CYCLE) */
4317 if (iLocality == InteractionLocality::Local)
4319 nbs->search_count++;
4321 if (nbs->print_cycles &&
4322 (!nbs->DomDec || iLocality == InteractionLocality::NonLocal) &&
4323 nbs->search_count % 100 == 0)
4325 nbs_cycle_print(stderr, nbs);
4328 /* If we have more than one list, they either got rebalancing (CPU)
4329 * or combined (GPU), so we should dump the final result to debug.
4331 if (debug && nbl_list->nnbl > 1)
4333 if (nbl_list->bSimple)
4335 for (int t = 0; t < nbl_list->nnbl; t++)
4337 print_nblist_statistics(debug, nbl_list->nbl[t], nbs, rlist);
4342 print_nblist_statistics(debug, nbl_list->nblGpu[0], nbs, rlist);
4350 if (nbl_list->bSimple)
4352 for (int t = 0; t < nbl_list->nnbl; t++)
4354 print_nblist_ci_cj(debug, nbl_list->nbl[t]);
4359 print_nblist_sci_cj(debug, nbl_list->nblGpu[0]);
4363 if (nbat->bUseBufferFlags)
4365 print_reduction_cost(&nbat->buffer_flags, nbl_list->nnbl);
4369 if (params_.useDynamicPruning && nbl_list->bSimple)
4371 nbnxnPrepareListForDynamicPruning(nbl_list);
4376 nonbonded_verlet_t::constructPairlist(const Nbnxm::InteractionLocality iLocality,
4377 const t_blocka *excl,
4381 pairlistSets_->construct(iLocality, nbs.get(), nbat.get(), excl,
4382 kernelSetup_.kernelType,
4387 /* Launch the transfer of the pairlist to the GPU.
4389 * NOTE: The launch overhead is currently not timed separately
4391 Nbnxm::gpu_init_pairlist(gpu_nbv,
4392 pairlistSets().pairlistSet(iLocality).nblGpu[0],
4397 void nbnxnPrepareListForDynamicPruning(nbnxn_pairlist_set_t *listSet)
4399 GMX_RELEASE_ASSERT(listSet->bSimple, "Should only be called for simple lists");
4401 /* TODO: Restructure the lists so we have actual outer and inner
4402 * list objects so we can set a single pointer instead of
4403 * swapping several pointers.
4406 for (int i = 0; i < listSet->nnbl; i++)
4408 NbnxnPairlistCpu &list = *listSet->nbl[i];
4410 /* The search produced a list in ci/cj.
4411 * Swap the list pointers so we get the outer list is ciOuter,cjOuter
4412 * and we can prune that to get an inner list in ci/cj.
4414 GMX_RELEASE_ASSERT(list.ciOuter.empty() && list.cjOuter.empty(),
4415 "The outer lists should be empty before preparation");
4417 std::swap(list.ci, list.ciOuter);
4418 std::swap(list.cj, list.cjOuter);