2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, 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"
49 #include "gromacs/domdec/domdec_struct.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/utilities.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/mdlib/gmx_omp_nthreads.h"
55 #include "gromacs/mdlib/nb_verlet.h"
56 #include "gromacs/mdlib/nbnxn_atomdata.h"
57 #include "gromacs/mdlib/nbnxn_consts.h"
58 #include "gromacs/mdlib/nbnxn_grid.h"
59 #include "gromacs/mdlib/nbnxn_internal.h"
60 #include "gromacs/mdlib/nbnxn_simd.h"
61 #include "gromacs/mdlib/nbnxn_util.h"
62 #include "gromacs/mdlib/ns.h"
63 #include "gromacs/mdtypes/group.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/pbcutil/ishift.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/simd/simd.h"
68 #include "gromacs/simd/vector_operations.h"
69 #include "gromacs/topology/block.h"
70 #include "gromacs/utility/exceptions.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/gmxomp.h"
73 #include "gromacs/utility/smalloc.h"
75 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
78 /* We shift the i-particles backward for PBC.
79 * This leads to more conditionals than shifting forward.
80 * We do this to get more balanced pair lists.
82 constexpr bool c_pbcShiftBackward = true;
85 static void nbs_cycle_clear(nbnxn_cycle_t *cc)
87 for (int i = 0; i < enbsCCnr; i++)
94 static double Mcyc_av(const nbnxn_cycle_t *cc)
96 return (double)cc->c*1e-6/cc->count;
99 static void nbs_cycle_print(FILE *fp, const nbnxn_search *nbs)
102 fprintf(fp, "ns %4d grid %4.1f search %4.1f red.f %5.3f",
103 nbs->cc[enbsCCgrid].count,
104 Mcyc_av(&nbs->cc[enbsCCgrid]),
105 Mcyc_av(&nbs->cc[enbsCCsearch]),
106 Mcyc_av(&nbs->cc[enbsCCreducef]));
108 if (nbs->work.size() > 1)
110 if (nbs->cc[enbsCCcombine].count > 0)
112 fprintf(fp, " comb %5.2f",
113 Mcyc_av(&nbs->cc[enbsCCcombine]));
115 fprintf(fp, " s. th");
116 for (const nbnxn_search_work_t &work : nbs->work)
118 fprintf(fp, " %4.1f",
119 Mcyc_av(&work.cc[enbsCCsearch]));
125 /* Layout for the nonbonded NxN pair lists */
126 enum class NbnxnLayout
128 NoSimd4x4, // i-cluster size 4, j-cluster size 4
129 Simd4xN, // i-cluster size 4, j-cluster size SIMD width
130 Simd2xNN, // i-cluster size 4, j-cluster size half SIMD width
131 Gpu8x8x8 // i-cluster size 8, j-cluster size 8 + super-clustering
135 /* Returns the j-cluster size */
136 template <NbnxnLayout layout>
137 static constexpr int jClusterSize()
139 static_assert(layout == NbnxnLayout::NoSimd4x4 || layout == NbnxnLayout::Simd4xN || layout == NbnxnLayout::Simd2xNN, "Currently jClusterSize only supports CPU layouts");
141 return layout == NbnxnLayout::Simd4xN ? GMX_SIMD_REAL_WIDTH : (layout == NbnxnLayout::Simd2xNN ? GMX_SIMD_REAL_WIDTH/2 : NBNXN_CPU_CLUSTER_I_SIZE);
144 /* Returns the j-cluster index given the i-cluster index */
145 template <int jClusterSize>
146 static inline int cjFromCi(int ci)
148 static_assert(jClusterSize == NBNXN_CPU_CLUSTER_I_SIZE/2 || jClusterSize == NBNXN_CPU_CLUSTER_I_SIZE || jClusterSize == NBNXN_CPU_CLUSTER_I_SIZE*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
150 if (jClusterSize == NBNXN_CPU_CLUSTER_I_SIZE/2)
154 else if (jClusterSize == NBNXN_CPU_CLUSTER_I_SIZE)
164 /* Returns the j-cluster index given the i-cluster index */
165 template <NbnxnLayout layout>
166 static inline int cjFromCi(int ci)
168 constexpr int clusterSize = jClusterSize<layout>();
170 return cjFromCi<clusterSize>(ci);
173 /* Returns the nbnxn coordinate data index given the i-cluster index */
174 template <NbnxnLayout layout>
175 static inline int xIndexFromCi(int ci)
177 constexpr int clusterSize = jClusterSize<layout>();
179 static_assert(clusterSize == NBNXN_CPU_CLUSTER_I_SIZE/2 || clusterSize == NBNXN_CPU_CLUSTER_I_SIZE || clusterSize == NBNXN_CPU_CLUSTER_I_SIZE*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
181 if (clusterSize <= NBNXN_CPU_CLUSTER_I_SIZE)
183 /* Coordinates are stored packed in groups of 4 */
188 /* Coordinates packed in 8, i-cluster size is half the packing width */
189 return (ci >> 1)*STRIDE_P8 + (ci & 1)*(c_packX8 >> 1);
193 /* Returns the nbnxn coordinate data index given the j-cluster index */
194 template <NbnxnLayout layout>
195 static inline int xIndexFromCj(int cj)
197 constexpr int clusterSize = jClusterSize<layout>();
199 static_assert(clusterSize == NBNXN_CPU_CLUSTER_I_SIZE/2 || clusterSize == NBNXN_CPU_CLUSTER_I_SIZE || clusterSize == NBNXN_CPU_CLUSTER_I_SIZE*2, "Only j-cluster sizes 2, 4 and 8 are currently implemented");
201 if (clusterSize == NBNXN_CPU_CLUSTER_I_SIZE/2)
203 /* Coordinates are stored packed in groups of 4 */
204 return (cj >> 1)*STRIDE_P4 + (cj & 1)*(c_packX4 >> 1);
206 else if (clusterSize == NBNXN_CPU_CLUSTER_I_SIZE)
208 /* Coordinates are stored packed in groups of 4 */
213 /* Coordinates are stored packed in groups of 8 */
219 gmx_bool nbnxn_kernel_pairlist_simple(int nb_kernel_type)
221 if (nb_kernel_type == nbnxnkNotSet)
223 gmx_fatal(FARGS, "Non-bonded kernel type not set for Verlet-style pair-list.");
226 switch (nb_kernel_type)
228 case nbnxnk8x8x8_GPU:
229 case nbnxnk8x8x8_PlainC:
232 case nbnxnk4x4_PlainC:
233 case nbnxnk4xN_SIMD_4xN:
234 case nbnxnk4xN_SIMD_2xNN:
238 gmx_incons("Invalid nonbonded kernel type passed!");
243 /* Initializes a single nbnxn_pairlist_t data structure */
244 static void nbnxn_init_pairlist_fep(t_nblist *nl)
246 nl->type = GMX_NBLIST_INTERACTION_FREE_ENERGY;
247 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
248 /* The interaction functions are set in the free energy kernel fuction */
261 nl->jindex = nullptr;
263 nl->excl_fep = nullptr;
267 static void free_nblist(t_nblist *nl)
277 nbnxn_search_work_t::nbnxn_search_work_t() :
280 buffer_flags({0, nullptr, 0}),
282 nbl_fep(new t_nblist),
285 nbnxn_init_pairlist_fep(nbl_fep.get());
290 nbnxn_search_work_t::~nbnxn_search_work_t()
292 sfree(buffer_flags.flag);
294 free_nblist(nbl_fep.get());
297 nbnxn_search::nbnxn_search(const ivec *n_dd_cells,
298 const gmx_domdec_zones_t *zones,
302 ePBC(epbcNONE), // The correct value will be set during the gridding
309 // The correct value will be set during the gridding
313 DomDec = n_dd_cells != nullptr;
316 for (int d = 0; d < DIM; d++)
318 if ((*n_dd_cells)[d] > 1)
321 /* Each grid matches a DD zone */
327 grid.resize(numGrids);
329 /* Initialize detailed nbsearch cycle counting */
330 print_cycles = (getenv("GMX_NBNXN_CYCLE") != nullptr);
334 nbnxn_search *nbnxn_init_search(const ivec *n_dd_cells,
335 const gmx_domdec_zones_t *zones,
339 return new nbnxn_search(n_dd_cells, zones, bFEP, nthread_max);
342 static void init_buffer_flags(nbnxn_buffer_flags_t *flags,
345 flags->nflag = (natoms + NBNXN_BUFFERFLAG_SIZE - 1)/NBNXN_BUFFERFLAG_SIZE;
346 if (flags->nflag > flags->flag_nalloc)
348 flags->flag_nalloc = over_alloc_large(flags->nflag);
349 srenew(flags->flag, flags->flag_nalloc);
351 for (int b = 0; b < flags->nflag; b++)
353 bitmask_clear(&(flags->flag[b]));
357 /* Determines the cell range along one dimension that
358 * the bounding box b0 - b1 sees.
361 static void get_cell_range(real b0, real b1,
362 const nbnxn_grid_t &gridj,
363 real d2, real r2, int *cf, int *cl)
365 real distanceInCells = (b0 - gridj.c0[dim])*gridj.invCellSize[dim];
366 *cf = std::max(static_cast<int>(distanceInCells), 0);
369 d2 + gmx::square((b0 - gridj.c0[dim]) - (*cf - 1 + 1)*gridj.cellSize[dim]) < r2)
374 *cl = std::min(static_cast<int>((b1 - gridj.c0[dim])*gridj.invCellSize[dim]), gridj.numCells[dim] - 1);
375 while (*cl < gridj.numCells[dim] - 1 &&
376 d2 + gmx::square((*cl + 1)*gridj.cellSize[dim] - (b1 - gridj.c0[dim])) < r2)
382 /* Reference code calculating the distance^2 between two bounding boxes */
384 static float box_dist2(float bx0, float bx1, float by0,
385 float by1, float bz0, float bz1,
386 const nbnxn_bb_t *bb)
389 float dl, dh, dm, dm0;
393 dl = bx0 - bb->upper[BB_X];
394 dh = bb->lower[BB_X] - bx1;
395 dm = std::max(dl, dh);
396 dm0 = std::max(dm, 0.0f);
399 dl = by0 - bb->upper[BB_Y];
400 dh = bb->lower[BB_Y] - by1;
401 dm = std::max(dl, dh);
402 dm0 = std::max(dm, 0.0f);
405 dl = bz0 - bb->upper[BB_Z];
406 dh = bb->lower[BB_Z] - bz1;
407 dm = std::max(dl, dh);
408 dm0 = std::max(dm, 0.0f);
415 /* Plain C code calculating the distance^2 between two bounding boxes */
416 static float subc_bb_dist2(int si,
417 const nbnxn_bb_t *bb_i_ci,
419 gmx::ArrayRef<const nbnxn_bb_t> bb_j_all)
421 const nbnxn_bb_t *bb_i = bb_i_ci + si;
422 const nbnxn_bb_t *bb_j = bb_j_all.data() + csj;
425 float dl, dh, dm, dm0;
427 dl = bb_i->lower[BB_X] - bb_j->upper[BB_X];
428 dh = bb_j->lower[BB_X] - bb_i->upper[BB_X];
429 dm = std::max(dl, dh);
430 dm0 = std::max(dm, 0.0f);
433 dl = bb_i->lower[BB_Y] - bb_j->upper[BB_Y];
434 dh = bb_j->lower[BB_Y] - bb_i->upper[BB_Y];
435 dm = std::max(dl, dh);
436 dm0 = std::max(dm, 0.0f);
439 dl = bb_i->lower[BB_Z] - bb_j->upper[BB_Z];
440 dh = bb_j->lower[BB_Z] - bb_i->upper[BB_Z];
441 dm = std::max(dl, dh);
442 dm0 = std::max(dm, 0.0f);
448 #if NBNXN_SEARCH_BB_SIMD4
450 /* 4-wide SIMD code for bb distance for bb format xyz0 */
451 static float subc_bb_dist2_simd4(int si,
452 const nbnxn_bb_t *bb_i_ci,
454 gmx::ArrayRef<const nbnxn_bb_t> bb_j_all)
456 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
459 Simd4Float bb_i_S0, bb_i_S1;
460 Simd4Float bb_j_S0, bb_j_S1;
466 bb_i_S0 = load4(&bb_i_ci[si].lower[0]);
467 bb_i_S1 = load4(&bb_i_ci[si].upper[0]);
468 bb_j_S0 = load4(&bb_j_all[csj].lower[0]);
469 bb_j_S1 = load4(&bb_j_all[csj].upper[0]);
471 dl_S = bb_i_S0 - bb_j_S1;
472 dh_S = bb_j_S0 - bb_i_S1;
474 dm_S = max(dl_S, dh_S);
475 dm0_S = max(dm_S, simd4SetZeroF());
477 return dotProduct(dm0_S, dm0_S);
480 /* Calculate bb bounding distances of bb_i[si,...,si+3] and store them in d2 */
481 #define SUBC_BB_DIST2_SIMD4_XXXX_INNER(si, bb_i, d2) \
485 Simd4Float dx_0, dy_0, dz_0; \
486 Simd4Float dx_1, dy_1, dz_1; \
488 Simd4Float mx, my, mz; \
489 Simd4Float m0x, m0y, m0z; \
491 Simd4Float d2x, d2y, d2z; \
492 Simd4Float d2s, d2t; \
494 shi = (si)*NNBSBB_D*DIM; \
496 xi_l = load4((bb_i)+shi+0*STRIDE_PBB); \
497 yi_l = load4((bb_i)+shi+1*STRIDE_PBB); \
498 zi_l = load4((bb_i)+shi+2*STRIDE_PBB); \
499 xi_h = load4((bb_i)+shi+3*STRIDE_PBB); \
500 yi_h = load4((bb_i)+shi+4*STRIDE_PBB); \
501 zi_h = load4((bb_i)+shi+5*STRIDE_PBB); \
503 dx_0 = xi_l - xj_h; \
504 dy_0 = yi_l - yj_h; \
505 dz_0 = zi_l - zj_h; \
507 dx_1 = xj_l - xi_h; \
508 dy_1 = yj_l - yi_h; \
509 dz_1 = zj_l - zi_h; \
511 mx = max(dx_0, dx_1); \
512 my = max(dy_0, dy_1); \
513 mz = max(dz_0, dz_1); \
515 m0x = max(mx, zero); \
516 m0y = max(my, zero); \
517 m0z = max(mz, zero); \
526 store4((d2)+(si), d2t); \
529 /* 4-wide SIMD code for nsi bb distances for bb format xxxxyyyyzzzz */
530 static void subc_bb_dist2_simd4_xxxx(const float *bb_j,
531 int nsi, const float *bb_i,
534 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
537 Simd4Float xj_l, yj_l, zj_l;
538 Simd4Float xj_h, yj_h, zj_h;
539 Simd4Float xi_l, yi_l, zi_l;
540 Simd4Float xi_h, yi_h, zi_h;
546 xj_l = Simd4Float(bb_j[0*STRIDE_PBB]);
547 yj_l = Simd4Float(bb_j[1*STRIDE_PBB]);
548 zj_l = Simd4Float(bb_j[2*STRIDE_PBB]);
549 xj_h = Simd4Float(bb_j[3*STRIDE_PBB]);
550 yj_h = Simd4Float(bb_j[4*STRIDE_PBB]);
551 zj_h = Simd4Float(bb_j[5*STRIDE_PBB]);
553 /* Here we "loop" over si (0,STRIDE_PBB) from 0 to nsi with step STRIDE_PBB.
554 * But as we know the number of iterations is 1 or 2, we unroll manually.
556 SUBC_BB_DIST2_SIMD4_XXXX_INNER(0, bb_i, d2);
557 if (STRIDE_PBB < nsi)
559 SUBC_BB_DIST2_SIMD4_XXXX_INNER(STRIDE_PBB, bb_i, d2);
563 #endif /* NBNXN_SEARCH_BB_SIMD4 */
566 /* Returns if any atom pair from two clusters is within distance sqrt(rlist2) */
567 static inline gmx_bool
568 clusterpair_in_range(const nbnxn_list_work_t *work,
570 int csj, int stride, const real *x_j,
573 #if !GMX_SIMD4_HAVE_REAL
576 * All coordinates are stored as xyzxyz...
579 const real *x_i = work->x_ci;
581 for (int i = 0; i < c_nbnxnGpuClusterSize; i++)
583 int i0 = (si*c_nbnxnGpuClusterSize + i)*DIM;
584 for (int j = 0; j < c_nbnxnGpuClusterSize; j++)
586 int j0 = (csj*c_nbnxnGpuClusterSize + j)*stride;
588 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]);
599 #else /* !GMX_SIMD4_HAVE_REAL */
601 /* 4-wide SIMD version.
602 * The coordinates x_i are stored as xxxxyyyy..., x_j is stored xyzxyz...
603 * Using 8-wide AVX(2) is not faster on Intel Sandy Bridge and Haswell.
605 static_assert(c_nbnxnGpuClusterSize == 8 || c_nbnxnGpuClusterSize == 4,
606 "A cluster is hard-coded to 4/8 atoms.");
608 Simd4Real rc2_S = Simd4Real(rlist2);
610 const real *x_i = work->x_ci_simd;
612 int dim_stride = c_nbnxnGpuClusterSize*DIM;
613 Simd4Real ix_S0 = load4(x_i + si*dim_stride + 0*GMX_SIMD4_WIDTH);
614 Simd4Real iy_S0 = load4(x_i + si*dim_stride + 1*GMX_SIMD4_WIDTH);
615 Simd4Real iz_S0 = load4(x_i + si*dim_stride + 2*GMX_SIMD4_WIDTH);
617 Simd4Real ix_S1, iy_S1, iz_S1;
618 if (c_nbnxnGpuClusterSize == 8)
620 ix_S1 = load4(x_i + si*dim_stride + 3*GMX_SIMD4_WIDTH);
621 iy_S1 = load4(x_i + si*dim_stride + 4*GMX_SIMD4_WIDTH);
622 iz_S1 = load4(x_i + si*dim_stride + 5*GMX_SIMD4_WIDTH);
624 /* We loop from the outer to the inner particles to maximize
625 * the chance that we find a pair in range quickly and return.
627 int j0 = csj*c_nbnxnGpuClusterSize;
628 int j1 = j0 + c_nbnxnGpuClusterSize - 1;
631 Simd4Real jx0_S, jy0_S, jz0_S;
632 Simd4Real jx1_S, jy1_S, jz1_S;
634 Simd4Real dx_S0, dy_S0, dz_S0;
635 Simd4Real dx_S1, dy_S1, dz_S1;
636 Simd4Real dx_S2, dy_S2, dz_S2;
637 Simd4Real dx_S3, dy_S3, dz_S3;
648 Simd4Bool wco_any_S01, wco_any_S23, wco_any_S;
650 jx0_S = Simd4Real(x_j[j0*stride+0]);
651 jy0_S = Simd4Real(x_j[j0*stride+1]);
652 jz0_S = Simd4Real(x_j[j0*stride+2]);
654 jx1_S = Simd4Real(x_j[j1*stride+0]);
655 jy1_S = Simd4Real(x_j[j1*stride+1]);
656 jz1_S = Simd4Real(x_j[j1*stride+2]);
658 /* Calculate distance */
659 dx_S0 = ix_S0 - jx0_S;
660 dy_S0 = iy_S0 - jy0_S;
661 dz_S0 = iz_S0 - jz0_S;
662 dx_S2 = ix_S0 - jx1_S;
663 dy_S2 = iy_S0 - jy1_S;
664 dz_S2 = iz_S0 - jz1_S;
665 if (c_nbnxnGpuClusterSize == 8)
667 dx_S1 = ix_S1 - jx0_S;
668 dy_S1 = iy_S1 - jy0_S;
669 dz_S1 = iz_S1 - jz0_S;
670 dx_S3 = ix_S1 - jx1_S;
671 dy_S3 = iy_S1 - jy1_S;
672 dz_S3 = iz_S1 - jz1_S;
675 /* rsq = dx*dx+dy*dy+dz*dz */
676 rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
677 rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
678 if (c_nbnxnGpuClusterSize == 8)
680 rsq_S1 = norm2(dx_S1, dy_S1, dz_S1);
681 rsq_S3 = norm2(dx_S3, dy_S3, dz_S3);
684 wco_S0 = (rsq_S0 < rc2_S);
685 wco_S2 = (rsq_S2 < rc2_S);
686 if (c_nbnxnGpuClusterSize == 8)
688 wco_S1 = (rsq_S1 < rc2_S);
689 wco_S3 = (rsq_S3 < rc2_S);
691 if (c_nbnxnGpuClusterSize == 8)
693 wco_any_S01 = wco_S0 || wco_S1;
694 wco_any_S23 = wco_S2 || wco_S3;
695 wco_any_S = wco_any_S01 || wco_any_S23;
699 wco_any_S = wco_S0 || wco_S2;
702 if (anyTrue(wco_any_S))
713 #endif /* !GMX_SIMD4_HAVE_REAL */
716 /* Returns the j-cluster index for index cjIndex in a cj list */
717 static inline int nblCj(const nbnxn_cj_t *cjList, int cjIndex)
719 return cjList[cjIndex].cj;
722 /* Returns the j-cluster index for index cjIndex in a cj4 list */
723 static inline int nblCj(const nbnxn_cj4_t *cj4List, int cjIndex)
725 return cj4List[cjIndex/c_nbnxnGpuJgroupSize].cj[cjIndex & (c_nbnxnGpuJgroupSize - 1)];
728 /* Returns the i-interaction mask of the j sub-cell for index cj_ind */
729 static unsigned int nbl_imask0(const nbnxn_pairlist_t *nbl, int cj_ind)
731 return nbl->cj4[cj_ind/c_nbnxnGpuJgroupSize].imei[0].imask;
734 /* Ensures there is enough space for extra extra exclusion masks */
735 static void check_excl_space(nbnxn_pairlist_t *nbl, int extra)
737 if (nbl->nexcl+extra > nbl->excl_nalloc)
739 nbl->excl_nalloc = over_alloc_small(nbl->nexcl+extra);
740 nbnxn_realloc_void((void **)&nbl->excl,
741 nbl->nexcl*sizeof(*nbl->excl),
742 nbl->excl_nalloc*sizeof(*nbl->excl),
743 nbl->alloc, nbl->free);
747 /* Ensures there is enough space for maxNumExtraClusters extra j-clusters in the list */
748 static void check_cell_list_space_simple(nbnxn_pairlist_t *nbl,
749 int maxNumExtraClusters)
753 cj_max = nbl->ncj + maxNumExtraClusters;
755 if (cj_max > nbl->cj_nalloc)
757 nbl->cj_nalloc = over_alloc_small(cj_max);
758 nbnxn_realloc_void((void **)&nbl->cj,
759 nbl->ncj*sizeof(*nbl->cj),
760 nbl->cj_nalloc*sizeof(*nbl->cj),
761 nbl->alloc, nbl->free);
763 nbnxn_realloc_void((void **)&nbl->cjOuter,
764 nbl->ncj*sizeof(*nbl->cjOuter),
765 nbl->cj_nalloc*sizeof(*nbl->cjOuter),
766 nbl->alloc, nbl->free);
770 /* Ensures there is enough space for ncell extra j-clusters in the list */
771 static void check_cell_list_space_supersub(nbnxn_pairlist_t *nbl,
776 /* We can have maximally nsupercell*c_gpuNumClusterPerCell sj lists */
777 /* We can store 4 j-subcell - i-supercell pairs in one struct.
778 * since we round down, we need one extra entry.
780 ncj4_max = ((nbl->work->cj_ind + ncell*c_gpuNumClusterPerCell + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize);
782 if (ncj4_max > nbl->cj4_nalloc)
784 nbl->cj4_nalloc = over_alloc_small(ncj4_max);
785 nbnxn_realloc_void((void **)&nbl->cj4,
786 nbl->work->cj4_init*sizeof(*nbl->cj4),
787 nbl->cj4_nalloc*sizeof(*nbl->cj4),
788 nbl->alloc, nbl->free);
791 if (ncj4_max > nbl->work->cj4_init)
793 for (int j4 = nbl->work->cj4_init; j4 < ncj4_max; j4++)
795 /* No i-subcells and no excl's in the list initially */
796 for (w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
798 nbl->cj4[j4].imei[w].imask = 0U;
799 nbl->cj4[j4].imei[w].excl_ind = 0;
803 nbl->work->cj4_init = ncj4_max;
807 /* Set all excl masks for one GPU warp no exclusions */
808 static void set_no_excls(nbnxn_excl_t *excl)
810 for (int t = 0; t < c_nbnxnGpuExclSize; t++)
812 /* Turn all interaction bits on */
813 excl->pair[t] = NBNXN_INTERACTION_MASK_ALL;
817 /* Initializes a single nbnxn_pairlist_t data structure */
818 static void nbnxn_init_pairlist(nbnxn_pairlist_t *nbl,
820 nbnxn_alloc_t *alloc,
823 if (alloc == nullptr)
825 nbl->alloc = nbnxn_alloc_aligned;
833 nbl->free = nbnxn_free_aligned;
840 nbl->bSimple = bSimple;
855 /* We need one element extra in sj, so alloc initially with 1 */
862 GMX_ASSERT(c_nbnxnGpuNumClusterPerSupercluster == c_gpuNumClusterPerCell, "The search code assumes that the a super-cluster matches a search grid cell");
864 GMX_ASSERT(sizeof(nbl->cj4[0].imei[0].imask)*8 >= c_nbnxnGpuJgroupSize*c_gpuNumClusterPerCell, "The i super-cluster cluster interaction mask does not contain a sufficient number of bits");
865 GMX_ASSERT(sizeof(nbl->excl[0])*8 >= c_nbnxnGpuJgroupSize*c_gpuNumClusterPerCell, "The GPU exclusion mask does not contain a sufficient number of bits");
868 nbl->excl_nalloc = 0;
870 check_excl_space(nbl, 1);
872 set_no_excls(&nbl->excl[0]);
878 snew_aligned(nbl->work->bb_ci, 1, NBNXN_SEARCH_BB_MEM_ALIGN);
883 snew_aligned(nbl->work->pbb_ci, c_gpuNumClusterPerCell/STRIDE_PBB*NNBSBB_XXXX, NBNXN_SEARCH_BB_MEM_ALIGN);
885 snew_aligned(nbl->work->bb_ci, c_gpuNumClusterPerCell, NBNXN_SEARCH_BB_MEM_ALIGN);
888 int gpu_clusterpair_nc = c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize*DIM;
889 snew(nbl->work->x_ci, gpu_clusterpair_nc);
891 snew_aligned(nbl->work->x_ci_simd,
892 std::max(NBNXN_CPU_CLUSTER_I_SIZE*DIM*GMX_SIMD_REAL_WIDTH,
894 GMX_SIMD_REAL_WIDTH);
896 snew_aligned(nbl->work->d2, c_gpuNumClusterPerCell, NBNXN_SEARCH_BB_MEM_ALIGN);
898 nbl->work->sort = nullptr;
899 nbl->work->sort_nalloc = 0;
900 nbl->work->sci_sort = nullptr;
901 nbl->work->sci_sort_nalloc = 0;
904 void nbnxn_init_pairlist_set(nbnxn_pairlist_set_t *nbl_list,
905 gmx_bool bSimple, gmx_bool bCombined,
906 nbnxn_alloc_t *alloc,
909 nbl_list->bSimple = bSimple;
910 nbl_list->bCombined = bCombined;
912 nbl_list->nnbl = gmx_omp_nthreads_get(emntNonbonded);
914 if (!nbl_list->bCombined &&
915 nbl_list->nnbl > NBNXN_BUFFERFLAG_MAX_THREADS)
917 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.",
918 nbl_list->nnbl, NBNXN_BUFFERFLAG_MAX_THREADS, NBNXN_BUFFERFLAG_MAX_THREADS);
921 snew(nbl_list->nbl, nbl_list->nnbl);
922 if (bSimple && nbl_list->nnbl > 1)
924 snew(nbl_list->nbl_work, nbl_list->nnbl);
926 snew(nbl_list->nbl_fep, nbl_list->nnbl);
927 /* Execute in order to avoid memory interleaving between threads */
928 #pragma omp parallel for num_threads(nbl_list->nnbl) schedule(static)
929 for (int i = 0; i < nbl_list->nnbl; i++)
933 /* Allocate the nblist data structure locally on each thread
934 * to optimize memory access for NUMA architectures.
936 snew(nbl_list->nbl[i], 1);
938 /* Only list 0 is used on the GPU, use normal allocation for i>0 */
939 if (!bSimple && i == 0)
941 nbnxn_init_pairlist(nbl_list->nbl[i], nbl_list->bSimple, alloc, free);
945 nbnxn_init_pairlist(nbl_list->nbl[i], nbl_list->bSimple, nullptr, nullptr);
946 if (bSimple && nbl_list->nnbl > 1)
948 snew(nbl_list->nbl_work[i], 1);
949 nbnxn_init_pairlist(nbl_list->nbl_work[i], nbl_list->bSimple, nullptr, nullptr);
953 snew(nbl_list->nbl_fep[i], 1);
954 nbnxn_init_pairlist_fep(nbl_list->nbl_fep[i]);
956 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
960 /* Print statistics of a pair list, used for debug output */
961 static void print_nblist_statistics_simple(FILE *fp, const nbnxn_pairlist_t *nbl,
962 const nbnxn_search *nbs, real rl)
964 const nbnxn_grid_t *grid;
968 grid = &nbs->grid[0];
970 fprintf(fp, "nbl nci %d ncj %d\n",
971 nbl->nci, nbl->ncjInUse);
972 fprintf(fp, "nbl na_sc %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
973 nbl->na_sc, rl, nbl->ncjInUse, nbl->ncjInUse/(double)grid->nc,
974 nbl->ncjInUse/(double)grid->nc*grid->na_sc,
975 nbl->ncjInUse/(double)grid->nc*grid->na_sc/(0.5*4.0/3.0*M_PI*rl*rl*rl*grid->nc*grid->na_sc/(grid->size[XX]*grid->size[YY]*grid->size[ZZ])));
977 fprintf(fp, "nbl average j cell list length %.1f\n",
978 0.25*nbl->ncjInUse/(double)std::max(nbl->nci, 1));
980 for (int s = 0; s < SHIFTS; s++)
985 for (int i = 0; i < nbl->nci; i++)
987 cs[nbl->ci[i].shift & NBNXN_CI_SHIFT] +=
988 nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start;
990 int j = nbl->ci[i].cj_ind_start;
991 while (j < nbl->ci[i].cj_ind_end &&
992 nbl->cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
998 fprintf(fp, "nbl cell pairs, total: %d excl: %d %.1f%%\n",
999 nbl->ncj, npexcl, 100*npexcl/(double)std::max(nbl->ncj, 1));
1000 for (int s = 0; s < SHIFTS; s++)
1004 fprintf(fp, "nbl shift %2d ncj %3d\n", s, cs[s]);
1009 /* Print statistics of a pair lists, used for debug output */
1010 static void print_nblist_statistics_supersub(FILE *fp, const nbnxn_pairlist_t *nbl,
1011 const nbnxn_search *nbs, real rl)
1013 const nbnxn_grid_t *grid;
1015 int c[c_gpuNumClusterPerCell + 1];
1016 double sum_nsp, sum_nsp2;
1019 /* This code only produces correct statistics with domain decomposition */
1020 grid = &nbs->grid[0];
1022 fprintf(fp, "nbl nsci %d ncj4 %d nsi %d excl4 %d\n",
1023 nbl->nsci, nbl->ncj4, nbl->nci_tot, nbl->nexcl);
1024 fprintf(fp, "nbl na_c %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
1025 nbl->na_ci, rl, nbl->nci_tot, nbl->nci_tot/(double)grid->nsubc_tot,
1026 nbl->nci_tot/(double)grid->nsubc_tot*grid->na_c,
1027 nbl->nci_tot/(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])));
1032 for (int si = 0; si <= c_gpuNumClusterPerCell; si++)
1036 for (int i = 0; i < nbl->nsci; i++)
1041 for (int j4 = nbl->sci[i].cj4_ind_start; j4 < nbl->sci[i].cj4_ind_end; j4++)
1043 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
1046 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
1048 if (nbl->cj4[j4].imei[0].imask & (1U << (j*c_gpuNumClusterPerCell + si)))
1058 sum_nsp2 += nsp*nsp;
1059 nsp_max = std::max(nsp_max, nsp);
1063 sum_nsp /= nbl->nsci;
1064 sum_nsp2 /= nbl->nsci;
1066 fprintf(fp, "nbl #cluster-pairs: av %.1f stddev %.1f max %d\n",
1067 sum_nsp, std::sqrt(sum_nsp2 - sum_nsp*sum_nsp), nsp_max);
1071 for (b = 0; b <= c_gpuNumClusterPerCell; b++)
1073 fprintf(fp, "nbl j-list #i-subcell %d %7d %4.1f\n",
1075 100.0*c[b]/(double)(nbl->ncj4*c_nbnxnGpuJgroupSize));
1080 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp */
1081 static void low_get_nbl_exclusions(nbnxn_pairlist_t *nbl, int cj4,
1082 int warp, nbnxn_excl_t **excl)
1084 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
1086 /* No exclusions set, make a new list entry */
1087 nbl->cj4[cj4].imei[warp].excl_ind = nbl->nexcl;
1089 *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
1090 set_no_excls(*excl);
1094 /* We already have some exclusions, new ones can be added to the list */
1095 *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
1099 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp,
1100 * generates a new element and allocates extra memory, if necessary.
1102 static void get_nbl_exclusions_1(nbnxn_pairlist_t *nbl, int cj4,
1103 int warp, nbnxn_excl_t **excl)
1105 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
1107 /* We need to make a new list entry, check if we have space */
1108 check_excl_space(nbl, 1);
1110 low_get_nbl_exclusions(nbl, cj4, warp, excl);
1113 /* Returns pointers to the exclusion masks for cj4-unit cj4 for both warps,
1114 * generates a new element and allocates extra memory, if necessary.
1116 static void get_nbl_exclusions_2(nbnxn_pairlist_t *nbl, int cj4,
1117 nbnxn_excl_t **excl_w0,
1118 nbnxn_excl_t **excl_w1)
1120 /* Check for space we might need */
1121 check_excl_space(nbl, 2);
1123 low_get_nbl_exclusions(nbl, cj4, 0, excl_w0);
1124 low_get_nbl_exclusions(nbl, cj4, 1, excl_w1);
1127 /* Sets the self exclusions i=j and pair exclusions i>j */
1128 static void set_self_and_newton_excls_supersub(nbnxn_pairlist_t *nbl,
1129 int cj4_ind, int sj_offset,
1130 int i_cluster_in_cell)
1132 nbnxn_excl_t *excl[c_nbnxnGpuClusterpairSplit];
1134 /* Here we only set the set self and double pair exclusions */
1136 static_assert(c_nbnxnGpuClusterpairSplit == 2, "");
1138 get_nbl_exclusions_2(nbl, cj4_ind, &excl[0], &excl[1]);
1140 /* Only minor < major bits set */
1141 for (int ej = 0; ej < nbl->na_ci; ej++)
1144 for (int ei = ej; ei < nbl->na_ci; ei++)
1146 excl[w]->pair[(ej & (c_nbnxnGpuJgroupSize-1))*nbl->na_ci + ei] &=
1147 ~(1U << (sj_offset*c_gpuNumClusterPerCell + i_cluster_in_cell));
1152 /* Returns a diagonal or off-diagonal interaction mask for plain C lists */
1153 static unsigned int get_imask(gmx_bool rdiag, int ci, int cj)
1155 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
1158 /* Returns a diagonal or off-diagonal interaction mask for cj-size=2 */
1159 gmx_unused static unsigned int get_imask_simd_j2(gmx_bool rdiag, int ci, int cj)
1161 return (rdiag && ci*2 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_0 :
1162 (rdiag && ci*2+1 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_1 :
1163 NBNXN_INTERACTION_MASK_ALL));
1166 /* Returns a diagonal or off-diagonal interaction mask for cj-size=4 */
1167 gmx_unused static unsigned int get_imask_simd_j4(gmx_bool rdiag, int ci, int cj)
1169 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
1172 /* Returns a diagonal or off-diagonal interaction mask for cj-size=8 */
1173 gmx_unused static unsigned int get_imask_simd_j8(gmx_bool rdiag, int ci, int cj)
1175 return (rdiag && ci == cj*2 ? NBNXN_INTERACTION_MASK_DIAG_J8_0 :
1176 (rdiag && ci == cj*2+1 ? NBNXN_INTERACTION_MASK_DIAG_J8_1 :
1177 NBNXN_INTERACTION_MASK_ALL));
1181 #if GMX_SIMD_REAL_WIDTH == 2
1182 #define get_imask_simd_4xn get_imask_simd_j2
1184 #if GMX_SIMD_REAL_WIDTH == 4
1185 #define get_imask_simd_4xn get_imask_simd_j4
1187 #if GMX_SIMD_REAL_WIDTH == 8
1188 #define get_imask_simd_4xn get_imask_simd_j8
1189 #define get_imask_simd_2xnn get_imask_simd_j4
1191 #if GMX_SIMD_REAL_WIDTH == 16
1192 #define get_imask_simd_2xnn get_imask_simd_j8
1196 /* Plain C code for checking and adding cluster-pairs to the list.
1198 * \param[in] gridj The j-grid
1199 * \param[in,out] nbl The pair-list to store the cluster pairs in
1200 * \param[in] icluster The index of the i-cluster
1201 * \param[in] jclusterFirst The first cluster in the j-range
1202 * \param[in] jclusterLast The last cluster in the j-range
1203 * \param[in] excludeSubDiagonal Exclude atom pairs with i-index > j-index
1204 * \param[in] x_j Coordinates for the j-atom, in xyz format
1205 * \param[in] rlist2 The squared list cut-off
1206 * \param[in] rbb2 The squared cut-off for putting cluster-pairs in the list based on bounding box distance only
1207 * \param[in,out] numDistanceChecks The number of distance checks performed
1210 makeClusterListSimple(const nbnxn_grid_t * gridj,
1211 nbnxn_pairlist_t * nbl,
1215 bool excludeSubDiagonal,
1216 const real * gmx_restrict x_j,
1219 int * gmx_restrict numDistanceChecks)
1221 const nbnxn_bb_t * gmx_restrict bb_ci = nbl->work->bb_ci;
1222 const real * gmx_restrict x_ci = nbl->work->x_ci;
1227 while (!InRange && jclusterFirst <= jclusterLast)
1229 real d2 = subc_bb_dist2(0, bb_ci, jclusterFirst, gridj->bb);
1230 *numDistanceChecks += 2;
1232 /* Check if the distance is within the distance where
1233 * we use only the bounding box distance rbb,
1234 * or within the cut-off and there is at least one atom pair
1235 * within the cut-off.
1241 else if (d2 < rlist2)
1243 int cjf_gl = gridj->cell0 + jclusterFirst;
1244 for (int i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
1246 for (int j = 0; j < NBNXN_CPU_CLUSTER_I_SIZE; j++)
1248 InRange = InRange ||
1249 (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
1250 gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
1251 gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rlist2);
1254 *numDistanceChecks += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
1267 while (!InRange && jclusterLast > jclusterFirst)
1269 real d2 = subc_bb_dist2(0, bb_ci, jclusterLast, gridj->bb);
1270 *numDistanceChecks += 2;
1272 /* Check if the distance is within the distance where
1273 * we use only the bounding box distance rbb,
1274 * or within the cut-off and there is at least one atom pair
1275 * within the cut-off.
1281 else if (d2 < rlist2)
1283 int cjl_gl = gridj->cell0 + jclusterLast;
1284 for (int i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
1286 for (int j = 0; j < NBNXN_CPU_CLUSTER_I_SIZE; j++)
1288 InRange = InRange ||
1289 (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
1290 gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
1291 gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rlist2);
1294 *numDistanceChecks += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
1302 if (jclusterFirst <= jclusterLast)
1304 for (int jcluster = jclusterFirst; jcluster <= jclusterLast; jcluster++)
1306 /* Store cj and the interaction mask */
1307 nbl->cj[nbl->ncj].cj = gridj->cell0 + jcluster;
1308 nbl->cj[nbl->ncj].excl = get_imask(excludeSubDiagonal, icluster, jcluster);
1311 /* Increase the closing index in i super-cell list */
1312 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
1316 #ifdef GMX_NBNXN_SIMD_4XN
1317 #include "gromacs/mdlib/nbnxn_search_simd_4xn.h"
1319 #ifdef GMX_NBNXN_SIMD_2XNN
1320 #include "gromacs/mdlib/nbnxn_search_simd_2xnn.h"
1323 /* Plain C or SIMD4 code for making a pair list of super-cell sci vs scj.
1324 * Checks bounding box distances and possibly atom pair distances.
1326 static void make_cluster_list_supersub(const nbnxn_grid_t *gridi,
1327 const nbnxn_grid_t *gridj,
1328 nbnxn_pairlist_t *nbl,
1330 gmx_bool sci_equals_scj,
1331 int stride, const real *x,
1332 real rlist2, float rbb2,
1333 int *numDistanceChecks)
1335 nbnxn_list_work_t *work = nbl->work;
1338 const float *pbb_ci = work->pbb_ci;
1340 const nbnxn_bb_t *bb_ci = work->bb_ci;
1343 assert(c_nbnxnGpuClusterSize == gridi->na_c);
1344 assert(c_nbnxnGpuClusterSize == gridj->na_c);
1346 /* We generate the pairlist mainly based on bounding-box distances
1347 * and do atom pair distance based pruning on the GPU.
1348 * Only if a j-group contains a single cluster-pair, we try to prune
1349 * that pair based on atom distances on the CPU to avoid empty j-groups.
1351 #define PRUNE_LIST_CPU_ONE 1
1352 #define PRUNE_LIST_CPU_ALL 0
1354 #if PRUNE_LIST_CPU_ONE
1358 float *d2l = work->d2;
1360 for (int subc = 0; subc < gridj->nsubc[scj]; subc++)
1362 int cj4_ind = nbl->work->cj_ind/c_nbnxnGpuJgroupSize;
1363 int cj_offset = nbl->work->cj_ind - cj4_ind*c_nbnxnGpuJgroupSize;
1364 nbnxn_cj4_t *cj4 = &nbl->cj4[cj4_ind];
1366 int cj = scj*c_gpuNumClusterPerCell + subc;
1368 int cj_gl = gridj->cell0*c_gpuNumClusterPerCell + cj;
1370 /* Initialize this j-subcell i-subcell list */
1371 cj4->cj[cj_offset] = cj_gl;
1380 ci1 = gridi->nsubc[sci];
1384 /* Determine all ci1 bb distances in one call with SIMD4 */
1385 subc_bb_dist2_simd4_xxxx(gridj->pbb.data() + (cj >> STRIDE_PBB_2LOG)*NNBSBB_XXXX + (cj & (STRIDE_PBB-1)),
1387 *numDistanceChecks += c_nbnxnGpuClusterSize*2;
1391 unsigned int imask = 0;
1392 /* We use a fixed upper-bound instead of ci1 to help optimization */
1393 for (int ci = 0; ci < c_gpuNumClusterPerCell; ci++)
1401 /* Determine the bb distance between ci and cj */
1402 d2l[ci] = subc_bb_dist2(ci, bb_ci, cj, gridj->bb);
1403 *numDistanceChecks += 2;
1407 #if PRUNE_LIST_CPU_ALL
1408 /* Check if the distance is within the distance where
1409 * we use only the bounding box distance rbb,
1410 * or within the cut-off and there is at least one atom pair
1411 * within the cut-off. This check is very costly.
1413 *numDistanceChecks += c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize;
1416 clusterpair_in_range(work, ci, cj_gl, stride, x, rlist2)))
1418 /* Check if the distance between the two bounding boxes
1419 * in within the pair-list cut-off.
1424 /* Flag this i-subcell to be taken into account */
1425 imask |= (1U << (cj_offset*c_gpuNumClusterPerCell + ci));
1427 #if PRUNE_LIST_CPU_ONE
1435 #if PRUNE_LIST_CPU_ONE
1436 /* If we only found 1 pair, check if any atoms are actually
1437 * within the cut-off, so we could get rid of it.
1439 if (npair == 1 && d2l[ci_last] >= rbb2 &&
1440 !clusterpair_in_range(work, ci_last, cj_gl, stride, x, rlist2))
1442 imask &= ~(1U << (cj_offset*c_gpuNumClusterPerCell + ci_last));
1449 /* We have a useful sj entry, close it now */
1451 /* Set the exclusions for the ci==sj entry.
1452 * Here we don't bother to check if this entry is actually flagged,
1453 * as it will nearly always be in the list.
1457 set_self_and_newton_excls_supersub(nbl, cj4_ind, cj_offset, subc);
1460 /* Copy the cluster interaction mask to the list */
1461 for (int w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
1463 cj4->imei[w].imask |= imask;
1466 nbl->work->cj_ind++;
1468 /* Keep the count */
1469 nbl->nci_tot += npair;
1471 /* Increase the closing index in i super-cell list */
1472 nbl->sci[nbl->nsci].cj4_ind_end =
1473 (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize;
1478 /* Returns how many contiguous j-clusters we have starting in the i-list */
1479 template <typename CjListType>
1480 static int numContiguousJClusters(const int cjIndexStart,
1481 const int cjIndexEnd,
1482 const CjListType &cjList)
1484 const int firstJCluster = nblCj(cjList, cjIndexStart);
1486 int numContiguous = 0;
1488 while (cjIndexStart + numContiguous < cjIndexEnd &&
1489 nblCj(cjList, cjIndexStart + numContiguous) == firstJCluster + numContiguous)
1494 return numContiguous;
1497 /* Helper struct for efficient searching for excluded atoms in a j-list */
1501 template <typename CjListType>
1502 JListRanges(int cjIndexStart,
1504 const CjListType &cjList);
1506 int cjIndexStart; // The start index in the j-list
1507 int cjIndexEnd; // The end index in the j-list
1508 int cjFirst; // The j-cluster with index cjIndexStart
1509 int cjLast; // The j-cluster with index cjIndexEnd-1
1510 int numDirect; // Up to cjIndexStart+numDirect the j-clusters are cjFirst + the index offset
1513 template <typename CjListType>
1514 JListRanges::JListRanges(int cjIndexStart,
1516 const CjListType &cjList) :
1517 cjIndexStart(cjIndexStart),
1518 cjIndexEnd(cjIndexEnd)
1520 GMX_ASSERT(cjIndexEnd > cjIndexStart, "JListRanges should only be called with non-empty lists");
1522 cjFirst = nblCj(cjList, cjIndexStart);
1523 cjLast = nblCj(cjList, cjIndexEnd - 1);
1525 /* Determine how many contiguous j-cells we have starting
1526 * from the first i-cell. This number can be used to directly
1527 * calculate j-cell indices for excluded atoms.
1529 numDirect = numContiguousJClusters(cjIndexStart, cjIndexEnd, cjList);
1532 /* Return the index of \p jCluster in the given range or -1 when not present
1534 * Note: This code is executed very often and therefore performance is
1535 * important. It should be inlined and fully optimized.
1537 template <typename CjListType>
1538 static inline int findJClusterInJList(int jCluster,
1539 const JListRanges &ranges,
1540 const CjListType &cjList)
1544 if (jCluster < ranges.cjFirst + ranges.numDirect)
1546 /* We can calculate the index directly using the offset */
1547 index = ranges.cjIndexStart + jCluster - ranges.cjFirst;
1551 /* Search for jCluster using bisection */
1553 int rangeStart = ranges.cjIndexStart + ranges.numDirect;
1554 int rangeEnd = ranges.cjIndexEnd;
1556 while (index == -1 && rangeStart < rangeEnd)
1558 rangeMiddle = (rangeStart + rangeEnd) >> 1;
1560 const int clusterMiddle = nblCj(cjList, rangeMiddle);
1562 if (jCluster == clusterMiddle)
1564 index = rangeMiddle;
1566 else if (jCluster < clusterMiddle)
1568 rangeEnd = rangeMiddle;
1572 rangeStart = rangeMiddle + 1;
1580 /* Set all atom-pair exclusions for a simple type list i-entry
1582 * Set all atom-pair exclusions from the topology stored in exclusions
1583 * as masks in the pair-list for simple list entry iEntry.
1586 setExclusionsForSimpleIentry(const nbnxn_search *nbs,
1587 nbnxn_pairlist_t *nbl,
1588 gmx_bool diagRemoved,
1590 const nbnxn_ci_t &iEntry,
1591 const t_blocka &exclusions)
1593 if (iEntry.cj_ind_end == iEntry.cj_ind_start)
1595 /* Empty list: no exclusions */
1599 const JListRanges ranges(iEntry.cj_ind_start, iEntry.cj_ind_end, nbl->cj);
1601 const int iCluster = iEntry.ci;
1603 gmx::ArrayRef<const int> cell = nbs->cell;
1605 /* Loop over the atoms in the i-cluster */
1606 for (int i = 0; i < nbl->na_sc; i++)
1608 const int iIndex = iCluster*nbl->na_sc + i;
1609 const int iAtom = nbs->a[iIndex];
1612 /* Loop over the topology-based exclusions for this i-atom */
1613 for (int exclIndex = exclusions.index[iAtom]; exclIndex < exclusions.index[iAtom + 1]; exclIndex++)
1615 const int jAtom = exclusions.a[exclIndex];
1619 /* The self exclusion are already set, save some time */
1623 /* Get the index of the j-atom in the nbnxn atom data */
1624 const int jIndex = cell[jAtom];
1626 /* Without shifts we only calculate interactions j>i
1627 * for one-way pair-lists.
1629 if (diagRemoved && jIndex <= iIndex)
1634 const int jCluster = (jIndex >> na_cj_2log);
1636 /* Could the cluster se be in our list? */
1637 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
1640 findJClusterInJList(jCluster, ranges, nbl->cj);
1644 /* We found an exclusion, clear the corresponding
1647 const int innerJ = jIndex - (jCluster << na_cj_2log);
1649 nbl->cj[index].excl &= ~(1U << ((i << na_cj_2log) + innerJ));
1657 /* Add a new i-entry to the FEP list and copy the i-properties */
1658 static inline void fep_list_new_nri_copy(t_nblist *nlist)
1660 /* Add a new i-entry */
1663 assert(nlist->nri < nlist->maxnri);
1665 /* Duplicate the last i-entry, except for jindex, which continues */
1666 nlist->iinr[nlist->nri] = nlist->iinr[nlist->nri-1];
1667 nlist->shift[nlist->nri] = nlist->shift[nlist->nri-1];
1668 nlist->gid[nlist->nri] = nlist->gid[nlist->nri-1];
1669 nlist->jindex[nlist->nri] = nlist->nrj;
1672 /* For load balancing of the free-energy lists over threads, we set
1673 * the maximum nrj size of an i-entry to 40. This leads to good
1674 * load balancing in the worst case scenario of a single perturbed
1675 * particle on 16 threads, while not introducing significant overhead.
1676 * Note that half of the perturbed pairs will anyhow end up in very small lists,
1677 * since non perturbed i-particles will see few perturbed j-particles).
1679 const int max_nrj_fep = 40;
1681 /* Exclude the perturbed pairs from the Verlet list. This is only done to avoid
1682 * singularities for overlapping particles (0/0), since the charges and
1683 * LJ parameters have been zeroed in the nbnxn data structure.
1684 * Simultaneously make a group pair list for the perturbed pairs.
1686 static void make_fep_list(const nbnxn_search *nbs,
1687 const nbnxn_atomdata_t *nbat,
1688 nbnxn_pairlist_t *nbl,
1689 gmx_bool bDiagRemoved,
1691 const nbnxn_grid_t *gridi,
1692 const nbnxn_grid_t *gridj,
1695 int ci, cj_ind_start, cj_ind_end, cja, cjr;
1697 int ngid, gid_i = 0, gid_j, gid;
1698 int egp_shift, egp_mask;
1700 int ind_i, ind_j, ai, aj;
1702 gmx_bool bFEP_i, bFEP_i_all;
1704 if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
1712 cj_ind_start = nbl_ci->cj_ind_start;
1713 cj_ind_end = nbl_ci->cj_ind_end;
1715 /* In worst case we have alternating energy groups
1716 * and create #atom-pair lists, which means we need the size
1717 * of a cluster pair (na_ci*na_cj) times the number of cj's.
1719 nri_max = nbl->na_ci*nbl->na_cj*(cj_ind_end - cj_ind_start);
1720 if (nlist->nri + nri_max > nlist->maxnri)
1722 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1723 reallocate_nblist(nlist);
1726 ngid = nbat->nenergrp;
1728 if (ngid*gridj->na_cj > gmx::index(sizeof(gid_cj)*8))
1730 gmx_fatal(FARGS, "The Verlet scheme with %dx%d kernels and free-energy only supports up to %lu energy groups",
1731 gridi->na_c, gridj->na_cj, (sizeof(gid_cj)*8)/gridj->na_cj);
1734 egp_shift = nbat->neg_2log;
1735 egp_mask = (1<<nbat->neg_2log) - 1;
1737 /* Loop over the atoms in the i sub-cell */
1739 for (int i = 0; i < nbl->na_ci; i++)
1741 ind_i = ci*nbl->na_ci + i;
1746 nlist->jindex[nri+1] = nlist->jindex[nri];
1747 nlist->iinr[nri] = ai;
1748 /* The actual energy group pair index is set later */
1749 nlist->gid[nri] = 0;
1750 nlist->shift[nri] = nbl_ci->shift & NBNXN_CI_SHIFT;
1752 bFEP_i = gridi->fep[ci - gridi->cell0] & (1 << i);
1754 bFEP_i_all = bFEP_i_all && bFEP_i;
1756 if (nlist->nrj + (cj_ind_end - cj_ind_start)*nbl->na_cj > nlist->maxnrj)
1758 nlist->maxnrj = over_alloc_small(nlist->nrj + (cj_ind_end - cj_ind_start)*nbl->na_cj);
1759 srenew(nlist->jjnr, nlist->maxnrj);
1760 srenew(nlist->excl_fep, nlist->maxnrj);
1765 gid_i = (nbat->energrp[ci] >> (egp_shift*i)) & egp_mask;
1768 for (int cj_ind = cj_ind_start; cj_ind < cj_ind_end; cj_ind++)
1770 unsigned int fep_cj;
1772 cja = nbl->cj[cj_ind].cj;
1774 if (gridj->na_cj == gridj->na_c)
1776 cjr = cja - gridj->cell0;
1777 fep_cj = gridj->fep[cjr];
1780 gid_cj = nbat->energrp[cja];
1783 else if (2*gridj->na_cj == gridj->na_c)
1785 cjr = cja - gridj->cell0*2;
1786 /* Extract half of the ci fep/energrp mask */
1787 fep_cj = (gridj->fep[cjr>>1] >> ((cjr&1)*gridj->na_cj)) & ((1<<gridj->na_cj) - 1);
1790 gid_cj = nbat->energrp[cja>>1] >> ((cja&1)*gridj->na_cj*egp_shift) & ((1<<(gridj->na_cj*egp_shift)) - 1);
1795 cjr = cja - (gridj->cell0>>1);
1796 /* Combine two ci fep masks/energrp */
1797 fep_cj = gridj->fep[cjr*2] + (gridj->fep[cjr*2+1] << gridj->na_c);
1800 gid_cj = nbat->energrp[cja*2] + (nbat->energrp[cja*2+1] << (gridj->na_c*egp_shift));
1804 if (bFEP_i || fep_cj != 0)
1806 for (int j = 0; j < nbl->na_cj; j++)
1808 /* Is this interaction perturbed and not excluded? */
1809 ind_j = cja*nbl->na_cj + j;
1812 (bFEP_i || (fep_cj & (1 << j))) &&
1813 (!bDiagRemoved || ind_j >= ind_i))
1817 gid_j = (gid_cj >> (j*egp_shift)) & egp_mask;
1818 gid = GID(gid_i, gid_j, ngid);
1820 if (nlist->nrj > nlist->jindex[nri] &&
1821 nlist->gid[nri] != gid)
1823 /* Energy group pair changed: new list */
1824 fep_list_new_nri_copy(nlist);
1827 nlist->gid[nri] = gid;
1830 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1832 fep_list_new_nri_copy(nlist);
1836 /* Add it to the FEP list */
1837 nlist->jjnr[nlist->nrj] = aj;
1838 nlist->excl_fep[nlist->nrj] = (nbl->cj[cj_ind].excl >> (i*nbl->na_cj + j)) & 1;
1841 /* Exclude it from the normal list.
1842 * Note that the charge has been set to zero,
1843 * but we need to avoid 0/0, as perturbed atoms
1844 * can be on top of each other.
1846 nbl->cj[cj_ind].excl &= ~(1U << (i*nbl->na_cj + j));
1852 if (nlist->nrj > nlist->jindex[nri])
1854 /* Actually add this new, non-empty, list */
1856 nlist->jindex[nlist->nri] = nlist->nrj;
1863 /* All interactions are perturbed, we can skip this entry */
1864 nbl_ci->cj_ind_end = cj_ind_start;
1865 nbl->ncjInUse -= cj_ind_end - cj_ind_start;
1869 /* Return the index of atom a within a cluster */
1870 static inline int cj_mod_cj4(int cj)
1872 return cj & (c_nbnxnGpuJgroupSize - 1);
1875 /* Convert a j-cluster to a cj4 group */
1876 static inline int cj_to_cj4(int cj)
1878 return cj/c_nbnxnGpuJgroupSize;
1881 /* Return the index of an j-atom within a warp */
1882 static inline int a_mod_wj(int a)
1884 return a & (c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit - 1);
1887 /* As make_fep_list above, but for super/sub lists. */
1888 static void make_fep_list_supersub(const nbnxn_search *nbs,
1889 const nbnxn_atomdata_t *nbat,
1890 nbnxn_pairlist_t *nbl,
1891 gmx_bool bDiagRemoved,
1892 const nbnxn_sci_t *nbl_sci,
1897 const nbnxn_grid_t *gridi,
1898 const nbnxn_grid_t *gridj,
1901 int sci, cj4_ind_start, cj4_ind_end, cjr;
1904 int ind_i, ind_j, ai, aj;
1908 const nbnxn_cj4_t *cj4;
1910 if (nbl_sci->cj4_ind_end == nbl_sci->cj4_ind_start)
1918 cj4_ind_start = nbl_sci->cj4_ind_start;
1919 cj4_ind_end = nbl_sci->cj4_ind_end;
1921 /* Here we process one super-cell, max #atoms na_sc, versus a list
1922 * cj4 entries, each with max c_nbnxnGpuJgroupSize cj's, each
1923 * of size na_cj atoms.
1924 * On the GPU we don't support energy groups (yet).
1925 * So for each of the na_sc i-atoms, we need max one FEP list
1926 * for each max_nrj_fep j-atoms.
1928 nri_max = nbl->na_sc*nbl->na_cj*(1 + ((cj4_ind_end - cj4_ind_start)*c_nbnxnGpuJgroupSize)/max_nrj_fep);
1929 if (nlist->nri + nri_max > nlist->maxnri)
1931 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1932 reallocate_nblist(nlist);
1935 /* Loop over the atoms in the i super-cluster */
1936 for (int c = 0; c < c_gpuNumClusterPerCell; c++)
1938 c_abs = sci*c_gpuNumClusterPerCell + c;
1940 for (int i = 0; i < nbl->na_ci; i++)
1942 ind_i = c_abs*nbl->na_ci + i;
1947 nlist->jindex[nri+1] = nlist->jindex[nri];
1948 nlist->iinr[nri] = ai;
1949 /* With GPUs, energy groups are not supported */
1950 nlist->gid[nri] = 0;
1951 nlist->shift[nri] = nbl_sci->shift & NBNXN_CI_SHIFT;
1953 bFEP_i = (gridi->fep[c_abs - gridi->cell0*c_gpuNumClusterPerCell] & (1 << i));
1955 xi = nbat->x[ind_i*nbat->xstride+XX] + shx;
1956 yi = nbat->x[ind_i*nbat->xstride+YY] + shy;
1957 zi = nbat->x[ind_i*nbat->xstride+ZZ] + shz;
1959 if ((nlist->nrj + cj4_ind_end - cj4_ind_start)*c_nbnxnGpuJgroupSize*nbl->na_cj > nlist->maxnrj)
1961 nlist->maxnrj = over_alloc_small((nlist->nrj + cj4_ind_end - cj4_ind_start)*c_nbnxnGpuJgroupSize*nbl->na_cj);
1962 srenew(nlist->jjnr, nlist->maxnrj);
1963 srenew(nlist->excl_fep, nlist->maxnrj);
1966 for (int cj4_ind = cj4_ind_start; cj4_ind < cj4_ind_end; cj4_ind++)
1968 cj4 = &nbl->cj4[cj4_ind];
1970 for (int gcj = 0; gcj < c_nbnxnGpuJgroupSize; gcj++)
1972 unsigned int fep_cj;
1974 if ((cj4->imei[0].imask & (1U << (gcj*c_gpuNumClusterPerCell + c))) == 0)
1976 /* Skip this ci for this cj */
1980 cjr = cj4->cj[gcj] - gridj->cell0*c_gpuNumClusterPerCell;
1982 fep_cj = gridj->fep[cjr];
1984 if (bFEP_i || fep_cj != 0)
1986 for (int j = 0; j < nbl->na_cj; j++)
1988 /* Is this interaction perturbed and not excluded? */
1989 ind_j = (gridj->cell0*c_gpuNumClusterPerCell + cjr)*nbl->na_cj + j;
1992 (bFEP_i || (fep_cj & (1 << j))) &&
1993 (!bDiagRemoved || ind_j >= ind_i))
1997 unsigned int excl_bit;
2000 const int jHalf = j/(c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit);
2001 get_nbl_exclusions_1(nbl, cj4_ind, jHalf, &excl);
2003 excl_pair = a_mod_wj(j)*nbl->na_ci + i;
2004 excl_bit = (1U << (gcj*c_gpuNumClusterPerCell + c));
2006 dx = nbat->x[ind_j*nbat->xstride+XX] - xi;
2007 dy = nbat->x[ind_j*nbat->xstride+YY] - yi;
2008 dz = nbat->x[ind_j*nbat->xstride+ZZ] - zi;
2010 /* The unpruned GPU list has more than 2/3
2011 * of the atom pairs beyond rlist. Using
2012 * this list will cause a lot of overhead
2013 * in the CPU FEP kernels, especially
2014 * relative to the fast GPU kernels.
2015 * So we prune the FEP list here.
2017 if (dx*dx + dy*dy + dz*dz < rlist_fep2)
2019 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
2021 fep_list_new_nri_copy(nlist);
2025 /* Add it to the FEP list */
2026 nlist->jjnr[nlist->nrj] = aj;
2027 nlist->excl_fep[nlist->nrj] = (excl->pair[excl_pair] & excl_bit) ? 1 : 0;
2031 /* Exclude it from the normal list.
2032 * Note that the charge and LJ parameters have
2033 * been set to zero, but we need to avoid 0/0,
2034 * as perturbed atoms can be on top of each other.
2036 excl->pair[excl_pair] &= ~excl_bit;
2040 /* Note that we could mask out this pair in imask
2041 * if all i- and/or all j-particles are perturbed.
2042 * But since the perturbed pairs on the CPU will
2043 * take an order of magnitude more time, the GPU
2044 * will finish before the CPU and there is no gain.
2050 if (nlist->nrj > nlist->jindex[nri])
2052 /* Actually add this new, non-empty, list */
2054 nlist->jindex[nlist->nri] = nlist->nrj;
2061 /* Set all atom-pair exclusions for a GPU type list i-entry
2063 * Sets all atom-pair exclusions from the topology stored in exclusions
2064 * as masks in the pair-list for i-super-cluster list entry iEntry.
2067 setExclusionsForGpuIentry(const nbnxn_search *nbs,
2068 nbnxn_pairlist_t *nbl,
2069 gmx_bool diagRemoved,
2070 const nbnxn_sci_t &iEntry,
2071 const t_blocka &exclusions)
2073 if (iEntry.cj4_ind_end == iEntry.cj4_ind_start)
2079 /* Set the search ranges using start and end j-cluster indices.
2080 * Note that here we can not use cj4_ind_end, since the last cj4
2081 * can be only partially filled, so we use cj_ind.
2083 const JListRanges ranges(iEntry.cj4_ind_start*c_nbnxnGpuJgroupSize,
2087 GMX_ASSERT(nbl->na_ci == c_nbnxnGpuClusterSize, "na_ci should match the GPU cluster size");
2088 constexpr int c_clusterSize = c_nbnxnGpuClusterSize;
2089 constexpr int c_superClusterSize = c_nbnxnGpuNumClusterPerSupercluster*c_nbnxnGpuClusterSize;
2091 const int iSuperCluster = iEntry.sci;
2093 gmx::ArrayRef<const int> cell = nbs->cell;
2095 /* Loop over the atoms in the i super-cluster */
2096 for (int i = 0; i < c_superClusterSize; i++)
2098 const int iIndex = iSuperCluster*c_superClusterSize + i;
2099 const int iAtom = nbs->a[iIndex];
2102 const int iCluster = i/c_clusterSize;
2104 /* Loop over the topology-based exclusions for this i-atom */
2105 for (int exclIndex = exclusions.index[iAtom]; exclIndex < exclusions.index[iAtom + 1]; exclIndex++)
2107 const int jAtom = exclusions.a[exclIndex];
2111 /* The self exclusions are already set, save some time */
2115 /* Get the index of the j-atom in the nbnxn atom data */
2116 const int jIndex = cell[jAtom];
2118 /* Without shifts we only calculate interactions j>i
2119 * for one-way pair-lists.
2121 /* NOTE: We would like to use iIndex on the right hand side,
2122 * but that makes this routine 25% slower with gcc6/7.
2123 * Even using c_superClusterSize makes it slower.
2124 * Either of these changes triggers peeling of the exclIndex
2125 * loop, which apparently leads to far less efficient code.
2127 if (diagRemoved && jIndex <= iSuperCluster*nbl->na_sc + i)
2132 const int jCluster = jIndex/c_clusterSize;
2134 /* Check whether the cluster is in our list? */
2135 if (jCluster >= ranges.cjFirst && jCluster <= ranges.cjLast)
2138 findJClusterInJList(jCluster, ranges, nbl->cj4);
2142 /* We found an exclusion, clear the corresponding
2145 const unsigned int pairMask = (1U << (cj_mod_cj4(index)*c_gpuNumClusterPerCell + iCluster));
2146 /* Check if the i-cluster interacts with the j-cluster */
2147 if (nbl_imask0(nbl, index) & pairMask)
2149 const int innerI = (i & (c_clusterSize - 1));
2150 const int innerJ = (jIndex & (c_clusterSize - 1));
2152 /* Determine which j-half (CUDA warp) we are in */
2153 const int jHalf = innerJ/(c_clusterSize/c_nbnxnGpuClusterpairSplit);
2155 nbnxn_excl_t *interactionMask;
2156 get_nbl_exclusions_1(nbl, cj_to_cj4(index), jHalf, &interactionMask);
2158 interactionMask->pair[a_mod_wj(innerJ)*c_clusterSize + innerI] &= ~pairMask;
2167 /* Reallocate the simple ci list for at least n entries */
2168 static void nb_realloc_ci(nbnxn_pairlist_t *nbl, int n)
2170 nbl->ci_nalloc = over_alloc_small(n);
2171 nbnxn_realloc_void((void **)&nbl->ci,
2172 nbl->nci*sizeof(*nbl->ci),
2173 nbl->ci_nalloc*sizeof(*nbl->ci),
2174 nbl->alloc, nbl->free);
2176 nbnxn_realloc_void((void **)&nbl->ciOuter,
2177 nbl->nci*sizeof(*nbl->ciOuter),
2178 nbl->ci_nalloc*sizeof(*nbl->ciOuter),
2179 nbl->alloc, nbl->free);
2182 /* Reallocate the super-cell sci list for at least n entries */
2183 static void nb_realloc_sci(nbnxn_pairlist_t *nbl, int n)
2185 nbl->sci_nalloc = over_alloc_small(n);
2186 nbnxn_realloc_void((void **)&nbl->sci,
2187 nbl->nsci*sizeof(*nbl->sci),
2188 nbl->sci_nalloc*sizeof(*nbl->sci),
2189 nbl->alloc, nbl->free);
2192 /* Make a new ci entry at index nbl->nci */
2193 static void new_ci_entry(nbnxn_pairlist_t *nbl, int ci, int shift, int flags)
2195 if (nbl->nci + 1 > nbl->ci_nalloc)
2197 nb_realloc_ci(nbl, nbl->nci+1);
2199 nbl->ci[nbl->nci].ci = ci;
2200 nbl->ci[nbl->nci].shift = shift;
2201 /* Store the interaction flags along with the shift */
2202 nbl->ci[nbl->nci].shift |= flags;
2203 nbl->ci[nbl->nci].cj_ind_start = nbl->ncj;
2204 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
2207 /* Make a new sci entry at index nbl->nsci */
2208 static void new_sci_entry(nbnxn_pairlist_t *nbl, int sci, int shift)
2210 if (nbl->nsci + 1 > nbl->sci_nalloc)
2212 nb_realloc_sci(nbl, nbl->nsci+1);
2214 nbl->sci[nbl->nsci].sci = sci;
2215 nbl->sci[nbl->nsci].shift = shift;
2216 nbl->sci[nbl->nsci].cj4_ind_start = nbl->ncj4;
2217 nbl->sci[nbl->nsci].cj4_ind_end = nbl->ncj4;
2220 /* Sort the simple j-list cj on exclusions.
2221 * Entries with exclusions will all be sorted to the beginning of the list.
2223 static void sort_cj_excl(nbnxn_cj_t *cj, int ncj,
2224 nbnxn_list_work_t *work)
2228 if (ncj > work->cj_nalloc)
2230 work->cj_nalloc = over_alloc_large(ncj);
2231 srenew(work->cj, work->cj_nalloc);
2234 /* Make a list of the j-cells involving exclusions */
2236 for (int j = 0; j < ncj; j++)
2238 if (cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
2240 work->cj[jnew++] = cj[j];
2243 /* Check if there are exclusions at all or not just the first entry */
2244 if (!((jnew == 0) ||
2245 (jnew == 1 && cj[0].excl != NBNXN_INTERACTION_MASK_ALL)))
2247 for (int j = 0; j < ncj; j++)
2249 if (cj[j].excl == NBNXN_INTERACTION_MASK_ALL)
2251 work->cj[jnew++] = cj[j];
2254 for (int j = 0; j < ncj; j++)
2256 cj[j] = work->cj[j];
2261 /* Close this simple list i entry */
2262 static void close_ci_entry_simple(nbnxn_pairlist_t *nbl)
2266 /* All content of the new ci entry have already been filled correctly,
2267 * we only need to increase the count here (for non empty lists).
2269 jlen = nbl->ci[nbl->nci].cj_ind_end - nbl->ci[nbl->nci].cj_ind_start;
2272 sort_cj_excl(nbl->cj+nbl->ci[nbl->nci].cj_ind_start, jlen, nbl->work);
2274 /* The counts below are used for non-bonded pair/flop counts
2275 * and should therefore match the available kernel setups.
2277 if (!(nbl->ci[nbl->nci].shift & NBNXN_CI_DO_COUL(0)))
2279 nbl->work->ncj_noq += jlen;
2281 else if ((nbl->ci[nbl->nci].shift & NBNXN_CI_HALF_LJ(0)) ||
2282 !(nbl->ci[nbl->nci].shift & NBNXN_CI_DO_LJ(0)))
2284 nbl->work->ncj_hlj += jlen;
2291 /* Split sci entry for load balancing on the GPU.
2292 * Splitting ensures we have enough lists to fully utilize the whole GPU.
2293 * With progBal we generate progressively smaller lists, which improves
2294 * load balancing. As we only know the current count on our own thread,
2295 * we will need to estimate the current total amount of i-entries.
2296 * As the lists get concatenated later, this estimate depends
2297 * both on nthread and our own thread index.
2299 static void split_sci_entry(nbnxn_pairlist_t *nbl,
2301 gmx_bool progBal, float nsp_tot_est,
2302 int thread, int nthread)
2305 int cj4_start, cj4_end, j4len;
2307 int nsp, nsp_sci, nsp_cj4, nsp_cj4_e, nsp_cj4_p;
2313 /* Estimate the total numbers of ci's of the nblist combined
2314 * over all threads using the target number of ci's.
2316 nsp_est = (nsp_tot_est*thread)/nthread + nbl->nci_tot;
2318 /* The first ci blocks should be larger, to avoid overhead.
2319 * The last ci blocks should be smaller, to improve load balancing.
2320 * The factor 3/2 makes the first block 3/2 times the target average
2321 * and ensures that the total number of blocks end up equal to
2322 * that of equally sized blocks of size nsp_target_av.
2324 nsp_max = static_cast<int>(nsp_target_av*(nsp_tot_est*1.5/(nsp_est + nsp_tot_est)));
2328 nsp_max = nsp_target_av;
2331 cj4_start = nbl->sci[nbl->nsci-1].cj4_ind_start;
2332 cj4_end = nbl->sci[nbl->nsci-1].cj4_ind_end;
2333 j4len = cj4_end - cj4_start;
2335 if (j4len > 1 && j4len*c_gpuNumClusterPerCell*c_nbnxnGpuJgroupSize > nsp_max)
2337 /* Remove the last ci entry and process the cj4's again */
2345 for (int cj4 = cj4_start; cj4 < cj4_end; cj4++)
2347 nsp_cj4_p = nsp_cj4;
2348 /* Count the number of cluster pairs in this cj4 group */
2350 for (int p = 0; p < c_gpuNumClusterPerCell*c_nbnxnGpuJgroupSize; p++)
2352 nsp_cj4 += (nbl->cj4[cj4].imei[0].imask >> p) & 1;
2355 /* If adding the current cj4 with nsp_cj4 pairs get us further
2356 * away from our target nsp_max, split the list before this cj4.
2358 if (nsp > 0 && nsp_max - nsp < nsp + nsp_cj4 - nsp_max)
2360 /* Split the list at cj4 */
2361 nbl->sci[sci].cj4_ind_end = cj4;
2362 /* Create a new sci entry */
2365 if (nbl->nsci+1 > nbl->sci_nalloc)
2367 nb_realloc_sci(nbl, nbl->nsci+1);
2369 nbl->sci[sci].sci = nbl->sci[nbl->nsci-1].sci;
2370 nbl->sci[sci].shift = nbl->sci[nbl->nsci-1].shift;
2371 nbl->sci[sci].cj4_ind_start = cj4;
2373 nsp_cj4_e = nsp_cj4_p;
2379 /* Put the remaining cj4's in the last sci entry */
2380 nbl->sci[sci].cj4_ind_end = cj4_end;
2382 /* Possibly balance out the last two sci's
2383 * by moving the last cj4 of the second last sci.
2385 if (nsp_sci - nsp_cj4_e >= nsp + nsp_cj4_e)
2387 nbl->sci[sci-1].cj4_ind_end--;
2388 nbl->sci[sci].cj4_ind_start--;
2395 /* Clost this super/sub list i entry */
2396 static void close_ci_entry_supersub(nbnxn_pairlist_t *nbl,
2398 gmx_bool progBal, float nsp_tot_est,
2399 int thread, int nthread)
2401 /* All content of the new ci entry have already been filled correctly,
2402 * we only need to increase the count here (for non empty lists).
2404 int j4len = nbl->sci[nbl->nsci].cj4_ind_end - nbl->sci[nbl->nsci].cj4_ind_start;
2407 /* We can only have complete blocks of 4 j-entries in a list,
2408 * so round the count up before closing.
2410 nbl->ncj4 = (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize;
2411 nbl->work->cj_ind = nbl->ncj4*c_nbnxnGpuJgroupSize;
2417 /* Measure the size of the new entry and potentially split it */
2418 split_sci_entry(nbl, nsp_max_av, progBal, nsp_tot_est,
2424 /* Syncs the working array before adding another grid pair to the list */
2425 static void sync_work(nbnxn_pairlist_t *nbl)
2429 nbl->work->cj_ind = nbl->ncj4*c_nbnxnGpuJgroupSize;
2430 nbl->work->cj4_init = nbl->ncj4;
2434 /* Clears an nbnxn_pairlist_t data structure */
2435 static void clear_pairlist(nbnxn_pairlist_t *nbl)
2446 nbl->work->ncj_noq = 0;
2447 nbl->work->ncj_hlj = 0;
2450 /* Clears a group scheme pair list */
2451 static void clear_pairlist_fep(t_nblist *nl)
2455 if (nl->jindex == nullptr)
2457 snew(nl->jindex, 1);
2462 /* Sets a simple list i-cell bounding box, including PBC shift */
2463 static inline void set_icell_bb_simple(gmx::ArrayRef<const nbnxn_bb_t> bb,
2465 real shx, real shy, real shz,
2468 bb_ci->lower[BB_X] = bb[ci].lower[BB_X] + shx;
2469 bb_ci->lower[BB_Y] = bb[ci].lower[BB_Y] + shy;
2470 bb_ci->lower[BB_Z] = bb[ci].lower[BB_Z] + shz;
2471 bb_ci->upper[BB_X] = bb[ci].upper[BB_X] + shx;
2472 bb_ci->upper[BB_Y] = bb[ci].upper[BB_Y] + shy;
2473 bb_ci->upper[BB_Z] = bb[ci].upper[BB_Z] + shz;
2477 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2478 static void set_icell_bbxxxx_supersub(gmx::ArrayRef<const float> bb,
2480 real shx, real shy, real shz,
2483 int ia = ci*(c_gpuNumClusterPerCell >> STRIDE_PBB_2LOG)*NNBSBB_XXXX;
2484 for (int m = 0; m < (c_gpuNumClusterPerCell >> STRIDE_PBB_2LOG)*NNBSBB_XXXX; m += NNBSBB_XXXX)
2486 for (int i = 0; i < STRIDE_PBB; i++)
2488 bb_ci[m+0*STRIDE_PBB+i] = bb[ia+m+0*STRIDE_PBB+i] + shx;
2489 bb_ci[m+1*STRIDE_PBB+i] = bb[ia+m+1*STRIDE_PBB+i] + shy;
2490 bb_ci[m+2*STRIDE_PBB+i] = bb[ia+m+2*STRIDE_PBB+i] + shz;
2491 bb_ci[m+3*STRIDE_PBB+i] = bb[ia+m+3*STRIDE_PBB+i] + shx;
2492 bb_ci[m+4*STRIDE_PBB+i] = bb[ia+m+4*STRIDE_PBB+i] + shy;
2493 bb_ci[m+5*STRIDE_PBB+i] = bb[ia+m+5*STRIDE_PBB+i] + shz;
2499 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2500 gmx_unused static void set_icell_bb_supersub(gmx::ArrayRef<const nbnxn_bb_t> bb,
2502 real shx, real shy, real shz,
2505 for (int i = 0; i < c_gpuNumClusterPerCell; i++)
2507 set_icell_bb_simple(bb, ci*c_gpuNumClusterPerCell+i,
2513 /* Copies PBC shifted i-cell atom coordinates x,y,z to working array */
2514 static void icell_set_x_simple(int ci,
2515 real shx, real shy, real shz,
2516 int stride, const real *x,
2517 nbnxn_list_work_t *work)
2519 int ia = ci*NBNXN_CPU_CLUSTER_I_SIZE;
2521 for (int i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE; i++)
2523 work->x_ci[i*STRIDE_XYZ+XX] = x[(ia+i)*stride+XX] + shx;
2524 work->x_ci[i*STRIDE_XYZ+YY] = x[(ia+i)*stride+YY] + shy;
2525 work->x_ci[i*STRIDE_XYZ+ZZ] = x[(ia+i)*stride+ZZ] + shz;
2529 /* Copies PBC shifted super-cell atom coordinates x,y,z to working array */
2530 static void icell_set_x_supersub(int ci,
2531 real shx, real shy, real shz,
2532 int stride, const real *x,
2533 nbnxn_list_work_t *work)
2535 #if !GMX_SIMD4_HAVE_REAL
2537 real * x_ci = work->x_ci;
2539 int ia = ci*c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize;
2540 for (int i = 0; i < c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize; i++)
2542 x_ci[i*DIM + XX] = x[(ia+i)*stride + XX] + shx;
2543 x_ci[i*DIM + YY] = x[(ia+i)*stride + YY] + shy;
2544 x_ci[i*DIM + ZZ] = x[(ia+i)*stride + ZZ] + shz;
2547 #else /* !GMX_SIMD4_HAVE_REAL */
2549 real * x_ci = work->x_ci_simd;
2551 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2553 for (int i = 0; i < c_nbnxnGpuClusterSize; i += GMX_SIMD4_WIDTH)
2555 int io = si*c_nbnxnGpuClusterSize + i;
2556 int ia = ci*c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize + io;
2557 for (int j = 0; j < GMX_SIMD4_WIDTH; j++)
2559 x_ci[io*DIM + j + XX*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + XX] + shx;
2560 x_ci[io*DIM + j + YY*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + YY] + shy;
2561 x_ci[io*DIM + j + ZZ*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + ZZ] + shz;
2566 #endif /* !GMX_SIMD4_HAVE_REAL */
2569 static real minimum_subgrid_size_xy(const nbnxn_grid_t *grid)
2573 return std::min(grid->cellSize[XX], grid->cellSize[YY]);
2577 return std::min(grid->cellSize[XX]/c_gpuNumClusterPerCellX,
2578 grid->cellSize[YY]/c_gpuNumClusterPerCellY);
2582 static real effective_buffer_1x1_vs_MxN(const nbnxn_grid_t *gridi,
2583 const nbnxn_grid_t *gridj)
2585 const real eff_1x1_buffer_fac_overest = 0.1;
2587 /* Determine an atom-pair list cut-off buffer size for atom pairs,
2588 * to be added to rlist (including buffer) used for MxN.
2589 * This is for converting an MxN list to a 1x1 list. This means we can't
2590 * use the normal buffer estimate, as we have an MxN list in which
2591 * some atom pairs beyond rlist are missing. We want to capture
2592 * the beneficial effect of buffering by extra pairs just outside rlist,
2593 * while removing the useless pairs that are further away from rlist.
2594 * (Also the buffer could have been set manually not using the estimate.)
2595 * This buffer size is an overestimate.
2596 * We add 10% of the smallest grid sub-cell dimensions.
2597 * Note that the z-size differs per cell and we don't use this,
2598 * so we overestimate.
2599 * With PME, the 10% value gives a buffer that is somewhat larger
2600 * than the effective buffer with a tolerance of 0.005 kJ/mol/ps.
2601 * Smaller tolerances or using RF lead to a smaller effective buffer,
2602 * so 10% gives a safe overestimate.
2604 return eff_1x1_buffer_fac_overest*(minimum_subgrid_size_xy(gridi) +
2605 minimum_subgrid_size_xy(gridj));
2608 /* Clusters at the cut-off only increase rlist by 60% of their size */
2609 static real nbnxn_rlist_inc_outside_fac = 0.6;
2611 /* Due to the cluster size the effective pair-list is longer than
2612 * that of a simple atom pair-list. This function gives the extra distance.
2614 real nbnxn_get_rlist_effective_inc(int cluster_size_j, real atom_density)
2617 real vol_inc_i, vol_inc_j;
2619 /* We should get this from the setup, but currently it's the same for
2620 * all setups, including GPUs.
2622 cluster_size_i = NBNXN_CPU_CLUSTER_I_SIZE;
2624 vol_inc_i = (cluster_size_i - 1)/atom_density;
2625 vol_inc_j = (cluster_size_j - 1)/atom_density;
2627 return nbnxn_rlist_inc_outside_fac*std::cbrt(vol_inc_i + vol_inc_j);
2630 /* Estimates the interaction volume^2 for non-local interactions */
2631 static real nonlocal_vol2(const struct gmx_domdec_zones_t *zones, const rvec ls, real r)
2639 /* Here we simply add up the volumes of 1, 2 or 3 1D decomposition
2640 * not home interaction volume^2. As these volumes are not additive,
2641 * this is an overestimate, but it would only be significant in the limit
2642 * of small cells, where we anyhow need to split the lists into
2643 * as small parts as possible.
2646 for (int z = 0; z < zones->n; z++)
2648 if (zones->shift[z][XX] + zones->shift[z][YY] + zones->shift[z][ZZ] == 1)
2653 for (int d = 0; d < DIM; d++)
2655 if (zones->shift[z][d] == 0)
2659 za *= zones->size[z].x1[d] - zones->size[z].x0[d];
2663 /* 4 octants of a sphere */
2664 vold_est = 0.25*M_PI*r*r*r*r;
2665 /* 4 quarter pie slices on the edges */
2666 vold_est += 4*cl*M_PI/6.0*r*r*r;
2667 /* One rectangular volume on a face */
2668 vold_est += ca*0.5*r*r;
2670 vol2_est_tot += vold_est*za;
2674 return vol2_est_tot;
2677 /* Estimates the average size of a full j-list for super/sub setup */
2678 static void get_nsubpair_target(const nbnxn_search *nbs,
2681 int min_ci_balanced,
2682 int *nsubpair_target,
2683 float *nsubpair_tot_est)
2685 /* The target value of 36 seems to be the optimum for Kepler.
2686 * Maxwell is less sensitive to the exact value.
2688 const int nsubpair_target_min = 36;
2689 const nbnxn_grid_t *grid;
2691 real r_eff_sup, vol_est, nsp_est, nsp_est_nl;
2693 grid = &nbs->grid[0];
2695 /* We don't need to balance list sizes if:
2696 * - We didn't request balancing.
2697 * - The number of grid cells >= the number of lists requested,
2698 * since we will always generate at least #cells lists.
2699 * - We don't have any cells, since then there won't be any lists.
2701 if (min_ci_balanced <= 0 || grid->nc >= min_ci_balanced || grid->nc == 0)
2703 /* nsubpair_target==0 signals no balancing */
2704 *nsubpair_target = 0;
2705 *nsubpair_tot_est = 0;
2710 ls[XX] = (grid->c1[XX] - grid->c0[XX])/(grid->numCells[XX]*c_gpuNumClusterPerCellX);
2711 ls[YY] = (grid->c1[YY] - grid->c0[YY])/(grid->numCells[YY]*c_gpuNumClusterPerCellY);
2712 ls[ZZ] = grid->na_c/(grid->atom_density*ls[XX]*ls[YY]);
2714 /* The average length of the diagonal of a sub cell */
2715 real diagonal = std::sqrt(ls[XX]*ls[XX] + ls[YY]*ls[YY] + ls[ZZ]*ls[ZZ]);
2717 /* The formulas below are a heuristic estimate of the average nsj per si*/
2718 r_eff_sup = rlist + nbnxn_rlist_inc_outside_fac*gmx::square((grid->na_c - 1.0)/grid->na_c)*0.5*diagonal;
2720 if (!nbs->DomDec || nbs->zones->n == 1)
2727 gmx::square(grid->atom_density/grid->na_c)*
2728 nonlocal_vol2(nbs->zones, ls, r_eff_sup);
2733 /* Sub-cell interacts with itself */
2734 vol_est = ls[XX]*ls[YY]*ls[ZZ];
2735 /* 6/2 rectangular volume on the faces */
2736 vol_est += (ls[XX]*ls[YY] + ls[XX]*ls[ZZ] + ls[YY]*ls[ZZ])*r_eff_sup;
2737 /* 12/2 quarter pie slices on the edges */
2738 vol_est += 2*(ls[XX] + ls[YY] + ls[ZZ])*0.25*M_PI*gmx::square(r_eff_sup);
2739 /* 4 octants of a sphere */
2740 vol_est += 0.5*4.0/3.0*M_PI*gmx::power3(r_eff_sup);
2742 /* Estimate the number of cluster pairs as the local number of
2743 * clusters times the volume they interact with times the density.
2745 nsp_est = grid->nsubc_tot*vol_est*grid->atom_density/grid->na_c;
2747 /* Subtract the non-local pair count */
2748 nsp_est -= nsp_est_nl;
2750 /* For small cut-offs nsp_est will be an underesimate.
2751 * With DD nsp_est_nl is an overestimate so nsp_est can get negative.
2752 * So to avoid too small or negative nsp_est we set a minimum of
2753 * all cells interacting with all 3^3 direct neighbors (3^3-1)/2+1=14.
2754 * This might be a slight overestimate for small non-periodic groups of
2755 * atoms as will occur for a local domain with DD, but for small
2756 * groups of atoms we'll anyhow be limited by nsubpair_target_min,
2757 * so this overestimation will not matter.
2759 nsp_est = std::max(nsp_est, grid->nsubc_tot*static_cast<real>(14));
2763 fprintf(debug, "nsp_est local %5.1f non-local %5.1f\n",
2764 nsp_est, nsp_est_nl);
2769 nsp_est = nsp_est_nl;
2772 /* Thus the (average) maximum j-list size should be as follows.
2773 * Since there is overhead, we shouldn't make the lists too small
2774 * (and we can't chop up j-groups) so we use a minimum target size of 36.
2776 *nsubpair_target = std::max(nsubpair_target_min,
2777 static_cast<int>(nsp_est/min_ci_balanced + 0.5));
2778 *nsubpair_tot_est = static_cast<int>(nsp_est);
2782 fprintf(debug, "nbl nsp estimate %.1f, nsubpair_target %d\n",
2783 nsp_est, *nsubpair_target);
2787 /* Debug list print function */
2788 static void print_nblist_ci_cj(FILE *fp, const nbnxn_pairlist_t *nbl)
2790 for (int i = 0; i < nbl->nci; i++)
2792 fprintf(fp, "ci %4d shift %2d ncj %3d\n",
2793 nbl->ci[i].ci, nbl->ci[i].shift,
2794 nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start);
2796 for (int j = nbl->ci[i].cj_ind_start; j < nbl->ci[i].cj_ind_end; j++)
2798 fprintf(fp, " cj %5d imask %x\n",
2805 /* Debug list print function */
2806 static void print_nblist_sci_cj(FILE *fp, const nbnxn_pairlist_t *nbl)
2808 for (int i = 0; i < nbl->nsci; i++)
2810 fprintf(fp, "ci %4d shift %2d ncj4 %2d\n",
2811 nbl->sci[i].sci, nbl->sci[i].shift,
2812 nbl->sci[i].cj4_ind_end - nbl->sci[i].cj4_ind_start);
2815 for (int j4 = nbl->sci[i].cj4_ind_start; j4 < nbl->sci[i].cj4_ind_end; j4++)
2817 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
2819 fprintf(fp, " sj %5d imask %x\n",
2821 nbl->cj4[j4].imei[0].imask);
2822 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2824 if (nbl->cj4[j4].imei[0].imask & (1U << (j*c_gpuNumClusterPerCell + si)))
2831 fprintf(fp, "ci %4d shift %2d ncj4 %2d ncp %3d\n",
2832 nbl->sci[i].sci, nbl->sci[i].shift,
2833 nbl->sci[i].cj4_ind_end - nbl->sci[i].cj4_ind_start,
2838 /* Combine pair lists *nbl generated on multiple threads nblc */
2839 static void combine_nblists(int nnbl, nbnxn_pairlist_t **nbl,
2840 nbnxn_pairlist_t *nblc)
2842 int nsci, ncj4, nexcl;
2846 gmx_incons("combine_nblists does not support simple lists");
2851 nexcl = nblc->nexcl;
2852 for (int i = 0; i < nnbl; i++)
2854 nsci += nbl[i]->nsci;
2855 ncj4 += nbl[i]->ncj4;
2856 nexcl += nbl[i]->nexcl;
2859 if (nsci > nblc->sci_nalloc)
2861 nb_realloc_sci(nblc, nsci);
2863 if (ncj4 > nblc->cj4_nalloc)
2865 nblc->cj4_nalloc = over_alloc_small(ncj4);
2866 nbnxn_realloc_void((void **)&nblc->cj4,
2867 nblc->ncj4*sizeof(*nblc->cj4),
2868 nblc->cj4_nalloc*sizeof(*nblc->cj4),
2869 nblc->alloc, nblc->free);
2871 if (nexcl > nblc->excl_nalloc)
2873 nblc->excl_nalloc = over_alloc_small(nexcl);
2874 nbnxn_realloc_void((void **)&nblc->excl,
2875 nblc->nexcl*sizeof(*nblc->excl),
2876 nblc->excl_nalloc*sizeof(*nblc->excl),
2877 nblc->alloc, nblc->free);
2880 /* Each thread should copy its own data to the combined arrays,
2881 * as otherwise data will go back and forth between different caches.
2883 #if GMX_OPENMP && !(defined __clang_analyzer__)
2884 // cppcheck-suppress unreadVariable
2885 int nthreads = gmx_omp_nthreads_get(emntPairsearch);
2888 #pragma omp parallel for num_threads(nthreads) schedule(static)
2889 for (int n = 0; n < nnbl; n++)
2896 const nbnxn_pairlist_t *nbli;
2898 /* Determine the offset in the combined data for our thread */
2899 sci_offset = nblc->nsci;
2900 cj4_offset = nblc->ncj4;
2901 excl_offset = nblc->nexcl;
2903 for (int i = 0; i < n; i++)
2905 sci_offset += nbl[i]->nsci;
2906 cj4_offset += nbl[i]->ncj4;
2907 excl_offset += nbl[i]->nexcl;
2912 for (int i = 0; i < nbli->nsci; i++)
2914 nblc->sci[sci_offset+i] = nbli->sci[i];
2915 nblc->sci[sci_offset+i].cj4_ind_start += cj4_offset;
2916 nblc->sci[sci_offset+i].cj4_ind_end += cj4_offset;
2919 for (int j4 = 0; j4 < nbli->ncj4; j4++)
2921 nblc->cj4[cj4_offset+j4] = nbli->cj4[j4];
2922 nblc->cj4[cj4_offset+j4].imei[0].excl_ind += excl_offset;
2923 nblc->cj4[cj4_offset+j4].imei[1].excl_ind += excl_offset;
2926 for (int j4 = 0; j4 < nbli->nexcl; j4++)
2928 nblc->excl[excl_offset+j4] = nbli->excl[j4];
2931 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2934 for (int n = 0; n < nnbl; n++)
2936 nblc->nsci += nbl[n]->nsci;
2937 nblc->ncj4 += nbl[n]->ncj4;
2938 nblc->nci_tot += nbl[n]->nci_tot;
2939 nblc->nexcl += nbl[n]->nexcl;
2943 static void balance_fep_lists(const nbnxn_search *nbs,
2944 nbnxn_pairlist_set_t *nbl_lists)
2947 int nri_tot, nrj_tot, nrj_target;
2951 nnbl = nbl_lists->nnbl;
2955 /* Nothing to balance */
2959 /* Count the total i-lists and pairs */
2962 for (int th = 0; th < nnbl; th++)
2964 nri_tot += nbl_lists->nbl_fep[th]->nri;
2965 nrj_tot += nbl_lists->nbl_fep[th]->nrj;
2968 nrj_target = (nrj_tot + nnbl - 1)/nnbl;
2970 assert(gmx_omp_nthreads_get(emntNonbonded) == nnbl);
2972 #pragma omp parallel for schedule(static) num_threads(nnbl)
2973 for (int th = 0; th < nnbl; th++)
2977 t_nblist *nbl = nbs->work[th].nbl_fep.get();
2979 /* Note that here we allocate for the total size, instead of
2980 * a per-thread esimate (which is hard to obtain).
2982 if (nri_tot > nbl->maxnri)
2984 nbl->maxnri = over_alloc_large(nri_tot);
2985 reallocate_nblist(nbl);
2987 if (nri_tot > nbl->maxnri || nrj_tot > nbl->maxnrj)
2989 nbl->maxnrj = over_alloc_small(nrj_tot);
2990 srenew(nbl->jjnr, nbl->maxnrj);
2991 srenew(nbl->excl_fep, nbl->maxnrj);
2994 clear_pairlist_fep(nbl);
2996 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2999 /* Loop over the source lists and assign and copy i-entries */
3001 nbld = nbs->work[th_dest].nbl_fep.get();
3002 for (int th = 0; th < nnbl; th++)
3006 nbls = nbl_lists->nbl_fep[th];
3008 for (int i = 0; i < nbls->nri; i++)
3012 /* The number of pairs in this i-entry */
3013 nrj = nbls->jindex[i+1] - nbls->jindex[i];
3015 /* Decide if list th_dest is too large and we should procede
3016 * to the next destination list.
3018 if (th_dest+1 < nnbl && nbld->nrj > 0 &&
3019 nbld->nrj + nrj - nrj_target > nrj_target - nbld->nrj)
3022 nbld = nbs->work[th_dest].nbl_fep.get();
3025 nbld->iinr[nbld->nri] = nbls->iinr[i];
3026 nbld->gid[nbld->nri] = nbls->gid[i];
3027 nbld->shift[nbld->nri] = nbls->shift[i];
3029 for (int j = nbls->jindex[i]; j < nbls->jindex[i+1]; j++)
3031 nbld->jjnr[nbld->nrj] = nbls->jjnr[j];
3032 nbld->excl_fep[nbld->nrj] = nbls->excl_fep[j];
3036 nbld->jindex[nbld->nri] = nbld->nrj;
3040 /* Swap the list pointers */
3041 for (int th = 0; th < nnbl; th++)
3043 t_nblist *nbl_tmp = nbs->work[th].nbl_fep.release();
3044 nbs->work[th].nbl_fep.reset(nbl_lists->nbl_fep[th]);
3045 nbl_lists->nbl_fep[th] = nbl_tmp;
3049 fprintf(debug, "nbl_fep[%d] nri %4d nrj %4d\n",
3051 nbl_lists->nbl_fep[th]->nri,
3052 nbl_lists->nbl_fep[th]->nrj);
3057 /* Returns the next ci to be processes by our thread */
3058 static gmx_bool next_ci(const nbnxn_grid_t *grid,
3059 int nth, int ci_block,
3060 int *ci_x, int *ci_y,
3066 if (*ci_b == ci_block)
3068 /* Jump to the next block assigned to this task */
3069 *ci += (nth - 1)*ci_block;
3073 if (*ci >= grid->nc)
3078 while (*ci >= grid->cxy_ind[*ci_x*grid->numCells[YY] + *ci_y + 1])
3081 if (*ci_y == grid->numCells[YY])
3091 /* Returns the distance^2 for which we put cell pairs in the list
3092 * without checking atom pair distances. This is usually < rlist^2.
3094 static float boundingbox_only_distance2(const nbnxn_grid_t *gridi,
3095 const nbnxn_grid_t *gridj,
3099 /* If the distance between two sub-cell bounding boxes is less
3100 * than this distance, do not check the distance between
3101 * all particle pairs in the sub-cell, since then it is likely
3102 * that the box pair has atom pairs within the cut-off.
3103 * We use the nblist cut-off minus 0.5 times the average x/y diagonal
3104 * spacing of the sub-cells. Around 40% of the checked pairs are pruned.
3105 * Using more than 0.5 gains at most 0.5%.
3106 * If forces are calculated more than twice, the performance gain
3107 * in the force calculation outweighs the cost of checking.
3108 * Note that with subcell lists, the atom-pair distance check
3109 * is only performed when only 1 out of 8 sub-cells in within range,
3110 * this is because the GPU is much faster than the cpu.
3115 bbx = 0.5*(gridi->cellSize[XX] + gridj->cellSize[XX]);
3116 bby = 0.5*(gridi->cellSize[YY] + gridj->cellSize[YY]);
3119 bbx /= c_gpuNumClusterPerCellX;
3120 bby /= c_gpuNumClusterPerCellY;
3123 rbb2 = std::max(0.0, rlist - 0.5*std::sqrt(bbx*bbx + bby*bby));
3129 return (float)((1+GMX_FLOAT_EPS)*rbb2);
3133 static int get_ci_block_size(const nbnxn_grid_t *gridi,
3134 gmx_bool bDomDec, int nth)
3136 const int ci_block_enum = 5;
3137 const int ci_block_denom = 11;
3138 const int ci_block_min_atoms = 16;
3141 /* Here we decide how to distribute the blocks over the threads.
3142 * We use prime numbers to try to avoid that the grid size becomes
3143 * a multiple of the number of threads, which would lead to some
3144 * threads getting "inner" pairs and others getting boundary pairs,
3145 * which in turns will lead to load imbalance between threads.
3146 * Set the block size as 5/11/ntask times the average number of cells
3147 * in a y,z slab. This should ensure a quite uniform distribution
3148 * of the grid parts of the different thread along all three grid
3149 * zone boundaries with 3D domain decomposition. At the same time
3150 * the blocks will not become too small.
3152 ci_block = (gridi->nc*ci_block_enum)/(ci_block_denom*gridi->numCells[XX]*nth);
3154 /* Ensure the blocks are not too small: avoids cache invalidation */
3155 if (ci_block*gridi->na_sc < ci_block_min_atoms)
3157 ci_block = (ci_block_min_atoms + gridi->na_sc - 1)/gridi->na_sc;
3160 /* Without domain decomposition
3161 * or with less than 3 blocks per task, divide in nth blocks.
3163 if (!bDomDec || nth*3*ci_block > gridi->nc)
3165 ci_block = (gridi->nc + nth - 1)/nth;
3168 if (ci_block > 1 && (nth - 1)*ci_block >= gridi->nc)
3170 /* Some threads have no work. Although reducing the block size
3171 * does not decrease the block count on the first few threads,
3172 * with GPUs better mixing of "upper" cells that have more empty
3173 * clusters results in a somewhat lower max load over all threads.
3174 * Without GPUs the regime of so few atoms per thread is less
3175 * performance relevant, but with 8-wide SIMD the same reasoning
3176 * applies, since the pair list uses 4 i-atom "sub-clusters".
3184 /* Returns the number of bits to right-shift a cluster index to obtain
3185 * the corresponding force buffer flag index.
3187 static int getBufferFlagShift(int numAtomsPerCluster)
3189 int bufferFlagShift = 0;
3190 while ((numAtomsPerCluster << bufferFlagShift) < NBNXN_BUFFERFLAG_SIZE)
3195 return bufferFlagShift;
3198 /* Generates the part of pair-list nbl assigned to our thread */
3199 static void nbnxn_make_pairlist_part(const nbnxn_search *nbs,
3200 const nbnxn_grid_t *gridi,
3201 const nbnxn_grid_t *gridj,
3202 nbnxn_search_work_t *work,
3203 const nbnxn_atomdata_t *nbat,
3204 const t_blocka &exclusions,
3208 gmx_bool bFBufferFlag,
3211 float nsubpair_tot_est,
3213 nbnxn_pairlist_t *nbl,
3218 real rlist2, rl_fep2 = 0;
3220 int ci_b, ci, ci_x, ci_y, ci_xy, cj;
3224 real bx0, bx1, by0, by1, bz0, bz1;
3226 real d2cx, d2z, d2z_cx, d2z_cy, d2zx, d2zxy, d2xy;
3227 int cxf, cxl, cyf, cyf_x, cyl;
3228 int numDistanceChecks;
3229 int gridi_flag_shift = 0, gridj_flag_shift = 0;
3230 gmx_bitmask_t *gridj_flag = nullptr;
3231 int ncj_old_i, ncj_old_j;
3233 nbs_cycle_start(&work->cc[enbsCCsearch]);
3235 if (gridj->bSimple != nbl->bSimple || gridi->bSimple != nbl->bSimple)
3237 gmx_incons("Grid incompatible with pair-list");
3241 nbl->na_sc = gridj->na_sc;
3242 nbl->na_ci = gridj->na_c;
3243 nbl->na_cj = nbnxn_kernel_to_cluster_j_size(nb_kernel_type);
3244 na_cj_2log = get_2log(nbl->na_cj);
3250 /* Determine conversion of clusters to flag blocks */
3251 gridi_flag_shift = getBufferFlagShift(nbl->na_ci);
3252 gridj_flag_shift = getBufferFlagShift(nbl->na_cj);
3254 gridj_flag = work->buffer_flags.flag;
3257 copy_mat(nbs->box, box);
3259 rlist2 = nbl->rlist*nbl->rlist;
3261 if (nbs->bFEP && !nbl->bSimple)
3263 /* Determine an atom-pair list cut-off distance for FEP atom pairs.
3264 * We should not simply use rlist, since then we would not have
3265 * the small, effective buffering of the NxN lists.
3266 * The buffer is on overestimate, but the resulting cost for pairs
3267 * beyond rlist is neglible compared to the FEP pairs within rlist.
3269 rl_fep2 = nbl->rlist + effective_buffer_1x1_vs_MxN(gridi, gridj);
3273 fprintf(debug, "nbl_fep atom-pair rlist %f\n", rl_fep2);
3275 rl_fep2 = rl_fep2*rl_fep2;
3278 rbb2 = boundingbox_only_distance2(gridi, gridj, nbl->rlist, nbl->bSimple);
3282 fprintf(debug, "nbl bounding box only distance %f\n", std::sqrt(rbb2));
3285 /* Set the shift range */
3286 for (int d = 0; d < DIM; d++)
3288 /* Check if we need periodicity shifts.
3289 * Without PBC or with domain decomposition we don't need them.
3291 if (d >= ePBC2npbcdim(nbs->ePBC) || nbs->dd_dim[d])
3298 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < std::sqrt(rlist2))
3308 const bool bSimple = nbl->bSimple;
3309 gmx::ArrayRef<const nbnxn_bb_t> bb_i;
3311 gmx::ArrayRef<const float> pbb_i;
3321 /* We use the normal bounding box format for both grid types */
3324 gmx::ArrayRef<const float> bbcz_i = gridi->bbcz;
3325 gmx::ArrayRef<const int> flags_i = gridi->flags;
3326 gmx::ArrayRef<const float> bbcz_j = gridj->bbcz;
3327 int cell0_i = gridi->cell0;
3331 fprintf(debug, "nbl nc_i %d col.av. %.1f ci_block %d\n",
3332 gridi->nc, gridi->nc/(double)(gridi->numCells[XX]*gridi->numCells[YY]), ci_block);
3335 numDistanceChecks = 0;
3337 /* Initially ci_b and ci to 1 before where we want them to start,
3338 * as they will both be incremented in next_ci.
3341 ci = th*ci_block - 1;
3344 while (next_ci(gridi, nth, ci_block, &ci_x, &ci_y, &ci_b, &ci))
3346 if (bSimple && flags_i[ci] == 0)
3351 ncj_old_i = nbl->ncj;
3354 if (gridj != gridi && shp[XX] == 0)
3358 bx1 = bb_i[ci].upper[BB_X];
3362 bx1 = gridi->c0[XX] + (ci_x+1)*gridi->cellSize[XX];
3364 if (bx1 < gridj->c0[XX])
3366 d2cx = gmx::square(gridj->c0[XX] - bx1);
3375 ci_xy = ci_x*gridi->numCells[YY] + ci_y;
3377 /* Loop over shift vectors in three dimensions */
3378 for (int tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
3380 shz = tz*box[ZZ][ZZ];
3382 bz0 = bbcz_i[ci*NNBSBB_D ] + shz;
3383 bz1 = bbcz_i[ci*NNBSBB_D+1] + shz;
3391 d2z = gmx::square(bz1);
3395 d2z = gmx::square(bz0 - box[ZZ][ZZ]);
3398 d2z_cx = d2z + d2cx;
3400 if (d2z_cx >= rlist2)
3405 bz1_frac = bz1/(gridi->cxy_ind[ci_xy+1] - gridi->cxy_ind[ci_xy]);
3410 /* The check with bz1_frac close to or larger than 1 comes later */
3412 for (int ty = -shp[YY]; ty <= shp[YY]; ty++)
3414 shy = ty*box[YY][YY] + tz*box[ZZ][YY];
3418 by0 = bb_i[ci].lower[BB_Y] + shy;
3419 by1 = bb_i[ci].upper[BB_Y] + shy;
3423 by0 = gridi->c0[YY] + (ci_y )*gridi->cellSize[YY] + shy;
3424 by1 = gridi->c0[YY] + (ci_y+1)*gridi->cellSize[YY] + shy;
3427 get_cell_range<YY>(by0, by1,
3438 if (by1 < gridj->c0[YY])
3440 d2z_cy += gmx::square(gridj->c0[YY] - by1);
3442 else if (by0 > gridj->c1[YY])
3444 d2z_cy += gmx::square(by0 - gridj->c1[YY]);
3447 for (int tx = -shp[XX]; tx <= shp[XX]; tx++)
3449 shift = XYZ2IS(tx, ty, tz);
3451 if (c_pbcShiftBackward && gridi == gridj && shift > CENTRAL)
3456 shx = tx*box[XX][XX] + ty*box[YY][XX] + tz*box[ZZ][XX];
3460 bx0 = bb_i[ci].lower[BB_X] + shx;
3461 bx1 = bb_i[ci].upper[BB_X] + shx;
3465 bx0 = gridi->c0[XX] + (ci_x )*gridi->cellSize[XX] + shx;
3466 bx1 = gridi->c0[XX] + (ci_x+1)*gridi->cellSize[XX] + shx;
3469 get_cell_range<XX>(bx0, bx1,
3481 new_ci_entry(nbl, cell0_i+ci, shift, flags_i[ci]);
3485 new_sci_entry(nbl, cell0_i+ci, shift);
3488 if ((!c_pbcShiftBackward || (shift == CENTRAL &&
3492 /* Leave the pairs with i > j.
3493 * x is the major index, so skip half of it.
3500 set_icell_bb_simple(bb_i, ci, shx, shy, shz,
3506 set_icell_bbxxxx_supersub(pbb_i, ci, shx, shy, shz,
3509 set_icell_bb_supersub(bb_i, ci, shx, shy, shz,
3514 nbs->icell_set_x(cell0_i+ci, shx, shy, shz,
3515 nbat->xstride, nbat->x,
3518 for (int cx = cxf; cx <= cxl; cx++)
3521 if (gridj->c0[XX] + cx*gridj->cellSize[XX] > bx1)
3523 d2zx += gmx::square(gridj->c0[XX] + cx*gridj->cellSize[XX] - bx1);
3525 else if (gridj->c0[XX] + (cx+1)*gridj->cellSize[XX] < bx0)
3527 d2zx += gmx::square(gridj->c0[XX] + (cx+1)*gridj->cellSize[XX] - bx0);
3530 if (gridi == gridj &&
3532 (!c_pbcShiftBackward || shift == CENTRAL) &&
3535 /* Leave the pairs with i > j.
3536 * Skip half of y when i and j have the same x.
3545 for (int cy = cyf_x; cy <= cyl; cy++)
3547 const int columnStart = gridj->cxy_ind[cx*gridj->numCells[YY] + cy];
3548 const int columnEnd = gridj->cxy_ind[cx*gridj->numCells[YY] + cy + 1];
3551 if (gridj->c0[YY] + cy*gridj->cellSize[YY] > by1)
3553 d2zxy += gmx::square(gridj->c0[YY] + cy*gridj->cellSize[YY] - by1);
3555 else if (gridj->c0[YY] + (cy+1)*gridj->cellSize[YY] < by0)
3557 d2zxy += gmx::square(gridj->c0[YY] + (cy+1)*gridj->cellSize[YY] - by0);
3559 if (columnStart < columnEnd && d2zxy < rlist2)
3561 /* To improve efficiency in the common case
3562 * of a homogeneous particle distribution,
3563 * we estimate the index of the middle cell
3564 * in range (midCell). We search down and up
3565 * starting from this index.
3567 * Note that the bbcz_j array contains bounds
3568 * for i-clusters, thus for clusters of 4 atoms.
3569 * For the common case where the j-cluster size
3570 * is 8, we could step with a stride of 2,
3571 * but we do not do this because it would
3572 * complicate this code even more.
3574 int midCell = columnStart + static_cast<int>(bz1_frac*(columnEnd - columnStart));
3575 if (midCell >= columnEnd)
3577 midCell = columnEnd - 1;
3582 /* Find the lowest cell that can possibly
3584 * Check if we hit the bottom of the grid,
3585 * if the j-cell is below the i-cell and if so,
3586 * if it is within range.
3588 int downTestCell = midCell;
3589 while (downTestCell >= columnStart &&
3590 (bbcz_j[downTestCell*NNBSBB_D + 1] >= bz0 ||
3591 d2xy + gmx::square(bbcz_j[downTestCell*NNBSBB_D + 1] - bz0) < rlist2))
3595 int firstCell = downTestCell + 1;
3597 /* Find the highest cell that can possibly
3599 * Check if we hit the top of the grid,
3600 * if the j-cell is above the i-cell and if so,
3601 * if it is within range.
3603 int upTestCell = midCell + 1;
3604 while (upTestCell < columnEnd &&
3605 (bbcz_j[upTestCell*NNBSBB_D] <= bz1 ||
3606 d2xy + gmx::square(bbcz_j[upTestCell*NNBSBB_D] - bz1) < rlist2))
3610 int lastCell = upTestCell - 1;
3612 #define NBNXN_REFCODE 0
3615 /* Simple reference code, for debugging,
3616 * overrides the more complex code above.
3618 firstCell = columnEnd;
3620 for (int k = columnStart; k < columnEnd; k++)
3622 if (d2xy + gmx::square(bbcz_j[k*NNBSBB_D + 1] - bz0) < rlist2 &&
3627 if (d2xy + gmx::square(bbcz_j[k*NNBSBB_D] - bz1) < rlist2 &&
3638 /* We want each atom/cell pair only once,
3639 * only use cj >= ci.
3641 if (!c_pbcShiftBackward || shift == CENTRAL)
3643 firstCell = std::max(firstCell, ci);
3647 if (firstCell <= lastCell)
3649 GMX_ASSERT(firstCell >= columnStart && lastCell < columnEnd, "The range should reside within the current grid column");
3651 /* For f buffer flags with simple lists */
3652 ncj_old_j = nbl->ncj;
3656 /* We have a maximum of 2 j-clusters
3657 * per i-cluster sized cell.
3659 check_cell_list_space_simple(nbl, 2*(lastCell - firstCell + 1));
3663 check_cell_list_space_supersub(nbl, lastCell - firstCell + 1);
3666 switch (nb_kernel_type)
3668 case nbnxnk4x4_PlainC:
3669 makeClusterListSimple(gridj,
3670 nbl, ci, firstCell, lastCell,
3671 (gridi == gridj && shift == CENTRAL),
3674 &numDistanceChecks);
3676 #ifdef GMX_NBNXN_SIMD_4XN
3677 case nbnxnk4xN_SIMD_4xN:
3678 makeClusterListSimd4xn(gridj,
3679 nbl, ci, firstCell, lastCell,
3680 (gridi == gridj && shift == CENTRAL),
3683 &numDistanceChecks);
3686 #ifdef GMX_NBNXN_SIMD_2XNN
3687 case nbnxnk4xN_SIMD_2xNN:
3688 makeClusterListSimd2xnn(gridj,
3689 nbl, ci, firstCell, lastCell,
3690 (gridi == gridj && shift == CENTRAL),
3693 &numDistanceChecks);
3696 case nbnxnk8x8x8_PlainC:
3697 case nbnxnk8x8x8_GPU:
3698 for (cj = firstCell; cj <= lastCell; cj++)
3700 make_cluster_list_supersub(gridi, gridj,
3702 (gridi == gridj && shift == CENTRAL && ci == cj),
3703 nbat->xstride, nbat->x,
3705 &numDistanceChecks);
3710 if (bFBufferFlag && nbl->ncj > ncj_old_j)
3712 int cbf = nbl->cj[ncj_old_j].cj >> gridj_flag_shift;
3713 int cbl = nbl->cj[nbl->ncj-1].cj >> gridj_flag_shift;
3714 for (int cb = cbf; cb <= cbl; cb++)
3716 bitmask_init_bit(&gridj_flag[cb], th);
3720 nbl->ncjInUse += nbl->ncj - ncj_old_j;
3726 /* Set the exclusions for this ci list */
3729 setExclusionsForSimpleIentry(nbs,
3731 shift == CENTRAL && gridi == gridj,
3738 make_fep_list(nbs, nbat, nbl,
3739 shift == CENTRAL && gridi == gridj,
3740 &(nbl->ci[nbl->nci]),
3741 gridi, gridj, nbl_fep);
3746 setExclusionsForGpuIentry(nbs,
3748 shift == CENTRAL && gridi == gridj,
3749 nbl->sci[nbl->nsci],
3754 make_fep_list_supersub(nbs, nbat, nbl,
3755 shift == CENTRAL && gridi == gridj,
3756 &(nbl->sci[nbl->nsci]),
3759 gridi, gridj, nbl_fep);
3763 /* Close this ci list */
3766 close_ci_entry_simple(nbl);
3770 close_ci_entry_supersub(nbl,
3772 progBal, nsubpair_tot_est,
3779 if (bFBufferFlag && nbl->ncj > ncj_old_i)
3781 bitmask_init_bit(&(work->buffer_flags.flag[(gridi->cell0+ci)>>gridi_flag_shift]), th);
3785 work->ndistc = numDistanceChecks;
3787 nbs_cycle_stop(&work->cc[enbsCCsearch]);
3789 GMX_ASSERT(nbl->ncjInUse == nbl->ncj || nbs->bFEP, "Without free-energy all cj pair-list entries should be in use. Note that subsequent code does not make use of the equality, this check is only here to catch bugs");
3793 fprintf(debug, "number of distance checks %d\n", numDistanceChecks);
3797 print_nblist_statistics_simple(debug, nbl, nbs, rlist);
3801 print_nblist_statistics_supersub(debug, nbl, nbs, rlist);
3806 fprintf(debug, "nbl FEP list pairs: %d\n", nbl_fep->nrj);
3811 static void reduce_buffer_flags(const nbnxn_search *nbs,
3813 const nbnxn_buffer_flags_t *dest)
3815 for (int s = 0; s < nsrc; s++)
3817 gmx_bitmask_t * flag = nbs->work[s].buffer_flags.flag;
3819 for (int b = 0; b < dest->nflag; b++)
3821 bitmask_union(&(dest->flag[b]), flag[b]);
3826 static void print_reduction_cost(const nbnxn_buffer_flags_t *flags, int nout)
3828 int nelem, nkeep, ncopy, nred, out;
3829 gmx_bitmask_t mask_0;
3835 bitmask_init_bit(&mask_0, 0);
3836 for (int b = 0; b < flags->nflag; b++)
3838 if (bitmask_is_equal(flags->flag[b], mask_0))
3840 /* Only flag 0 is set, no copy of reduction required */
3844 else if (!bitmask_is_zero(flags->flag[b]))
3847 for (out = 0; out < nout; out++)
3849 if (bitmask_is_set(flags->flag[b], out))
3866 fprintf(debug, "nbnxn reduction: #flag %d #list %d elem %4.2f, keep %4.2f copy %4.2f red %4.2f\n",
3868 nelem/(double)(flags->nflag),
3869 nkeep/(double)(flags->nflag),
3870 ncopy/(double)(flags->nflag),
3871 nred/(double)(flags->nflag));
3874 /* Copies the list entries from src to dest when cjStart <= *cjGlobal < cjEnd.
3875 * *cjGlobal is updated with the cj count in src.
3876 * When setFlags==true, flag bit t is set in flag for all i and j clusters.
3878 template<bool setFlags>
3879 static void copySelectedListRange(const nbnxn_ci_t * gmx_restrict srcCi,
3880 const nbnxn_pairlist_t * gmx_restrict src,
3881 nbnxn_pairlist_t * gmx_restrict dest,
3882 gmx_bitmask_t *flag,
3883 int iFlagShift, int jFlagShift, int t)
3885 int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3887 if (dest->nci + 1 >= dest->ci_nalloc)
3889 nb_realloc_ci(dest, dest->nci + 1);
3891 check_cell_list_space_simple(dest, ncj);
3893 dest->ci[dest->nci] = *srcCi;
3894 dest->ci[dest->nci].cj_ind_start = dest->ncj;
3895 dest->ci[dest->nci].cj_ind_end = dest->ncj + ncj;
3899 bitmask_init_bit(&flag[srcCi->ci >> iFlagShift], t);
3902 for (int j = srcCi->cj_ind_start; j < srcCi->cj_ind_end; j++)
3904 dest->cj[dest->ncj++] = src->cj[j];
3908 /* NOTE: This is relatively expensive, since this
3909 * operation is done for all elements in the list,
3910 * whereas at list generation this is done only
3911 * once for each flag entry.
3913 bitmask_init_bit(&flag[src->cj[j].cj >> jFlagShift], t);
3920 /* This routine re-balances the pairlists such that all are nearly equally
3921 * sized. Only whole i-entries are moved between lists. These are moved
3922 * between the ends of the lists, such that the buffer reduction cost should
3923 * not change significantly.
3924 * Note that all original reduction flags are currently kept. This can lead
3925 * to reduction of parts of the force buffer that could be avoided. But since
3926 * the original lists are quite balanced, this will only give minor overhead.
3928 static void rebalanceSimpleLists(int numLists,
3929 nbnxn_pairlist_t * const * const srcSet,
3930 nbnxn_pairlist_t **destSet,
3931 gmx::ArrayRef<nbnxn_search_work_t> searchWork)
3934 for (int s = 0; s < numLists; s++)
3936 ncjTotal += srcSet[s]->ncjInUse;
3938 int ncjTarget = (ncjTotal + numLists - 1)/numLists;
3940 #pragma omp parallel num_threads(numLists)
3942 int t = gmx_omp_get_thread_num();
3944 int cjStart = ncjTarget* t;
3945 int cjEnd = ncjTarget*(t + 1);
3947 /* The destination pair-list for task/thread t */
3948 nbnxn_pairlist_t *dest = destSet[t];
3950 clear_pairlist(dest);
3951 dest->bSimple = srcSet[0]->bSimple;
3952 dest->na_ci = srcSet[0]->na_ci;
3953 dest->na_cj = srcSet[0]->na_cj;
3955 /* Note that the flags in the work struct (still) contain flags
3956 * for all entries that are present in srcSet->nbl[t].
3958 gmx_bitmask_t *flag = searchWork[t].buffer_flags.flag;
3960 int iFlagShift = getBufferFlagShift(dest->na_ci);
3961 int jFlagShift = getBufferFlagShift(dest->na_cj);
3964 for (int s = 0; s < numLists && cjGlobal < cjEnd; s++)
3966 const nbnxn_pairlist_t *src = srcSet[s];
3968 if (cjGlobal + src->ncjInUse > cjStart)
3970 for (int i = 0; i < src->nci && cjGlobal < cjEnd; i++)
3972 const nbnxn_ci_t *srcCi = &src->ci[i];
3973 int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3974 if (cjGlobal >= cjStart)
3976 /* If the source list is not our own, we need to set
3977 * extra flags (the template bool parameter).
3981 copySelectedListRange
3984 flag, iFlagShift, jFlagShift, t);
3988 copySelectedListRange
3991 dest, flag, iFlagShift, jFlagShift, t);
3999 cjGlobal += src->ncjInUse;
4003 dest->ncjInUse = dest->ncj;
4007 int ncjTotalNew = 0;
4008 for (int s = 0; s < numLists; s++)
4010 ncjTotalNew += destSet[s]->ncjInUse;
4012 GMX_RELEASE_ASSERT(ncjTotalNew == ncjTotal, "The total size of the lists before and after rebalancing should match");
4016 /* Returns if the pairlists are so imbalanced that it is worth rebalancing. */
4017 static bool checkRebalanceSimpleLists(const nbnxn_pairlist_set_t *listSet)
4019 int numLists = listSet->nnbl;
4022 for (int s = 0; s < numLists; s++)
4024 ncjMax = std::max(ncjMax, listSet->nbl[s]->ncjInUse);
4025 ncjTotal += listSet->nbl[s]->ncjInUse;
4029 fprintf(debug, "Pair-list ncjMax %d ncjTotal %d\n", ncjMax, ncjTotal);
4031 /* The rebalancing adds 3% extra time to the search. Heuristically we
4032 * determined that under common conditions the non-bonded kernel balance
4033 * improvement will outweigh this when the imbalance is more than 3%.
4034 * But this will, obviously, depend on search vs kernel time and nstlist.
4036 const real rebalanceTolerance = 1.03;
4038 return numLists*ncjMax > ncjTotal*rebalanceTolerance;
4041 /* Perform a count (linear) sort to sort the smaller lists to the end.
4042 * This avoids load imbalance on the GPU, as large lists will be
4043 * scheduled and executed first and the smaller lists later.
4044 * Load balancing between multi-processors only happens at the end
4045 * and there smaller lists lead to more effective load balancing.
4046 * The sorting is done on the cj4 count, not on the actual pair counts.
4047 * Not only does this make the sort faster, but it also results in
4048 * better load balancing than using a list sorted on exact load.
4049 * This function swaps the pointer in the pair list to avoid a copy operation.
4051 static void sort_sci(nbnxn_pairlist_t *nbl)
4053 nbnxn_list_work_t *work;
4055 nbnxn_sci_t *sci_sort;
4057 if (nbl->ncj4 <= nbl->nsci)
4059 /* nsci = 0 or all sci have size 1, sorting won't change the order */
4065 /* We will distinguish differences up to double the average */
4066 m = (2*nbl->ncj4)/nbl->nsci;
4068 if (m + 1 > work->sort_nalloc)
4070 work->sort_nalloc = over_alloc_large(m + 1);
4071 srenew(work->sort, work->sort_nalloc);
4074 if (work->sci_sort_nalloc != nbl->sci_nalloc)
4076 work->sci_sort_nalloc = nbl->sci_nalloc;
4077 nbnxn_realloc_void((void **)&work->sci_sort,
4079 work->sci_sort_nalloc*sizeof(*work->sci_sort),
4080 nbl->alloc, nbl->free);
4083 /* Count the entries of each size */
4084 for (int i = 0; i <= m; i++)
4088 for (int s = 0; s < nbl->nsci; s++)
4090 int i = std::min(m, nbl->sci[s].cj4_ind_end - nbl->sci[s].cj4_ind_start);
4093 /* Calculate the offset for each count */
4096 for (int i = m - 1; i >= 0; i--)
4099 work->sort[i] = work->sort[i + 1] + s0;
4103 /* Sort entries directly into place */
4104 sci_sort = work->sci_sort;
4105 for (int s = 0; s < nbl->nsci; s++)
4107 int i = std::min(m, nbl->sci[s].cj4_ind_end - nbl->sci[s].cj4_ind_start);
4108 sci_sort[work->sort[i]++] = nbl->sci[s];
4111 /* Swap the sci pointers so we use the new, sorted list */
4112 work->sci_sort = nbl->sci;
4113 nbl->sci = sci_sort;
4116 /* Make a local or non-local pair-list, depending on iloc */
4117 void nbnxn_make_pairlist(nbnxn_search *nbs,
4118 nbnxn_atomdata_t *nbat,
4119 const t_blocka *excl,
4121 int min_ci_balanced,
4122 nbnxn_pairlist_set_t *nbl_list,
4127 nbnxn_grid_t *gridi, *gridj;
4129 int nsubpair_target;
4130 float nsubpair_tot_est;
4132 nbnxn_pairlist_t **nbl;
4134 gmx_bool CombineNBLists;
4136 int np_tot, np_noq, np_hlj, nap;
4138 nnbl = nbl_list->nnbl;
4139 nbl = nbl_list->nbl;
4140 CombineNBLists = nbl_list->bCombined;
4144 fprintf(debug, "ns making %d nblists\n", nnbl);
4147 nbat->bUseBufferFlags = (nbat->nout > 1);
4148 /* We should re-init the flags before making the first list */
4149 if (nbat->bUseBufferFlags && LOCAL_I(iloc))
4151 init_buffer_flags(&nbat->buffer_flags, nbat->natoms);
4154 if (nbl_list->bSimple)
4157 switch (nb_kernel_type)
4159 #ifdef GMX_NBNXN_SIMD_4XN
4160 case nbnxnk4xN_SIMD_4xN:
4161 nbs->icell_set_x = icell_set_x_simd_4xn;
4164 #ifdef GMX_NBNXN_SIMD_2XNN
4165 case nbnxnk4xN_SIMD_2xNN:
4166 nbs->icell_set_x = icell_set_x_simd_2xnn;
4170 nbs->icell_set_x = icell_set_x_simple;
4174 /* MSVC 2013 complains about switch statements without case */
4175 nbs->icell_set_x = icell_set_x_simple;
4180 nbs->icell_set_x = icell_set_x_supersub;
4185 /* Only zone (grid) 0 vs 0 */
4192 nzi = nbs->zones->nizone;
4195 if (!nbl_list->bSimple && min_ci_balanced > 0)
4197 get_nsubpair_target(nbs, iloc, rlist, min_ci_balanced,
4198 &nsubpair_target, &nsubpair_tot_est);
4202 nsubpair_target = 0;
4203 nsubpair_tot_est = 0;
4206 /* Clear all pair-lists */
4207 for (int th = 0; th < nnbl; th++)
4209 clear_pairlist(nbl[th]);
4213 clear_pairlist_fep(nbl_list->nbl_fep[th]);
4217 for (int zi = 0; zi < nzi; zi++)
4219 gridi = &nbs->grid[zi];
4221 if (NONLOCAL_I(iloc))
4223 zj0 = nbs->zones->izone[zi].j0;
4224 zj1 = nbs->zones->izone[zi].j1;
4230 for (int zj = zj0; zj < zj1; zj++)
4232 gridj = &nbs->grid[zj];
4236 fprintf(debug, "ns search grid %d vs %d\n", zi, zj);
4239 nbs_cycle_start(&nbs->cc[enbsCCsearch]);
4241 ci_block = get_ci_block_size(gridi, nbs->DomDec, nnbl);
4243 /* With GPU: generate progressively smaller lists for
4244 * load balancing for local only or non-local with 2 zones.
4246 progBal = (LOCAL_I(iloc) || nbs->zones->n <= 2);
4248 #pragma omp parallel for num_threads(nnbl) schedule(static)
4249 for (int th = 0; th < nnbl; th++)
4253 /* Re-init the thread-local work flag data before making
4254 * the first list (not an elegant conditional).
4256 if (nbat->bUseBufferFlags && ((zi == 0 && zj == 0)))
4258 init_buffer_flags(&nbs->work[th].buffer_flags, nbat->natoms);
4261 if (CombineNBLists && th > 0)
4263 clear_pairlist(nbl[th]);
4266 /* Divide the i super cell equally over the nblists */
4267 nbnxn_make_pairlist_part(nbs, gridi, gridj,
4268 &nbs->work[th], nbat, *excl,
4272 nbat->bUseBufferFlags,
4274 progBal, nsubpair_tot_est,
4277 nbl_list->nbl_fep[th]);
4279 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4281 nbs_cycle_stop(&nbs->cc[enbsCCsearch]);
4286 for (int th = 0; th < nnbl; th++)
4288 inc_nrnb(nrnb, eNR_NBNXN_DIST2, nbs->work[th].ndistc);
4290 if (nbl_list->bSimple)
4292 np_tot += nbl[th]->ncj;
4293 np_noq += nbl[th]->work->ncj_noq;
4294 np_hlj += nbl[th]->work->ncj_hlj;
4298 /* This count ignores potential subsequent pair pruning */
4299 np_tot += nbl[th]->nci_tot;
4302 nap = nbl[0]->na_ci*nbl[0]->na_cj;
4303 nbl_list->natpair_ljq = (np_tot - np_noq)*nap - np_hlj*nap/2;
4304 nbl_list->natpair_lj = np_noq*nap;
4305 nbl_list->natpair_q = np_hlj*nap/2;
4307 if (CombineNBLists && nnbl > 1)
4309 nbs_cycle_start(&nbs->cc[enbsCCcombine]);
4311 combine_nblists(nnbl-1, nbl+1, nbl[0]);
4313 nbs_cycle_stop(&nbs->cc[enbsCCcombine]);
4318 if (nbl_list->bSimple)
4320 if (nnbl > 1 && checkRebalanceSimpleLists(nbl_list))
4322 rebalanceSimpleLists(nbl_list->nnbl, nbl_list->nbl, nbl_list->nbl_work, nbs->work);
4324 /* Swap the pointer of the sets of pair lists */
4325 nbnxn_pairlist_t **tmp = nbl_list->nbl;
4326 nbl_list->nbl = nbl_list->nbl_work;
4327 nbl_list->nbl_work = tmp;
4332 /* Sort the entries on size, large ones first */
4333 if (CombineNBLists || nnbl == 1)
4339 #pragma omp parallel for num_threads(nnbl) schedule(static)
4340 for (int th = 0; th < nnbl; th++)
4346 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4351 if (nbat->bUseBufferFlags)
4353 reduce_buffer_flags(nbs, nbl_list->nnbl, &nbat->buffer_flags);
4358 /* Balance the free-energy lists over all the threads */
4359 balance_fep_lists(nbs, nbl_list);
4362 /* This is a fresh list, so not pruned, stored using ci and nci.
4363 * ciOuter and nciOuter are invalid at this point.
4365 GMX_ASSERT(nbl_list->nbl[0]->nciOuter == -1, "nciOuter should have been set to -1 to signal that it is invalid");
4367 /* Special performance logging stuff (env.var. GMX_NBNXN_CYCLE) */
4370 nbs->search_count++;
4372 if (nbs->print_cycles &&
4373 (!nbs->DomDec || !LOCAL_I(iloc)) &&
4374 nbs->search_count % 100 == 0)
4376 nbs_cycle_print(stderr, nbs);
4379 /* If we have more than one list, they either got rebalancing (CPU)
4380 * or combined (GPU), so we should dump the final result to debug.
4382 if (debug && nbl_list->nnbl > 1)
4384 if (nbl_list->bSimple)
4386 for (int t = 0; t < nbl_list->nnbl; t++)
4388 print_nblist_statistics_simple(debug, nbl_list->nbl[t], nbs, rlist);
4393 print_nblist_statistics_supersub(debug, nbl_list->nbl[0], nbs, rlist);
4401 if (nbl_list->bSimple)
4403 for (int t = 0; t < nbl_list->nnbl; t++)
4405 print_nblist_ci_cj(debug, nbl_list->nbl[t]);
4410 print_nblist_sci_cj(debug, nbl_list->nbl[0]);
4414 if (nbat->bUseBufferFlags)
4416 print_reduction_cost(&nbat->buffer_flags, nbl_list->nnbl);
4421 void nbnxnPrepareListForDynamicPruning(nbnxn_pairlist_set_t *listSet)
4423 /* TODO: Restructure the lists so we have actual outer and inner
4424 * list objects so we can set a single pointer instead of
4425 * swapping several pointers.
4428 for (int i = 0; i < listSet->nnbl; i++)
4430 /* The search produced a list in ci/cj.
4431 * Swap the list pointers so we get the outer list is ciOuter,cjOuter
4432 * and we can prune that to get an inner list in ci/cj.
4434 nbnxn_pairlist_t *list = listSet->nbl[i];
4435 list->nciOuter = list->nci;
4437 nbnxn_ci_t *ciTmp = list->ciOuter;
4438 list->ciOuter = list->ci;
4441 nbnxn_cj_t *cjTmp = list->cjOuter;
4442 list->cjOuter = list->cj;
4445 /* Signal that this inner list is currently invalid */