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.
38 #include "nbnxn_search.h"
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/nb_verlet.h"
55 #include "gromacs/mdlib/nbnxn_atomdata.h"
56 #include "gromacs/mdlib/nbnxn_consts.h"
57 #include "gromacs/mdlib/nbnxn_grid.h"
58 #include "gromacs/mdlib/nbnxn_internal.h"
59 #include "gromacs/mdlib/nbnxn_simd.h"
60 #include "gromacs/mdlib/nbnxn_util.h"
61 #include "gromacs/mdlib/ns.h"
62 #include "gromacs/mdtypes/group.h"
63 #include "gromacs/mdtypes/md_enums.h"
64 #include "gromacs/pbcutil/ishift.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/simd/simd.h"
67 #include "gromacs/simd/vector_operations.h"
68 #include "gromacs/topology/block.h"
69 #include "gromacs/utility/exceptions.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/gmxomp.h"
72 #include "gromacs/utility/smalloc.h"
74 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
77 /* We shift the i-particles backward for PBC.
78 * This leads to more conditionals than shifting forward.
79 * We do this to get more balanced pair lists.
81 constexpr bool c_pbcShiftBackward = true;
84 static void nbs_cycle_clear(nbnxn_cycle_t *cc)
86 for (int i = 0; i < enbsCCnr; i++)
93 static double Mcyc_av(const nbnxn_cycle_t *cc)
95 return static_cast<double>(cc->c)*1e-6/cc->count;
98 static void nbs_cycle_print(FILE *fp, const nbnxn_search *nbs)
101 fprintf(fp, "ns %4d grid %4.1f search %4.1f red.f %5.3f",
102 nbs->cc[enbsCCgrid].count,
103 Mcyc_av(&nbs->cc[enbsCCgrid]),
104 Mcyc_av(&nbs->cc[enbsCCsearch]),
105 Mcyc_av(&nbs->cc[enbsCCreducef]));
107 if (nbs->work.size() > 1)
109 if (nbs->cc[enbsCCcombine].count > 0)
111 fprintf(fp, " comb %5.2f",
112 Mcyc_av(&nbs->cc[enbsCCcombine]));
114 fprintf(fp, " s. th");
115 for (const nbnxn_search_work_t &work : nbs->work)
117 fprintf(fp, " %4.1f",
118 Mcyc_av(&work.cc[enbsCCsearch]));
124 /* Layout for the nonbonded NxN pair lists */
125 enum class NbnxnLayout
127 NoSimd4x4, // i-cluster size 4, j-cluster size 4
128 Simd4xN, // i-cluster size 4, j-cluster size SIMD width
129 Simd2xNN, // i-cluster size 4, j-cluster size half SIMD width
130 Gpu8x8x8 // i-cluster size 8, j-cluster size 8 + super-clustering
134 /* Returns the j-cluster size */
135 template <NbnxnLayout layout>
136 static constexpr int jClusterSize()
138 static_assert(layout == NbnxnLayout::NoSimd4x4 || layout == NbnxnLayout::Simd4xN || layout == NbnxnLayout::Simd2xNN, "Currently jClusterSize only supports CPU layouts");
140 return layout == NbnxnLayout::Simd4xN ? GMX_SIMD_REAL_WIDTH : (layout == NbnxnLayout::Simd2xNN ? GMX_SIMD_REAL_WIDTH/2 : c_nbnxnCpuIClusterSize);
143 /*! \brief Returns the j-cluster index given the i-cluster index.
145 * \tparam jClusterSize The number of atoms in a j-cluster
146 * \tparam jSubClusterIndex The j-sub-cluster index (0/1), used when size(j-cluster) < size(i-cluster)
147 * \param[in] ci The i-cluster index
149 template <int jClusterSize, int jSubClusterIndex>
150 static inline int cjFromCi(int ci)
152 static_assert(jClusterSize == c_nbnxnCpuIClusterSize/2 || jClusterSize == c_nbnxnCpuIClusterSize || jClusterSize == c_nbnxnCpuIClusterSize*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
154 static_assert(jSubClusterIndex == 0 || jSubClusterIndex == 1,
155 "Only sub-cluster indices 0 and 1 are supported");
157 if (jClusterSize == c_nbnxnCpuIClusterSize/2)
159 if (jSubClusterIndex == 0)
165 return ((ci + 1) << 1) - 1;
168 else if (jClusterSize == c_nbnxnCpuIClusterSize)
178 /*! \brief Returns the j-cluster index given the i-cluster index.
180 * \tparam layout The pair-list layout
181 * \tparam jSubClusterIndex The j-sub-cluster index (0/1), used when size(j-cluster) < size(i-cluster)
182 * \param[in] ci The i-cluster index
184 template <NbnxnLayout layout, int jSubClusterIndex>
185 static inline int cjFromCi(int ci)
187 constexpr int clusterSize = jClusterSize<layout>();
189 return cjFromCi<clusterSize, jSubClusterIndex>(ci);
192 /* Returns the nbnxn coordinate data index given the i-cluster index */
193 template <NbnxnLayout layout>
194 static inline int xIndexFromCi(int ci)
196 constexpr int clusterSize = jClusterSize<layout>();
198 static_assert(clusterSize == c_nbnxnCpuIClusterSize/2 || clusterSize == c_nbnxnCpuIClusterSize || clusterSize == c_nbnxnCpuIClusterSize*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
200 if (clusterSize <= c_nbnxnCpuIClusterSize)
202 /* Coordinates are stored packed in groups of 4 */
207 /* Coordinates packed in 8, i-cluster size is half the packing width */
208 return (ci >> 1)*STRIDE_P8 + (ci & 1)*(c_packX8 >> 1);
212 /* Returns the nbnxn coordinate data index given the j-cluster index */
213 template <NbnxnLayout layout>
214 static inline int xIndexFromCj(int cj)
216 constexpr int clusterSize = jClusterSize<layout>();
218 static_assert(clusterSize == c_nbnxnCpuIClusterSize/2 || clusterSize == c_nbnxnCpuIClusterSize || clusterSize == c_nbnxnCpuIClusterSize*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
220 if (clusterSize == c_nbnxnCpuIClusterSize/2)
222 /* Coordinates are stored packed in groups of 4 */
223 return (cj >> 1)*STRIDE_P4 + (cj & 1)*(c_packX4 >> 1);
225 else if (clusterSize == c_nbnxnCpuIClusterSize)
227 /* Coordinates are stored packed in groups of 4 */
232 /* Coordinates are stored packed in groups of 8 */
238 gmx_bool nbnxn_kernel_pairlist_simple(int nb_kernel_type)
240 if (nb_kernel_type == nbnxnkNotSet)
242 gmx_fatal(FARGS, "Non-bonded kernel type not set for Verlet-style pair-list.");
245 switch (nb_kernel_type)
247 case nbnxnk8x8x8_GPU:
248 case nbnxnk8x8x8_PlainC:
251 case nbnxnk4x4_PlainC:
252 case nbnxnk4xN_SIMD_4xN:
253 case nbnxnk4xN_SIMD_2xNN:
257 gmx_incons("Invalid nonbonded kernel type passed!");
262 /* Initializes a single nbnxn_pairlist_t data structure */
263 static void nbnxn_init_pairlist_fep(t_nblist *nl)
265 nl->type = GMX_NBLIST_INTERACTION_FREE_ENERGY;
266 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
267 /* The interaction functions are set in the free energy kernel fuction */
280 nl->jindex = nullptr;
282 nl->excl_fep = nullptr;
286 static void free_nblist(t_nblist *nl)
296 nbnxn_search_work_t::nbnxn_search_work_t() :
299 buffer_flags({0, nullptr, 0}),
301 nbl_fep(new t_nblist),
304 nbnxn_init_pairlist_fep(nbl_fep.get());
309 nbnxn_search_work_t::~nbnxn_search_work_t()
311 sfree(buffer_flags.flag);
313 free_nblist(nbl_fep.get());
316 nbnxn_search::nbnxn_search(const ivec *n_dd_cells,
317 const gmx_domdec_zones_t *zones,
321 ePBC(epbcNONE), // The correct value will be set during the gridding
328 // The correct value will be set during the gridding
332 DomDec = n_dd_cells != nullptr;
335 for (int d = 0; d < DIM; d++)
337 if ((*n_dd_cells)[d] > 1)
340 /* Each grid matches a DD zone */
346 grid.resize(numGrids);
348 /* Initialize detailed nbsearch cycle counting */
349 print_cycles = (getenv("GMX_NBNXN_CYCLE") != nullptr);
353 nbnxn_search *nbnxn_init_search(const ivec *n_dd_cells,
354 const gmx_domdec_zones_t *zones,
358 return new nbnxn_search(n_dd_cells, zones, bFEP, nthread_max);
361 static void init_buffer_flags(nbnxn_buffer_flags_t *flags,
364 flags->nflag = (natoms + NBNXN_BUFFERFLAG_SIZE - 1)/NBNXN_BUFFERFLAG_SIZE;
365 if (flags->nflag > flags->flag_nalloc)
367 flags->flag_nalloc = over_alloc_large(flags->nflag);
368 srenew(flags->flag, flags->flag_nalloc);
370 for (int b = 0; b < flags->nflag; b++)
372 bitmask_clear(&(flags->flag[b]));
376 /* Returns the pair-list cutoff between a bounding box and a grid cell given an atom-to-atom pair-list cutoff
378 * Given a cutoff distance between atoms, this functions returns the cutoff
379 * distance2 between a bounding box of a group of atoms and a grid cell.
380 * Since atoms can be geometrically outside of the cell they have been
381 * assigned to (when atom groups instead of individual atoms are assigned
382 * to cells), this distance returned can be larger than the input.
384 static real listRangeForBoundingBoxToGridCell(real rlist,
385 const nbnxn_grid_t &grid)
387 return rlist + grid.maxAtomGroupRadius;
390 /* Returns the pair-list cutoff between a grid cells given an atom-to-atom pair-list cutoff
392 * Given a cutoff distance between atoms, this functions returns the cutoff
393 * distance2 between two grid cells.
394 * Since atoms can be geometrically outside of the cell they have been
395 * assigned to (when atom groups instead of individual atoms are assigned
396 * to cells), this distance returned can be larger than the input.
398 static real listRangeForGridCellToGridCell(real rlist,
399 const nbnxn_grid_t &iGrid,
400 const nbnxn_grid_t &jGrid)
402 return rlist + iGrid.maxAtomGroupRadius + jGrid.maxAtomGroupRadius;
405 /* Determines the cell range along one dimension that
406 * the bounding box b0 - b1 sees.
409 static void get_cell_range(real b0, real b1,
410 const nbnxn_grid_t &jGrid,
411 real d2, real rlist, int *cf, int *cl)
413 real listRangeBBToCell2 = gmx::square(listRangeForBoundingBoxToGridCell(rlist, jGrid));
414 real distanceInCells = (b0 - jGrid.c0[dim])*jGrid.invCellSize[dim];
415 *cf = std::max(static_cast<int>(distanceInCells), 0);
418 d2 + gmx::square((b0 - jGrid.c0[dim]) - (*cf - 1 + 1)*jGrid.cellSize[dim]) < listRangeBBToCell2)
423 *cl = std::min(static_cast<int>((b1 - jGrid.c0[dim])*jGrid.invCellSize[dim]), jGrid.numCells[dim] - 1);
424 while (*cl < jGrid.numCells[dim] - 1 &&
425 d2 + gmx::square((*cl + 1)*jGrid.cellSize[dim] - (b1 - jGrid.c0[dim])) < listRangeBBToCell2)
431 /* Reference code calculating the distance^2 between two bounding boxes */
433 static float box_dist2(float bx0, float bx1, float by0,
434 float by1, float bz0, float bz1,
435 const nbnxn_bb_t *bb)
438 float dl, dh, dm, dm0;
442 dl = bx0 - bb->upper[BB_X];
443 dh = bb->lower[BB_X] - bx1;
444 dm = std::max(dl, dh);
445 dm0 = std::max(dm, 0.0f);
448 dl = by0 - bb->upper[BB_Y];
449 dh = bb->lower[BB_Y] - by1;
450 dm = std::max(dl, dh);
451 dm0 = std::max(dm, 0.0f);
454 dl = bz0 - bb->upper[BB_Z];
455 dh = bb->lower[BB_Z] - bz1;
456 dm = std::max(dl, dh);
457 dm0 = std::max(dm, 0.0f);
464 /* Plain C code calculating the distance^2 between two bounding boxes */
465 static float subc_bb_dist2(int si,
466 const nbnxn_bb_t *bb_i_ci,
468 gmx::ArrayRef<const nbnxn_bb_t> bb_j_all)
470 const nbnxn_bb_t *bb_i = bb_i_ci + si;
471 const nbnxn_bb_t *bb_j = bb_j_all.data() + csj;
474 float dl, dh, dm, dm0;
476 dl = bb_i->lower[BB_X] - bb_j->upper[BB_X];
477 dh = bb_j->lower[BB_X] - bb_i->upper[BB_X];
478 dm = std::max(dl, dh);
479 dm0 = std::max(dm, 0.0f);
482 dl = bb_i->lower[BB_Y] - bb_j->upper[BB_Y];
483 dh = bb_j->lower[BB_Y] - bb_i->upper[BB_Y];
484 dm = std::max(dl, dh);
485 dm0 = std::max(dm, 0.0f);
488 dl = bb_i->lower[BB_Z] - bb_j->upper[BB_Z];
489 dh = bb_j->lower[BB_Z] - bb_i->upper[BB_Z];
490 dm = std::max(dl, dh);
491 dm0 = std::max(dm, 0.0f);
497 #if NBNXN_SEARCH_BB_SIMD4
499 /* 4-wide SIMD code for bb distance for bb format xyz0 */
500 static float subc_bb_dist2_simd4(int si,
501 const nbnxn_bb_t *bb_i_ci,
503 gmx::ArrayRef<const nbnxn_bb_t> bb_j_all)
505 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
508 Simd4Float bb_i_S0, bb_i_S1;
509 Simd4Float bb_j_S0, bb_j_S1;
515 bb_i_S0 = load4(&bb_i_ci[si].lower[0]);
516 bb_i_S1 = load4(&bb_i_ci[si].upper[0]);
517 bb_j_S0 = load4(&bb_j_all[csj].lower[0]);
518 bb_j_S1 = load4(&bb_j_all[csj].upper[0]);
520 dl_S = bb_i_S0 - bb_j_S1;
521 dh_S = bb_j_S0 - bb_i_S1;
523 dm_S = max(dl_S, dh_S);
524 dm0_S = max(dm_S, simd4SetZeroF());
526 return dotProduct(dm0_S, dm0_S);
529 /* Calculate bb bounding distances of bb_i[si,...,si+3] and store them in d2 */
530 #define SUBC_BB_DIST2_SIMD4_XXXX_INNER(si, bb_i, d2) \
534 Simd4Float dx_0, dy_0, dz_0; \
535 Simd4Float dx_1, dy_1, dz_1; \
537 Simd4Float mx, my, mz; \
538 Simd4Float m0x, m0y, m0z; \
540 Simd4Float d2x, d2y, d2z; \
541 Simd4Float d2s, d2t; \
543 shi = (si)*NNBSBB_D*DIM; \
545 xi_l = load4((bb_i)+shi+0*STRIDE_PBB); \
546 yi_l = load4((bb_i)+shi+1*STRIDE_PBB); \
547 zi_l = load4((bb_i)+shi+2*STRIDE_PBB); \
548 xi_h = load4((bb_i)+shi+3*STRIDE_PBB); \
549 yi_h = load4((bb_i)+shi+4*STRIDE_PBB); \
550 zi_h = load4((bb_i)+shi+5*STRIDE_PBB); \
552 dx_0 = xi_l - xj_h; \
553 dy_0 = yi_l - yj_h; \
554 dz_0 = zi_l - zj_h; \
556 dx_1 = xj_l - xi_h; \
557 dy_1 = yj_l - yi_h; \
558 dz_1 = zj_l - zi_h; \
560 mx = max(dx_0, dx_1); \
561 my = max(dy_0, dy_1); \
562 mz = max(dz_0, dz_1); \
564 m0x = max(mx, zero); \
565 m0y = max(my, zero); \
566 m0z = max(mz, zero); \
575 store4((d2)+(si), d2t); \
578 /* 4-wide SIMD code for nsi bb distances for bb format xxxxyyyyzzzz */
579 static void subc_bb_dist2_simd4_xxxx(const float *bb_j,
580 int nsi, const float *bb_i,
583 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
586 Simd4Float xj_l, yj_l, zj_l;
587 Simd4Float xj_h, yj_h, zj_h;
588 Simd4Float xi_l, yi_l, zi_l;
589 Simd4Float xi_h, yi_h, zi_h;
595 xj_l = Simd4Float(bb_j[0*STRIDE_PBB]);
596 yj_l = Simd4Float(bb_j[1*STRIDE_PBB]);
597 zj_l = Simd4Float(bb_j[2*STRIDE_PBB]);
598 xj_h = Simd4Float(bb_j[3*STRIDE_PBB]);
599 yj_h = Simd4Float(bb_j[4*STRIDE_PBB]);
600 zj_h = Simd4Float(bb_j[5*STRIDE_PBB]);
602 /* Here we "loop" over si (0,STRIDE_PBB) from 0 to nsi with step STRIDE_PBB.
603 * But as we know the number of iterations is 1 or 2, we unroll manually.
605 SUBC_BB_DIST2_SIMD4_XXXX_INNER(0, bb_i, d2);
606 if (STRIDE_PBB < nsi)
608 SUBC_BB_DIST2_SIMD4_XXXX_INNER(STRIDE_PBB, bb_i, d2);
612 #endif /* NBNXN_SEARCH_BB_SIMD4 */
615 /* Returns if any atom pair from two clusters is within distance sqrt(rlist2) */
616 static inline gmx_bool
617 clusterpair_in_range(const NbnxnPairlistGpuWork &work,
619 int csj, int stride, const real *x_j,
622 #if !GMX_SIMD4_HAVE_REAL
625 * All coordinates are stored as xyzxyz...
628 const real *x_i = work.iSuperClusterData.x.data();
630 for (int i = 0; i < c_nbnxnGpuClusterSize; i++)
632 int i0 = (si*c_nbnxnGpuClusterSize + i)*DIM;
633 for (int j = 0; j < c_nbnxnGpuClusterSize; j++)
635 int j0 = (csj*c_nbnxnGpuClusterSize + j)*stride;
637 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]);
648 #else /* !GMX_SIMD4_HAVE_REAL */
650 /* 4-wide SIMD version.
651 * The coordinates x_i are stored as xxxxyyyy..., x_j is stored xyzxyz...
652 * Using 8-wide AVX(2) is not faster on Intel Sandy Bridge and Haswell.
654 static_assert(c_nbnxnGpuClusterSize == 8 || c_nbnxnGpuClusterSize == 4,
655 "A cluster is hard-coded to 4/8 atoms.");
657 Simd4Real rc2_S = Simd4Real(rlist2);
659 const real *x_i = work.iSuperClusterData.xSimd.data();
661 int dim_stride = c_nbnxnGpuClusterSize*DIM;
662 Simd4Real ix_S0 = load4(x_i + si*dim_stride + 0*GMX_SIMD4_WIDTH);
663 Simd4Real iy_S0 = load4(x_i + si*dim_stride + 1*GMX_SIMD4_WIDTH);
664 Simd4Real iz_S0 = load4(x_i + si*dim_stride + 2*GMX_SIMD4_WIDTH);
666 Simd4Real ix_S1, iy_S1, iz_S1;
667 if (c_nbnxnGpuClusterSize == 8)
669 ix_S1 = load4(x_i + si*dim_stride + 3*GMX_SIMD4_WIDTH);
670 iy_S1 = load4(x_i + si*dim_stride + 4*GMX_SIMD4_WIDTH);
671 iz_S1 = load4(x_i + si*dim_stride + 5*GMX_SIMD4_WIDTH);
673 /* We loop from the outer to the inner particles to maximize
674 * the chance that we find a pair in range quickly and return.
676 int j0 = csj*c_nbnxnGpuClusterSize;
677 int j1 = j0 + c_nbnxnGpuClusterSize - 1;
680 Simd4Real jx0_S, jy0_S, jz0_S;
681 Simd4Real jx1_S, jy1_S, jz1_S;
683 Simd4Real dx_S0, dy_S0, dz_S0;
684 Simd4Real dx_S1, dy_S1, dz_S1;
685 Simd4Real dx_S2, dy_S2, dz_S2;
686 Simd4Real dx_S3, dy_S3, dz_S3;
697 Simd4Bool wco_any_S01, wco_any_S23, wco_any_S;
699 jx0_S = Simd4Real(x_j[j0*stride+0]);
700 jy0_S = Simd4Real(x_j[j0*stride+1]);
701 jz0_S = Simd4Real(x_j[j0*stride+2]);
703 jx1_S = Simd4Real(x_j[j1*stride+0]);
704 jy1_S = Simd4Real(x_j[j1*stride+1]);
705 jz1_S = Simd4Real(x_j[j1*stride+2]);
707 /* Calculate distance */
708 dx_S0 = ix_S0 - jx0_S;
709 dy_S0 = iy_S0 - jy0_S;
710 dz_S0 = iz_S0 - jz0_S;
711 dx_S2 = ix_S0 - jx1_S;
712 dy_S2 = iy_S0 - jy1_S;
713 dz_S2 = iz_S0 - jz1_S;
714 if (c_nbnxnGpuClusterSize == 8)
716 dx_S1 = ix_S1 - jx0_S;
717 dy_S1 = iy_S1 - jy0_S;
718 dz_S1 = iz_S1 - jz0_S;
719 dx_S3 = ix_S1 - jx1_S;
720 dy_S3 = iy_S1 - jy1_S;
721 dz_S3 = iz_S1 - jz1_S;
724 /* rsq = dx*dx+dy*dy+dz*dz */
725 rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
726 rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
727 if (c_nbnxnGpuClusterSize == 8)
729 rsq_S1 = norm2(dx_S1, dy_S1, dz_S1);
730 rsq_S3 = norm2(dx_S3, dy_S3, dz_S3);
733 wco_S0 = (rsq_S0 < rc2_S);
734 wco_S2 = (rsq_S2 < rc2_S);
735 if (c_nbnxnGpuClusterSize == 8)
737 wco_S1 = (rsq_S1 < rc2_S);
738 wco_S3 = (rsq_S3 < rc2_S);
740 if (c_nbnxnGpuClusterSize == 8)
742 wco_any_S01 = wco_S0 || wco_S1;
743 wco_any_S23 = wco_S2 || wco_S3;
744 wco_any_S = wco_any_S01 || wco_any_S23;
748 wco_any_S = wco_S0 || wco_S2;
751 if (anyTrue(wco_any_S))
762 #endif /* !GMX_SIMD4_HAVE_REAL */
765 /* Returns the j-cluster index for index cjIndex in a cj list */
766 static inline int nblCj(gmx::ArrayRef<const nbnxn_cj_t> cjList,
769 return cjList[cjIndex].cj;
772 /* Returns the j-cluster index for index cjIndex in a cj4 list */
773 static inline int nblCj(gmx::ArrayRef<const nbnxn_cj4_t> cj4List,
776 return cj4List[cjIndex/c_nbnxnGpuJgroupSize].cj[cjIndex & (c_nbnxnGpuJgroupSize - 1)];
779 /* Returns the i-interaction mask of the j sub-cell for index cj_ind */
780 static unsigned int nbl_imask0(const NbnxnPairlistGpu *nbl, int cj_ind)
782 return nbl->cj4[cj_ind/c_nbnxnGpuJgroupSize].imei[0].imask;
785 /* Initializes a single NbnxnPairlistCpu data structure */
786 static void nbnxn_init_pairlist(NbnxnPairlistCpu *nbl)
788 nbl->na_ci = c_nbnxnCpuIClusterSize;
791 nbl->ciOuter.clear();
794 nbl->cjOuter.clear();
797 nbl->work = new NbnxnPairlistCpuWork();
800 NbnxnPairlistGpu::NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy) :
801 na_ci(c_nbnxnGpuClusterSize),
802 na_cj(c_nbnxnGpuClusterSize),
803 na_sc(c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize),
805 sci({}, {pinningPolicy}),
806 cj4({}, {pinningPolicy}),
807 excl({}, {pinningPolicy}),
810 static_assert(c_nbnxnGpuNumClusterPerSupercluster == c_gpuNumClusterPerCell,
811 "The search code assumes that the a super-cluster matches a search grid cell");
813 static_assert(sizeof(cj4[0].imei[0].imask)*8 >= c_nbnxnGpuJgroupSize*c_gpuNumClusterPerCell,
814 "The i super-cluster cluster interaction mask does not contain a sufficient number of bits");
816 static_assert(sizeof(excl[0])*8 >= c_nbnxnGpuJgroupSize*c_gpuNumClusterPerCell, "The GPU exclusion mask does not contain a sufficient number of bits");
818 // We always want a first entry without any exclusions
821 work = new NbnxnPairlistGpuWork();
824 void nbnxn_init_pairlist_set(nbnxn_pairlist_set_t *nbl_list,
825 gmx_bool bSimple, gmx_bool bCombined)
827 GMX_RELEASE_ASSERT(!bSimple || !bCombined, "Can only combine non-simple lists");
829 nbl_list->bSimple = bSimple;
830 nbl_list->bCombined = bCombined;
832 nbl_list->nnbl = gmx_omp_nthreads_get(emntNonbonded);
834 if (!nbl_list->bCombined &&
835 nbl_list->nnbl > NBNXN_BUFFERFLAG_MAX_THREADS)
837 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.",
838 nbl_list->nnbl, NBNXN_BUFFERFLAG_MAX_THREADS, NBNXN_BUFFERFLAG_MAX_THREADS);
843 snew(nbl_list->nbl, nbl_list->nnbl);
844 if (nbl_list->nnbl > 1)
846 snew(nbl_list->nbl_work, nbl_list->nnbl);
851 snew(nbl_list->nblGpu, nbl_list->nnbl);
853 snew(nbl_list->nbl_fep, nbl_list->nnbl);
854 /* Execute in order to avoid memory interleaving between threads */
855 #pragma omp parallel for num_threads(nbl_list->nnbl) schedule(static)
856 for (int i = 0; i < nbl_list->nnbl; i++)
860 /* Allocate the nblist data structure locally on each thread
861 * to optimize memory access for NUMA architectures.
865 nbl_list->nbl[i] = new NbnxnPairlistCpu();
867 nbnxn_init_pairlist(nbl_list->nbl[i]);
868 if (nbl_list->nnbl > 1)
870 nbl_list->nbl_work[i] = new NbnxnPairlistCpu();
871 nbnxn_init_pairlist(nbl_list->nbl_work[i]);
876 /* Only list 0 is used on the GPU, use normal allocation for i>0 */
877 auto pinningPolicy = (i == 0 ? gmx::PinningPolicy::PinnedIfSupported : gmx::PinningPolicy::CannotBePinned);
879 nbl_list->nblGpu[i] = new NbnxnPairlistGpu(pinningPolicy);
882 snew(nbl_list->nbl_fep[i], 1);
883 nbnxn_init_pairlist_fep(nbl_list->nbl_fep[i]);
885 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
889 /* Print statistics of a pair list, used for debug output */
890 static void print_nblist_statistics(FILE *fp, const NbnxnPairlistCpu *nbl,
891 const nbnxn_search *nbs, real rl)
893 const nbnxn_grid_t *grid;
897 grid = &nbs->grid[0];
899 fprintf(fp, "nbl nci %zu ncj %d\n",
900 nbl->ci.size(), nbl->ncjInUse);
901 fprintf(fp, "nbl na_cj %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
902 nbl->na_cj, rl, nbl->ncjInUse, nbl->ncjInUse/static_cast<double>(grid->nc),
903 nbl->ncjInUse/static_cast<double>(grid->nc)*grid->na_cj,
904 nbl->ncjInUse/static_cast<double>(grid->nc)*grid->na_cj/(0.5*4.0/3.0*M_PI*rl*rl*rl*grid->nc*grid->na_cj/(grid->size[XX]*grid->size[YY]*grid->size[ZZ])));
906 fprintf(fp, "nbl average j cell list length %.1f\n",
907 0.25*nbl->ncjInUse/std::max(static_cast<double>(nbl->ci.size()), 1.0));
909 for (int s = 0; s < SHIFTS; s++)
914 for (const nbnxn_ci_t &ciEntry : nbl->ci)
916 cs[ciEntry.shift & NBNXN_CI_SHIFT] +=
917 ciEntry.cj_ind_end - ciEntry.cj_ind_start;
919 int j = ciEntry.cj_ind_start;
920 while (j < ciEntry.cj_ind_end &&
921 nbl->cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
927 fprintf(fp, "nbl cell pairs, total: %zu excl: %d %.1f%%\n",
928 nbl->cj.size(), npexcl, 100*npexcl/std::max(static_cast<double>(nbl->cj.size()), 1.0));
929 for (int s = 0; s < SHIFTS; s++)
933 fprintf(fp, "nbl shift %2d ncj %3d\n", s, cs[s]);
938 /* Print statistics of a pair lists, used for debug output */
939 static void print_nblist_statistics(FILE *fp, const NbnxnPairlistGpu *nbl,
940 const nbnxn_search *nbs, real rl)
942 const nbnxn_grid_t *grid;
944 int c[c_gpuNumClusterPerCell + 1];
945 double sum_nsp, sum_nsp2;
948 /* This code only produces correct statistics with domain decomposition */
949 grid = &nbs->grid[0];
951 fprintf(fp, "nbl nsci %zu ncj4 %zu nsi %d excl4 %zu\n",
952 nbl->sci.size(), nbl->cj4.size(), nbl->nci_tot, nbl->excl.size());
953 fprintf(fp, "nbl na_c %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
954 nbl->na_ci, rl, nbl->nci_tot, nbl->nci_tot/static_cast<double>(grid->nsubc_tot),
955 nbl->nci_tot/static_cast<double>(grid->nsubc_tot)*grid->na_c,
956 nbl->nci_tot/static_cast<double>(grid->nsubc_tot)*grid->na_c/(0.5*4.0/3.0*M_PI*rl*rl*rl*grid->nsubc_tot*grid->na_c/(grid->size[XX]*grid->size[YY]*grid->size[ZZ])));
961 for (int si = 0; si <= c_gpuNumClusterPerCell; si++)
965 for (const nbnxn_sci_t &sci : nbl->sci)
968 for (int j4 = sci.cj4_ind_start; j4 < sci.cj4_ind_end; j4++)
970 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
973 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
975 if (nbl->cj4[j4].imei[0].imask & (1U << (j*c_gpuNumClusterPerCell + si)))
986 nsp_max = std::max(nsp_max, nsp);
988 if (!nbl->sci.empty())
990 sum_nsp /= nbl->sci.size();
991 sum_nsp2 /= nbl->sci.size();
993 fprintf(fp, "nbl #cluster-pairs: av %.1f stddev %.1f max %d\n",
994 sum_nsp, std::sqrt(sum_nsp2 - sum_nsp*sum_nsp), nsp_max);
996 if (!nbl->cj4.empty())
998 for (b = 0; b <= c_gpuNumClusterPerCell; b++)
1000 fprintf(fp, "nbl j-list #i-subcell %d %7d %4.1f\n",
1001 b, c[b], 100.0*c[b]/size_t {nbl->cj4.size()*c_nbnxnGpuJgroupSize});
1006 /* Returns a pointer to the exclusion mask for j-cluster-group \p cj4 and warp \p warp
1007 * Generates a new exclusion entry when the j-cluster-group uses
1008 * the default all-interaction mask at call time, so the returned mask
1009 * can be modified when needed.
1011 static nbnxn_excl_t *get_exclusion_mask(NbnxnPairlistGpu *nbl,
1015 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
1017 /* No exclusions set, make a new list entry */
1018 const size_t oldSize = nbl->excl.size();
1019 GMX_ASSERT(oldSize >= 1, "We should always have entry [0]");
1020 /* Add entry with default values: no exclusions */
1021 nbl->excl.resize(oldSize + 1);
1022 nbl->cj4[cj4].imei[warp].excl_ind = oldSize;
1025 return &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
1028 static void set_self_and_newton_excls_supersub(NbnxnPairlistGpu *nbl,
1029 int cj4_ind, int sj_offset,
1030 int i_cluster_in_cell)
1032 nbnxn_excl_t *excl[c_nbnxnGpuClusterpairSplit];
1034 /* Here we only set the set self and double pair exclusions */
1036 /* Reserve extra elements, so the resize() in get_exclusion_mask()
1037 * will not invalidate excl entries in the loop below
1039 nbl->excl.reserve(nbl->excl.size() + c_nbnxnGpuClusterpairSplit);
1040 for (int w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
1042 excl[w] = get_exclusion_mask(nbl, cj4_ind, w);
1045 /* Only minor < major bits set */
1046 for (int ej = 0; ej < nbl->na_ci; ej++)
1049 for (int ei = ej; ei < nbl->na_ci; ei++)
1051 excl[w]->pair[(ej & (c_nbnxnGpuJgroupSize-1))*nbl->na_ci + ei] &=
1052 ~(1U << (sj_offset*c_gpuNumClusterPerCell + i_cluster_in_cell));
1057 /* Returns a diagonal or off-diagonal interaction mask for plain C lists */
1058 static unsigned int get_imask(gmx_bool rdiag, int ci, int cj)
1060 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
1063 /* Returns a diagonal or off-diagonal interaction mask for cj-size=2 */
1064 gmx_unused static unsigned int get_imask_simd_j2(gmx_bool rdiag, int ci, int cj)
1066 return (rdiag && ci*2 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_0 :
1067 (rdiag && ci*2+1 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_1 :
1068 NBNXN_INTERACTION_MASK_ALL));
1071 /* Returns a diagonal or off-diagonal interaction mask for cj-size=4 */
1072 gmx_unused static unsigned int get_imask_simd_j4(gmx_bool rdiag, int ci, int cj)
1074 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
1077 /* Returns a diagonal or off-diagonal interaction mask for cj-size=8 */
1078 gmx_unused static unsigned int get_imask_simd_j8(gmx_bool rdiag, int ci, int cj)
1080 return (rdiag && ci == cj*2 ? NBNXN_INTERACTION_MASK_DIAG_J8_0 :
1081 (rdiag && ci == cj*2+1 ? NBNXN_INTERACTION_MASK_DIAG_J8_1 :
1082 NBNXN_INTERACTION_MASK_ALL));
1086 #if GMX_SIMD_REAL_WIDTH == 2
1087 #define get_imask_simd_4xn get_imask_simd_j2
1089 #if GMX_SIMD_REAL_WIDTH == 4
1090 #define get_imask_simd_4xn get_imask_simd_j4
1092 #if GMX_SIMD_REAL_WIDTH == 8
1093 #define get_imask_simd_4xn get_imask_simd_j8
1094 #define get_imask_simd_2xnn get_imask_simd_j4
1096 #if GMX_SIMD_REAL_WIDTH == 16
1097 #define get_imask_simd_2xnn get_imask_simd_j8
1101 /* Plain C code for checking and adding cluster-pairs to the list.
1103 * \param[in] gridj The j-grid
1104 * \param[in,out] nbl The pair-list to store the cluster pairs in
1105 * \param[in] icluster The index of the i-cluster
1106 * \param[in] jclusterFirst The first cluster in the j-range
1107 * \param[in] jclusterLast The last cluster in the j-range
1108 * \param[in] excludeSubDiagonal Exclude atom pairs with i-index > j-index
1109 * \param[in] x_j Coordinates for the j-atom, in xyz format
1110 * \param[in] rlist2 The squared list cut-off
1111 * \param[in] rbb2 The squared cut-off for putting cluster-pairs in the list based on bounding box distance only
1112 * \param[in,out] numDistanceChecks The number of distance checks performed
1115 makeClusterListSimple(const nbnxn_grid_t &jGrid,
1116 NbnxnPairlistCpu * nbl,
1120 bool excludeSubDiagonal,
1121 const real * gmx_restrict x_j,
1124 int * gmx_restrict numDistanceChecks)
1126 const nbnxn_bb_t * gmx_restrict bb_ci = nbl->work->iClusterData.bb.data();
1127 const real * gmx_restrict x_ci = nbl->work->iClusterData.x.data();
1132 while (!InRange && jclusterFirst <= jclusterLast)
1134 real d2 = subc_bb_dist2(0, bb_ci, jclusterFirst, jGrid.bb);
1135 *numDistanceChecks += 2;
1137 /* Check if the distance is within the distance where
1138 * we use only the bounding box distance rbb,
1139 * or within the cut-off and there is at least one atom pair
1140 * within the cut-off.
1146 else if (d2 < rlist2)
1148 int cjf_gl = jGrid.cell0 + jclusterFirst;
1149 for (int i = 0; i < c_nbnxnCpuIClusterSize && !InRange; i++)
1151 for (int j = 0; j < c_nbnxnCpuIClusterSize; j++)
1153 InRange = InRange ||
1154 (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjf_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+XX]) +
1155 gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjf_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+YY]) +
1156 gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjf_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+ZZ]) < rlist2);
1159 *numDistanceChecks += c_nbnxnCpuIClusterSize*c_nbnxnCpuIClusterSize;
1172 while (!InRange && jclusterLast > jclusterFirst)
1174 real d2 = subc_bb_dist2(0, bb_ci, jclusterLast, jGrid.bb);
1175 *numDistanceChecks += 2;
1177 /* Check if the distance is within the distance where
1178 * we use only the bounding box distance rbb,
1179 * or within the cut-off and there is at least one atom pair
1180 * within the cut-off.
1186 else if (d2 < rlist2)
1188 int cjl_gl = jGrid.cell0 + jclusterLast;
1189 for (int i = 0; i < c_nbnxnCpuIClusterSize && !InRange; i++)
1191 for (int j = 0; j < c_nbnxnCpuIClusterSize; j++)
1193 InRange = InRange ||
1194 (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjl_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+XX]) +
1195 gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjl_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+YY]) +
1196 gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjl_gl*c_nbnxnCpuIClusterSize+j)*STRIDE_XYZ+ZZ]) < rlist2);
1199 *numDistanceChecks += c_nbnxnCpuIClusterSize*c_nbnxnCpuIClusterSize;
1207 if (jclusterFirst <= jclusterLast)
1209 for (int jcluster = jclusterFirst; jcluster <= jclusterLast; jcluster++)
1211 /* Store cj and the interaction mask */
1213 cjEntry.cj = jGrid.cell0 + jcluster;
1214 cjEntry.excl = get_imask(excludeSubDiagonal, icluster, jcluster);
1215 nbl->cj.push_back(cjEntry);
1217 /* Increase the closing index in the i list */
1218 nbl->ci.back().cj_ind_end = nbl->cj.size();
1222 #ifdef GMX_NBNXN_SIMD_4XN
1223 #include "gromacs/mdlib/nbnxn_search_simd_4xn.h"
1225 #ifdef GMX_NBNXN_SIMD_2XNN
1226 #include "gromacs/mdlib/nbnxn_search_simd_2xnn.h"
1229 /* Plain C or SIMD4 code for making a pair list of super-cell sci vs scj.
1230 * Checks bounding box distances and possibly atom pair distances.
1232 static void make_cluster_list_supersub(const nbnxn_grid_t &iGrid,
1233 const nbnxn_grid_t &jGrid,
1234 NbnxnPairlistGpu *nbl,
1237 const bool excludeSubDiagonal,
1242 int *numDistanceChecks)
1244 NbnxnPairlistGpuWork &work = *nbl->work;
1247 const float *pbb_ci = work.iSuperClusterData.bbPacked.data();
1249 const nbnxn_bb_t *bb_ci = work.iSuperClusterData.bb.data();
1252 assert(c_nbnxnGpuClusterSize == iGrid.na_c);
1253 assert(c_nbnxnGpuClusterSize == jGrid.na_c);
1255 /* We generate the pairlist mainly based on bounding-box distances
1256 * and do atom pair distance based pruning on the GPU.
1257 * Only if a j-group contains a single cluster-pair, we try to prune
1258 * that pair based on atom distances on the CPU to avoid empty j-groups.
1260 #define PRUNE_LIST_CPU_ONE 1
1261 #define PRUNE_LIST_CPU_ALL 0
1263 #if PRUNE_LIST_CPU_ONE
1267 float *d2l = work.distanceBuffer.data();
1269 for (int subc = 0; subc < jGrid.nsubc[scj]; subc++)
1271 const int cj4_ind = work.cj_ind/c_nbnxnGpuJgroupSize;
1272 const int cj_offset = work.cj_ind - cj4_ind*c_nbnxnGpuJgroupSize;
1273 const int cj = scj*c_gpuNumClusterPerCell + subc;
1275 const int cj_gl = jGrid.cell0*c_gpuNumClusterPerCell + cj;
1278 if (excludeSubDiagonal && sci == scj)
1284 ci1 = iGrid.nsubc[sci];
1288 /* Determine all ci1 bb distances in one call with SIMD4 */
1289 subc_bb_dist2_simd4_xxxx(jGrid.pbb.data() + (cj >> STRIDE_PBB_2LOG)*NNBSBB_XXXX + (cj & (STRIDE_PBB-1)),
1291 *numDistanceChecks += c_nbnxnGpuClusterSize*2;
1295 unsigned int imask = 0;
1296 /* We use a fixed upper-bound instead of ci1 to help optimization */
1297 for (int ci = 0; ci < c_gpuNumClusterPerCell; ci++)
1305 /* Determine the bb distance between ci and cj */
1306 d2l[ci] = subc_bb_dist2(ci, bb_ci, cj, jGrid.bb);
1307 *numDistanceChecks += 2;
1311 #if PRUNE_LIST_CPU_ALL
1312 /* Check if the distance is within the distance where
1313 * we use only the bounding box distance rbb,
1314 * or within the cut-off and there is at least one atom pair
1315 * within the cut-off. This check is very costly.
1317 *numDistanceChecks += c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize;
1320 clusterpair_in_range(work, ci, cj_gl, stride, x, rlist2)))
1322 /* Check if the distance between the two bounding boxes
1323 * in within the pair-list cut-off.
1328 /* Flag this i-subcell to be taken into account */
1329 imask |= (1U << (cj_offset*c_gpuNumClusterPerCell + ci));
1331 #if PRUNE_LIST_CPU_ONE
1339 #if PRUNE_LIST_CPU_ONE
1340 /* If we only found 1 pair, check if any atoms are actually
1341 * within the cut-off, so we could get rid of it.
1343 if (npair == 1 && d2l[ci_last] >= rbb2 &&
1344 !clusterpair_in_range(work, ci_last, cj_gl, stride, x, rlist2))
1346 imask &= ~(1U << (cj_offset*c_gpuNumClusterPerCell + ci_last));
1353 /* We have at least one cluster pair: add a j-entry */
1354 if (static_cast<size_t>(cj4_ind) == nbl->cj4.size())
1356 nbl->cj4.resize(nbl->cj4.size() + 1);
1358 nbnxn_cj4_t *cj4 = &nbl->cj4[cj4_ind];
1360 cj4->cj[cj_offset] = cj_gl;
1362 /* Set the exclusions for the ci==sj entry.
1363 * Here we don't bother to check if this entry is actually flagged,
1364 * as it will nearly always be in the list.
1366 if (excludeSubDiagonal && sci == scj)
1368 set_self_and_newton_excls_supersub(nbl, cj4_ind, cj_offset, subc);
1371 /* Copy the cluster interaction mask to the list */
1372 for (int w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
1374 cj4->imei[w].imask |= imask;
1377 nbl->work->cj_ind++;
1379 /* Keep the count */
1380 nbl->nci_tot += npair;
1382 /* Increase the closing index in i super-cell list */
1383 nbl->sci.back().cj4_ind_end =
1384 (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize;
1389 /* Returns how many contiguous j-clusters we have starting in the i-list */
1390 template <typename CjListType>
1391 static int numContiguousJClusters(const int cjIndexStart,
1392 const int cjIndexEnd,
1393 gmx::ArrayRef<const CjListType> cjList)
1395 const int firstJCluster = nblCj(cjList, cjIndexStart);
1397 int numContiguous = 0;
1399 while (cjIndexStart + numContiguous < cjIndexEnd &&
1400 nblCj(cjList, cjIndexStart + numContiguous) == firstJCluster + numContiguous)
1405 return numContiguous;
1409 * \brief Helper struct for efficient searching for excluded atoms in a j-list
1413 /*! \brief Constructs a j-list range from \p cjList with the given index range */
1414 template <typename CjListType>
1415 JListRanges(int cjIndexStart,
1417 gmx::ArrayRef<const CjListType> cjList);
1419 int cjIndexStart; //!< The start index in the j-list
1420 int cjIndexEnd; //!< The end index in the j-list
1421 int cjFirst; //!< The j-cluster with index cjIndexStart
1422 int cjLast; //!< The j-cluster with index cjIndexEnd-1
1423 int numDirect; //!< Up to cjIndexStart+numDirect the j-clusters are cjFirst + the index offset
1427 template <typename CjListType>
1428 JListRanges::JListRanges(int cjIndexStart,
1430 gmx::ArrayRef<const CjListType> cjList) :
1431 cjIndexStart(cjIndexStart),
1432 cjIndexEnd(cjIndexEnd)
1434 GMX_ASSERT(cjIndexEnd > cjIndexStart, "JListRanges should only be called with non-empty lists");
1436 cjFirst = nblCj(cjList, cjIndexStart);
1437 cjLast = nblCj(cjList, cjIndexEnd - 1);
1439 /* Determine how many contiguous j-cells we have starting
1440 * from the first i-cell. This number can be used to directly
1441 * calculate j-cell indices for excluded atoms.
1443 numDirect = numContiguousJClusters(cjIndexStart, cjIndexEnd, cjList);
1447 /* Return the index of \p jCluster in the given range or -1 when not present
1449 * Note: This code is executed very often and therefore performance is
1450 * important. It should be inlined and fully optimized.
1452 template <typename CjListType>
1454 findJClusterInJList(int jCluster,
1455 const JListRanges &ranges,
1456 gmx::ArrayRef<const CjListType> cjList)
1460 if (jCluster < ranges.cjFirst + ranges.numDirect)
1462 /* We can calculate the index directly using the offset */
1463 index = ranges.cjIndexStart + jCluster - ranges.cjFirst;
1467 /* Search for jCluster using bisection */
1469 int rangeStart = ranges.cjIndexStart + ranges.numDirect;
1470 int rangeEnd = ranges.cjIndexEnd;
1472 while (index == -1 && rangeStart < rangeEnd)
1474 rangeMiddle = (rangeStart + rangeEnd) >> 1;
1476 const int clusterMiddle = nblCj(cjList, rangeMiddle);
1478 if (jCluster == clusterMiddle)
1480 index = rangeMiddle;
1482 else if (jCluster < clusterMiddle)
1484 rangeEnd = rangeMiddle;
1488 rangeStart = rangeMiddle + 1;
1496 // TODO: Get rid of the two functions below by renaming sci to ci (or something better)
1498 /* Return the i-entry in the list we are currently operating on */
1499 static nbnxn_ci_t *getOpenIEntry(NbnxnPairlistCpu *nbl)
1501 return &nbl->ci.back();
1504 /* Return the i-entry in the list we are currently operating on */
1505 static nbnxn_sci_t *getOpenIEntry(NbnxnPairlistGpu *nbl)
1507 return &nbl->sci.back();
1510 /* Set all atom-pair exclusions for a simple type list i-entry
1512 * Set all atom-pair exclusions from the topology stored in exclusions
1513 * as masks in the pair-list for simple list entry iEntry.
1516 setExclusionsForIEntry(const nbnxn_search *nbs,
1517 NbnxnPairlistCpu *nbl,
1518 gmx_bool diagRemoved,
1520 const nbnxn_ci_t &iEntry,
1521 const t_blocka &exclusions)
1523 if (iEntry.cj_ind_end == iEntry.cj_ind_start)
1525 /* Empty list: no exclusions */
1529 const JListRanges ranges(iEntry.cj_ind_start, iEntry.cj_ind_end, gmx::makeConstArrayRef(nbl->cj));
1531 const int iCluster = iEntry.ci;
1533 gmx::ArrayRef<const int> cell = nbs->cell;
1535 /* Loop over the atoms in the i-cluster */
1536 for (int i = 0; i < nbl->na_ci; i++)
1538 const int iIndex = iCluster*nbl->na_ci + i;
1539 const int iAtom = nbs->a[iIndex];
1542 /* Loop over the topology-based exclusions for this i-atom */
1543 for (int exclIndex = exclusions.index[iAtom]; exclIndex < exclusions.index[iAtom + 1]; exclIndex++)
1545 const int jAtom = exclusions.a[exclIndex];
1549 /* The self exclusion are already set, save some time */
1553 /* Get the index of the j-atom in the nbnxn atom data */
1554 const int jIndex = cell[jAtom];
1556 /* Without shifts we only calculate interactions j>i
1557 * for one-way pair-lists.
1559 if (diagRemoved && jIndex <= iIndex)
1564 const int jCluster = (jIndex >> na_cj_2log);
1566 /* Could the cluster se be in our list? */
1567 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
1570 findJClusterInJList(jCluster, ranges,
1571 gmx::makeConstArrayRef(nbl->cj));
1575 /* We found an exclusion, clear the corresponding
1578 const int innerJ = jIndex - (jCluster << na_cj_2log);
1580 nbl->cj[index].excl &= ~(1U << ((i << na_cj_2log) + innerJ));
1588 /* Add a new i-entry to the FEP list and copy the i-properties */
1589 static inline void fep_list_new_nri_copy(t_nblist *nlist)
1591 /* Add a new i-entry */
1594 assert(nlist->nri < nlist->maxnri);
1596 /* Duplicate the last i-entry, except for jindex, which continues */
1597 nlist->iinr[nlist->nri] = nlist->iinr[nlist->nri-1];
1598 nlist->shift[nlist->nri] = nlist->shift[nlist->nri-1];
1599 nlist->gid[nlist->nri] = nlist->gid[nlist->nri-1];
1600 nlist->jindex[nlist->nri] = nlist->nrj;
1603 /* For load balancing of the free-energy lists over threads, we set
1604 * the maximum nrj size of an i-entry to 40. This leads to good
1605 * load balancing in the worst case scenario of a single perturbed
1606 * particle on 16 threads, while not introducing significant overhead.
1607 * Note that half of the perturbed pairs will anyhow end up in very small lists,
1608 * since non perturbed i-particles will see few perturbed j-particles).
1610 const int max_nrj_fep = 40;
1612 /* Exclude the perturbed pairs from the Verlet list. This is only done to avoid
1613 * singularities for overlapping particles (0/0), since the charges and
1614 * LJ parameters have been zeroed in the nbnxn data structure.
1615 * Simultaneously make a group pair list for the perturbed pairs.
1617 static void make_fep_list(const nbnxn_search *nbs,
1618 const nbnxn_atomdata_t *nbat,
1619 NbnxnPairlistCpu *nbl,
1620 gmx_bool bDiagRemoved,
1622 real gmx_unused shx,
1623 real gmx_unused shy,
1624 real gmx_unused shz,
1625 real gmx_unused rlist_fep2,
1626 const nbnxn_grid_t &iGrid,
1627 const nbnxn_grid_t &jGrid,
1630 int ci, cj_ind_start, cj_ind_end, cja, cjr;
1632 int ngid, gid_i = 0, gid_j, gid;
1633 int egp_shift, egp_mask;
1635 int ind_i, ind_j, ai, aj;
1637 gmx_bool bFEP_i, bFEP_i_all;
1639 if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
1647 cj_ind_start = nbl_ci->cj_ind_start;
1648 cj_ind_end = nbl_ci->cj_ind_end;
1650 /* In worst case we have alternating energy groups
1651 * and create #atom-pair lists, which means we need the size
1652 * of a cluster pair (na_ci*na_cj) times the number of cj's.
1654 nri_max = nbl->na_ci*nbl->na_cj*(cj_ind_end - cj_ind_start);
1655 if (nlist->nri + nri_max > nlist->maxnri)
1657 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1658 reallocate_nblist(nlist);
1661 const nbnxn_atomdata_t::Params &nbatParams = nbat->params();
1663 ngid = nbatParams.nenergrp;
1665 if (ngid*jGrid.na_cj > gmx::index(sizeof(gid_cj)*8))
1667 gmx_fatal(FARGS, "The Verlet scheme with %dx%d kernels and free-energy only supports up to %zu energy groups",
1668 iGrid.na_c, jGrid.na_cj, (sizeof(gid_cj)*8)/jGrid.na_cj);
1671 egp_shift = nbatParams.neg_2log;
1672 egp_mask = (1 << egp_shift) - 1;
1674 /* Loop over the atoms in the i sub-cell */
1676 for (int i = 0; i < nbl->na_ci; i++)
1678 ind_i = ci*nbl->na_ci + i;
1683 nlist->jindex[nri+1] = nlist->jindex[nri];
1684 nlist->iinr[nri] = ai;
1685 /* The actual energy group pair index is set later */
1686 nlist->gid[nri] = 0;
1687 nlist->shift[nri] = nbl_ci->shift & NBNXN_CI_SHIFT;
1689 bFEP_i = ((iGrid.fep[ci - iGrid.cell0] & (1 << i)) != 0u);
1691 bFEP_i_all = bFEP_i_all && bFEP_i;
1693 if (nlist->nrj + (cj_ind_end - cj_ind_start)*nbl->na_cj > nlist->maxnrj)
1695 nlist->maxnrj = over_alloc_small(nlist->nrj + (cj_ind_end - cj_ind_start)*nbl->na_cj);
1696 srenew(nlist->jjnr, nlist->maxnrj);
1697 srenew(nlist->excl_fep, nlist->maxnrj);
1702 gid_i = (nbatParams.energrp[ci] >> (egp_shift*i)) & egp_mask;
1705 for (int cj_ind = cj_ind_start; cj_ind < cj_ind_end; cj_ind++)
1707 unsigned int fep_cj;
1709 cja = nbl->cj[cj_ind].cj;
1711 if (jGrid.na_cj == jGrid.na_c)
1713 cjr = cja - jGrid.cell0;
1714 fep_cj = jGrid.fep[cjr];
1717 gid_cj = nbatParams.energrp[cja];
1720 else if (2*jGrid.na_cj == jGrid.na_c)
1722 cjr = cja - jGrid.cell0*2;
1723 /* Extract half of the ci fep/energrp mask */
1724 fep_cj = (jGrid.fep[cjr>>1] >> ((cjr&1)*jGrid.na_cj)) & ((1<<jGrid.na_cj) - 1);
1727 gid_cj = nbatParams.energrp[cja>>1] >> ((cja&1)*jGrid.na_cj*egp_shift) & ((1<<(jGrid.na_cj*egp_shift)) - 1);
1732 cjr = cja - (jGrid.cell0>>1);
1733 /* Combine two ci fep masks/energrp */
1734 fep_cj = jGrid.fep[cjr*2] + (jGrid.fep[cjr*2+1] << jGrid.na_c);
1737 gid_cj = nbatParams.energrp[cja*2] + (nbatParams.energrp[cja*2+1] << (jGrid.na_c*egp_shift));
1741 if (bFEP_i || fep_cj != 0)
1743 for (int j = 0; j < nbl->na_cj; j++)
1745 /* Is this interaction perturbed and not excluded? */
1746 ind_j = cja*nbl->na_cj + j;
1749 (bFEP_i || (fep_cj & (1 << j))) &&
1750 (!bDiagRemoved || ind_j >= ind_i))
1754 gid_j = (gid_cj >> (j*egp_shift)) & egp_mask;
1755 gid = GID(gid_i, gid_j, ngid);
1757 if (nlist->nrj > nlist->jindex[nri] &&
1758 nlist->gid[nri] != gid)
1760 /* Energy group pair changed: new list */
1761 fep_list_new_nri_copy(nlist);
1764 nlist->gid[nri] = gid;
1767 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1769 fep_list_new_nri_copy(nlist);
1773 /* Add it to the FEP list */
1774 nlist->jjnr[nlist->nrj] = aj;
1775 nlist->excl_fep[nlist->nrj] = (nbl->cj[cj_ind].excl >> (i*nbl->na_cj + j)) & 1;
1778 /* Exclude it from the normal list.
1779 * Note that the charge has been set to zero,
1780 * but we need to avoid 0/0, as perturbed atoms
1781 * can be on top of each other.
1783 nbl->cj[cj_ind].excl &= ~(1U << (i*nbl->na_cj + j));
1789 if (nlist->nrj > nlist->jindex[nri])
1791 /* Actually add this new, non-empty, list */
1793 nlist->jindex[nlist->nri] = nlist->nrj;
1800 /* All interactions are perturbed, we can skip this entry */
1801 nbl_ci->cj_ind_end = cj_ind_start;
1802 nbl->ncjInUse -= cj_ind_end - cj_ind_start;
1806 /* Return the index of atom a within a cluster */
1807 static inline int cj_mod_cj4(int cj)
1809 return cj & (c_nbnxnGpuJgroupSize - 1);
1812 /* Convert a j-cluster to a cj4 group */
1813 static inline int cj_to_cj4(int cj)
1815 return cj/c_nbnxnGpuJgroupSize;
1818 /* Return the index of an j-atom within a warp */
1819 static inline int a_mod_wj(int a)
1821 return a & (c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit - 1);
1824 /* As make_fep_list above, but for super/sub lists. */
1825 static void make_fep_list(const nbnxn_search *nbs,
1826 const nbnxn_atomdata_t *nbat,
1827 NbnxnPairlistGpu *nbl,
1828 gmx_bool bDiagRemoved,
1829 const nbnxn_sci_t *nbl_sci,
1834 const nbnxn_grid_t &iGrid,
1835 const nbnxn_grid_t &jGrid,
1840 int ind_i, ind_j, ai, aj;
1844 const nbnxn_cj4_t *cj4;
1846 const int numJClusterGroups = nbl_sci->numJClusterGroups();
1847 if (numJClusterGroups == 0)
1853 const int sci = nbl_sci->sci;
1855 const int cj4_ind_start = nbl_sci->cj4_ind_start;
1856 const int cj4_ind_end = nbl_sci->cj4_ind_end;
1858 /* Here we process one super-cell, max #atoms na_sc, versus a list
1859 * cj4 entries, each with max c_nbnxnGpuJgroupSize cj's, each
1860 * of size na_cj atoms.
1861 * On the GPU we don't support energy groups (yet).
1862 * So for each of the na_sc i-atoms, we need max one FEP list
1863 * for each max_nrj_fep j-atoms.
1865 nri_max = nbl->na_sc*nbl->na_cj*(1 + (numJClusterGroups*c_nbnxnGpuJgroupSize)/max_nrj_fep);
1866 if (nlist->nri + nri_max > nlist->maxnri)
1868 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1869 reallocate_nblist(nlist);
1872 /* Loop over the atoms in the i super-cluster */
1873 for (int c = 0; c < c_gpuNumClusterPerCell; c++)
1875 c_abs = sci*c_gpuNumClusterPerCell + c;
1877 for (int i = 0; i < nbl->na_ci; i++)
1879 ind_i = c_abs*nbl->na_ci + i;
1884 nlist->jindex[nri+1] = nlist->jindex[nri];
1885 nlist->iinr[nri] = ai;
1886 /* With GPUs, energy groups are not supported */
1887 nlist->gid[nri] = 0;
1888 nlist->shift[nri] = nbl_sci->shift & NBNXN_CI_SHIFT;
1890 bFEP_i = ((iGrid.fep[c_abs - iGrid.cell0*c_gpuNumClusterPerCell] & (1 << i)) != 0u);
1892 xi = nbat->x()[ind_i*nbat->xstride+XX] + shx;
1893 yi = nbat->x()[ind_i*nbat->xstride+YY] + shy;
1894 zi = nbat->x()[ind_i*nbat->xstride+ZZ] + shz;
1896 const int nrjMax = nlist->nrj + numJClusterGroups*c_nbnxnGpuJgroupSize*nbl->na_cj;
1897 if (nrjMax > nlist->maxnrj)
1899 nlist->maxnrj = over_alloc_small(nrjMax);
1900 srenew(nlist->jjnr, nlist->maxnrj);
1901 srenew(nlist->excl_fep, nlist->maxnrj);
1904 for (int cj4_ind = cj4_ind_start; cj4_ind < cj4_ind_end; cj4_ind++)
1906 cj4 = &nbl->cj4[cj4_ind];
1908 for (int gcj = 0; gcj < c_nbnxnGpuJgroupSize; gcj++)
1910 unsigned int fep_cj;
1912 if ((cj4->imei[0].imask & (1U << (gcj*c_gpuNumClusterPerCell + c))) == 0)
1914 /* Skip this ci for this cj */
1919 cj4->cj[gcj] - jGrid.cell0*c_gpuNumClusterPerCell;
1921 fep_cj = jGrid.fep[cjr];
1923 if (bFEP_i || fep_cj != 0)
1925 for (int j = 0; j < nbl->na_cj; j++)
1927 /* Is this interaction perturbed and not excluded? */
1928 ind_j = (jGrid.cell0*c_gpuNumClusterPerCell + cjr)*nbl->na_cj + j;
1931 (bFEP_i || (fep_cj & (1 << j))) &&
1932 (!bDiagRemoved || ind_j >= ind_i))
1935 unsigned int excl_bit;
1938 const int jHalf = j/(c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit);
1939 nbnxn_excl_t *excl =
1940 get_exclusion_mask(nbl, cj4_ind, jHalf);
1942 excl_pair = a_mod_wj(j)*nbl->na_ci + i;
1943 excl_bit = (1U << (gcj*c_gpuNumClusterPerCell + c));
1945 dx = nbat->x()[ind_j*nbat->xstride+XX] - xi;
1946 dy = nbat->x()[ind_j*nbat->xstride+YY] - yi;
1947 dz = nbat->x()[ind_j*nbat->xstride+ZZ] - zi;
1949 /* The unpruned GPU list has more than 2/3
1950 * of the atom pairs beyond rlist. Using
1951 * this list will cause a lot of overhead
1952 * in the CPU FEP kernels, especially
1953 * relative to the fast GPU kernels.
1954 * So we prune the FEP list here.
1956 if (dx*dx + dy*dy + dz*dz < rlist_fep2)
1958 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1960 fep_list_new_nri_copy(nlist);
1964 /* Add it to the FEP list */
1965 nlist->jjnr[nlist->nrj] = aj;
1966 nlist->excl_fep[nlist->nrj] = (excl->pair[excl_pair] & excl_bit) ? 1 : 0;
1970 /* Exclude it from the normal list.
1971 * Note that the charge and LJ parameters have
1972 * been set to zero, but we need to avoid 0/0,
1973 * as perturbed atoms can be on top of each other.
1975 excl->pair[excl_pair] &= ~excl_bit;
1979 /* Note that we could mask out this pair in imask
1980 * if all i- and/or all j-particles are perturbed.
1981 * But since the perturbed pairs on the CPU will
1982 * take an order of magnitude more time, the GPU
1983 * will finish before the CPU and there is no gain.
1989 if (nlist->nrj > nlist->jindex[nri])
1991 /* Actually add this new, non-empty, list */
1993 nlist->jindex[nlist->nri] = nlist->nrj;
2000 /* Set all atom-pair exclusions for a GPU type list i-entry
2002 * Sets all atom-pair exclusions from the topology stored in exclusions
2003 * as masks in the pair-list for i-super-cluster list entry iEntry.
2006 setExclusionsForIEntry(const nbnxn_search *nbs,
2007 NbnxnPairlistGpu *nbl,
2008 gmx_bool diagRemoved,
2009 int gmx_unused na_cj_2log,
2010 const nbnxn_sci_t &iEntry,
2011 const t_blocka &exclusions)
2013 if (iEntry.numJClusterGroups() == 0)
2019 /* Set the search ranges using start and end j-cluster indices.
2020 * Note that here we can not use cj4_ind_end, since the last cj4
2021 * can be only partially filled, so we use cj_ind.
2023 const JListRanges ranges(iEntry.cj4_ind_start*c_nbnxnGpuJgroupSize,
2025 gmx::makeConstArrayRef(nbl->cj4));
2027 GMX_ASSERT(nbl->na_ci == c_nbnxnGpuClusterSize, "na_ci should match the GPU cluster size");
2028 constexpr int c_clusterSize = c_nbnxnGpuClusterSize;
2029 constexpr int c_superClusterSize = c_nbnxnGpuNumClusterPerSupercluster*c_nbnxnGpuClusterSize;
2031 const int iSuperCluster = iEntry.sci;
2033 gmx::ArrayRef<const int> cell = nbs->cell;
2035 /* Loop over the atoms in the i super-cluster */
2036 for (int i = 0; i < c_superClusterSize; i++)
2038 const int iIndex = iSuperCluster*c_superClusterSize + i;
2039 const int iAtom = nbs->a[iIndex];
2042 const int iCluster = i/c_clusterSize;
2044 /* Loop over the topology-based exclusions for this i-atom */
2045 for (int exclIndex = exclusions.index[iAtom]; exclIndex < exclusions.index[iAtom + 1]; exclIndex++)
2047 const int jAtom = exclusions.a[exclIndex];
2051 /* The self exclusions are already set, save some time */
2055 /* Get the index of the j-atom in the nbnxn atom data */
2056 const int jIndex = cell[jAtom];
2058 /* Without shifts we only calculate interactions j>i
2059 * for one-way pair-lists.
2061 /* NOTE: We would like to use iIndex on the right hand side,
2062 * but that makes this routine 25% slower with gcc6/7.
2063 * Even using c_superClusterSize makes it slower.
2064 * Either of these changes triggers peeling of the exclIndex
2065 * loop, which apparently leads to far less efficient code.
2067 if (diagRemoved && jIndex <= iSuperCluster*nbl->na_sc + i)
2072 const int jCluster = jIndex/c_clusterSize;
2074 /* Check whether the cluster is in our list? */
2075 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
2078 findJClusterInJList(jCluster, ranges,
2079 gmx::makeConstArrayRef(nbl->cj4));
2083 /* We found an exclusion, clear the corresponding
2086 const unsigned int pairMask = (1U << (cj_mod_cj4(index)*c_gpuNumClusterPerCell + iCluster));
2087 /* Check if the i-cluster interacts with the j-cluster */
2088 if (nbl_imask0(nbl, index) & pairMask)
2090 const int innerI = (i & (c_clusterSize - 1));
2091 const int innerJ = (jIndex & (c_clusterSize - 1));
2093 /* Determine which j-half (CUDA warp) we are in */
2094 const int jHalf = innerJ/(c_clusterSize/c_nbnxnGpuClusterpairSplit);
2096 nbnxn_excl_t *interactionMask =
2097 get_exclusion_mask(nbl, cj_to_cj4(index), jHalf);
2099 interactionMask->pair[a_mod_wj(innerJ)*c_clusterSize + innerI] &= ~pairMask;
2108 /* Make a new ci entry at the back of nbl->ci */
2109 static void addNewIEntry(NbnxnPairlistCpu *nbl, int ci, int shift, int flags)
2113 ciEntry.shift = shift;
2114 /* Store the interaction flags along with the shift */
2115 ciEntry.shift |= flags;
2116 ciEntry.cj_ind_start = nbl->cj.size();
2117 ciEntry.cj_ind_end = nbl->cj.size();
2118 nbl->ci.push_back(ciEntry);
2121 /* Make a new sci entry at index nbl->nsci */
2122 static void addNewIEntry(NbnxnPairlistGpu *nbl, int sci, int shift, int gmx_unused flags)
2124 nbnxn_sci_t sciEntry;
2126 sciEntry.shift = shift;
2127 sciEntry.cj4_ind_start = nbl->cj4.size();
2128 sciEntry.cj4_ind_end = nbl->cj4.size();
2130 nbl->sci.push_back(sciEntry);
2133 /* Sort the simple j-list cj on exclusions.
2134 * Entries with exclusions will all be sorted to the beginning of the list.
2136 static void sort_cj_excl(nbnxn_cj_t *cj, int ncj,
2137 NbnxnPairlistCpuWork *work)
2139 work->cj.resize(ncj);
2141 /* Make a list of the j-cells involving exclusions */
2143 for (int j = 0; j < ncj; j++)
2145 if (cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
2147 work->cj[jnew++] = cj[j];
2150 /* Check if there are exclusions at all or not just the first entry */
2151 if (!((jnew == 0) ||
2152 (jnew == 1 && cj[0].excl != NBNXN_INTERACTION_MASK_ALL)))
2154 for (int j = 0; j < ncj; j++)
2156 if (cj[j].excl == NBNXN_INTERACTION_MASK_ALL)
2158 work->cj[jnew++] = cj[j];
2161 for (int j = 0; j < ncj; j++)
2163 cj[j] = work->cj[j];
2168 /* Close this simple list i entry */
2169 static void closeIEntry(NbnxnPairlistCpu *nbl,
2170 int gmx_unused sp_max_av,
2171 gmx_bool gmx_unused progBal,
2172 float gmx_unused nsp_tot_est,
2173 int gmx_unused thread,
2174 int gmx_unused nthread)
2176 nbnxn_ci_t &ciEntry = nbl->ci.back();
2178 /* All content of the new ci entry have already been filled correctly,
2179 * we only need to sort and increase counts or remove the entry when empty.
2181 const int jlen = ciEntry.cj_ind_end - ciEntry.cj_ind_start;
2184 sort_cj_excl(nbl->cj.data() + ciEntry.cj_ind_start, jlen, nbl->work);
2186 /* The counts below are used for non-bonded pair/flop counts
2187 * and should therefore match the available kernel setups.
2189 if (!(ciEntry.shift & NBNXN_CI_DO_COUL(0)))
2191 nbl->work->ncj_noq += jlen;
2193 else if ((ciEntry.shift & NBNXN_CI_HALF_LJ(0)) ||
2194 !(ciEntry.shift & NBNXN_CI_DO_LJ(0)))
2196 nbl->work->ncj_hlj += jlen;
2201 /* Entry is empty: remove it */
2206 /* Split sci entry for load balancing on the GPU.
2207 * Splitting ensures we have enough lists to fully utilize the whole GPU.
2208 * With progBal we generate progressively smaller lists, which improves
2209 * load balancing. As we only know the current count on our own thread,
2210 * we will need to estimate the current total amount of i-entries.
2211 * As the lists get concatenated later, this estimate depends
2212 * both on nthread and our own thread index.
2214 static void split_sci_entry(NbnxnPairlistGpu *nbl,
2216 gmx_bool progBal, float nsp_tot_est,
2217 int thread, int nthread)
2225 /* Estimate the total numbers of ci's of the nblist combined
2226 * over all threads using the target number of ci's.
2228 nsp_est = (nsp_tot_est*thread)/nthread + nbl->nci_tot;
2230 /* The first ci blocks should be larger, to avoid overhead.
2231 * The last ci blocks should be smaller, to improve load balancing.
2232 * The factor 3/2 makes the first block 3/2 times the target average
2233 * and ensures that the total number of blocks end up equal to
2234 * that of equally sized blocks of size nsp_target_av.
2236 nsp_max = static_cast<int>(nsp_target_av*(nsp_tot_est*1.5/(nsp_est + nsp_tot_est)));
2240 nsp_max = nsp_target_av;
2243 const int cj4_start = nbl->sci.back().cj4_ind_start;
2244 const int cj4_end = nbl->sci.back().cj4_ind_end;
2245 const int j4len = cj4_end - cj4_start;
2247 if (j4len > 1 && j4len*c_gpuNumClusterPerCell*c_nbnxnGpuJgroupSize > nsp_max)
2249 /* Modify the last ci entry and process the cj4's again */
2255 for (int cj4 = cj4_start; cj4 < cj4_end; cj4++)
2257 int nsp_cj4_p = nsp_cj4;
2258 /* Count the number of cluster pairs in this cj4 group */
2260 for (int p = 0; p < c_gpuNumClusterPerCell*c_nbnxnGpuJgroupSize; p++)
2262 nsp_cj4 += (nbl->cj4[cj4].imei[0].imask >> p) & 1;
2265 /* If adding the current cj4 with nsp_cj4 pairs get us further
2266 * away from our target nsp_max, split the list before this cj4.
2268 if (nsp > 0 && nsp_max - nsp < nsp + nsp_cj4 - nsp_max)
2270 /* Split the list at cj4 */
2271 nbl->sci.back().cj4_ind_end = cj4;
2272 /* Create a new sci entry */
2274 sciNew.sci = nbl->sci.back().sci;
2275 sciNew.shift = nbl->sci.back().shift;
2276 sciNew.cj4_ind_start = cj4;
2277 nbl->sci.push_back(sciNew);
2280 nsp_cj4_e = nsp_cj4_p;
2286 /* Put the remaining cj4's in the last sci entry */
2287 nbl->sci.back().cj4_ind_end = cj4_end;
2289 /* Possibly balance out the last two sci's
2290 * by moving the last cj4 of the second last sci.
2292 if (nsp_sci - nsp_cj4_e >= nsp + nsp_cj4_e)
2294 GMX_ASSERT(nbl->sci.size() >= 2, "We expect at least two elements");
2295 nbl->sci[nbl->sci.size() - 2].cj4_ind_end--;
2296 nbl->sci[nbl->sci.size() - 1].cj4_ind_start--;
2301 /* Clost this super/sub list i entry */
2302 static void closeIEntry(NbnxnPairlistGpu *nbl,
2304 gmx_bool progBal, float nsp_tot_est,
2305 int thread, int nthread)
2307 nbnxn_sci_t &sciEntry = *getOpenIEntry(nbl);
2309 /* All content of the new ci entry have already been filled correctly,
2310 * we only need to, potentially, split or remove the entry when empty.
2312 int j4len = sciEntry.numJClusterGroups();
2315 /* We can only have complete blocks of 4 j-entries in a list,
2316 * so round the count up before closing.
2318 int ncj4 = (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize;
2319 nbl->work->cj_ind = ncj4*c_nbnxnGpuJgroupSize;
2323 /* Measure the size of the new entry and potentially split it */
2324 split_sci_entry(nbl, nsp_max_av, progBal, nsp_tot_est,
2330 /* Entry is empty: remove it */
2331 nbl->sci.pop_back();
2335 /* Syncs the working array before adding another grid pair to the GPU list */
2336 static void sync_work(NbnxnPairlistCpu gmx_unused *nbl)
2340 /* Syncs the working array before adding another grid pair to the GPU list */
2341 static void sync_work(NbnxnPairlistGpu *nbl)
2343 nbl->work->cj_ind = nbl->cj4.size()*c_nbnxnGpuJgroupSize;
2346 /* Clears an NbnxnPairlistCpu data structure */
2347 static void clear_pairlist(NbnxnPairlistCpu *nbl)
2353 nbl->ciOuter.clear();
2354 nbl->cjOuter.clear();
2356 nbl->work->ncj_noq = 0;
2357 nbl->work->ncj_hlj = 0;
2360 /* Clears an NbnxnPairlistGpu data structure */
2361 static void clear_pairlist(NbnxnPairlistGpu *nbl)
2365 nbl->excl.resize(1);
2369 /* Clears a group scheme pair list */
2370 static void clear_pairlist_fep(t_nblist *nl)
2374 if (nl->jindex == nullptr)
2376 snew(nl->jindex, 1);
2381 /* Sets a simple list i-cell bounding box, including PBC shift */
2382 static inline void set_icell_bb_simple(gmx::ArrayRef<const nbnxn_bb_t> bb,
2384 real shx, real shy, real shz,
2387 bb_ci->lower[BB_X] = bb[ci].lower[BB_X] + shx;
2388 bb_ci->lower[BB_Y] = bb[ci].lower[BB_Y] + shy;
2389 bb_ci->lower[BB_Z] = bb[ci].lower[BB_Z] + shz;
2390 bb_ci->upper[BB_X] = bb[ci].upper[BB_X] + shx;
2391 bb_ci->upper[BB_Y] = bb[ci].upper[BB_Y] + shy;
2392 bb_ci->upper[BB_Z] = bb[ci].upper[BB_Z] + shz;
2395 /* Sets a simple list i-cell bounding box, including PBC shift */
2396 static inline void set_icell_bb(const nbnxn_grid_t &iGrid,
2398 real shx, real shy, real shz,
2399 NbnxnPairlistCpuWork *work)
2401 set_icell_bb_simple(iGrid.bb, ci, shx, shy, shz, &work->iClusterData.bb[0]);
2405 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2406 static void set_icell_bbxxxx_supersub(gmx::ArrayRef<const float> bb,
2408 real shx, real shy, real shz,
2411 int ia = ci*(c_gpuNumClusterPerCell >> STRIDE_PBB_2LOG)*NNBSBB_XXXX;
2412 for (int m = 0; m < (c_gpuNumClusterPerCell >> STRIDE_PBB_2LOG)*NNBSBB_XXXX; m += NNBSBB_XXXX)
2414 for (int i = 0; i < STRIDE_PBB; i++)
2416 bb_ci[m+0*STRIDE_PBB+i] = bb[ia+m+0*STRIDE_PBB+i] + shx;
2417 bb_ci[m+1*STRIDE_PBB+i] = bb[ia+m+1*STRIDE_PBB+i] + shy;
2418 bb_ci[m+2*STRIDE_PBB+i] = bb[ia+m+2*STRIDE_PBB+i] + shz;
2419 bb_ci[m+3*STRIDE_PBB+i] = bb[ia+m+3*STRIDE_PBB+i] + shx;
2420 bb_ci[m+4*STRIDE_PBB+i] = bb[ia+m+4*STRIDE_PBB+i] + shy;
2421 bb_ci[m+5*STRIDE_PBB+i] = bb[ia+m+5*STRIDE_PBB+i] + shz;
2427 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2428 gmx_unused static void set_icell_bb_supersub(gmx::ArrayRef<const nbnxn_bb_t> bb,
2430 real shx, real shy, real shz,
2433 for (int i = 0; i < c_gpuNumClusterPerCell; i++)
2435 set_icell_bb_simple(bb, ci*c_gpuNumClusterPerCell+i,
2441 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2442 gmx_unused static void set_icell_bb(const nbnxn_grid_t &iGrid,
2444 real shx, real shy, real shz,
2445 NbnxnPairlistGpuWork *work)
2448 set_icell_bbxxxx_supersub(iGrid.pbb, ci, shx, shy, shz,
2449 work->iSuperClusterData.bbPacked.data());
2451 set_icell_bb_supersub(iGrid.bb, ci, shx, shy, shz,
2452 work->iSuperClusterData.bb.data());
2456 /* Copies PBC shifted i-cell atom coordinates x,y,z to working array */
2457 static void icell_set_x_simple(int ci,
2458 real shx, real shy, real shz,
2459 int stride, const real *x,
2460 NbnxnPairlistCpuWork::IClusterData *iClusterData)
2462 const int ia = ci*c_nbnxnCpuIClusterSize;
2464 for (int i = 0; i < c_nbnxnCpuIClusterSize; i++)
2466 iClusterData->x[i*STRIDE_XYZ+XX] = x[(ia+i)*stride+XX] + shx;
2467 iClusterData->x[i*STRIDE_XYZ+YY] = x[(ia+i)*stride+YY] + shy;
2468 iClusterData->x[i*STRIDE_XYZ+ZZ] = x[(ia+i)*stride+ZZ] + shz;
2472 static void icell_set_x(int ci,
2473 real shx, real shy, real shz,
2474 int stride, const real *x,
2476 NbnxnPairlistCpuWork *work)
2478 switch (nb_kernel_type)
2481 #ifdef GMX_NBNXN_SIMD_4XN
2482 case nbnxnk4xN_SIMD_4xN:
2483 icell_set_x_simd_4xn(ci, shx, shy, shz, stride, x, work);
2486 #ifdef GMX_NBNXN_SIMD_2XNN
2487 case nbnxnk4xN_SIMD_2xNN:
2488 icell_set_x_simd_2xnn(ci, shx, shy, shz, stride, x, work);
2492 case nbnxnk4x4_PlainC:
2493 icell_set_x_simple(ci, shx, shy, shz, stride, x, &work->iClusterData);
2496 GMX_ASSERT(false, "Unhandled case");
2501 /* Copies PBC shifted super-cell atom coordinates x,y,z to working array */
2502 static void icell_set_x(int ci,
2503 real shx, real shy, real shz,
2504 int stride, const real *x,
2505 int gmx_unused nb_kernel_type,
2506 NbnxnPairlistGpuWork *work)
2508 #if !GMX_SIMD4_HAVE_REAL
2510 real * x_ci = work->iSuperClusterData.x.data();
2512 int ia = ci*c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize;
2513 for (int i = 0; i < c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize; i++)
2515 x_ci[i*DIM + XX] = x[(ia+i)*stride + XX] + shx;
2516 x_ci[i*DIM + YY] = x[(ia+i)*stride + YY] + shy;
2517 x_ci[i*DIM + ZZ] = x[(ia+i)*stride + ZZ] + shz;
2520 #else /* !GMX_SIMD4_HAVE_REAL */
2522 real * x_ci = work->iSuperClusterData.xSimd.data();
2524 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2526 for (int i = 0; i < c_nbnxnGpuClusterSize; i += GMX_SIMD4_WIDTH)
2528 int io = si*c_nbnxnGpuClusterSize + i;
2529 int ia = ci*c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize + io;
2530 for (int j = 0; j < GMX_SIMD4_WIDTH; j++)
2532 x_ci[io*DIM + j + XX*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + XX] + shx;
2533 x_ci[io*DIM + j + YY*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + YY] + shy;
2534 x_ci[io*DIM + j + ZZ*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + ZZ] + shz;
2539 #endif /* !GMX_SIMD4_HAVE_REAL */
2542 static real minimum_subgrid_size_xy(const nbnxn_grid_t &grid)
2546 return std::min(grid.cellSize[XX], grid.cellSize[YY]);
2550 return std::min(grid.cellSize[XX]/c_gpuNumClusterPerCellX,
2551 grid.cellSize[YY]/c_gpuNumClusterPerCellY);
2555 static real effective_buffer_1x1_vs_MxN(const nbnxn_grid_t &iGrid,
2556 const nbnxn_grid_t &jGrid)
2558 const real eff_1x1_buffer_fac_overest = 0.1;
2560 /* Determine an atom-pair list cut-off buffer size for atom pairs,
2561 * to be added to rlist (including buffer) used for MxN.
2562 * This is for converting an MxN list to a 1x1 list. This means we can't
2563 * use the normal buffer estimate, as we have an MxN list in which
2564 * some atom pairs beyond rlist are missing. We want to capture
2565 * the beneficial effect of buffering by extra pairs just outside rlist,
2566 * while removing the useless pairs that are further away from rlist.
2567 * (Also the buffer could have been set manually not using the estimate.)
2568 * This buffer size is an overestimate.
2569 * We add 10% of the smallest grid sub-cell dimensions.
2570 * Note that the z-size differs per cell and we don't use this,
2571 * so we overestimate.
2572 * With PME, the 10% value gives a buffer that is somewhat larger
2573 * than the effective buffer with a tolerance of 0.005 kJ/mol/ps.
2574 * Smaller tolerances or using RF lead to a smaller effective buffer,
2575 * so 10% gives a safe overestimate.
2577 return eff_1x1_buffer_fac_overest*(minimum_subgrid_size_xy(iGrid) +
2578 minimum_subgrid_size_xy(jGrid));
2581 /* Clusters at the cut-off only increase rlist by 60% of their size */
2582 static real nbnxn_rlist_inc_outside_fac = 0.6;
2584 /* Due to the cluster size the effective pair-list is longer than
2585 * that of a simple atom pair-list. This function gives the extra distance.
2587 real nbnxn_get_rlist_effective_inc(int cluster_size_j, real atom_density)
2590 real vol_inc_i, vol_inc_j;
2592 /* We should get this from the setup, but currently it's the same for
2593 * all setups, including GPUs.
2595 cluster_size_i = c_nbnxnCpuIClusterSize;
2597 vol_inc_i = (cluster_size_i - 1)/atom_density;
2598 vol_inc_j = (cluster_size_j - 1)/atom_density;
2600 return nbnxn_rlist_inc_outside_fac*std::cbrt(vol_inc_i + vol_inc_j);
2603 /* Estimates the interaction volume^2 for non-local interactions */
2604 static real nonlocal_vol2(const struct gmx_domdec_zones_t *zones, const rvec ls, real r)
2612 /* Here we simply add up the volumes of 1, 2 or 3 1D decomposition
2613 * not home interaction volume^2. As these volumes are not additive,
2614 * this is an overestimate, but it would only be significant in the limit
2615 * of small cells, where we anyhow need to split the lists into
2616 * as small parts as possible.
2619 for (int z = 0; z < zones->n; z++)
2621 if (zones->shift[z][XX] + zones->shift[z][YY] + zones->shift[z][ZZ] == 1)
2626 for (int d = 0; d < DIM; d++)
2628 if (zones->shift[z][d] == 0)
2632 za *= zones->size[z].x1[d] - zones->size[z].x0[d];
2636 /* 4 octants of a sphere */
2637 vold_est = 0.25*M_PI*r*r*r*r;
2638 /* 4 quarter pie slices on the edges */
2639 vold_est += 4*cl*M_PI/6.0*r*r*r;
2640 /* One rectangular volume on a face */
2641 vold_est += ca*0.5*r*r;
2643 vol2_est_tot += vold_est*za;
2647 return vol2_est_tot;
2650 /* Estimates the average size of a full j-list for super/sub setup */
2651 static void get_nsubpair_target(const nbnxn_search *nbs,
2654 int min_ci_balanced,
2655 int *nsubpair_target,
2656 float *nsubpair_tot_est)
2658 /* The target value of 36 seems to be the optimum for Kepler.
2659 * Maxwell is less sensitive to the exact value.
2661 const int nsubpair_target_min = 36;
2663 real r_eff_sup, vol_est, nsp_est, nsp_est_nl;
2665 const nbnxn_grid_t &grid = nbs->grid[0];
2667 /* We don't need to balance list sizes if:
2668 * - We didn't request balancing.
2669 * - The number of grid cells >= the number of lists requested,
2670 * since we will always generate at least #cells lists.
2671 * - We don't have any cells, since then there won't be any lists.
2673 if (min_ci_balanced <= 0 || grid.nc >= min_ci_balanced || grid.nc == 0)
2675 /* nsubpair_target==0 signals no balancing */
2676 *nsubpair_target = 0;
2677 *nsubpair_tot_est = 0;
2682 ls[XX] = (grid.c1[XX] - grid.c0[XX])/(grid.numCells[XX]*c_gpuNumClusterPerCellX);
2683 ls[YY] = (grid.c1[YY] - grid.c0[YY])/(grid.numCells[YY]*c_gpuNumClusterPerCellY);
2684 ls[ZZ] = grid.na_c/(grid.atom_density*ls[XX]*ls[YY]);
2686 /* The average length of the diagonal of a sub cell */
2687 real diagonal = std::sqrt(ls[XX]*ls[XX] + ls[YY]*ls[YY] + ls[ZZ]*ls[ZZ]);
2689 /* The formulas below are a heuristic estimate of the average nsj per si*/
2690 r_eff_sup = rlist + nbnxn_rlist_inc_outside_fac*gmx::square((grid.na_c - 1.0)/grid.na_c)*0.5*diagonal;
2692 if (!nbs->DomDec || nbs->zones->n == 1)
2699 gmx::square(grid.atom_density/grid.na_c)*
2700 nonlocal_vol2(nbs->zones, ls, r_eff_sup);
2705 /* Sub-cell interacts with itself */
2706 vol_est = ls[XX]*ls[YY]*ls[ZZ];
2707 /* 6/2 rectangular volume on the faces */
2708 vol_est += (ls[XX]*ls[YY] + ls[XX]*ls[ZZ] + ls[YY]*ls[ZZ])*r_eff_sup;
2709 /* 12/2 quarter pie slices on the edges */
2710 vol_est += 2*(ls[XX] + ls[YY] + ls[ZZ])*0.25*M_PI*gmx::square(r_eff_sup);
2711 /* 4 octants of a sphere */
2712 vol_est += 0.5*4.0/3.0*M_PI*gmx::power3(r_eff_sup);
2714 /* Estimate the number of cluster pairs as the local number of
2715 * clusters times the volume they interact with times the density.
2717 nsp_est = grid.nsubc_tot*vol_est*grid.atom_density/grid.na_c;
2719 /* Subtract the non-local pair count */
2720 nsp_est -= nsp_est_nl;
2722 /* For small cut-offs nsp_est will be an underesimate.
2723 * With DD nsp_est_nl is an overestimate so nsp_est can get negative.
2724 * So to avoid too small or negative nsp_est we set a minimum of
2725 * all cells interacting with all 3^3 direct neighbors (3^3-1)/2+1=14.
2726 * This might be a slight overestimate for small non-periodic groups of
2727 * atoms as will occur for a local domain with DD, but for small
2728 * groups of atoms we'll anyhow be limited by nsubpair_target_min,
2729 * so this overestimation will not matter.
2731 nsp_est = std::max(nsp_est, grid.nsubc_tot*14._real);
2735 fprintf(debug, "nsp_est local %5.1f non-local %5.1f\n",
2736 nsp_est, nsp_est_nl);
2741 nsp_est = nsp_est_nl;
2744 /* Thus the (average) maximum j-list size should be as follows.
2745 * Since there is overhead, we shouldn't make the lists too small
2746 * (and we can't chop up j-groups) so we use a minimum target size of 36.
2748 *nsubpair_target = std::max(nsubpair_target_min,
2749 roundToInt(nsp_est/min_ci_balanced));
2750 *nsubpair_tot_est = static_cast<int>(nsp_est);
2754 fprintf(debug, "nbl nsp estimate %.1f, nsubpair_target %d\n",
2755 nsp_est, *nsubpair_target);
2759 /* Debug list print function */
2760 static void print_nblist_ci_cj(FILE *fp, const NbnxnPairlistCpu *nbl)
2762 for (const nbnxn_ci_t &ciEntry : nbl->ci)
2764 fprintf(fp, "ci %4d shift %2d ncj %3d\n",
2765 ciEntry.ci, ciEntry.shift,
2766 ciEntry.cj_ind_end - ciEntry.cj_ind_start);
2768 for (int j = ciEntry.cj_ind_start; j < ciEntry.cj_ind_end; j++)
2770 fprintf(fp, " cj %5d imask %x\n",
2777 /* Debug list print function */
2778 static void print_nblist_sci_cj(FILE *fp, const NbnxnPairlistGpu *nbl)
2780 for (const nbnxn_sci_t &sci : nbl->sci)
2782 fprintf(fp, "ci %4d shift %2d ncj4 %2d\n",
2784 sci.numJClusterGroups());
2787 for (int j4 = sci.cj4_ind_start; j4 < sci.cj4_ind_end; j4++)
2789 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
2791 fprintf(fp, " sj %5d imask %x\n",
2793 nbl->cj4[j4].imei[0].imask);
2794 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2796 if (nbl->cj4[j4].imei[0].imask & (1U << (j*c_gpuNumClusterPerCell + si)))
2803 fprintf(fp, "ci %4d shift %2d ncj4 %2d ncp %3d\n",
2805 sci.numJClusterGroups(),
2810 /* Combine pair lists *nbl generated on multiple threads nblc */
2811 static void combine_nblists(int nnbl, NbnxnPairlistGpu **nbl,
2812 NbnxnPairlistGpu *nblc)
2814 int nsci = nblc->sci.size();
2815 int ncj4 = nblc->cj4.size();
2816 int nexcl = nblc->excl.size();
2817 for (int i = 0; i < nnbl; i++)
2819 nsci += nbl[i]->sci.size();
2820 ncj4 += nbl[i]->cj4.size();
2821 nexcl += nbl[i]->excl.size();
2824 /* Resize with the final, combined size, so we can fill in parallel */
2825 /* NOTE: For better performance we should use default initialization */
2826 nblc->sci.resize(nsci);
2827 nblc->cj4.resize(ncj4);
2828 nblc->excl.resize(nexcl);
2830 /* Each thread should copy its own data to the combined arrays,
2831 * as otherwise data will go back and forth between different caches.
2833 #if GMX_OPENMP && !(defined __clang_analyzer__)
2834 int nthreads = gmx_omp_nthreads_get(emntPairsearch);
2837 #pragma omp parallel for num_threads(nthreads) schedule(static)
2838 for (int n = 0; n < nnbl; n++)
2842 /* Determine the offset in the combined data for our thread.
2843 * Note that the original sizes in nblc are lost.
2845 int sci_offset = nsci;
2846 int cj4_offset = ncj4;
2847 int excl_offset = nexcl;
2849 for (int i = n; i < nnbl; i++)
2851 sci_offset -= nbl[i]->sci.size();
2852 cj4_offset -= nbl[i]->cj4.size();
2853 excl_offset -= nbl[i]->excl.size();
2856 const NbnxnPairlistGpu &nbli = *nbl[n];
2858 for (size_t i = 0; i < nbli.sci.size(); i++)
2860 nblc->sci[sci_offset + i] = nbli.sci[i];
2861 nblc->sci[sci_offset + i].cj4_ind_start += cj4_offset;
2862 nblc->sci[sci_offset + i].cj4_ind_end += cj4_offset;
2865 for (size_t j4 = 0; j4 < nbli.cj4.size(); j4++)
2867 nblc->cj4[cj4_offset + j4] = nbli.cj4[j4];
2868 nblc->cj4[cj4_offset + j4].imei[0].excl_ind += excl_offset;
2869 nblc->cj4[cj4_offset + j4].imei[1].excl_ind += excl_offset;
2872 for (size_t j4 = 0; j4 < nbli.excl.size(); j4++)
2874 nblc->excl[excl_offset + j4] = nbli.excl[j4];
2877 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2880 for (int n = 0; n < nnbl; n++)
2882 nblc->nci_tot += nbl[n]->nci_tot;
2886 static void balance_fep_lists(const nbnxn_search *nbs,
2887 nbnxn_pairlist_set_t *nbl_lists)
2890 int nri_tot, nrj_tot, nrj_target;
2894 nnbl = nbl_lists->nnbl;
2898 /* Nothing to balance */
2902 /* Count the total i-lists and pairs */
2905 for (int th = 0; th < nnbl; th++)
2907 nri_tot += nbl_lists->nbl_fep[th]->nri;
2908 nrj_tot += nbl_lists->nbl_fep[th]->nrj;
2911 nrj_target = (nrj_tot + nnbl - 1)/nnbl;
2913 assert(gmx_omp_nthreads_get(emntNonbonded) == nnbl);
2915 #pragma omp parallel for schedule(static) num_threads(nnbl)
2916 for (int th = 0; th < nnbl; th++)
2920 t_nblist *nbl = nbs->work[th].nbl_fep.get();
2922 /* Note that here we allocate for the total size, instead of
2923 * a per-thread esimate (which is hard to obtain).
2925 if (nri_tot > nbl->maxnri)
2927 nbl->maxnri = over_alloc_large(nri_tot);
2928 reallocate_nblist(nbl);
2930 if (nri_tot > nbl->maxnri || nrj_tot > nbl->maxnrj)
2932 nbl->maxnrj = over_alloc_small(nrj_tot);
2933 srenew(nbl->jjnr, nbl->maxnrj);
2934 srenew(nbl->excl_fep, nbl->maxnrj);
2937 clear_pairlist_fep(nbl);
2939 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2942 /* Loop over the source lists and assign and copy i-entries */
2944 nbld = nbs->work[th_dest].nbl_fep.get();
2945 for (int th = 0; th < nnbl; th++)
2949 nbls = nbl_lists->nbl_fep[th];
2951 for (int i = 0; i < nbls->nri; i++)
2955 /* The number of pairs in this i-entry */
2956 nrj = nbls->jindex[i+1] - nbls->jindex[i];
2958 /* Decide if list th_dest is too large and we should procede
2959 * to the next destination list.
2961 if (th_dest+1 < nnbl && nbld->nrj > 0 &&
2962 nbld->nrj + nrj - nrj_target > nrj_target - nbld->nrj)
2965 nbld = nbs->work[th_dest].nbl_fep.get();
2968 nbld->iinr[nbld->nri] = nbls->iinr[i];
2969 nbld->gid[nbld->nri] = nbls->gid[i];
2970 nbld->shift[nbld->nri] = nbls->shift[i];
2972 for (int j = nbls->jindex[i]; j < nbls->jindex[i+1]; j++)
2974 nbld->jjnr[nbld->nrj] = nbls->jjnr[j];
2975 nbld->excl_fep[nbld->nrj] = nbls->excl_fep[j];
2979 nbld->jindex[nbld->nri] = nbld->nrj;
2983 /* Swap the list pointers */
2984 for (int th = 0; th < nnbl; th++)
2986 t_nblist *nbl_tmp = nbs->work[th].nbl_fep.release();
2987 nbs->work[th].nbl_fep.reset(nbl_lists->nbl_fep[th]);
2988 nbl_lists->nbl_fep[th] = nbl_tmp;
2992 fprintf(debug, "nbl_fep[%d] nri %4d nrj %4d\n",
2994 nbl_lists->nbl_fep[th]->nri,
2995 nbl_lists->nbl_fep[th]->nrj);
3000 /* Returns the next ci to be processes by our thread */
3001 static gmx_bool next_ci(const nbnxn_grid_t &grid,
3002 int nth, int ci_block,
3003 int *ci_x, int *ci_y,
3009 if (*ci_b == ci_block)
3011 /* Jump to the next block assigned to this task */
3012 *ci += (nth - 1)*ci_block;
3021 while (*ci >= grid.cxy_ind[*ci_x*grid.numCells[YY] + *ci_y + 1])
3024 if (*ci_y == grid.numCells[YY])
3034 /* Returns the distance^2 for which we put cell pairs in the list
3035 * without checking atom pair distances. This is usually < rlist^2.
3037 static float boundingbox_only_distance2(const nbnxn_grid_t &iGrid,
3038 const nbnxn_grid_t &jGrid,
3042 /* If the distance between two sub-cell bounding boxes is less
3043 * than this distance, do not check the distance between
3044 * all particle pairs in the sub-cell, since then it is likely
3045 * that the box pair has atom pairs within the cut-off.
3046 * We use the nblist cut-off minus 0.5 times the average x/y diagonal
3047 * spacing of the sub-cells. Around 40% of the checked pairs are pruned.
3048 * Using more than 0.5 gains at most 0.5%.
3049 * If forces are calculated more than twice, the performance gain
3050 * in the force calculation outweighs the cost of checking.
3051 * Note that with subcell lists, the atom-pair distance check
3052 * is only performed when only 1 out of 8 sub-cells in within range,
3053 * this is because the GPU is much faster than the cpu.
3058 bbx = 0.5*(iGrid.cellSize[XX] + jGrid.cellSize[XX]);
3059 bby = 0.5*(iGrid.cellSize[YY] + jGrid.cellSize[YY]);
3062 bbx /= c_gpuNumClusterPerCellX;
3063 bby /= c_gpuNumClusterPerCellY;
3066 rbb2 = std::max(0.0, rlist - 0.5*std::sqrt(bbx*bbx + bby*bby));
3072 return (float)((1+GMX_FLOAT_EPS)*rbb2);
3076 static int get_ci_block_size(const nbnxn_grid_t &iGrid,
3077 gmx_bool bDomDec, int nth)
3079 const int ci_block_enum = 5;
3080 const int ci_block_denom = 11;
3081 const int ci_block_min_atoms = 16;
3084 /* Here we decide how to distribute the blocks over the threads.
3085 * We use prime numbers to try to avoid that the grid size becomes
3086 * a multiple of the number of threads, which would lead to some
3087 * threads getting "inner" pairs and others getting boundary pairs,
3088 * which in turns will lead to load imbalance between threads.
3089 * Set the block size as 5/11/ntask times the average number of cells
3090 * in a y,z slab. This should ensure a quite uniform distribution
3091 * of the grid parts of the different thread along all three grid
3092 * zone boundaries with 3D domain decomposition. At the same time
3093 * the blocks will not become too small.
3095 ci_block = (iGrid.nc*ci_block_enum)/(ci_block_denom*iGrid.numCells[XX]*nth);
3097 /* Ensure the blocks are not too small: avoids cache invalidation */
3098 if (ci_block*iGrid.na_sc < ci_block_min_atoms)
3100 ci_block = (ci_block_min_atoms + iGrid.na_sc - 1)/iGrid.na_sc;
3103 /* Without domain decomposition
3104 * or with less than 3 blocks per task, divide in nth blocks.
3106 if (!bDomDec || nth*3*ci_block > iGrid.nc)
3108 ci_block = (iGrid.nc + nth - 1)/nth;
3111 if (ci_block > 1 && (nth - 1)*ci_block >= iGrid.nc)
3113 /* Some threads have no work. Although reducing the block size
3114 * does not decrease the block count on the first few threads,
3115 * with GPUs better mixing of "upper" cells that have more empty
3116 * clusters results in a somewhat lower max load over all threads.
3117 * Without GPUs the regime of so few atoms per thread is less
3118 * performance relevant, but with 8-wide SIMD the same reasoning
3119 * applies, since the pair list uses 4 i-atom "sub-clusters".
3127 /* Returns the number of bits to right-shift a cluster index to obtain
3128 * the corresponding force buffer flag index.
3130 static int getBufferFlagShift(int numAtomsPerCluster)
3132 int bufferFlagShift = 0;
3133 while ((numAtomsPerCluster << bufferFlagShift) < NBNXN_BUFFERFLAG_SIZE)
3138 return bufferFlagShift;
3141 static bool pairlistIsSimple(const NbnxnPairlistCpu gmx_unused &pairlist)
3146 static bool pairlistIsSimple(const NbnxnPairlistGpu gmx_unused &pairlist)
3151 static void makeClusterListWrapper(NbnxnPairlistCpu *nbl,
3152 const nbnxn_grid_t gmx_unused &iGrid,
3154 const nbnxn_grid_t &jGrid,
3155 const int firstCell,
3157 const bool excludeSubDiagonal,
3158 const nbnxn_atomdata_t *nbat,
3161 const int nb_kernel_type,
3162 int *numDistanceChecks)
3164 switch (nb_kernel_type)
3166 case nbnxnk4x4_PlainC:
3167 makeClusterListSimple(jGrid,
3168 nbl, ci, firstCell, lastCell,
3174 #ifdef GMX_NBNXN_SIMD_4XN
3175 case nbnxnk4xN_SIMD_4xN:
3176 makeClusterListSimd4xn(jGrid,
3177 nbl, ci, firstCell, lastCell,
3184 #ifdef GMX_NBNXN_SIMD_2XNN
3185 case nbnxnk4xN_SIMD_2xNN:
3186 makeClusterListSimd2xnn(jGrid,
3187 nbl, ci, firstCell, lastCell,
3197 static void makeClusterListWrapper(NbnxnPairlistGpu *nbl,
3198 const nbnxn_grid_t &gmx_unused iGrid,
3200 const nbnxn_grid_t &jGrid,
3201 const int firstCell,
3203 const bool excludeSubDiagonal,
3204 const nbnxn_atomdata_t *nbat,
3207 const int gmx_unused nb_kernel_type,
3208 int *numDistanceChecks)
3210 for (int cj = firstCell; cj <= lastCell; cj++)
3212 make_cluster_list_supersub(iGrid, jGrid,
3215 nbat->xstride, nbat->x().data(),
3221 static int getNumSimpleJClustersInList(const NbnxnPairlistCpu &nbl)
3223 return nbl.cj.size();
3226 static int getNumSimpleJClustersInList(const gmx_unused NbnxnPairlistGpu &nbl)
3231 static void incrementNumSimpleJClustersInList(NbnxnPairlistCpu *nbl,
3234 nbl->ncjInUse += nbl->cj.size() - ncj_old_j;
3237 static void incrementNumSimpleJClustersInList(NbnxnPairlistGpu gmx_unused *nbl,
3238 int gmx_unused ncj_old_j)
3242 static void checkListSizeConsistency(const NbnxnPairlistCpu &nbl,
3243 const bool haveFreeEnergy)
3245 GMX_RELEASE_ASSERT(static_cast<size_t>(nbl.ncjInUse) == nbl.cj.size() || haveFreeEnergy,
3246 "Without free-energy all cj pair-list entries should be in use. "
3247 "Note that subsequent code does not make use of the equality, "
3248 "this check is only here to catch bugs");
3251 static void checkListSizeConsistency(const NbnxnPairlistGpu gmx_unused &nbl,
3252 bool gmx_unused haveFreeEnergy)
3254 /* We currently can not check consistency here */
3257 /* Set the buffer flags for newly added entries in the list */
3258 static void setBufferFlags(const NbnxnPairlistCpu &nbl,
3259 const int ncj_old_j,
3260 const int gridj_flag_shift,
3261 gmx_bitmask_t *gridj_flag,
3264 if (gmx::ssize(nbl.cj) > ncj_old_j)
3266 int cbFirst = nbl.cj[ncj_old_j].cj >> gridj_flag_shift;
3267 int cbLast = nbl.cj.back().cj >> gridj_flag_shift;
3268 for (int cb = cbFirst; cb <= cbLast; cb++)
3270 bitmask_init_bit(&gridj_flag[cb], th);
3275 static void setBufferFlags(const NbnxnPairlistGpu gmx_unused &nbl,
3276 int gmx_unused ncj_old_j,
3277 int gmx_unused gridj_flag_shift,
3278 gmx_bitmask_t gmx_unused *gridj_flag,
3281 GMX_ASSERT(false, "This function should never be called");
3284 /* Generates the part of pair-list nbl assigned to our thread */
3285 template <typename T>
3286 static void nbnxn_make_pairlist_part(const nbnxn_search *nbs,
3287 const nbnxn_grid_t &iGrid,
3288 const nbnxn_grid_t &jGrid,
3289 nbnxn_search_work_t *work,
3290 const nbnxn_atomdata_t *nbat,
3291 const t_blocka &exclusions,
3295 gmx_bool bFBufferFlag,
3298 float nsubpair_tot_est,
3305 real rlist2, rl_fep2 = 0;
3307 int ci_b, ci, ci_x, ci_y, ci_xy;
3309 real bx0, bx1, by0, by1, bz0, bz1;
3311 real d2cx, d2z, d2z_cx, d2z_cy, d2zx, d2zxy, d2xy;
3312 int cxf, cxl, cyf, cyf_x, cyl;
3313 int numDistanceChecks;
3314 int gridi_flag_shift = 0, gridj_flag_shift = 0;
3315 gmx_bitmask_t *gridj_flag = nullptr;
3316 int ncj_old_i, ncj_old_j;
3318 nbs_cycle_start(&work->cc[enbsCCsearch]);
3320 if (jGrid.bSimple != pairlistIsSimple(*nbl) ||
3321 iGrid.bSimple != pairlistIsSimple(*nbl))
3323 gmx_incons("Grid incompatible with pair-list");
3327 GMX_ASSERT(nbl->na_ci == jGrid.na_c, "The cluster sizes in the list and grid should match");
3328 nbl->na_cj = nbnxn_kernel_to_cluster_j_size(nb_kernel_type);
3329 na_cj_2log = get_2log(nbl->na_cj);
3335 /* Determine conversion of clusters to flag blocks */
3336 gridi_flag_shift = getBufferFlagShift(nbl->na_ci);
3337 gridj_flag_shift = getBufferFlagShift(nbl->na_cj);
3339 gridj_flag = work->buffer_flags.flag;
3342 copy_mat(nbs->box, box);
3344 rlist2 = nbl->rlist*nbl->rlist;
3346 if (nbs->bFEP && !pairlistIsSimple(*nbl))
3348 /* Determine an atom-pair list cut-off distance for FEP atom pairs.
3349 * We should not simply use rlist, since then we would not have
3350 * the small, effective buffering of the NxN lists.
3351 * The buffer is on overestimate, but the resulting cost for pairs
3352 * beyond rlist is neglible compared to the FEP pairs within rlist.
3354 rl_fep2 = nbl->rlist + effective_buffer_1x1_vs_MxN(iGrid, jGrid);
3358 fprintf(debug, "nbl_fep atom-pair rlist %f\n", rl_fep2);
3360 rl_fep2 = rl_fep2*rl_fep2;
3363 rbb2 = boundingbox_only_distance2(iGrid, jGrid, nbl->rlist, pairlistIsSimple(*nbl));
3367 fprintf(debug, "nbl bounding box only distance %f\n", std::sqrt(rbb2));
3370 const bool isIntraGridList = (&iGrid == &jGrid);
3372 /* Set the shift range */
3373 for (int d = 0; d < DIM; d++)
3375 /* Check if we need periodicity shifts.
3376 * Without PBC or with domain decomposition we don't need them.
3378 if (d >= ePBC2npbcdim(nbs->ePBC) || nbs->dd_dim[d])
3384 const real listRangeCellToCell = listRangeForGridCellToGridCell(rlist, iGrid, jGrid);
3386 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < listRangeCellToCell)
3396 const bool bSimple = pairlistIsSimple(*nbl);
3397 gmx::ArrayRef<const nbnxn_bb_t> bb_i;
3399 gmx::ArrayRef<const float> pbb_i;
3409 /* We use the normal bounding box format for both grid types */
3412 gmx::ArrayRef<const float> bbcz_i = iGrid.bbcz;
3413 gmx::ArrayRef<const int> flags_i = iGrid.flags;
3414 gmx::ArrayRef<const float> bbcz_j = jGrid.bbcz;
3415 int cell0_i = iGrid.cell0;
3419 fprintf(debug, "nbl nc_i %d col.av. %.1f ci_block %d\n",
3420 iGrid.nc, iGrid.nc/static_cast<double>(iGrid.numCells[XX]*iGrid.numCells[YY]), ci_block);
3423 numDistanceChecks = 0;
3425 const real listRangeBBToJCell2 = gmx::square(listRangeForBoundingBoxToGridCell(rlist, jGrid));
3427 /* Initially ci_b and ci to 1 before where we want them to start,
3428 * as they will both be incremented in next_ci.
3431 ci = th*ci_block - 1;
3434 while (next_ci(iGrid, nth, ci_block, &ci_x, &ci_y, &ci_b, &ci))
3436 if (bSimple && flags_i[ci] == 0)
3441 ncj_old_i = getNumSimpleJClustersInList(*nbl);
3444 if (!isIntraGridList && shp[XX] == 0)
3448 bx1 = bb_i[ci].upper[BB_X];
3452 bx1 = iGrid.c0[XX] + (ci_x+1)*iGrid.cellSize[XX];
3454 if (bx1 < jGrid.c0[XX])
3456 d2cx = gmx::square(jGrid.c0[XX] - bx1);
3458 if (d2cx >= listRangeBBToJCell2)
3465 ci_xy = ci_x*iGrid.numCells[YY] + ci_y;
3467 /* Loop over shift vectors in three dimensions */
3468 for (int tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
3470 const real shz = tz*box[ZZ][ZZ];
3472 bz0 = bbcz_i[ci*NNBSBB_D ] + shz;
3473 bz1 = bbcz_i[ci*NNBSBB_D+1] + shz;
3481 d2z = gmx::square(bz1);
3485 d2z = gmx::square(bz0 - box[ZZ][ZZ]);
3488 d2z_cx = d2z + d2cx;
3490 if (d2z_cx >= rlist2)
3495 bz1_frac = bz1/(iGrid.cxy_ind[ci_xy+1] - iGrid.cxy_ind[ci_xy]);
3500 /* The check with bz1_frac close to or larger than 1 comes later */
3502 for (int ty = -shp[YY]; ty <= shp[YY]; ty++)
3504 const real shy = ty*box[YY][YY] + tz*box[ZZ][YY];
3508 by0 = bb_i[ci].lower[BB_Y] + shy;
3509 by1 = bb_i[ci].upper[BB_Y] + shy;
3513 by0 = iGrid.c0[YY] + (ci_y )*iGrid.cellSize[YY] + shy;
3514 by1 = iGrid.c0[YY] + (ci_y+1)*iGrid.cellSize[YY] + shy;
3517 get_cell_range<YY>(by0, by1,
3528 if (by1 < jGrid.c0[YY])
3530 d2z_cy += gmx::square(jGrid.c0[YY] - by1);
3532 else if (by0 > jGrid.c1[YY])
3534 d2z_cy += gmx::square(by0 - jGrid.c1[YY]);
3537 for (int tx = -shp[XX]; tx <= shp[XX]; tx++)
3539 const int shift = XYZ2IS(tx, ty, tz);
3541 const bool excludeSubDiagonal = (isIntraGridList && shift == CENTRAL);
3543 if (c_pbcShiftBackward && isIntraGridList && shift > CENTRAL)
3548 const real shx = tx*box[XX][XX] + ty*box[YY][XX] + tz*box[ZZ][XX];
3552 bx0 = bb_i[ci].lower[BB_X] + shx;
3553 bx1 = bb_i[ci].upper[BB_X] + shx;
3557 bx0 = iGrid.c0[XX] + (ci_x )*iGrid.cellSize[XX] + shx;
3558 bx1 = iGrid.c0[XX] + (ci_x+1)*iGrid.cellSize[XX] + shx;
3561 get_cell_range<XX>(bx0, bx1,
3571 addNewIEntry(nbl, cell0_i+ci, shift, flags_i[ci]);
3573 if ((!c_pbcShiftBackward || excludeSubDiagonal) &&
3576 /* Leave the pairs with i > j.
3577 * x is the major index, so skip half of it.
3582 set_icell_bb(iGrid, ci, shx, shy, shz,
3585 icell_set_x(cell0_i+ci, shx, shy, shz,
3586 nbat->xstride, nbat->x().data(),
3590 for (int cx = cxf; cx <= cxl; cx++)
3593 if (jGrid.c0[XX] + cx*jGrid.cellSize[XX] > bx1)
3595 d2zx += gmx::square(jGrid.c0[XX] + cx*jGrid.cellSize[XX] - bx1);
3597 else if (jGrid.c0[XX] + (cx+1)*jGrid.cellSize[XX] < bx0)
3599 d2zx += gmx::square(jGrid.c0[XX] + (cx+1)*jGrid.cellSize[XX] - bx0);
3602 if (isIntraGridList &&
3604 (!c_pbcShiftBackward || shift == CENTRAL) &&
3607 /* Leave the pairs with i > j.
3608 * Skip half of y when i and j have the same x.
3617 for (int cy = cyf_x; cy <= cyl; cy++)
3619 const int columnStart = jGrid.cxy_ind[cx*jGrid.numCells[YY] + cy];
3620 const int columnEnd = jGrid.cxy_ind[cx*jGrid.numCells[YY] + cy + 1];
3623 if (jGrid.c0[YY] + cy*jGrid.cellSize[YY] > by1)
3625 d2zxy += gmx::square(jGrid.c0[YY] + cy*jGrid.cellSize[YY] - by1);
3627 else if (jGrid.c0[YY] + (cy+1)*jGrid.cellSize[YY] < by0)
3629 d2zxy += gmx::square(jGrid.c0[YY] + (cy+1)*jGrid.cellSize[YY] - by0);
3631 if (columnStart < columnEnd && d2zxy < listRangeBBToJCell2)
3633 /* To improve efficiency in the common case
3634 * of a homogeneous particle distribution,
3635 * we estimate the index of the middle cell
3636 * in range (midCell). We search down and up
3637 * starting from this index.
3639 * Note that the bbcz_j array contains bounds
3640 * for i-clusters, thus for clusters of 4 atoms.
3641 * For the common case where the j-cluster size
3642 * is 8, we could step with a stride of 2,
3643 * but we do not do this because it would
3644 * complicate this code even more.
3646 int midCell = columnStart + static_cast<int>(bz1_frac*(columnEnd - columnStart));
3647 if (midCell >= columnEnd)
3649 midCell = columnEnd - 1;
3654 /* Find the lowest cell that can possibly
3656 * Check if we hit the bottom of the grid,
3657 * if the j-cell is below the i-cell and if so,
3658 * if it is within range.
3660 int downTestCell = midCell;
3661 while (downTestCell >= columnStart &&
3662 (bbcz_j[downTestCell*NNBSBB_D + 1] >= bz0 ||
3663 d2xy + gmx::square(bbcz_j[downTestCell*NNBSBB_D + 1] - bz0) < rlist2))
3667 int firstCell = downTestCell + 1;
3669 /* Find the highest cell that can possibly
3671 * Check if we hit the top of the grid,
3672 * if the j-cell is above the i-cell and if so,
3673 * if it is within range.
3675 int upTestCell = midCell + 1;
3676 while (upTestCell < columnEnd &&
3677 (bbcz_j[upTestCell*NNBSBB_D] <= bz1 ||
3678 d2xy + gmx::square(bbcz_j[upTestCell*NNBSBB_D] - bz1) < rlist2))
3682 int lastCell = upTestCell - 1;
3684 #define NBNXN_REFCODE 0
3687 /* Simple reference code, for debugging,
3688 * overrides the more complex code above.
3690 firstCell = columnEnd;
3692 for (int k = columnStart; k < columnEnd; k++)
3694 if (d2xy + gmx::square(bbcz_j[k*NNBSBB_D + 1] - bz0) < rlist2 &&
3699 if (d2xy + gmx::square(bbcz_j[k*NNBSBB_D] - bz1) < rlist2 &&
3708 if (isIntraGridList)
3710 /* We want each atom/cell pair only once,
3711 * only use cj >= ci.
3713 if (!c_pbcShiftBackward || shift == CENTRAL)
3715 firstCell = std::max(firstCell, ci);
3719 if (firstCell <= lastCell)
3721 GMX_ASSERT(firstCell >= columnStart && lastCell < columnEnd, "The range should reside within the current grid column");
3723 /* For f buffer flags with simple lists */
3724 ncj_old_j = getNumSimpleJClustersInList(*nbl);
3726 makeClusterListWrapper(nbl,
3728 jGrid, firstCell, lastCell,
3733 &numDistanceChecks);
3737 setBufferFlags(*nbl, ncj_old_j, gridj_flag_shift,
3741 incrementNumSimpleJClustersInList(nbl, ncj_old_j);
3747 /* Set the exclusions for this ci list */
3748 setExclusionsForIEntry(nbs,
3752 *getOpenIEntry(nbl),
3757 make_fep_list(nbs, nbat, nbl,
3762 iGrid, jGrid, nbl_fep);
3765 /* Close this ci list */
3768 progBal, nsubpair_tot_est,
3774 if (bFBufferFlag && getNumSimpleJClustersInList(*nbl) > ncj_old_i)
3776 bitmask_init_bit(&(work->buffer_flags.flag[(iGrid.cell0+ci) >> gridi_flag_shift]), th);
3780 work->ndistc = numDistanceChecks;
3782 nbs_cycle_stop(&work->cc[enbsCCsearch]);
3784 checkListSizeConsistency(*nbl, nbs->bFEP);
3788 fprintf(debug, "number of distance checks %d\n", numDistanceChecks);
3790 print_nblist_statistics(debug, nbl, nbs, rlist);
3794 fprintf(debug, "nbl FEP list pairs: %d\n", nbl_fep->nrj);
3799 static void reduce_buffer_flags(const nbnxn_search *nbs,
3801 const nbnxn_buffer_flags_t *dest)
3803 for (int s = 0; s < nsrc; s++)
3805 gmx_bitmask_t * flag = nbs->work[s].buffer_flags.flag;
3807 for (int b = 0; b < dest->nflag; b++)
3809 bitmask_union(&(dest->flag[b]), flag[b]);
3814 static void print_reduction_cost(const nbnxn_buffer_flags_t *flags, int nout)
3816 int nelem, nkeep, ncopy, nred, out;
3817 gmx_bitmask_t mask_0;
3823 bitmask_init_bit(&mask_0, 0);
3824 for (int b = 0; b < flags->nflag; b++)
3826 if (bitmask_is_equal(flags->flag[b], mask_0))
3828 /* Only flag 0 is set, no copy of reduction required */
3832 else if (!bitmask_is_zero(flags->flag[b]))
3835 for (out = 0; out < nout; out++)
3837 if (bitmask_is_set(flags->flag[b], out))
3854 fprintf(debug, "nbnxn reduction: #flag %d #list %d elem %4.2f, keep %4.2f copy %4.2f red %4.2f\n",
3856 nelem/static_cast<double>(flags->nflag),
3857 nkeep/static_cast<double>(flags->nflag),
3858 ncopy/static_cast<double>(flags->nflag),
3859 nred/static_cast<double>(flags->nflag));
3862 /* Copies the list entries from src to dest when cjStart <= *cjGlobal < cjEnd.
3863 * *cjGlobal is updated with the cj count in src.
3864 * When setFlags==true, flag bit t is set in flag for all i and j clusters.
3866 template<bool setFlags>
3867 static void copySelectedListRange(const nbnxn_ci_t * gmx_restrict srcCi,
3868 const NbnxnPairlistCpu * gmx_restrict src,
3869 NbnxnPairlistCpu * gmx_restrict dest,
3870 gmx_bitmask_t *flag,
3871 int iFlagShift, int jFlagShift, int t)
3873 const int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3875 dest->ci.push_back(*srcCi);
3876 dest->ci.back().cj_ind_start = dest->cj.size();
3877 dest->ci.back().cj_ind_end = dest->cj.size() + ncj;
3881 bitmask_init_bit(&flag[srcCi->ci >> iFlagShift], t);
3884 for (int j = srcCi->cj_ind_start; j < srcCi->cj_ind_end; j++)
3886 dest->cj.push_back(src->cj[j]);
3890 /* NOTE: This is relatively expensive, since this
3891 * operation is done for all elements in the list,
3892 * whereas at list generation this is done only
3893 * once for each flag entry.
3895 bitmask_init_bit(&flag[src->cj[j].cj >> jFlagShift], t);
3900 /* This routine re-balances the pairlists such that all are nearly equally
3901 * sized. Only whole i-entries are moved between lists. These are moved
3902 * between the ends of the lists, such that the buffer reduction cost should
3903 * not change significantly.
3904 * Note that all original reduction flags are currently kept. This can lead
3905 * to reduction of parts of the force buffer that could be avoided. But since
3906 * the original lists are quite balanced, this will only give minor overhead.
3908 static void rebalanceSimpleLists(int numLists,
3909 NbnxnPairlistCpu * const * const srcSet,
3910 NbnxnPairlistCpu **destSet,
3911 gmx::ArrayRef<nbnxn_search_work_t> searchWork)
3914 for (int s = 0; s < numLists; s++)
3916 ncjTotal += srcSet[s]->ncjInUse;
3918 int ncjTarget = (ncjTotal + numLists - 1)/numLists;
3920 #pragma omp parallel num_threads(numLists)
3922 int t = gmx_omp_get_thread_num();
3924 int cjStart = ncjTarget* t;
3925 int cjEnd = ncjTarget*(t + 1);
3927 /* The destination pair-list for task/thread t */
3928 NbnxnPairlistCpu *dest = destSet[t];
3930 clear_pairlist(dest);
3931 dest->na_cj = srcSet[0]->na_cj;
3933 /* Note that the flags in the work struct (still) contain flags
3934 * for all entries that are present in srcSet->nbl[t].
3936 gmx_bitmask_t *flag = searchWork[t].buffer_flags.flag;
3938 int iFlagShift = getBufferFlagShift(dest->na_ci);
3939 int jFlagShift = getBufferFlagShift(dest->na_cj);
3942 for (int s = 0; s < numLists && cjGlobal < cjEnd; s++)
3944 const NbnxnPairlistCpu *src = srcSet[s];
3946 if (cjGlobal + src->ncjInUse > cjStart)
3948 for (gmx::index i = 0; i < gmx::ssize(src->ci) && cjGlobal < cjEnd; i++)
3950 const nbnxn_ci_t *srcCi = &src->ci[i];
3951 int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3952 if (cjGlobal >= cjStart)
3954 /* If the source list is not our own, we need to set
3955 * extra flags (the template bool parameter).
3959 copySelectedListRange
3962 flag, iFlagShift, jFlagShift, t);
3966 copySelectedListRange
3969 dest, flag, iFlagShift, jFlagShift, t);
3977 cjGlobal += src->ncjInUse;
3981 dest->ncjInUse = dest->cj.size();
3985 int ncjTotalNew = 0;
3986 for (int s = 0; s < numLists; s++)
3988 ncjTotalNew += destSet[s]->ncjInUse;
3990 GMX_RELEASE_ASSERT(ncjTotalNew == ncjTotal, "The total size of the lists before and after rebalancing should match");
3994 /* Returns if the pairlists are so imbalanced that it is worth rebalancing. */
3995 static bool checkRebalanceSimpleLists(const nbnxn_pairlist_set_t *listSet)
3997 int numLists = listSet->nnbl;
4000 for (int s = 0; s < numLists; s++)
4002 ncjMax = std::max(ncjMax, listSet->nbl[s]->ncjInUse);
4003 ncjTotal += listSet->nbl[s]->ncjInUse;
4007 fprintf(debug, "Pair-list ncjMax %d ncjTotal %d\n", ncjMax, ncjTotal);
4009 /* The rebalancing adds 3% extra time to the search. Heuristically we
4010 * determined that under common conditions the non-bonded kernel balance
4011 * improvement will outweigh this when the imbalance is more than 3%.
4012 * But this will, obviously, depend on search vs kernel time and nstlist.
4014 const real rebalanceTolerance = 1.03;
4016 return numLists*ncjMax > ncjTotal*rebalanceTolerance;
4019 /* Perform a count (linear) sort to sort the smaller lists to the end.
4020 * This avoids load imbalance on the GPU, as large lists will be
4021 * scheduled and executed first and the smaller lists later.
4022 * Load balancing between multi-processors only happens at the end
4023 * and there smaller lists lead to more effective load balancing.
4024 * The sorting is done on the cj4 count, not on the actual pair counts.
4025 * Not only does this make the sort faster, but it also results in
4026 * better load balancing than using a list sorted on exact load.
4027 * This function swaps the pointer in the pair list to avoid a copy operation.
4029 static void sort_sci(NbnxnPairlistGpu *nbl)
4031 if (nbl->cj4.size() <= nbl->sci.size())
4033 /* nsci = 0 or all sci have size 1, sorting won't change the order */
4037 NbnxnPairlistGpuWork &work = *nbl->work;
4039 /* We will distinguish differences up to double the average */
4040 const int m = (2*nbl->cj4.size())/nbl->sci.size();
4042 /* Resize work.sci_sort so we can sort into it */
4043 work.sci_sort.resize(nbl->sci.size());
4045 std::vector<int> &sort = work.sortBuffer;
4046 /* Set up m + 1 entries in sort, initialized at 0 */
4048 sort.resize(m + 1, 0);
4049 /* Count the entries of each size */
4050 for (const nbnxn_sci_t &sci : nbl->sci)
4052 int i = std::min(m, sci.numJClusterGroups());
4055 /* Calculate the offset for each count */
4058 for (int i = m - 1; i >= 0; i--)
4061 sort[i] = sort[i + 1] + s0;
4065 /* Sort entries directly into place */
4066 gmx::ArrayRef<nbnxn_sci_t> sci_sort = work.sci_sort;
4067 for (const nbnxn_sci_t &sci : nbl->sci)
4069 int i = std::min(m, sci.numJClusterGroups());
4070 sci_sort[sort[i]++] = sci;
4073 /* Swap the sci pointers so we use the new, sorted list */
4074 std::swap(nbl->sci, work.sci_sort);
4077 /* Make a local or non-local pair-list, depending on iloc */
4078 void nbnxn_make_pairlist(nbnxn_search *nbs,
4079 nbnxn_atomdata_t *nbat,
4080 const t_blocka *excl,
4082 int min_ci_balanced,
4083 nbnxn_pairlist_set_t *nbl_list,
4088 int nsubpair_target;
4089 float nsubpair_tot_est;
4092 gmx_bool CombineNBLists;
4094 int np_tot, np_noq, np_hlj, nap;
4096 nnbl = nbl_list->nnbl;
4097 CombineNBLists = nbl_list->bCombined;
4101 fprintf(debug, "ns making %d nblists\n", nnbl);
4104 nbat->bUseBufferFlags = (nbat->out.size() > 1);
4105 /* We should re-init the flags before making the first list */
4106 if (nbat->bUseBufferFlags && LOCAL_I(iloc))
4108 init_buffer_flags(&nbat->buffer_flags, nbat->numAtoms());
4114 /* Only zone (grid) 0 vs 0 */
4119 nzi = nbs->zones->nizone;
4122 if (!nbl_list->bSimple && min_ci_balanced > 0)
4124 get_nsubpair_target(nbs, iloc, rlist, min_ci_balanced,
4125 &nsubpair_target, &nsubpair_tot_est);
4129 nsubpair_target = 0;
4130 nsubpair_tot_est = 0;
4133 /* Clear all pair-lists */
4134 for (int th = 0; th < nnbl; th++)
4136 if (nbl_list->bSimple)
4138 clear_pairlist(nbl_list->nbl[th]);
4142 clear_pairlist(nbl_list->nblGpu[th]);
4147 clear_pairlist_fep(nbl_list->nbl_fep[th]);
4151 for (int zi = 0; zi < nzi; zi++)
4153 const nbnxn_grid_t &iGrid = nbs->grid[zi];
4164 zj0 = nbs->zones->izone[zi].j0;
4165 zj1 = nbs->zones->izone[zi].j1;
4171 for (int zj = zj0; zj < zj1; zj++)
4173 const nbnxn_grid_t &jGrid = nbs->grid[zj];
4177 fprintf(debug, "ns search grid %d vs %d\n", zi, zj);
4180 nbs_cycle_start(&nbs->cc[enbsCCsearch]);
4182 ci_block = get_ci_block_size(iGrid, nbs->DomDec, nnbl);
4184 /* With GPU: generate progressively smaller lists for
4185 * load balancing for local only or non-local with 2 zones.
4187 progBal = (LOCAL_I(iloc) || nbs->zones->n <= 2);
4189 #pragma omp parallel for num_threads(nnbl) schedule(static)
4190 for (int th = 0; th < nnbl; th++)
4194 /* Re-init the thread-local work flag data before making
4195 * the first list (not an elegant conditional).
4197 if (nbat->bUseBufferFlags && ((zi == 0 && zj == 0)))
4199 init_buffer_flags(&nbs->work[th].buffer_flags, nbat->numAtoms());
4202 if (CombineNBLists && th > 0)
4204 GMX_ASSERT(!nbl_list->bSimple, "Can only combine GPU lists");
4206 clear_pairlist(nbl_list->nblGpu[th]);
4209 /* Divide the i super cell equally over the nblists */
4210 if (nbl_list->bSimple)
4212 nbnxn_make_pairlist_part(nbs, iGrid, jGrid,
4213 &nbs->work[th], nbat, *excl,
4217 nbat->bUseBufferFlags,
4219 progBal, nsubpair_tot_est,
4222 nbl_list->nbl_fep[th]);
4226 nbnxn_make_pairlist_part(nbs, iGrid, jGrid,
4227 &nbs->work[th], nbat, *excl,
4231 nbat->bUseBufferFlags,
4233 progBal, nsubpair_tot_est,
4235 nbl_list->nblGpu[th],
4236 nbl_list->nbl_fep[th]);
4239 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4241 nbs_cycle_stop(&nbs->cc[enbsCCsearch]);
4246 for (int th = 0; th < nnbl; th++)
4248 inc_nrnb(nrnb, eNR_NBNXN_DIST2, nbs->work[th].ndistc);
4250 if (nbl_list->bSimple)
4252 NbnxnPairlistCpu *nbl = nbl_list->nbl[th];
4253 np_tot += nbl->cj.size();
4254 np_noq += nbl->work->ncj_noq;
4255 np_hlj += nbl->work->ncj_hlj;
4259 NbnxnPairlistGpu *nbl = nbl_list->nblGpu[th];
4260 /* This count ignores potential subsequent pair pruning */
4261 np_tot += nbl->nci_tot;
4264 if (nbl_list->bSimple)
4266 nap = nbl_list->nbl[0]->na_ci*nbl_list->nbl[0]->na_cj;
4270 nap = gmx::square(nbl_list->nblGpu[0]->na_ci);
4272 nbl_list->natpair_ljq = (np_tot - np_noq)*nap - np_hlj*nap/2;
4273 nbl_list->natpair_lj = np_noq*nap;
4274 nbl_list->natpair_q = np_hlj*nap/2;
4276 if (CombineNBLists && nnbl > 1)
4278 GMX_ASSERT(!nbl_list->bSimple, "Can only combine GPU lists");
4279 NbnxnPairlistGpu **nbl = nbl_list->nblGpu;
4281 nbs_cycle_start(&nbs->cc[enbsCCcombine]);
4283 combine_nblists(nnbl-1, nbl+1, nbl[0]);
4285 nbs_cycle_stop(&nbs->cc[enbsCCcombine]);
4290 if (nbl_list->bSimple)
4292 if (nnbl > 1 && checkRebalanceSimpleLists(nbl_list))
4294 rebalanceSimpleLists(nbl_list->nnbl, nbl_list->nbl, nbl_list->nbl_work, nbs->work);
4296 /* Swap the pointer of the sets of pair lists */
4297 NbnxnPairlistCpu **tmp = nbl_list->nbl;
4298 nbl_list->nbl = nbl_list->nbl_work;
4299 nbl_list->nbl_work = tmp;
4304 /* Sort the entries on size, large ones first */
4305 if (CombineNBLists || nnbl == 1)
4307 sort_sci(nbl_list->nblGpu[0]);
4311 #pragma omp parallel for num_threads(nnbl) schedule(static)
4312 for (int th = 0; th < nnbl; th++)
4316 sort_sci(nbl_list->nblGpu[th]);
4318 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4323 if (nbat->bUseBufferFlags)
4325 reduce_buffer_flags(nbs, nbl_list->nnbl, &nbat->buffer_flags);
4330 /* Balance the free-energy lists over all the threads */
4331 balance_fep_lists(nbs, nbl_list);
4334 if (nbl_list->bSimple)
4336 /* This is a fresh list, so not pruned, stored using ci.
4337 * ciOuter is invalid at this point.
4339 GMX_ASSERT(nbl_list->nbl[0]->ciOuter.empty(), "ciOuter is invalid so it should be empty");
4342 /* Special performance logging stuff (env.var. GMX_NBNXN_CYCLE) */
4345 nbs->search_count++;
4347 if (nbs->print_cycles &&
4348 (!nbs->DomDec || !LOCAL_I(iloc)) &&
4349 nbs->search_count % 100 == 0)
4351 nbs_cycle_print(stderr, nbs);
4354 /* If we have more than one list, they either got rebalancing (CPU)
4355 * or combined (GPU), so we should dump the final result to debug.
4357 if (debug && nbl_list->nnbl > 1)
4359 if (nbl_list->bSimple)
4361 for (int t = 0; t < nbl_list->nnbl; t++)
4363 print_nblist_statistics(debug, nbl_list->nbl[t], nbs, rlist);
4368 print_nblist_statistics(debug, nbl_list->nblGpu[0], nbs, rlist);
4376 if (nbl_list->bSimple)
4378 for (int t = 0; t < nbl_list->nnbl; t++)
4380 print_nblist_ci_cj(debug, nbl_list->nbl[t]);
4385 print_nblist_sci_cj(debug, nbl_list->nblGpu[0]);
4389 if (nbat->bUseBufferFlags)
4391 print_reduction_cost(&nbat->buffer_flags, nbl_list->nnbl);
4396 void nbnxnPrepareListForDynamicPruning(nbnxn_pairlist_set_t *listSet)
4398 GMX_RELEASE_ASSERT(listSet->bSimple, "Should only be called for simple lists");
4400 /* TODO: Restructure the lists so we have actual outer and inner
4401 * list objects so we can set a single pointer instead of
4402 * swapping several pointers.
4405 for (int i = 0; i < listSet->nnbl; i++)
4407 NbnxnPairlistCpu &list = *listSet->nbl[i];
4409 /* The search produced a list in ci/cj.
4410 * Swap the list pointers so we get the outer list is ciOuter,cjOuter
4411 * and we can prune that to get an inner list in ci/cj.
4413 GMX_RELEASE_ASSERT(list.ciOuter.empty() && list.cjOuter.empty(),
4414 "The outer lists should be empty before preparation");
4416 std::swap(list.ci, list.ciOuter);
4417 std::swap(list.cj, list.cjOuter);