2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017, 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/utility/exceptions.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/gmxomp.h"
72 #include "gromacs/utility/smalloc.h"
74 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
77 /* We shift the i-particles backward for PBC.
78 * This leads to more conditionals than shifting forward.
79 * We do this to get more balanced pair lists.
81 constexpr bool c_pbcShiftBackward = true;
84 static void nbs_cycle_clear(nbnxn_cycle_t *cc)
86 for (int i = 0; i < enbsCCnr; i++)
93 static double Mcyc_av(const nbnxn_cycle_t *cc)
95 return (double)cc->c*1e-6/cc->count;
98 static void nbs_cycle_print(FILE *fp, const nbnxn_search_t nbs)
101 fprintf(fp, "ns %4d grid %4.1f search %4.1f red.f %5.3f",
102 nbs->cc[enbsCCgrid].count,
103 Mcyc_av(&nbs->cc[enbsCCgrid]),
104 Mcyc_av(&nbs->cc[enbsCCsearch]),
105 Mcyc_av(&nbs->cc[enbsCCreducef]));
107 if (nbs->nthread_max > 1)
109 if (nbs->cc[enbsCCcombine].count > 0)
111 fprintf(fp, " comb %5.2f",
112 Mcyc_av(&nbs->cc[enbsCCcombine]));
114 fprintf(fp, " s. th");
115 for (int t = 0; t < nbs->nthread_max; t++)
117 fprintf(fp, " %4.1f",
118 Mcyc_av(&nbs->work[t].cc[enbsCCsearch]));
124 static gmx_inline int ci_to_cj(int ci, int na_cj_2log)
128 case 2: return ci; break;
129 case 3: return (ci>>1); break;
130 case 1: return (ci<<1); break;
138 /* Returns the j-cluster index corresponding to the i-cluster index */
139 template<int cj_size> static gmx_inline int ci_to_cj(int ci)
145 else if (cj_size == 4)
149 else if (cj_size == 8)
155 GMX_ASSERT(false, "Only j-cluster sizes 2, 4, 8 are implemented");
160 /* Returns the index in the coordinate array corresponding to the i-cluster index */
161 template<int cj_size> static gmx_inline int x_ind_ci(int ci)
165 /* Coordinates are stored packed in groups of 4 */
168 else if (cj_size == 8)
170 /* Coordinates packed in 8, i-cluster size is half the packing width */
171 return (ci >> 1)*STRIDE_P8 + (ci & 1)*(c_packX8 >> 1);
175 GMX_ASSERT(false, "Only j-cluster sizes 2, 4, 8 are implemented");
180 /* Returns the index in the coordinate array corresponding to the j-cluster index */
181 template<int cj_size> static gmx_inline int x_ind_cj(int cj)
185 /* Coordinates are stored packed in groups of 4 */
186 return (cj >> 1)*STRIDE_P4 + (cj & 1)*(c_packX4 >> 1);
188 else if (cj_size <= 4)
190 /* Coordinates are stored packed in groups of 4 */
193 else if (cj_size == 8)
195 /* Coordinates are stored packed in groups of 8 */
200 GMX_ASSERT(false, "Only j-cluster sizes 2, 4, 8 are implemented");
205 /* The 6 functions below are only introduced to make the code more readable */
207 static gmx_inline int ci_to_cj_simd_4xn(int ci)
209 return ci_to_cj<GMX_SIMD_REAL_WIDTH>(ci);
212 static gmx_inline int x_ind_ci_simd_4xn(int ci)
214 return x_ind_ci<GMX_SIMD_REAL_WIDTH>(ci);
217 static gmx_inline int x_ind_cj_simd_4xn(int cj)
219 return x_ind_cj<GMX_SIMD_REAL_WIDTH>(cj);
222 static gmx_inline int ci_to_cj_simd_2xnn(int ci)
224 return ci_to_cj<GMX_SIMD_REAL_WIDTH/2>(ci);
227 static gmx_inline int x_ind_ci_simd_2xnn(int ci)
229 return x_ind_ci<GMX_SIMD_REAL_WIDTH/2>(ci);
232 static gmx_inline int x_ind_cj_simd_2xnn(int cj)
234 return x_ind_cj<GMX_SIMD_REAL_WIDTH/2>(cj);
239 gmx_bool nbnxn_kernel_pairlist_simple(int nb_kernel_type)
241 if (nb_kernel_type == nbnxnkNotSet)
243 gmx_fatal(FARGS, "Non-bonded kernel type not set for Verlet-style pair-list.");
246 switch (nb_kernel_type)
248 case nbnxnk8x8x8_GPU:
249 case nbnxnk8x8x8_PlainC:
252 case nbnxnk4x4_PlainC:
253 case nbnxnk4xN_SIMD_4xN:
254 case nbnxnk4xN_SIMD_2xNN:
258 gmx_incons("Invalid nonbonded kernel type passed!");
263 /* Initializes a single nbnxn_pairlist_t data structure */
264 static void nbnxn_init_pairlist_fep(t_nblist *nl)
266 nl->type = GMX_NBLIST_INTERACTION_FREE_ENERGY;
267 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
268 /* The interaction functions are set in the free energy kernel fuction */
281 nl->jindex = nullptr;
283 nl->excl_fep = nullptr;
287 void nbnxn_init_search(nbnxn_search_t * nbs_ptr,
289 struct gmx_domdec_zones_t *zones,
301 nbs->DomDec = (n_dd_cells != nullptr);
303 clear_ivec(nbs->dd_dim);
309 for (int d = 0; d < DIM; d++)
311 if ((*n_dd_cells)[d] > 1)
314 /* Each grid matches a DD zone */
320 nbnxn_grids_init(nbs, ngrid);
323 nbs->cell_nalloc = 0;
327 nbs->nthread_max = nthread_max;
329 /* Initialize the work data structures for each thread */
330 snew(nbs->work, nbs->nthread_max);
331 for (int t = 0; t < nbs->nthread_max; t++)
333 nbs->work[t].cxy_na = nullptr;
334 nbs->work[t].cxy_na_nalloc = 0;
335 nbs->work[t].sort_work = nullptr;
336 nbs->work[t].sort_work_nalloc = 0;
338 snew(nbs->work[t].nbl_fep, 1);
339 nbnxn_init_pairlist_fep(nbs->work[t].nbl_fep);
342 /* Initialize detailed nbsearch cycle counting */
343 nbs->print_cycles = (getenv("GMX_NBNXN_CYCLE") != nullptr);
344 nbs->search_count = 0;
345 nbs_cycle_clear(nbs->cc);
346 for (int t = 0; t < nbs->nthread_max; t++)
348 nbs_cycle_clear(nbs->work[t].cc);
352 static void init_buffer_flags(nbnxn_buffer_flags_t *flags,
355 flags->nflag = (natoms + NBNXN_BUFFERFLAG_SIZE - 1)/NBNXN_BUFFERFLAG_SIZE;
356 if (flags->nflag > flags->flag_nalloc)
358 flags->flag_nalloc = over_alloc_large(flags->nflag);
359 srenew(flags->flag, flags->flag_nalloc);
361 for (int b = 0; b < flags->nflag; b++)
363 bitmask_clear(&(flags->flag[b]));
367 /* Determines the cell range along one dimension that
368 * the bounding box b0 - b1 sees.
370 static void get_cell_range(real b0, real b1,
371 int nc, real c0, real s, real invs,
372 real d2, real r2, int *cf, int *cl)
374 *cf = std::max(static_cast<int>((b0 - c0)*invs), 0);
376 while (*cf > 0 && d2 + gmx::square((b0 - c0) - (*cf-1+1)*s) < r2)
381 *cl = std::min(static_cast<int>((b1 - c0)*invs), nc-1);
382 while (*cl < nc-1 && d2 + gmx::square((*cl+1)*s - (b1 - c0)) < r2)
388 /* Reference code calculating the distance^2 between two bounding boxes */
389 static float box_dist2(float bx0, float bx1, float by0,
390 float by1, float bz0, float bz1,
391 const nbnxn_bb_t *bb)
394 float dl, dh, dm, dm0;
398 dl = bx0 - bb->upper[BB_X];
399 dh = bb->lower[BB_X] - bx1;
400 dm = std::max(dl, dh);
401 dm0 = std::max(dm, 0.0f);
404 dl = by0 - bb->upper[BB_Y];
405 dh = bb->lower[BB_Y] - by1;
406 dm = std::max(dl, dh);
407 dm0 = std::max(dm, 0.0f);
410 dl = bz0 - bb->upper[BB_Z];
411 dh = bb->lower[BB_Z] - bz1;
412 dm = std::max(dl, dh);
413 dm0 = std::max(dm, 0.0f);
419 /* Plain C code calculating the distance^2 between two bounding boxes */
420 static float subc_bb_dist2(int si, const nbnxn_bb_t *bb_i_ci,
421 int csj, const nbnxn_bb_t *bb_j_all)
423 const nbnxn_bb_t *bb_i, *bb_j;
425 float dl, dh, dm, dm0;
428 bb_j = bb_j_all + csj;
432 dl = bb_i->lower[BB_X] - bb_j->upper[BB_X];
433 dh = bb_j->lower[BB_X] - bb_i->upper[BB_X];
434 dm = std::max(dl, dh);
435 dm0 = std::max(dm, 0.0f);
438 dl = bb_i->lower[BB_Y] - bb_j->upper[BB_Y];
439 dh = bb_j->lower[BB_Y] - bb_i->upper[BB_Y];
440 dm = std::max(dl, dh);
441 dm0 = std::max(dm, 0.0f);
444 dl = bb_i->lower[BB_Z] - bb_j->upper[BB_Z];
445 dh = bb_j->lower[BB_Z] - bb_i->upper[BB_Z];
446 dm = std::max(dl, dh);
447 dm0 = std::max(dm, 0.0f);
453 #if NBNXN_SEARCH_BB_SIMD4
455 /* 4-wide SIMD code for bb distance for bb format xyz0 */
456 static float subc_bb_dist2_simd4(int si, const nbnxn_bb_t *bb_i_ci,
457 int csj, const nbnxn_bb_t *bb_j_all)
459 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
462 Simd4Float bb_i_S0, bb_i_S1;
463 Simd4Float bb_j_S0, bb_j_S1;
469 bb_i_S0 = load4(&bb_i_ci[si].lower[0]);
470 bb_i_S1 = load4(&bb_i_ci[si].upper[0]);
471 bb_j_S0 = load4(&bb_j_all[csj].lower[0]);
472 bb_j_S1 = load4(&bb_j_all[csj].upper[0]);
474 dl_S = bb_i_S0 - bb_j_S1;
475 dh_S = bb_j_S0 - bb_i_S1;
477 dm_S = max(dl_S, dh_S);
478 dm0_S = max(dm_S, simd4SetZeroF());
480 return dotProduct(dm0_S, dm0_S);
483 /* Calculate bb bounding distances of bb_i[si,...,si+3] and store them in d2 */
484 #define SUBC_BB_DIST2_SIMD4_XXXX_INNER(si, bb_i, d2) \
488 Simd4Float dx_0, dy_0, dz_0; \
489 Simd4Float dx_1, dy_1, dz_1; \
491 Simd4Float mx, my, mz; \
492 Simd4Float m0x, m0y, m0z; \
494 Simd4Float d2x, d2y, d2z; \
495 Simd4Float d2s, d2t; \
497 shi = si*NNBSBB_D*DIM; \
499 xi_l = load4(bb_i+shi+0*STRIDE_PBB); \
500 yi_l = load4(bb_i+shi+1*STRIDE_PBB); \
501 zi_l = load4(bb_i+shi+2*STRIDE_PBB); \
502 xi_h = load4(bb_i+shi+3*STRIDE_PBB); \
503 yi_h = load4(bb_i+shi+4*STRIDE_PBB); \
504 zi_h = load4(bb_i+shi+5*STRIDE_PBB); \
506 dx_0 = xi_l - xj_h; \
507 dy_0 = yi_l - yj_h; \
508 dz_0 = zi_l - zj_h; \
510 dx_1 = xj_l - xi_h; \
511 dy_1 = yj_l - yi_h; \
512 dz_1 = zj_l - zi_h; \
514 mx = max(dx_0, dx_1); \
515 my = max(dy_0, dy_1); \
516 mz = max(dz_0, dz_1); \
518 m0x = max(mx, zero); \
519 m0y = max(my, zero); \
520 m0z = max(mz, zero); \
529 store4(d2+si, d2t); \
532 /* 4-wide SIMD code for nsi bb distances for bb format xxxxyyyyzzzz */
533 static void subc_bb_dist2_simd4_xxxx(const float *bb_j,
534 int nsi, const float *bb_i,
537 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
540 Simd4Float xj_l, yj_l, zj_l;
541 Simd4Float xj_h, yj_h, zj_h;
542 Simd4Float xi_l, yi_l, zi_l;
543 Simd4Float xi_h, yi_h, zi_h;
549 xj_l = Simd4Float(bb_j[0*STRIDE_PBB]);
550 yj_l = Simd4Float(bb_j[1*STRIDE_PBB]);
551 zj_l = Simd4Float(bb_j[2*STRIDE_PBB]);
552 xj_h = Simd4Float(bb_j[3*STRIDE_PBB]);
553 yj_h = Simd4Float(bb_j[4*STRIDE_PBB]);
554 zj_h = Simd4Float(bb_j[5*STRIDE_PBB]);
556 /* Here we "loop" over si (0,STRIDE_PBB) from 0 to nsi with step STRIDE_PBB.
557 * But as we know the number of iterations is 1 or 2, we unroll manually.
559 SUBC_BB_DIST2_SIMD4_XXXX_INNER(0, bb_i, d2);
560 if (STRIDE_PBB < nsi)
562 SUBC_BB_DIST2_SIMD4_XXXX_INNER(STRIDE_PBB, bb_i, d2);
566 #endif /* NBNXN_SEARCH_BB_SIMD4 */
569 /* Returns if any atom pair from two clusters is within distance sqrt(rlist2) */
570 static gmx_inline gmx_bool
571 clusterpair_in_range(const nbnxn_list_work_t *work,
573 int csj, int stride, const real *x_j,
576 #if !GMX_SIMD4_HAVE_REAL
579 * All coordinates are stored as xyzxyz...
582 const real *x_i = work->x_ci;
584 for (int i = 0; i < c_nbnxnGpuClusterSize; i++)
586 int i0 = (si*c_nbnxnGpuClusterSize + i)*DIM;
587 for (int j = 0; j < c_nbnxnGpuClusterSize; j++)
589 int j0 = (csj*c_nbnxnGpuClusterSize + j)*stride;
591 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]);
602 #else /* !GMX_SIMD4_HAVE_REAL */
604 /* 4-wide SIMD version.
605 * A cluster is hard-coded to 8 atoms.
606 * The coordinates x_i are stored as xxxxyyyy..., x_j is stored xyzxyz...
607 * Using 8-wide AVX(2) is not faster on Intel Sandy Bridge and Haswell.
609 assert(c_nbnxnGpuClusterSize == 8);
611 Simd4Real rc2_S = Simd4Real(rlist2);
613 const real *x_i = work->x_ci_simd;
615 int dim_stride = c_nbnxnGpuClusterSize*DIM;
616 Simd4Real ix_S0 = load4(x_i + si*dim_stride + 0*GMX_SIMD4_WIDTH);
617 Simd4Real iy_S0 = load4(x_i + si*dim_stride + 1*GMX_SIMD4_WIDTH);
618 Simd4Real iz_S0 = load4(x_i + si*dim_stride + 2*GMX_SIMD4_WIDTH);
619 Simd4Real ix_S1 = load4(x_i + si*dim_stride + 3*GMX_SIMD4_WIDTH);
620 Simd4Real iy_S1 = load4(x_i + si*dim_stride + 4*GMX_SIMD4_WIDTH);
621 Simd4Real iz_S1 = load4(x_i + si*dim_stride + 5*GMX_SIMD4_WIDTH);
623 /* We loop from the outer to the inner particles to maximize
624 * the chance that we find a pair in range quickly and return.
626 int j0 = csj*c_nbnxnGpuClusterSize;
627 int j1 = j0 + c_nbnxnGpuClusterSize - 1;
630 Simd4Real jx0_S, jy0_S, jz0_S;
631 Simd4Real jx1_S, jy1_S, jz1_S;
633 Simd4Real dx_S0, dy_S0, dz_S0;
634 Simd4Real dx_S1, dy_S1, dz_S1;
635 Simd4Real dx_S2, dy_S2, dz_S2;
636 Simd4Real dx_S3, dy_S3, dz_S3;
647 Simd4Bool wco_any_S01, wco_any_S23, wco_any_S;
649 jx0_S = Simd4Real(x_j[j0*stride+0]);
650 jy0_S = Simd4Real(x_j[j0*stride+1]);
651 jz0_S = Simd4Real(x_j[j0*stride+2]);
653 jx1_S = Simd4Real(x_j[j1*stride+0]);
654 jy1_S = Simd4Real(x_j[j1*stride+1]);
655 jz1_S = Simd4Real(x_j[j1*stride+2]);
657 /* Calculate distance */
658 dx_S0 = ix_S0 - jx0_S;
659 dy_S0 = iy_S0 - jy0_S;
660 dz_S0 = iz_S0 - jz0_S;
661 dx_S1 = ix_S1 - jx0_S;
662 dy_S1 = iy_S1 - jy0_S;
663 dz_S1 = iz_S1 - jz0_S;
664 dx_S2 = ix_S0 - jx1_S;
665 dy_S2 = iy_S0 - jy1_S;
666 dz_S2 = iz_S0 - jz1_S;
667 dx_S3 = ix_S1 - jx1_S;
668 dy_S3 = iy_S1 - jy1_S;
669 dz_S3 = iz_S1 - jz1_S;
671 /* rsq = dx*dx+dy*dy+dz*dz */
672 rsq_S0 = norm2(dx_S0, dy_S0, dz_S0);
673 rsq_S1 = norm2(dx_S1, dy_S1, dz_S1);
674 rsq_S2 = norm2(dx_S2, dy_S2, dz_S2);
675 rsq_S3 = norm2(dx_S3, dy_S3, dz_S3);
677 wco_S0 = (rsq_S0 < rc2_S);
678 wco_S1 = (rsq_S1 < rc2_S);
679 wco_S2 = (rsq_S2 < rc2_S);
680 wco_S3 = (rsq_S3 < rc2_S);
682 wco_any_S01 = wco_S0 || wco_S1;
683 wco_any_S23 = wco_S2 || wco_S3;
684 wco_any_S = wco_any_S01 || wco_any_S23;
686 if (anyTrue(wco_any_S))
697 #endif /* !GMX_SIMD4_HAVE_REAL */
700 /* Returns the j sub-cell for index cj_ind */
701 static int nbl_cj(const nbnxn_pairlist_t *nbl, int cj_ind)
703 return nbl->cj4[cj_ind/c_nbnxnGpuJgroupSize].cj[cj_ind & (c_nbnxnGpuJgroupSize - 1)];
706 /* Returns the i-interaction mask of the j sub-cell for index cj_ind */
707 static unsigned int nbl_imask0(const nbnxn_pairlist_t *nbl, int cj_ind)
709 return nbl->cj4[cj_ind/c_nbnxnGpuJgroupSize].imei[0].imask;
712 /* Ensures there is enough space for extra extra exclusion masks */
713 static void check_excl_space(nbnxn_pairlist_t *nbl, int extra)
715 if (nbl->nexcl+extra > nbl->excl_nalloc)
717 nbl->excl_nalloc = over_alloc_small(nbl->nexcl+extra);
718 nbnxn_realloc_void((void **)&nbl->excl,
719 nbl->nexcl*sizeof(*nbl->excl),
720 nbl->excl_nalloc*sizeof(*nbl->excl),
721 nbl->alloc, nbl->free);
725 /* Ensures there is enough space for maxNumExtraClusters extra j-clusters in the list */
726 static void check_cell_list_space_simple(nbnxn_pairlist_t *nbl,
727 int maxNumExtraClusters)
731 cj_max = nbl->ncj + maxNumExtraClusters;
733 if (cj_max > nbl->cj_nalloc)
735 nbl->cj_nalloc = over_alloc_small(cj_max);
736 nbnxn_realloc_void((void **)&nbl->cj,
737 nbl->ncj*sizeof(*nbl->cj),
738 nbl->cj_nalloc*sizeof(*nbl->cj),
739 nbl->alloc, nbl->free);
743 /* Ensures there is enough space for ncell extra j-clusters in the list */
744 static void check_cell_list_space_supersub(nbnxn_pairlist_t *nbl,
749 /* We can have maximally nsupercell*c_gpuNumClusterPerCell sj lists */
750 /* We can store 4 j-subcell - i-supercell pairs in one struct.
751 * since we round down, we need one extra entry.
753 ncj4_max = ((nbl->work->cj_ind + ncell*c_gpuNumClusterPerCell + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize);
755 if (ncj4_max > nbl->cj4_nalloc)
757 nbl->cj4_nalloc = over_alloc_small(ncj4_max);
758 nbnxn_realloc_void((void **)&nbl->cj4,
759 nbl->work->cj4_init*sizeof(*nbl->cj4),
760 nbl->cj4_nalloc*sizeof(*nbl->cj4),
761 nbl->alloc, nbl->free);
764 if (ncj4_max > nbl->work->cj4_init)
766 for (int j4 = nbl->work->cj4_init; j4 < ncj4_max; j4++)
768 /* No i-subcells and no excl's in the list initially */
769 for (w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
771 nbl->cj4[j4].imei[w].imask = 0U;
772 nbl->cj4[j4].imei[w].excl_ind = 0;
776 nbl->work->cj4_init = ncj4_max;
780 /* Set all excl masks for one GPU warp no exclusions */
781 static void set_no_excls(nbnxn_excl_t *excl)
783 for (int t = 0; t < c_nbnxnGpuExclSize; t++)
785 /* Turn all interaction bits on */
786 excl->pair[t] = NBNXN_INTERACTION_MASK_ALL;
790 /* Initializes a single nbnxn_pairlist_t data structure */
791 static void nbnxn_init_pairlist(nbnxn_pairlist_t *nbl,
793 nbnxn_alloc_t *alloc,
796 if (alloc == nullptr)
798 nbl->alloc = nbnxn_alloc_aligned;
806 nbl->free = nbnxn_free_aligned;
813 nbl->bSimple = bSimple;
828 /* We need one element extra in sj, so alloc initially with 1 */
835 GMX_ASSERT(c_nbnxnGpuNumClusterPerSupercluster == c_gpuNumClusterPerCell, "The search code assumes that the a super-cluster matches a search grid cell");
837 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");
838 GMX_ASSERT(sizeof(nbl->excl[0])*8 >= c_nbnxnGpuJgroupSize*c_gpuNumClusterPerCell, "The GPU exclusion mask does not contain a sufficient number of bits");
841 nbl->excl_nalloc = 0;
843 check_excl_space(nbl, 1);
845 set_no_excls(&nbl->excl[0]);
851 snew_aligned(nbl->work->bb_ci, 1, NBNXN_SEARCH_BB_MEM_ALIGN);
856 snew_aligned(nbl->work->pbb_ci, c_gpuNumClusterPerCell/STRIDE_PBB*NNBSBB_XXXX, NBNXN_SEARCH_BB_MEM_ALIGN);
858 snew_aligned(nbl->work->bb_ci, c_gpuNumClusterPerCell, NBNXN_SEARCH_BB_MEM_ALIGN);
861 int gpu_clusterpair_nc = c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize*DIM;
862 snew(nbl->work->x_ci, gpu_clusterpair_nc);
864 snew_aligned(nbl->work->x_ci_simd,
865 std::max(NBNXN_CPU_CLUSTER_I_SIZE*DIM*GMX_SIMD_REAL_WIDTH,
867 GMX_SIMD_REAL_WIDTH);
869 snew_aligned(nbl->work->d2, c_gpuNumClusterPerCell, NBNXN_SEARCH_BB_MEM_ALIGN);
871 nbl->work->sort = nullptr;
872 nbl->work->sort_nalloc = 0;
873 nbl->work->sci_sort = nullptr;
874 nbl->work->sci_sort_nalloc = 0;
877 void nbnxn_init_pairlist_set(nbnxn_pairlist_set_t *nbl_list,
878 gmx_bool bSimple, gmx_bool bCombined,
879 nbnxn_alloc_t *alloc,
882 nbl_list->bSimple = bSimple;
883 nbl_list->bCombined = bCombined;
885 nbl_list->nnbl = gmx_omp_nthreads_get(emntNonbonded);
887 if (!nbl_list->bCombined &&
888 nbl_list->nnbl > NBNXN_BUFFERFLAG_MAX_THREADS)
890 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.",
891 nbl_list->nnbl, NBNXN_BUFFERFLAG_MAX_THREADS, NBNXN_BUFFERFLAG_MAX_THREADS);
894 snew(nbl_list->nbl, nbl_list->nnbl);
895 if (bSimple && nbl_list->nnbl > 1)
897 snew(nbl_list->nbl_work, nbl_list->nnbl);
899 snew(nbl_list->nbl_fep, nbl_list->nnbl);
900 /* Execute in order to avoid memory interleaving between threads */
901 #pragma omp parallel for num_threads(nbl_list->nnbl) schedule(static)
902 for (int i = 0; i < nbl_list->nnbl; i++)
906 /* Allocate the nblist data structure locally on each thread
907 * to optimize memory access for NUMA architectures.
909 snew(nbl_list->nbl[i], 1);
911 /* Only list 0 is used on the GPU, use normal allocation for i>0 */
912 if (!bSimple && i == 0)
914 nbnxn_init_pairlist(nbl_list->nbl[i], nbl_list->bSimple, alloc, free);
918 nbnxn_init_pairlist(nbl_list->nbl[i], nbl_list->bSimple, nullptr, nullptr);
919 if (bSimple && nbl_list->nnbl > 1)
921 snew(nbl_list->nbl_work[i], 1);
922 nbnxn_init_pairlist(nbl_list->nbl_work[i], nbl_list->bSimple, nullptr, nullptr);
926 snew(nbl_list->nbl_fep[i], 1);
927 nbnxn_init_pairlist_fep(nbl_list->nbl_fep[i]);
929 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
933 /* Print statistics of a pair list, used for debug output */
934 static void print_nblist_statistics_simple(FILE *fp, const nbnxn_pairlist_t *nbl,
935 const nbnxn_search_t nbs, real rl)
937 const nbnxn_grid_t *grid;
941 grid = &nbs->grid[0];
943 fprintf(fp, "nbl nci %d ncj %d\n",
944 nbl->nci, nbl->ncjInUse);
945 fprintf(fp, "nbl na_sc %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
946 nbl->na_sc, rl, nbl->ncjInUse, nbl->ncjInUse/(double)grid->nc,
947 nbl->ncjInUse/(double)grid->nc*grid->na_sc,
948 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])));
950 fprintf(fp, "nbl average j cell list length %.1f\n",
951 0.25*nbl->ncjInUse/(double)std::max(nbl->nci, 1));
953 for (int s = 0; s < SHIFTS; s++)
958 for (int i = 0; i < nbl->nci; i++)
960 cs[nbl->ci[i].shift & NBNXN_CI_SHIFT] +=
961 nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start;
963 int j = nbl->ci[i].cj_ind_start;
964 while (j < nbl->ci[i].cj_ind_end &&
965 nbl->cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
971 fprintf(fp, "nbl cell pairs, total: %d excl: %d %.1f%%\n",
972 nbl->ncj, npexcl, 100*npexcl/(double)std::max(nbl->ncj, 1));
973 for (int s = 0; s < SHIFTS; s++)
977 fprintf(fp, "nbl shift %2d ncj %3d\n", s, cs[s]);
982 /* Print statistics of a pair lists, used for debug output */
983 static void print_nblist_statistics_supersub(FILE *fp, const nbnxn_pairlist_t *nbl,
984 const nbnxn_search_t nbs, real rl)
986 const nbnxn_grid_t *grid;
988 int c[c_gpuNumClusterPerCell + 1];
989 double sum_nsp, sum_nsp2;
992 /* This code only produces correct statistics with domain decomposition */
993 grid = &nbs->grid[0];
995 fprintf(fp, "nbl nsci %d ncj4 %d nsi %d excl4 %d\n",
996 nbl->nsci, nbl->ncj4, nbl->nci_tot, nbl->nexcl);
997 fprintf(fp, "nbl na_c %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
998 nbl->na_ci, rl, nbl->nci_tot, nbl->nci_tot/(double)grid->nsubc_tot,
999 nbl->nci_tot/(double)grid->nsubc_tot*grid->na_c,
1000 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])));
1005 for (int si = 0; si <= c_gpuNumClusterPerCell; si++)
1009 for (int i = 0; i < nbl->nsci; i++)
1014 for (int j4 = nbl->sci[i].cj4_ind_start; j4 < nbl->sci[i].cj4_ind_end; j4++)
1016 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
1019 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
1021 if (nbl->cj4[j4].imei[0].imask & (1U << (j*c_gpuNumClusterPerCell + si)))
1031 sum_nsp2 += nsp*nsp;
1032 nsp_max = std::max(nsp_max, nsp);
1036 sum_nsp /= nbl->nsci;
1037 sum_nsp2 /= nbl->nsci;
1039 fprintf(fp, "nbl #cluster-pairs: av %.1f stddev %.1f max %d\n",
1040 sum_nsp, std::sqrt(sum_nsp2 - sum_nsp*sum_nsp), nsp_max);
1044 for (b = 0; b <= c_gpuNumClusterPerCell; b++)
1046 fprintf(fp, "nbl j-list #i-subcell %d %7d %4.1f\n",
1048 100.0*c[b]/(double)(nbl->ncj4*c_nbnxnGpuJgroupSize));
1053 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp */
1054 static void low_get_nbl_exclusions(nbnxn_pairlist_t *nbl, int cj4,
1055 int warp, nbnxn_excl_t **excl)
1057 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
1059 /* No exclusions set, make a new list entry */
1060 nbl->cj4[cj4].imei[warp].excl_ind = nbl->nexcl;
1062 *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
1063 set_no_excls(*excl);
1067 /* We already have some exclusions, new ones can be added to the list */
1068 *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
1072 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp,
1073 * generates a new element and allocates extra memory, if necessary.
1075 static void get_nbl_exclusions_1(nbnxn_pairlist_t *nbl, int cj4,
1076 int warp, nbnxn_excl_t **excl)
1078 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
1080 /* We need to make a new list entry, check if we have space */
1081 check_excl_space(nbl, 1);
1083 low_get_nbl_exclusions(nbl, cj4, warp, excl);
1086 /* Returns pointers to the exclusion masks for cj4-unit cj4 for both warps,
1087 * generates a new element and allocates extra memory, if necessary.
1089 static void get_nbl_exclusions_2(nbnxn_pairlist_t *nbl, int cj4,
1090 nbnxn_excl_t **excl_w0,
1091 nbnxn_excl_t **excl_w1)
1093 /* Check for space we might need */
1094 check_excl_space(nbl, 2);
1096 low_get_nbl_exclusions(nbl, cj4, 0, excl_w0);
1097 low_get_nbl_exclusions(nbl, cj4, 1, excl_w1);
1100 /* Sets the self exclusions i=j and pair exclusions i>j */
1101 static void set_self_and_newton_excls_supersub(nbnxn_pairlist_t *nbl,
1102 int cj4_ind, int sj_offset,
1103 int i_cluster_in_cell)
1105 nbnxn_excl_t *excl[c_nbnxnGpuClusterpairSplit];
1107 /* Here we only set the set self and double pair exclusions */
1109 assert(c_nbnxnGpuClusterpairSplit == 2);
1111 get_nbl_exclusions_2(nbl, cj4_ind, &excl[0], &excl[1]);
1113 /* Only minor < major bits set */
1114 for (int ej = 0; ej < nbl->na_ci; ej++)
1117 for (int ei = ej; ei < nbl->na_ci; ei++)
1119 excl[w]->pair[(ej & (c_nbnxnGpuJgroupSize-1))*nbl->na_ci + ei] &=
1120 ~(1U << (sj_offset*c_gpuNumClusterPerCell + i_cluster_in_cell));
1125 /* Returns a diagonal or off-diagonal interaction mask for plain C lists */
1126 static unsigned int get_imask(gmx_bool rdiag, int ci, int cj)
1128 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
1131 /* Returns a diagonal or off-diagonal interaction mask for cj-size=2 */
1132 static unsigned int get_imask_simd_j2(gmx_bool rdiag, int ci, int cj)
1134 return (rdiag && ci*2 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_0 :
1135 (rdiag && ci*2+1 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_1 :
1136 NBNXN_INTERACTION_MASK_ALL));
1139 /* Returns a diagonal or off-diagonal interaction mask for cj-size=4 */
1140 static unsigned int get_imask_simd_j4(gmx_bool rdiag, int ci, int cj)
1142 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
1145 /* Returns a diagonal or off-diagonal interaction mask for cj-size=8 */
1146 static unsigned int get_imask_simd_j8(gmx_bool rdiag, int ci, int cj)
1148 return (rdiag && ci == cj*2 ? NBNXN_INTERACTION_MASK_DIAG_J8_0 :
1149 (rdiag && ci == cj*2+1 ? NBNXN_INTERACTION_MASK_DIAG_J8_1 :
1150 NBNXN_INTERACTION_MASK_ALL));
1154 #if GMX_SIMD_REAL_WIDTH == 2
1155 #define get_imask_simd_4xn get_imask_simd_j2
1157 #if GMX_SIMD_REAL_WIDTH == 4
1158 #define get_imask_simd_4xn get_imask_simd_j4
1160 #if GMX_SIMD_REAL_WIDTH == 8
1161 #define get_imask_simd_4xn get_imask_simd_j8
1162 #define get_imask_simd_2xnn get_imask_simd_j4
1164 #if GMX_SIMD_REAL_WIDTH == 16
1165 #define get_imask_simd_2xnn get_imask_simd_j8
1169 /* Plain C code for checking and adding cluster-pairs to the list.
1171 * \param[in] gridj The j-grid
1172 * \param[in,out] nbl The pair-list to store the cluster pairs in
1173 * \param[in] icluster The index of the i-cluster
1174 * \param[in] jclusterFirst The first cluster in the j-range
1175 * \param[in] jclusterLast The last cluster in the j-range
1176 * \param[in] excludeSubDiagonal Exclude atom pairs with i-index > j-index
1177 * \param[in] x_j Coordinates for the j-atom, in xyz format
1178 * \param[in] rlist2 The squared list cut-off
1179 * \param[in] rbb2 The squared cut-off for putting cluster-pairs in the list based on bounding box distance only
1180 * \param[in,out] numDistanceChecks The number of distance checks performed
1183 makeClusterListSimple(const nbnxn_grid_t * gridj,
1184 nbnxn_pairlist_t * nbl,
1188 bool excludeSubDiagonal,
1189 const real * gmx_restrict x_j,
1192 int * gmx_restrict numDistanceChecks)
1194 const nbnxn_bb_t * gmx_restrict bb_ci = nbl->work->bb_ci;
1195 const real * gmx_restrict x_ci = nbl->work->x_ci;
1200 while (!InRange && jclusterFirst <= jclusterLast)
1202 real d2 = subc_bb_dist2(0, bb_ci, jclusterFirst, gridj->bb);
1203 *numDistanceChecks += 2;
1205 /* Check if the distance is within the distance where
1206 * we use only the bounding box distance rbb,
1207 * or within the cut-off and there is at least one atom pair
1208 * within the cut-off.
1214 else if (d2 < rlist2)
1216 int cjf_gl = gridj->cell0 + jclusterFirst;
1217 for (int i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
1219 for (int j = 0; j < NBNXN_CPU_CLUSTER_I_SIZE; j++)
1221 InRange = InRange ||
1222 (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
1223 gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
1224 gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rlist2);
1227 *numDistanceChecks += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
1240 while (!InRange && jclusterLast > jclusterFirst)
1242 real d2 = subc_bb_dist2(0, bb_ci, jclusterLast, gridj->bb);
1243 *numDistanceChecks += 2;
1245 /* Check if the distance is within the distance where
1246 * we use only the bounding box distance rbb,
1247 * or within the cut-off and there is at least one atom pair
1248 * within the cut-off.
1254 else if (d2 < rlist2)
1256 int cjl_gl = gridj->cell0 + jclusterLast;
1257 for (int i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
1259 for (int j = 0; j < NBNXN_CPU_CLUSTER_I_SIZE; j++)
1261 InRange = InRange ||
1262 (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
1263 gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
1264 gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rlist2);
1267 *numDistanceChecks += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
1275 if (jclusterFirst <= jclusterLast)
1277 for (int jcluster = jclusterFirst; jcluster <= jclusterLast; jcluster++)
1279 /* Store cj and the interaction mask */
1280 nbl->cj[nbl->ncj].cj = gridj->cell0 + jcluster;
1281 nbl->cj[nbl->ncj].excl = get_imask(excludeSubDiagonal, icluster, jcluster);
1284 /* Increase the closing index in i super-cell list */
1285 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
1289 #ifdef GMX_NBNXN_SIMD_4XN
1290 #include "gromacs/mdlib/nbnxn_search_simd_4xn.h"
1292 #ifdef GMX_NBNXN_SIMD_2XNN
1293 #include "gromacs/mdlib/nbnxn_search_simd_2xnn.h"
1296 /* Plain C or SIMD4 code for making a pair list of super-cell sci vs scj.
1297 * Checks bounding box distances and possibly atom pair distances.
1299 static void make_cluster_list_supersub(const nbnxn_grid_t *gridi,
1300 const nbnxn_grid_t *gridj,
1301 nbnxn_pairlist_t *nbl,
1303 gmx_bool sci_equals_scj,
1304 int stride, const real *x,
1305 real rlist2, float rbb2,
1306 int *numDistanceChecks)
1308 nbnxn_list_work_t *work = nbl->work;
1311 const float *pbb_ci = work->pbb_ci;
1313 const nbnxn_bb_t *bb_ci = work->bb_ci;
1316 assert(c_nbnxnGpuClusterSize == gridi->na_c);
1317 assert(c_nbnxnGpuClusterSize == gridj->na_c);
1319 /* We generate the pairlist mainly based on bounding-box distances
1320 * and do atom pair distance based pruning on the GPU.
1321 * Only if a j-group contains a single cluster-pair, we try to prune
1322 * that pair based on atom distances on the CPU to avoid empty j-groups.
1324 #define PRUNE_LIST_CPU_ONE 1
1325 #define PRUNE_LIST_CPU_ALL 0
1327 #if PRUNE_LIST_CPU_ONE
1331 float *d2l = work->d2;
1333 for (int subc = 0; subc < gridj->nsubc[scj]; subc++)
1335 int cj4_ind = nbl->work->cj_ind/c_nbnxnGpuJgroupSize;
1336 int cj_offset = nbl->work->cj_ind - cj4_ind*c_nbnxnGpuJgroupSize;
1337 nbnxn_cj4_t *cj4 = &nbl->cj4[cj4_ind];
1339 int cj = scj*c_gpuNumClusterPerCell + subc;
1341 int cj_gl = gridj->cell0*c_gpuNumClusterPerCell + cj;
1343 /* Initialize this j-subcell i-subcell list */
1344 cj4->cj[cj_offset] = cj_gl;
1353 ci1 = gridi->nsubc[sci];
1357 /* Determine all ci1 bb distances in one call with SIMD4 */
1358 subc_bb_dist2_simd4_xxxx(gridj->pbb+(cj>>STRIDE_PBB_2LOG)*NNBSBB_XXXX+(cj & (STRIDE_PBB-1)),
1360 *numDistanceChecks += c_nbnxnGpuClusterSize*2;
1364 unsigned int imask = 0;
1365 /* We use a fixed upper-bound instead of ci1 to help optimization */
1366 for (int ci = 0; ci < c_gpuNumClusterPerCell; ci++)
1374 /* Determine the bb distance between ci and cj */
1375 d2l[ci] = subc_bb_dist2(ci, bb_ci, cj, gridj->bb);
1376 *numDistanceChecks += 2;
1380 #if PRUNE_LIST_CPU_ALL
1381 /* Check if the distance is within the distance where
1382 * we use only the bounding box distance rbb,
1383 * or within the cut-off and there is at least one atom pair
1384 * within the cut-off. This check is very costly.
1386 *numDistanceChecks += c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize;
1389 clusterpair_in_range(work, ci, cj_gl, stride, x, rlist2)))
1391 /* Check if the distance between the two bounding boxes
1392 * in within the pair-list cut-off.
1397 /* Flag this i-subcell to be taken into account */
1398 imask |= (1U << (cj_offset*c_gpuNumClusterPerCell + ci));
1400 #if PRUNE_LIST_CPU_ONE
1408 #if PRUNE_LIST_CPU_ONE
1409 /* If we only found 1 pair, check if any atoms are actually
1410 * within the cut-off, so we could get rid of it.
1412 if (npair == 1 && d2l[ci_last] >= rbb2 &&
1413 !clusterpair_in_range(work, ci_last, cj_gl, stride, x, rlist2))
1415 imask &= ~(1U << (cj_offset*c_gpuNumClusterPerCell + ci_last));
1422 /* We have a useful sj entry, close it now */
1424 /* Set the exclusions for the ci==sj entry.
1425 * Here we don't bother to check if this entry is actually flagged,
1426 * as it will nearly always be in the list.
1430 set_self_and_newton_excls_supersub(nbl, cj4_ind, cj_offset, subc);
1433 /* Copy the cluster interaction mask to the list */
1434 for (int w = 0; w < c_nbnxnGpuClusterpairSplit; w++)
1436 cj4->imei[w].imask |= imask;
1439 nbl->work->cj_ind++;
1441 /* Keep the count */
1442 nbl->nci_tot += npair;
1444 /* Increase the closing index in i super-cell list */
1445 nbl->sci[nbl->nsci].cj4_ind_end =
1446 (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize;
1451 /* Set all atom-pair exclusions from the topology stored in excl
1452 * as masks in the pair-list for simple list i-entry nbl_ci
1454 static void set_ci_top_excls(const nbnxn_search_t nbs,
1455 nbnxn_pairlist_t *nbl,
1456 gmx_bool diagRemoved,
1459 const nbnxn_ci_t *nbl_ci,
1460 const t_blocka *excl)
1464 int cj_ind_first, cj_ind_last;
1465 int cj_first, cj_last;
1467 int ai, aj, si, ge, se;
1468 int found, cj_ind_0, cj_ind_1, cj_ind_m;
1470 int inner_i, inner_e;
1474 if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
1482 cj_ind_first = nbl_ci->cj_ind_start;
1483 cj_ind_last = nbl->ncj - 1;
1485 cj_first = nbl->cj[cj_ind_first].cj;
1486 cj_last = nbl->cj[cj_ind_last].cj;
1488 /* Determine how many contiguous j-cells we have starting
1489 * from the first i-cell. This number can be used to directly
1490 * calculate j-cell indices for excluded atoms.
1493 if (na_ci_2log == na_cj_2log)
1495 while (cj_ind_first + ndirect <= cj_ind_last &&
1496 nbl->cj[cj_ind_first+ndirect].cj == ci + ndirect)
1501 #if NBNXN_SEARCH_BB_SIMD4
1504 while (cj_ind_first + ndirect <= cj_ind_last &&
1505 nbl->cj[cj_ind_first+ndirect].cj == ci_to_cj(ci, na_cj_2log) + ndirect)
1512 /* Loop over the atoms in the i super-cell */
1513 for (int i = 0; i < nbl->na_sc; i++)
1515 ai = nbs->a[ci*nbl->na_sc+i];
1518 si = (i>>na_ci_2log);
1520 /* Loop over the topology-based exclusions for this i-atom */
1521 for (int eind = excl->index[ai]; eind < excl->index[ai+1]; eind++)
1527 /* The self exclusion are already set, save some time */
1533 /* Without shifts we only calculate interactions j>i
1534 * for one-way pair-lists.
1536 if (diagRemoved && ge <= ci*nbl->na_sc + i)
1541 se = (ge >> na_cj_2log);
1543 /* Could the cluster se be in our list? */
1544 if (se >= cj_first && se <= cj_last)
1546 if (se < cj_first + ndirect)
1548 /* We can calculate cj_ind directly from se */
1549 found = cj_ind_first + se - cj_first;
1553 /* Search for se using bisection */
1555 cj_ind_0 = cj_ind_first + ndirect;
1556 cj_ind_1 = cj_ind_last + 1;
1557 while (found == -1 && cj_ind_0 < cj_ind_1)
1559 cj_ind_m = (cj_ind_0 + cj_ind_1)>>1;
1561 cj_m = nbl->cj[cj_ind_m].cj;
1569 cj_ind_1 = cj_ind_m;
1573 cj_ind_0 = cj_ind_m + 1;
1580 inner_i = i - (si << na_ci_2log);
1581 inner_e = ge - (se << na_cj_2log);
1583 nbl->cj[found].excl &= ~(1U<<((inner_i<<na_cj_2log) + inner_e));
1591 /* Add a new i-entry to the FEP list and copy the i-properties */
1592 static gmx_inline void fep_list_new_nri_copy(t_nblist *nlist)
1594 /* Add a new i-entry */
1597 assert(nlist->nri < nlist->maxnri);
1599 /* Duplicate the last i-entry, except for jindex, which continues */
1600 nlist->iinr[nlist->nri] = nlist->iinr[nlist->nri-1];
1601 nlist->shift[nlist->nri] = nlist->shift[nlist->nri-1];
1602 nlist->gid[nlist->nri] = nlist->gid[nlist->nri-1];
1603 nlist->jindex[nlist->nri] = nlist->nrj;
1606 /* For load balancing of the free-energy lists over threads, we set
1607 * the maximum nrj size of an i-entry to 40. This leads to good
1608 * load balancing in the worst case scenario of a single perturbed
1609 * particle on 16 threads, while not introducing significant overhead.
1610 * Note that half of the perturbed pairs will anyhow end up in very small lists,
1611 * since non perturbed i-particles will see few perturbed j-particles).
1613 const int max_nrj_fep = 40;
1615 /* Exclude the perturbed pairs from the Verlet list. This is only done to avoid
1616 * singularities for overlapping particles (0/0), since the charges and
1617 * LJ parameters have been zeroed in the nbnxn data structure.
1618 * Simultaneously make a group pair list for the perturbed pairs.
1620 static void make_fep_list(const nbnxn_search_t nbs,
1621 const nbnxn_atomdata_t *nbat,
1622 nbnxn_pairlist_t *nbl,
1623 gmx_bool bDiagRemoved,
1625 const nbnxn_grid_t *gridi,
1626 const nbnxn_grid_t *gridj,
1629 int ci, cj_ind_start, cj_ind_end, cja, cjr;
1631 int ngid, gid_i = 0, gid_j, gid;
1632 int egp_shift, egp_mask;
1634 int ind_i, ind_j, ai, aj;
1636 gmx_bool bFEP_i, bFEP_i_all;
1638 if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
1646 cj_ind_start = nbl_ci->cj_ind_start;
1647 cj_ind_end = nbl_ci->cj_ind_end;
1649 /* In worst case we have alternating energy groups
1650 * and create #atom-pair lists, which means we need the size
1651 * of a cluster pair (na_ci*na_cj) times the number of cj's.
1653 nri_max = nbl->na_ci*nbl->na_cj*(cj_ind_end - cj_ind_start);
1654 if (nlist->nri + nri_max > nlist->maxnri)
1656 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1657 reallocate_nblist(nlist);
1660 ngid = nbat->nenergrp;
1662 if (static_cast<std::size_t>(ngid*gridj->na_cj) > sizeof(gid_cj)*8)
1664 gmx_fatal(FARGS, "The Verlet scheme with %dx%d kernels and free-energy only supports up to %d energy groups",
1665 gridi->na_c, gridj->na_cj, (sizeof(gid_cj)*8)/gridj->na_cj);
1668 egp_shift = nbat->neg_2log;
1669 egp_mask = (1<<nbat->neg_2log) - 1;
1671 /* Loop over the atoms in the i sub-cell */
1673 for (int i = 0; i < nbl->na_ci; i++)
1675 ind_i = ci*nbl->na_ci + i;
1680 nlist->jindex[nri+1] = nlist->jindex[nri];
1681 nlist->iinr[nri] = ai;
1682 /* The actual energy group pair index is set later */
1683 nlist->gid[nri] = 0;
1684 nlist->shift[nri] = nbl_ci->shift & NBNXN_CI_SHIFT;
1686 bFEP_i = gridi->fep[ci - gridi->cell0] & (1 << i);
1688 bFEP_i_all = bFEP_i_all && bFEP_i;
1690 if (nlist->nrj + (cj_ind_end - cj_ind_start)*nbl->na_cj > nlist->maxnrj)
1692 nlist->maxnrj = over_alloc_small(nlist->nrj + (cj_ind_end - cj_ind_start)*nbl->na_cj);
1693 srenew(nlist->jjnr, nlist->maxnrj);
1694 srenew(nlist->excl_fep, nlist->maxnrj);
1699 gid_i = (nbat->energrp[ci] >> (egp_shift*i)) & egp_mask;
1702 for (int cj_ind = cj_ind_start; cj_ind < cj_ind_end; cj_ind++)
1704 unsigned int fep_cj;
1706 cja = nbl->cj[cj_ind].cj;
1708 if (gridj->na_cj == gridj->na_c)
1710 cjr = cja - gridj->cell0;
1711 fep_cj = gridj->fep[cjr];
1714 gid_cj = nbat->energrp[cja];
1717 else if (2*gridj->na_cj == gridj->na_c)
1719 cjr = cja - gridj->cell0*2;
1720 /* Extract half of the ci fep/energrp mask */
1721 fep_cj = (gridj->fep[cjr>>1] >> ((cjr&1)*gridj->na_cj)) & ((1<<gridj->na_cj) - 1);
1724 gid_cj = nbat->energrp[cja>>1] >> ((cja&1)*gridj->na_cj*egp_shift) & ((1<<(gridj->na_cj*egp_shift)) - 1);
1729 cjr = cja - (gridj->cell0>>1);
1730 /* Combine two ci fep masks/energrp */
1731 fep_cj = gridj->fep[cjr*2] + (gridj->fep[cjr*2+1] << gridj->na_c);
1734 gid_cj = nbat->energrp[cja*2] + (nbat->energrp[cja*2+1] << (gridj->na_c*egp_shift));
1738 if (bFEP_i || fep_cj != 0)
1740 for (int j = 0; j < nbl->na_cj; j++)
1742 /* Is this interaction perturbed and not excluded? */
1743 ind_j = cja*nbl->na_cj + j;
1746 (bFEP_i || (fep_cj & (1 << j))) &&
1747 (!bDiagRemoved || ind_j >= ind_i))
1751 gid_j = (gid_cj >> (j*egp_shift)) & egp_mask;
1752 gid = GID(gid_i, gid_j, ngid);
1754 if (nlist->nrj > nlist->jindex[nri] &&
1755 nlist->gid[nri] != gid)
1757 /* Energy group pair changed: new list */
1758 fep_list_new_nri_copy(nlist);
1761 nlist->gid[nri] = gid;
1764 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1766 fep_list_new_nri_copy(nlist);
1770 /* Add it to the FEP list */
1771 nlist->jjnr[nlist->nrj] = aj;
1772 nlist->excl_fep[nlist->nrj] = (nbl->cj[cj_ind].excl >> (i*nbl->na_cj + j)) & 1;
1775 /* Exclude it from the normal list.
1776 * Note that the charge has been set to zero,
1777 * but we need to avoid 0/0, as perturbed atoms
1778 * can be on top of each other.
1780 nbl->cj[cj_ind].excl &= ~(1U << (i*nbl->na_cj + j));
1786 if (nlist->nrj > nlist->jindex[nri])
1788 /* Actually add this new, non-empty, list */
1790 nlist->jindex[nlist->nri] = nlist->nrj;
1797 /* All interactions are perturbed, we can skip this entry */
1798 nbl_ci->cj_ind_end = cj_ind_start;
1799 nbl->ncjInUse -= cj_ind_end - cj_ind_start;
1803 /* Return the index of atom a within a cluster */
1804 static gmx_inline int cj_mod_cj4(int cj)
1806 return cj & (c_nbnxnGpuJgroupSize - 1);
1809 /* Convert a j-cluster to a cj4 group */
1810 static gmx_inline int cj_to_cj4(int cj)
1812 return cj/c_nbnxnGpuJgroupSize;
1815 /* Return the index of an j-atom within a warp */
1816 static gmx_inline int a_mod_wj(int a)
1818 return a & (c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit - 1);
1821 /* As make_fep_list above, but for super/sub lists. */
1822 static void make_fep_list_supersub(const nbnxn_search_t nbs,
1823 const nbnxn_atomdata_t *nbat,
1824 nbnxn_pairlist_t *nbl,
1825 gmx_bool bDiagRemoved,
1826 const nbnxn_sci_t *nbl_sci,
1831 const nbnxn_grid_t *gridi,
1832 const nbnxn_grid_t *gridj,
1835 int sci, cj4_ind_start, cj4_ind_end, cjr;
1838 int ind_i, ind_j, ai, aj;
1842 const nbnxn_cj4_t *cj4;
1844 if (nbl_sci->cj4_ind_end == nbl_sci->cj4_ind_start)
1852 cj4_ind_start = nbl_sci->cj4_ind_start;
1853 cj4_ind_end = nbl_sci->cj4_ind_end;
1855 /* Here we process one super-cell, max #atoms na_sc, versus a list
1856 * cj4 entries, each with max c_nbnxnGpuJgroupSize cj's, each
1857 * of size na_cj atoms.
1858 * On the GPU we don't support energy groups (yet).
1859 * So for each of the na_sc i-atoms, we need max one FEP list
1860 * for each max_nrj_fep j-atoms.
1862 nri_max = nbl->na_sc*nbl->na_cj*(1 + ((cj4_ind_end - cj4_ind_start)*c_nbnxnGpuJgroupSize)/max_nrj_fep);
1863 if (nlist->nri + nri_max > nlist->maxnri)
1865 nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1866 reallocate_nblist(nlist);
1869 /* Loop over the atoms in the i super-cluster */
1870 for (int c = 0; c < c_gpuNumClusterPerCell; c++)
1872 c_abs = sci*c_gpuNumClusterPerCell + c;
1874 for (int i = 0; i < nbl->na_ci; i++)
1876 ind_i = c_abs*nbl->na_ci + i;
1881 nlist->jindex[nri+1] = nlist->jindex[nri];
1882 nlist->iinr[nri] = ai;
1883 /* With GPUs, energy groups are not supported */
1884 nlist->gid[nri] = 0;
1885 nlist->shift[nri] = nbl_sci->shift & NBNXN_CI_SHIFT;
1887 bFEP_i = (gridi->fep[c_abs - gridi->cell0*c_gpuNumClusterPerCell] & (1 << i));
1889 xi = nbat->x[ind_i*nbat->xstride+XX] + shx;
1890 yi = nbat->x[ind_i*nbat->xstride+YY] + shy;
1891 zi = nbat->x[ind_i*nbat->xstride+ZZ] + shz;
1893 if ((nlist->nrj + cj4_ind_end - cj4_ind_start)*c_nbnxnGpuJgroupSize*nbl->na_cj > nlist->maxnrj)
1895 nlist->maxnrj = over_alloc_small((nlist->nrj + cj4_ind_end - cj4_ind_start)*c_nbnxnGpuJgroupSize*nbl->na_cj);
1896 srenew(nlist->jjnr, nlist->maxnrj);
1897 srenew(nlist->excl_fep, nlist->maxnrj);
1900 for (int cj4_ind = cj4_ind_start; cj4_ind < cj4_ind_end; cj4_ind++)
1902 cj4 = &nbl->cj4[cj4_ind];
1904 for (int gcj = 0; gcj < c_nbnxnGpuJgroupSize; gcj++)
1906 unsigned int fep_cj;
1908 if ((cj4->imei[0].imask & (1U << (gcj*c_gpuNumClusterPerCell + c))) == 0)
1910 /* Skip this ci for this cj */
1914 cjr = cj4->cj[gcj] - gridj->cell0*c_gpuNumClusterPerCell;
1916 fep_cj = gridj->fep[cjr];
1918 if (bFEP_i || fep_cj != 0)
1920 for (int j = 0; j < nbl->na_cj; j++)
1922 /* Is this interaction perturbed and not excluded? */
1923 ind_j = (gridj->cell0*c_gpuNumClusterPerCell + cjr)*nbl->na_cj + j;
1926 (bFEP_i || (fep_cj & (1 << j))) &&
1927 (!bDiagRemoved || ind_j >= ind_i))
1931 unsigned int excl_bit;
1934 get_nbl_exclusions_1(nbl, cj4_ind, j>>2, &excl);
1936 excl_pair = a_mod_wj(j)*nbl->na_ci + i;
1937 excl_bit = (1U << (gcj*c_gpuNumClusterPerCell + c));
1939 dx = nbat->x[ind_j*nbat->xstride+XX] - xi;
1940 dy = nbat->x[ind_j*nbat->xstride+YY] - yi;
1941 dz = nbat->x[ind_j*nbat->xstride+ZZ] - zi;
1943 /* The unpruned GPU list has more than 2/3
1944 * of the atom pairs beyond rlist. Using
1945 * this list will cause a lot of overhead
1946 * in the CPU FEP kernels, especially
1947 * relative to the fast GPU kernels.
1948 * So we prune the FEP list here.
1950 if (dx*dx + dy*dy + dz*dz < rlist_fep2)
1952 if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1954 fep_list_new_nri_copy(nlist);
1958 /* Add it to the FEP list */
1959 nlist->jjnr[nlist->nrj] = aj;
1960 nlist->excl_fep[nlist->nrj] = (excl->pair[excl_pair] & excl_bit) ? 1 : 0;
1964 /* Exclude it from the normal list.
1965 * Note that the charge and LJ parameters have
1966 * been set to zero, but we need to avoid 0/0,
1967 * as perturbed atoms can be on top of each other.
1969 excl->pair[excl_pair] &= ~excl_bit;
1973 /* Note that we could mask out this pair in imask
1974 * if all i- and/or all j-particles are perturbed.
1975 * But since the perturbed pairs on the CPU will
1976 * take an order of magnitude more time, the GPU
1977 * will finish before the CPU and there is no gain.
1983 if (nlist->nrj > nlist->jindex[nri])
1985 /* Actually add this new, non-empty, list */
1987 nlist->jindex[nlist->nri] = nlist->nrj;
1994 /* Set all atom-pair exclusions from the topology stored in excl
1995 * as masks in the pair-list for i-super-cell entry nbl_sci
1997 static void set_sci_top_excls(const nbnxn_search_t nbs,
1998 nbnxn_pairlist_t *nbl,
1999 gmx_bool diagRemoved,
2001 const nbnxn_sci_t *nbl_sci,
2002 const t_blocka *excl)
2007 int cj_ind_first, cj_ind_last;
2008 int cj_first, cj_last;
2010 int ai, aj, si, ge, se;
2011 int found, cj_ind_0, cj_ind_1, cj_ind_m;
2013 nbnxn_excl_t *nbl_excl;
2014 int inner_i, inner_e, w;
2020 if (nbl_sci->cj4_ind_end == nbl_sci->cj4_ind_start)
2028 cj_ind_first = nbl_sci->cj4_ind_start*c_nbnxnGpuJgroupSize;
2029 cj_ind_last = nbl->work->cj_ind - 1;
2031 cj_first = nbl->cj4[nbl_sci->cj4_ind_start].cj[0];
2032 cj_last = nbl_cj(nbl, cj_ind_last);
2034 /* Determine how many contiguous j-clusters we have starting
2035 * from the first i-cluster. This number can be used to directly
2036 * calculate j-cluster indices for excluded atoms.
2039 while (cj_ind_first + ndirect <= cj_ind_last &&
2040 nbl_cj(nbl, cj_ind_first+ndirect) == sci*c_gpuNumClusterPerCell + ndirect)
2045 /* Loop over the atoms in the i super-cell */
2046 for (int i = 0; i < nbl->na_sc; i++)
2048 ai = nbs->a[sci*nbl->na_sc+i];
2051 si = (i>>na_c_2log);
2053 /* Loop over the topology-based exclusions for this i-atom */
2054 for (int eind = excl->index[ai]; eind < excl->index[ai+1]; eind++)
2060 /* The self exclusion are already set, save some time */
2066 /* Without shifts we only calculate interactions j>i
2067 * for one-way pair-lists.
2069 if (diagRemoved && ge <= sci*nbl->na_sc + i)
2075 /* Could the cluster se be in our list? */
2076 if (se >= cj_first && se <= cj_last)
2078 if (se < cj_first + ndirect)
2080 /* We can calculate cj_ind directly from se */
2081 found = cj_ind_first + se - cj_first;
2085 /* Search for se using bisection */
2087 cj_ind_0 = cj_ind_first + ndirect;
2088 cj_ind_1 = cj_ind_last + 1;
2089 while (found == -1 && cj_ind_0 < cj_ind_1)
2091 cj_ind_m = (cj_ind_0 + cj_ind_1)>>1;
2093 cj_m = nbl_cj(nbl, cj_ind_m);
2101 cj_ind_1 = cj_ind_m;
2105 cj_ind_0 = cj_ind_m + 1;
2112 inner_i = i - si*na_c;
2113 inner_e = ge - se*na_c;
2115 if (nbl_imask0(nbl, found) & (1U << (cj_mod_cj4(found)*c_gpuNumClusterPerCell + si)))
2119 get_nbl_exclusions_1(nbl, cj_to_cj4(found), w, &nbl_excl);
2121 nbl_excl->pair[a_mod_wj(inner_e)*nbl->na_ci+inner_i] &=
2122 ~(1U << (cj_mod_cj4(found)*c_gpuNumClusterPerCell + si));
2131 /* Reallocate the simple ci list for at least n entries */
2132 static void nb_realloc_ci(nbnxn_pairlist_t *nbl, int n)
2134 nbl->ci_nalloc = over_alloc_small(n);
2135 nbnxn_realloc_void((void **)&nbl->ci,
2136 nbl->nci*sizeof(*nbl->ci),
2137 nbl->ci_nalloc*sizeof(*nbl->ci),
2138 nbl->alloc, nbl->free);
2141 /* Reallocate the super-cell sci list for at least n entries */
2142 static void nb_realloc_sci(nbnxn_pairlist_t *nbl, int n)
2144 nbl->sci_nalloc = over_alloc_small(n);
2145 nbnxn_realloc_void((void **)&nbl->sci,
2146 nbl->nsci*sizeof(*nbl->sci),
2147 nbl->sci_nalloc*sizeof(*nbl->sci),
2148 nbl->alloc, nbl->free);
2151 /* Make a new ci entry at index nbl->nci */
2152 static void new_ci_entry(nbnxn_pairlist_t *nbl, int ci, int shift, int flags)
2154 if (nbl->nci + 1 > nbl->ci_nalloc)
2156 nb_realloc_ci(nbl, nbl->nci+1);
2158 nbl->ci[nbl->nci].ci = ci;
2159 nbl->ci[nbl->nci].shift = shift;
2160 /* Store the interaction flags along with the shift */
2161 nbl->ci[nbl->nci].shift |= flags;
2162 nbl->ci[nbl->nci].cj_ind_start = nbl->ncj;
2163 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
2166 /* Make a new sci entry at index nbl->nsci */
2167 static void new_sci_entry(nbnxn_pairlist_t *nbl, int sci, int shift)
2169 if (nbl->nsci + 1 > nbl->sci_nalloc)
2171 nb_realloc_sci(nbl, nbl->nsci+1);
2173 nbl->sci[nbl->nsci].sci = sci;
2174 nbl->sci[nbl->nsci].shift = shift;
2175 nbl->sci[nbl->nsci].cj4_ind_start = nbl->ncj4;
2176 nbl->sci[nbl->nsci].cj4_ind_end = nbl->ncj4;
2179 /* Sort the simple j-list cj on exclusions.
2180 * Entries with exclusions will all be sorted to the beginning of the list.
2182 static void sort_cj_excl(nbnxn_cj_t *cj, int ncj,
2183 nbnxn_list_work_t *work)
2187 if (ncj > work->cj_nalloc)
2189 work->cj_nalloc = over_alloc_large(ncj);
2190 srenew(work->cj, work->cj_nalloc);
2193 /* Make a list of the j-cells involving exclusions */
2195 for (int j = 0; j < ncj; j++)
2197 if (cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
2199 work->cj[jnew++] = cj[j];
2202 /* Check if there are exclusions at all or not just the first entry */
2203 if (!((jnew == 0) ||
2204 (jnew == 1 && cj[0].excl != NBNXN_INTERACTION_MASK_ALL)))
2206 for (int j = 0; j < ncj; j++)
2208 if (cj[j].excl == NBNXN_INTERACTION_MASK_ALL)
2210 work->cj[jnew++] = cj[j];
2213 for (int j = 0; j < ncj; j++)
2215 cj[j] = work->cj[j];
2220 /* Close this simple list i entry */
2221 static void close_ci_entry_simple(nbnxn_pairlist_t *nbl)
2225 /* All content of the new ci entry have already been filled correctly,
2226 * we only need to increase the count here (for non empty lists).
2228 jlen = nbl->ci[nbl->nci].cj_ind_end - nbl->ci[nbl->nci].cj_ind_start;
2231 sort_cj_excl(nbl->cj+nbl->ci[nbl->nci].cj_ind_start, jlen, nbl->work);
2233 /* The counts below are used for non-bonded pair/flop counts
2234 * and should therefore match the available kernel setups.
2236 if (!(nbl->ci[nbl->nci].shift & NBNXN_CI_DO_COUL(0)))
2238 nbl->work->ncj_noq += jlen;
2240 else if ((nbl->ci[nbl->nci].shift & NBNXN_CI_HALF_LJ(0)) ||
2241 !(nbl->ci[nbl->nci].shift & NBNXN_CI_DO_LJ(0)))
2243 nbl->work->ncj_hlj += jlen;
2250 /* Split sci entry for load balancing on the GPU.
2251 * Splitting ensures we have enough lists to fully utilize the whole GPU.
2252 * With progBal we generate progressively smaller lists, which improves
2253 * load balancing. As we only know the current count on our own thread,
2254 * we will need to estimate the current total amount of i-entries.
2255 * As the lists get concatenated later, this estimate depends
2256 * both on nthread and our own thread index.
2258 static void split_sci_entry(nbnxn_pairlist_t *nbl,
2260 gmx_bool progBal, float nsp_tot_est,
2261 int thread, int nthread)
2264 int cj4_start, cj4_end, j4len;
2266 int nsp, nsp_sci, nsp_cj4, nsp_cj4_e, nsp_cj4_p;
2272 /* Estimate the total numbers of ci's of the nblist combined
2273 * over all threads using the target number of ci's.
2275 nsp_est = (nsp_tot_est*thread)/nthread + nbl->nci_tot;
2277 /* The first ci blocks should be larger, to avoid overhead.
2278 * The last ci blocks should be smaller, to improve load balancing.
2279 * The factor 3/2 makes the first block 3/2 times the target average
2280 * and ensures that the total number of blocks end up equal to
2281 * that of equally sized blocks of size nsp_target_av.
2283 nsp_max = static_cast<int>(nsp_target_av*(nsp_tot_est*1.5/(nsp_est + nsp_tot_est)));
2287 nsp_max = nsp_target_av;
2290 cj4_start = nbl->sci[nbl->nsci-1].cj4_ind_start;
2291 cj4_end = nbl->sci[nbl->nsci-1].cj4_ind_end;
2292 j4len = cj4_end - cj4_start;
2294 if (j4len > 1 && j4len*c_gpuNumClusterPerCell*c_nbnxnGpuJgroupSize > nsp_max)
2296 /* Remove the last ci entry and process the cj4's again */
2304 for (int cj4 = cj4_start; cj4 < cj4_end; cj4++)
2306 nsp_cj4_p = nsp_cj4;
2307 /* Count the number of cluster pairs in this cj4 group */
2309 for (int p = 0; p < c_gpuNumClusterPerCell*c_nbnxnGpuJgroupSize; p++)
2311 nsp_cj4 += (nbl->cj4[cj4].imei[0].imask >> p) & 1;
2314 /* If adding the current cj4 with nsp_cj4 pairs get us further
2315 * away from our target nsp_max, split the list before this cj4.
2317 if (nsp > 0 && nsp_max - nsp < nsp + nsp_cj4 - nsp_max)
2319 /* Split the list at cj4 */
2320 nbl->sci[sci].cj4_ind_end = cj4;
2321 /* Create a new sci entry */
2324 if (nbl->nsci+1 > nbl->sci_nalloc)
2326 nb_realloc_sci(nbl, nbl->nsci+1);
2328 nbl->sci[sci].sci = nbl->sci[nbl->nsci-1].sci;
2329 nbl->sci[sci].shift = nbl->sci[nbl->nsci-1].shift;
2330 nbl->sci[sci].cj4_ind_start = cj4;
2332 nsp_cj4_e = nsp_cj4_p;
2338 /* Put the remaining cj4's in the last sci entry */
2339 nbl->sci[sci].cj4_ind_end = cj4_end;
2341 /* Possibly balance out the last two sci's
2342 * by moving the last cj4 of the second last sci.
2344 if (nsp_sci - nsp_cj4_e >= nsp + nsp_cj4_e)
2346 nbl->sci[sci-1].cj4_ind_end--;
2347 nbl->sci[sci].cj4_ind_start--;
2354 /* Clost this super/sub list i entry */
2355 static void close_ci_entry_supersub(nbnxn_pairlist_t *nbl,
2357 gmx_bool progBal, float nsp_tot_est,
2358 int thread, int nthread)
2360 /* All content of the new ci entry have already been filled correctly,
2361 * we only need to increase the count here (for non empty lists).
2363 int j4len = nbl->sci[nbl->nsci].cj4_ind_end - nbl->sci[nbl->nsci].cj4_ind_start;
2366 /* We can only have complete blocks of 4 j-entries in a list,
2367 * so round the count up before closing.
2369 nbl->ncj4 = (nbl->work->cj_ind + c_nbnxnGpuJgroupSize - 1)/c_nbnxnGpuJgroupSize;
2370 nbl->work->cj_ind = nbl->ncj4*c_nbnxnGpuJgroupSize;
2376 /* Measure the size of the new entry and potentially split it */
2377 split_sci_entry(nbl, nsp_max_av, progBal, nsp_tot_est,
2383 /* Syncs the working array before adding another grid pair to the list */
2384 static void sync_work(nbnxn_pairlist_t *nbl)
2388 nbl->work->cj_ind = nbl->ncj4*c_nbnxnGpuJgroupSize;
2389 nbl->work->cj4_init = nbl->ncj4;
2393 /* Clears an nbnxn_pairlist_t data structure */
2394 static void clear_pairlist(nbnxn_pairlist_t *nbl)
2404 nbl->work->ncj_noq = 0;
2405 nbl->work->ncj_hlj = 0;
2408 /* Clears a group scheme pair list */
2409 static void clear_pairlist_fep(t_nblist *nl)
2413 if (nl->jindex == nullptr)
2415 snew(nl->jindex, 1);
2420 /* Sets a simple list i-cell bounding box, including PBC shift */
2421 static gmx_inline void set_icell_bb_simple(const nbnxn_bb_t *bb, int ci,
2422 real shx, real shy, real shz,
2425 bb_ci->lower[BB_X] = bb[ci].lower[BB_X] + shx;
2426 bb_ci->lower[BB_Y] = bb[ci].lower[BB_Y] + shy;
2427 bb_ci->lower[BB_Z] = bb[ci].lower[BB_Z] + shz;
2428 bb_ci->upper[BB_X] = bb[ci].upper[BB_X] + shx;
2429 bb_ci->upper[BB_Y] = bb[ci].upper[BB_Y] + shy;
2430 bb_ci->upper[BB_Z] = bb[ci].upper[BB_Z] + shz;
2434 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2435 static void set_icell_bbxxxx_supersub(const float *bb, int ci,
2436 real shx, real shy, real shz,
2439 int ia = ci*(c_gpuNumClusterPerCell >> STRIDE_PBB_2LOG)*NNBSBB_XXXX;
2440 for (int m = 0; m < (c_gpuNumClusterPerCell >> STRIDE_PBB_2LOG)*NNBSBB_XXXX; m += NNBSBB_XXXX)
2442 for (int i = 0; i < STRIDE_PBB; i++)
2444 bb_ci[m+0*STRIDE_PBB+i] = bb[ia+m+0*STRIDE_PBB+i] + shx;
2445 bb_ci[m+1*STRIDE_PBB+i] = bb[ia+m+1*STRIDE_PBB+i] + shy;
2446 bb_ci[m+2*STRIDE_PBB+i] = bb[ia+m+2*STRIDE_PBB+i] + shz;
2447 bb_ci[m+3*STRIDE_PBB+i] = bb[ia+m+3*STRIDE_PBB+i] + shx;
2448 bb_ci[m+4*STRIDE_PBB+i] = bb[ia+m+4*STRIDE_PBB+i] + shy;
2449 bb_ci[m+5*STRIDE_PBB+i] = bb[ia+m+5*STRIDE_PBB+i] + shz;
2455 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2456 static void set_icell_bb_supersub(const nbnxn_bb_t *bb, int ci,
2457 real shx, real shy, real shz,
2460 for (int i = 0; i < c_gpuNumClusterPerCell; i++)
2462 set_icell_bb_simple(bb, ci*c_gpuNumClusterPerCell+i,
2468 /* Copies PBC shifted i-cell atom coordinates x,y,z to working array */
2469 static void icell_set_x_simple(int ci,
2470 real shx, real shy, real shz,
2471 int stride, const real *x,
2472 nbnxn_list_work_t *work)
2474 int ia = ci*NBNXN_CPU_CLUSTER_I_SIZE;
2476 for (int i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE; i++)
2478 work->x_ci[i*STRIDE_XYZ+XX] = x[(ia+i)*stride+XX] + shx;
2479 work->x_ci[i*STRIDE_XYZ+YY] = x[(ia+i)*stride+YY] + shy;
2480 work->x_ci[i*STRIDE_XYZ+ZZ] = x[(ia+i)*stride+ZZ] + shz;
2484 /* Copies PBC shifted super-cell atom coordinates x,y,z to working array */
2485 static void icell_set_x_supersub(int ci,
2486 real shx, real shy, real shz,
2487 int stride, const real *x,
2488 nbnxn_list_work_t *work)
2490 #if !GMX_SIMD4_HAVE_REAL
2492 real * x_ci = work->x_ci;
2494 int ia = ci*c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize;
2495 for (int i = 0; i < c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize; i++)
2497 x_ci[i*DIM + XX] = x[(ia+i)*stride + XX] + shx;
2498 x_ci[i*DIM + YY] = x[(ia+i)*stride + YY] + shy;
2499 x_ci[i*DIM + ZZ] = x[(ia+i)*stride + ZZ] + shz;
2502 #else /* !GMX_SIMD4_HAVE_REAL */
2504 real * x_ci = work->x_ci_simd;
2506 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2508 for (int i = 0; i < c_nbnxnGpuClusterSize; i += GMX_SIMD4_WIDTH)
2510 int io = si*c_nbnxnGpuClusterSize + i;
2511 int ia = ci*c_gpuNumClusterPerCell*c_nbnxnGpuClusterSize + io;
2512 for (int j = 0; j < GMX_SIMD4_WIDTH; j++)
2514 x_ci[io*DIM + j + XX*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + XX] + shx;
2515 x_ci[io*DIM + j + YY*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + YY] + shy;
2516 x_ci[io*DIM + j + ZZ*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + ZZ] + shz;
2521 #endif /* !GMX_SIMD4_HAVE_REAL */
2524 static real minimum_subgrid_size_xy(const nbnxn_grid_t *grid)
2528 return std::min(grid->sx, grid->sy);
2532 return std::min(grid->sx/c_gpuNumClusterPerCellX,
2533 grid->sy/c_gpuNumClusterPerCellY);
2537 static real effective_buffer_1x1_vs_MxN(const nbnxn_grid_t *gridi,
2538 const nbnxn_grid_t *gridj)
2540 const real eff_1x1_buffer_fac_overest = 0.1;
2542 /* Determine an atom-pair list cut-off buffer size for atom pairs,
2543 * to be added to rlist (including buffer) used for MxN.
2544 * This is for converting an MxN list to a 1x1 list. This means we can't
2545 * use the normal buffer estimate, as we have an MxN list in which
2546 * some atom pairs beyond rlist are missing. We want to capture
2547 * the beneficial effect of buffering by extra pairs just outside rlist,
2548 * while removing the useless pairs that are further away from rlist.
2549 * (Also the buffer could have been set manually not using the estimate.)
2550 * This buffer size is an overestimate.
2551 * We add 10% of the smallest grid sub-cell dimensions.
2552 * Note that the z-size differs per cell and we don't use this,
2553 * so we overestimate.
2554 * With PME, the 10% value gives a buffer that is somewhat larger
2555 * than the effective buffer with a tolerance of 0.005 kJ/mol/ps.
2556 * Smaller tolerances or using RF lead to a smaller effective buffer,
2557 * so 10% gives a safe overestimate.
2559 return eff_1x1_buffer_fac_overest*(minimum_subgrid_size_xy(gridi) +
2560 minimum_subgrid_size_xy(gridj));
2563 /* Clusters at the cut-off only increase rlist by 60% of their size */
2564 static real nbnxn_rlist_inc_outside_fac = 0.6;
2566 /* Due to the cluster size the effective pair-list is longer than
2567 * that of a simple atom pair-list. This function gives the extra distance.
2569 real nbnxn_get_rlist_effective_inc(int cluster_size_j, real atom_density)
2572 real vol_inc_i, vol_inc_j;
2574 /* We should get this from the setup, but currently it's the same for
2575 * all setups, including GPUs.
2577 cluster_size_i = NBNXN_CPU_CLUSTER_I_SIZE;
2579 vol_inc_i = (cluster_size_i - 1)/atom_density;
2580 vol_inc_j = (cluster_size_j - 1)/atom_density;
2582 return nbnxn_rlist_inc_outside_fac*std::cbrt(vol_inc_i + vol_inc_j);
2585 /* Estimates the interaction volume^2 for non-local interactions */
2586 static real nonlocal_vol2(const struct gmx_domdec_zones_t *zones, rvec ls, real r)
2594 /* Here we simply add up the volumes of 1, 2 or 3 1D decomposition
2595 * not home interaction volume^2. As these volumes are not additive,
2596 * this is an overestimate, but it would only be significant in the limit
2597 * of small cells, where we anyhow need to split the lists into
2598 * as small parts as possible.
2601 for (int z = 0; z < zones->n; z++)
2603 if (zones->shift[z][XX] + zones->shift[z][YY] + zones->shift[z][ZZ] == 1)
2608 for (int d = 0; d < DIM; d++)
2610 if (zones->shift[z][d] == 0)
2614 za *= zones->size[z].x1[d] - zones->size[z].x0[d];
2618 /* 4 octants of a sphere */
2619 vold_est = 0.25*M_PI*r*r*r*r;
2620 /* 4 quarter pie slices on the edges */
2621 vold_est += 4*cl*M_PI/6.0*r*r*r;
2622 /* One rectangular volume on a face */
2623 vold_est += ca*0.5*r*r;
2625 vol2_est_tot += vold_est*za;
2629 return vol2_est_tot;
2632 /* Estimates the average size of a full j-list for super/sub setup */
2633 static void get_nsubpair_target(const nbnxn_search_t nbs,
2636 int min_ci_balanced,
2637 int *nsubpair_target,
2638 float *nsubpair_tot_est)
2640 /* The target value of 36 seems to be the optimum for Kepler.
2641 * Maxwell is less sensitive to the exact value.
2643 const int nsubpair_target_min = 36;
2644 const nbnxn_grid_t *grid;
2646 real r_eff_sup, vol_est, nsp_est, nsp_est_nl;
2648 grid = &nbs->grid[0];
2650 /* We don't need to balance list sizes if:
2651 * - We didn't request balancing.
2652 * - The number of grid cells >= the number of lists requested,
2653 * since we will always generate at least #cells lists.
2654 * - We don't have any cells, since then there won't be any lists.
2656 if (min_ci_balanced <= 0 || grid->nc >= min_ci_balanced || grid->nc == 0)
2658 /* nsubpair_target==0 signals no balancing */
2659 *nsubpair_target = 0;
2660 *nsubpair_tot_est = 0;
2665 ls[XX] = (grid->c1[XX] - grid->c0[XX])/(grid->ncx*c_gpuNumClusterPerCellX);
2666 ls[YY] = (grid->c1[YY] - grid->c0[YY])/(grid->ncy*c_gpuNumClusterPerCellY);
2667 ls[ZZ] = grid->na_c/(grid->atom_density*ls[XX]*ls[YY]);
2669 /* The average length of the diagonal of a sub cell */
2670 real diagonal = std::sqrt(ls[XX]*ls[XX] + ls[YY]*ls[YY] + ls[ZZ]*ls[ZZ]);
2672 /* The formulas below are a heuristic estimate of the average nsj per si*/
2673 r_eff_sup = rlist + nbnxn_rlist_inc_outside_fac*gmx::square((grid->na_c - 1.0)/grid->na_c)*0.5*diagonal;
2675 if (!nbs->DomDec || nbs->zones->n == 1)
2682 gmx::square(grid->atom_density/grid->na_c)*
2683 nonlocal_vol2(nbs->zones, ls, r_eff_sup);
2688 /* Sub-cell interacts with itself */
2689 vol_est = ls[XX]*ls[YY]*ls[ZZ];
2690 /* 6/2 rectangular volume on the faces */
2691 vol_est += (ls[XX]*ls[YY] + ls[XX]*ls[ZZ] + ls[YY]*ls[ZZ])*r_eff_sup;
2692 /* 12/2 quarter pie slices on the edges */
2693 vol_est += 2*(ls[XX] + ls[YY] + ls[ZZ])*0.25*M_PI*gmx::square(r_eff_sup);
2694 /* 4 octants of a sphere */
2695 vol_est += 0.5*4.0/3.0*M_PI*gmx::power3(r_eff_sup);
2697 /* Estimate the number of cluster pairs as the local number of
2698 * clusters times the volume they interact with times the density.
2700 nsp_est = grid->nsubc_tot*vol_est*grid->atom_density/grid->na_c;
2702 /* Subtract the non-local pair count */
2703 nsp_est -= nsp_est_nl;
2705 /* For small cut-offs nsp_est will be an underesimate.
2706 * With DD nsp_est_nl is an overestimate so nsp_est can get negative.
2707 * So to avoid too small or negative nsp_est we set a minimum of
2708 * all cells interacting with all 3^3 direct neighbors (3^3-1)/2+1=14.
2709 * This might be a slight overestimate for small non-periodic groups of
2710 * atoms as will occur for a local domain with DD, but for small
2711 * groups of atoms we'll anyhow be limited by nsubpair_target_min,
2712 * so this overestimation will not matter.
2714 nsp_est = std::max(nsp_est, grid->nsubc_tot*static_cast<real>(14));
2718 fprintf(debug, "nsp_est local %5.1f non-local %5.1f\n",
2719 nsp_est, nsp_est_nl);
2724 nsp_est = nsp_est_nl;
2727 /* Thus the (average) maximum j-list size should be as follows.
2728 * Since there is overhead, we shouldn't make the lists too small
2729 * (and we can't chop up j-groups) so we use a minimum target size of 36.
2731 *nsubpair_target = std::max(nsubpair_target_min,
2732 static_cast<int>(nsp_est/min_ci_balanced + 0.5));
2733 *nsubpair_tot_est = static_cast<int>(nsp_est);
2737 fprintf(debug, "nbl nsp estimate %.1f, nsubpair_target %d\n",
2738 nsp_est, *nsubpair_target);
2742 /* Debug list print function */
2743 static void print_nblist_ci_cj(FILE *fp, const nbnxn_pairlist_t *nbl)
2745 for (int i = 0; i < nbl->nci; i++)
2747 fprintf(fp, "ci %4d shift %2d ncj %3d\n",
2748 nbl->ci[i].ci, nbl->ci[i].shift,
2749 nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start);
2751 for (int j = nbl->ci[i].cj_ind_start; j < nbl->ci[i].cj_ind_end; j++)
2753 fprintf(fp, " cj %5d imask %x\n",
2760 /* Debug list print function */
2761 static void print_nblist_sci_cj(FILE *fp, const nbnxn_pairlist_t *nbl)
2763 for (int i = 0; i < nbl->nsci; i++)
2765 fprintf(fp, "ci %4d shift %2d ncj4 %2d\n",
2766 nbl->sci[i].sci, nbl->sci[i].shift,
2767 nbl->sci[i].cj4_ind_end - nbl->sci[i].cj4_ind_start);
2770 for (int j4 = nbl->sci[i].cj4_ind_start; j4 < nbl->sci[i].cj4_ind_end; j4++)
2772 for (int j = 0; j < c_nbnxnGpuJgroupSize; j++)
2774 fprintf(fp, " sj %5d imask %x\n",
2776 nbl->cj4[j4].imei[0].imask);
2777 for (int si = 0; si < c_gpuNumClusterPerCell; si++)
2779 if (nbl->cj4[j4].imei[0].imask & (1U << (j*c_gpuNumClusterPerCell + si)))
2786 fprintf(fp, "ci %4d shift %2d ncj4 %2d ncp %3d\n",
2787 nbl->sci[i].sci, nbl->sci[i].shift,
2788 nbl->sci[i].cj4_ind_end - nbl->sci[i].cj4_ind_start,
2793 /* Combine pair lists *nbl generated on multiple threads nblc */
2794 static void combine_nblists(int nnbl, nbnxn_pairlist_t **nbl,
2795 nbnxn_pairlist_t *nblc)
2797 int nsci, ncj4, nexcl;
2801 gmx_incons("combine_nblists does not support simple lists");
2806 nexcl = nblc->nexcl;
2807 for (int i = 0; i < nnbl; i++)
2809 nsci += nbl[i]->nsci;
2810 ncj4 += nbl[i]->ncj4;
2811 nexcl += nbl[i]->nexcl;
2814 if (nsci > nblc->sci_nalloc)
2816 nb_realloc_sci(nblc, nsci);
2818 if (ncj4 > nblc->cj4_nalloc)
2820 nblc->cj4_nalloc = over_alloc_small(ncj4);
2821 nbnxn_realloc_void((void **)&nblc->cj4,
2822 nblc->ncj4*sizeof(*nblc->cj4),
2823 nblc->cj4_nalloc*sizeof(*nblc->cj4),
2824 nblc->alloc, nblc->free);
2826 if (nexcl > nblc->excl_nalloc)
2828 nblc->excl_nalloc = over_alloc_small(nexcl);
2829 nbnxn_realloc_void((void **)&nblc->excl,
2830 nblc->nexcl*sizeof(*nblc->excl),
2831 nblc->excl_nalloc*sizeof(*nblc->excl),
2832 nblc->alloc, nblc->free);
2835 /* Each thread should copy its own data to the combined arrays,
2836 * as otherwise data will go back and forth between different caches.
2838 #if GMX_OPENMP && !(defined __clang_analyzer__)
2839 // cppcheck-suppress unreadVariable
2840 int nthreads = gmx_omp_nthreads_get(emntPairsearch);
2843 #pragma omp parallel for num_threads(nthreads) schedule(static)
2844 for (int n = 0; n < nnbl; n++)
2851 const nbnxn_pairlist_t *nbli;
2853 /* Determine the offset in the combined data for our thread */
2854 sci_offset = nblc->nsci;
2855 cj4_offset = nblc->ncj4;
2856 excl_offset = nblc->nexcl;
2858 for (int i = 0; i < n; i++)
2860 sci_offset += nbl[i]->nsci;
2861 cj4_offset += nbl[i]->ncj4;
2862 excl_offset += nbl[i]->nexcl;
2867 for (int i = 0; i < nbli->nsci; i++)
2869 nblc->sci[sci_offset+i] = nbli->sci[i];
2870 nblc->sci[sci_offset+i].cj4_ind_start += cj4_offset;
2871 nblc->sci[sci_offset+i].cj4_ind_end += cj4_offset;
2874 for (int j4 = 0; j4 < nbli->ncj4; j4++)
2876 nblc->cj4[cj4_offset+j4] = nbli->cj4[j4];
2877 nblc->cj4[cj4_offset+j4].imei[0].excl_ind += excl_offset;
2878 nblc->cj4[cj4_offset+j4].imei[1].excl_ind += excl_offset;
2881 for (int j4 = 0; j4 < nbli->nexcl; j4++)
2883 nblc->excl[excl_offset+j4] = nbli->excl[j4];
2886 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2889 for (int n = 0; n < nnbl; n++)
2891 nblc->nsci += nbl[n]->nsci;
2892 nblc->ncj4 += nbl[n]->ncj4;
2893 nblc->nci_tot += nbl[n]->nci_tot;
2894 nblc->nexcl += nbl[n]->nexcl;
2898 static void balance_fep_lists(const nbnxn_search_t nbs,
2899 nbnxn_pairlist_set_t *nbl_lists)
2902 int nri_tot, nrj_tot, nrj_target;
2906 nnbl = nbl_lists->nnbl;
2910 /* Nothing to balance */
2914 /* Count the total i-lists and pairs */
2917 for (int th = 0; th < nnbl; th++)
2919 nri_tot += nbl_lists->nbl_fep[th]->nri;
2920 nrj_tot += nbl_lists->nbl_fep[th]->nrj;
2923 nrj_target = (nrj_tot + nnbl - 1)/nnbl;
2925 assert(gmx_omp_nthreads_get(emntNonbonded) == nnbl);
2927 #pragma omp parallel for schedule(static) num_threads(nnbl)
2928 for (int th = 0; th < nnbl; th++)
2934 nbl = nbs->work[th].nbl_fep;
2936 /* Note that here we allocate for the total size, instead of
2937 * a per-thread esimate (which is hard to obtain).
2939 if (nri_tot > nbl->maxnri)
2941 nbl->maxnri = over_alloc_large(nri_tot);
2942 reallocate_nblist(nbl);
2944 if (nri_tot > nbl->maxnri || nrj_tot > nbl->maxnrj)
2946 nbl->maxnrj = over_alloc_small(nrj_tot);
2947 srenew(nbl->jjnr, nbl->maxnrj);
2948 srenew(nbl->excl_fep, nbl->maxnrj);
2951 clear_pairlist_fep(nbl);
2953 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2956 /* Loop over the source lists and assign and copy i-entries */
2958 nbld = nbs->work[th_dest].nbl_fep;
2959 for (int th = 0; th < nnbl; th++)
2963 nbls = nbl_lists->nbl_fep[th];
2965 for (int i = 0; i < nbls->nri; i++)
2969 /* The number of pairs in this i-entry */
2970 nrj = nbls->jindex[i+1] - nbls->jindex[i];
2972 /* Decide if list th_dest is too large and we should procede
2973 * to the next destination list.
2975 if (th_dest+1 < nnbl && nbld->nrj > 0 &&
2976 nbld->nrj + nrj - nrj_target > nrj_target - nbld->nrj)
2979 nbld = nbs->work[th_dest].nbl_fep;
2982 nbld->iinr[nbld->nri] = nbls->iinr[i];
2983 nbld->gid[nbld->nri] = nbls->gid[i];
2984 nbld->shift[nbld->nri] = nbls->shift[i];
2986 for (int j = nbls->jindex[i]; j < nbls->jindex[i+1]; j++)
2988 nbld->jjnr[nbld->nrj] = nbls->jjnr[j];
2989 nbld->excl_fep[nbld->nrj] = nbls->excl_fep[j];
2993 nbld->jindex[nbld->nri] = nbld->nrj;
2997 /* Swap the list pointers */
2998 for (int th = 0; th < nnbl; th++)
3002 nbl_tmp = nbl_lists->nbl_fep[th];
3003 nbl_lists->nbl_fep[th] = nbs->work[th].nbl_fep;
3004 nbs->work[th].nbl_fep = nbl_tmp;
3008 fprintf(debug, "nbl_fep[%d] nri %4d nrj %4d\n",
3010 nbl_lists->nbl_fep[th]->nri,
3011 nbl_lists->nbl_fep[th]->nrj);
3016 /* Returns the next ci to be processes by our thread */
3017 static gmx_bool next_ci(const nbnxn_grid_t *grid,
3019 int nth, int ci_block,
3020 int *ci_x, int *ci_y,
3026 if (*ci_b == ci_block)
3028 /* Jump to the next block assigned to this task */
3029 *ci += (nth - 1)*ci_block;
3033 if (*ci >= grid->nc*conv)
3038 while (*ci >= grid->cxy_ind[*ci_x*grid->ncy + *ci_y + 1]*conv)
3041 if (*ci_y == grid->ncy)
3051 /* Returns the distance^2 for which we put cell pairs in the list
3052 * without checking atom pair distances. This is usually < rlist^2.
3054 static float boundingbox_only_distance2(const nbnxn_grid_t *gridi,
3055 const nbnxn_grid_t *gridj,
3059 /* If the distance between two sub-cell bounding boxes is less
3060 * than this distance, do not check the distance between
3061 * all particle pairs in the sub-cell, since then it is likely
3062 * that the box pair has atom pairs within the cut-off.
3063 * We use the nblist cut-off minus 0.5 times the average x/y diagonal
3064 * spacing of the sub-cells. Around 40% of the checked pairs are pruned.
3065 * Using more than 0.5 gains at most 0.5%.
3066 * If forces are calculated more than twice, the performance gain
3067 * in the force calculation outweighs the cost of checking.
3068 * Note that with subcell lists, the atom-pair distance check
3069 * is only performed when only 1 out of 8 sub-cells in within range,
3070 * this is because the GPU is much faster than the cpu.
3075 bbx = 0.5*(gridi->sx + gridj->sx);
3076 bby = 0.5*(gridi->sy + gridj->sy);
3079 bbx /= c_gpuNumClusterPerCellX;
3080 bby /= c_gpuNumClusterPerCellY;
3083 rbb2 = std::max(0.0, rlist - 0.5*std::sqrt(bbx*bbx + bby*bby));
3089 return (float)((1+GMX_FLOAT_EPS)*rbb2);
3093 static int get_ci_block_size(const nbnxn_grid_t *gridi,
3094 gmx_bool bDomDec, int nth)
3096 const int ci_block_enum = 5;
3097 const int ci_block_denom = 11;
3098 const int ci_block_min_atoms = 16;
3101 /* Here we decide how to distribute the blocks over the threads.
3102 * We use prime numbers to try to avoid that the grid size becomes
3103 * a multiple of the number of threads, which would lead to some
3104 * threads getting "inner" pairs and others getting boundary pairs,
3105 * which in turns will lead to load imbalance between threads.
3106 * Set the block size as 5/11/ntask times the average number of cells
3107 * in a y,z slab. This should ensure a quite uniform distribution
3108 * of the grid parts of the different thread along all three grid
3109 * zone boundaries with 3D domain decomposition. At the same time
3110 * the blocks will not become too small.
3112 ci_block = (gridi->nc*ci_block_enum)/(ci_block_denom*gridi->ncx*nth);
3114 /* Ensure the blocks are not too small: avoids cache invalidation */
3115 if (ci_block*gridi->na_sc < ci_block_min_atoms)
3117 ci_block = (ci_block_min_atoms + gridi->na_sc - 1)/gridi->na_sc;
3120 /* Without domain decomposition
3121 * or with less than 3 blocks per task, divide in nth blocks.
3123 if (!bDomDec || nth*3*ci_block > gridi->nc)
3125 ci_block = (gridi->nc + nth - 1)/nth;
3128 if (ci_block > 1 && (nth - 1)*ci_block >= gridi->nc)
3130 /* Some threads have no work. Although reducing the block size
3131 * does not decrease the block count on the first few threads,
3132 * with GPUs better mixing of "upper" cells that have more empty
3133 * clusters results in a somewhat lower max load over all threads.
3134 * Without GPUs the regime of so few atoms per thread is less
3135 * performance relevant, but with 8-wide SIMD the same reasoning
3136 * applies, since the pair list uses 4 i-atom "sub-clusters".
3144 /* Returns the number of bits to right-shift a cluster index to obtain
3145 * the corresponding force buffer flag index.
3147 static int getBufferFlagShift(int numAtomsPerCluster)
3149 int bufferFlagShift = 0;
3150 while ((numAtomsPerCluster << bufferFlagShift) < NBNXN_BUFFERFLAG_SIZE)
3155 return bufferFlagShift;
3158 /* Generates the part of pair-list nbl assigned to our thread */
3159 static void nbnxn_make_pairlist_part(const nbnxn_search_t nbs,
3160 const nbnxn_grid_t *gridi,
3161 const nbnxn_grid_t *gridj,
3162 nbnxn_search_work_t *work,
3163 const nbnxn_atomdata_t *nbat,
3164 const t_blocka *excl,
3168 gmx_bool bFBufferFlag,
3171 float nsubpair_tot_est,
3173 nbnxn_pairlist_t *nbl,
3178 real rlist2, rl_fep2 = 0;
3180 int ci_b, ci, ci_x, ci_y, ci_xy, cj;
3184 int conv_i, cell0_i;
3185 const nbnxn_bb_t *bb_i = nullptr;
3187 const float *pbb_i = nullptr;
3189 const float *bbcz_i, *bbcz_j;
3191 real bx0, bx1, by0, by1, bz0, bz1;
3193 real d2cx, d2z, d2z_cx, d2z_cy, d2zx, d2zxy, d2xy;
3194 int cxf, cxl, cyf, cyf_x, cyl;
3195 int numDistanceChecks;
3196 int gridi_flag_shift = 0, gridj_flag_shift = 0;
3197 gmx_bitmask_t *gridj_flag = nullptr;
3198 int ncj_old_i, ncj_old_j;
3200 nbs_cycle_start(&work->cc[enbsCCsearch]);
3202 if (gridj->bSimple != nbl->bSimple)
3204 gmx_incons("Grid incompatible with pair-list");
3208 nbl->na_sc = gridj->na_sc;
3209 nbl->na_ci = gridj->na_c;
3210 nbl->na_cj = nbnxn_kernel_to_cluster_j_size(nb_kernel_type);
3211 na_cj_2log = get_2log(nbl->na_cj);
3217 /* Determine conversion of clusters to flag blocks */
3218 gridi_flag_shift = getBufferFlagShift(nbl->na_ci);
3219 gridj_flag_shift = getBufferFlagShift(nbl->na_cj);
3221 gridj_flag = work->buffer_flags.flag;
3224 copy_mat(nbs->box, box);
3226 rlist2 = nbl->rlist*nbl->rlist;
3228 if (nbs->bFEP && !nbl->bSimple)
3230 /* Determine an atom-pair list cut-off distance for FEP atom pairs.
3231 * We should not simply use rlist, since then we would not have
3232 * the small, effective buffering of the NxN lists.
3233 * The buffer is on overestimate, but the resulting cost for pairs
3234 * beyond rlist is neglible compared to the FEP pairs within rlist.
3236 rl_fep2 = nbl->rlist + effective_buffer_1x1_vs_MxN(gridi, gridj);
3240 fprintf(debug, "nbl_fep atom-pair rlist %f\n", rl_fep2);
3242 rl_fep2 = rl_fep2*rl_fep2;
3245 rbb2 = boundingbox_only_distance2(gridi, gridj, nbl->rlist, nbl->bSimple);
3249 fprintf(debug, "nbl bounding box only distance %f\n", std::sqrt(rbb2));
3252 /* Set the shift range */
3253 for (int d = 0; d < DIM; d++)
3255 /* Check if we need periodicity shifts.
3256 * Without PBC or with domain decomposition we don't need them.
3258 if (d >= ePBC2npbcdim(nbs->ePBC) || nbs->dd_dim[d])
3265 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < std::sqrt(rlist2))
3276 if (nbl->bSimple && !gridi->bSimple)
3278 conv_i = gridi->na_sc/gridj->na_sc;
3279 bb_i = gridi->bb_simple;
3280 bbcz_i = gridi->bbcz_simple;
3281 flags_i = gridi->flags_simple;
3296 /* We use the normal bounding box format for both grid types */
3299 bbcz_i = gridi->bbcz;
3300 flags_i = gridi->flags;
3302 cell0_i = gridi->cell0*conv_i;
3304 bbcz_j = gridj->bbcz;
3308 /* Blocks of the conversion factor - 1 give a large repeat count
3309 * combined with a small block size. This should result in good
3310 * load balancing for both small and large domains.
3312 ci_block = conv_i - 1;
3316 fprintf(debug, "nbl nc_i %d col.av. %.1f ci_block %d\n",
3317 gridi->nc, gridi->nc/(double)(gridi->ncx*gridi->ncy), ci_block);
3320 numDistanceChecks = 0;
3322 /* Initially ci_b and ci to 1 before where we want them to start,
3323 * as they will both be incremented in next_ci.
3326 ci = th*ci_block - 1;
3329 while (next_ci(gridi, conv_i, nth, ci_block, &ci_x, &ci_y, &ci_b, &ci))
3331 if (nbl->bSimple && flags_i[ci] == 0)
3336 ncj_old_i = nbl->ncj;
3339 if (gridj != gridi && shp[XX] == 0)
3343 bx1 = bb_i[ci].upper[BB_X];
3347 bx1 = gridi->c0[XX] + (ci_x+1)*gridi->sx;
3349 if (bx1 < gridj->c0[XX])
3351 d2cx = gmx::square(gridj->c0[XX] - bx1);
3360 ci_xy = ci_x*gridi->ncy + ci_y;
3362 /* Loop over shift vectors in three dimensions */
3363 for (int tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
3365 shz = tz*box[ZZ][ZZ];
3367 bz0 = bbcz_i[ci*NNBSBB_D ] + shz;
3368 bz1 = bbcz_i[ci*NNBSBB_D+1] + shz;
3376 d2z = gmx::square(bz1);
3380 d2z = gmx::square(bz0 - box[ZZ][ZZ]);
3383 d2z_cx = d2z + d2cx;
3385 if (d2z_cx >= rlist2)
3390 bz1_frac = bz1/(gridi->cxy_ind[ci_xy+1] - gridi->cxy_ind[ci_xy]);
3395 /* The check with bz1_frac close to or larger than 1 comes later */
3397 for (int ty = -shp[YY]; ty <= shp[YY]; ty++)
3399 shy = ty*box[YY][YY] + tz*box[ZZ][YY];
3403 by0 = bb_i[ci].lower[BB_Y] + shy;
3404 by1 = bb_i[ci].upper[BB_Y] + shy;
3408 by0 = gridi->c0[YY] + (ci_y )*gridi->sy + shy;
3409 by1 = gridi->c0[YY] + (ci_y+1)*gridi->sy + shy;
3412 get_cell_range(by0, by1,
3413 gridj->ncy, gridj->c0[YY], gridj->sy, gridj->inv_sy,
3423 if (by1 < gridj->c0[YY])
3425 d2z_cy += gmx::square(gridj->c0[YY] - by1);
3427 else if (by0 > gridj->c1[YY])
3429 d2z_cy += gmx::square(by0 - gridj->c1[YY]);
3432 for (int tx = -shp[XX]; tx <= shp[XX]; tx++)
3434 shift = XYZ2IS(tx, ty, tz);
3436 if (c_pbcShiftBackward && gridi == gridj && shift > CENTRAL)
3441 shx = tx*box[XX][XX] + ty*box[YY][XX] + tz*box[ZZ][XX];
3445 bx0 = bb_i[ci].lower[BB_X] + shx;
3446 bx1 = bb_i[ci].upper[BB_X] + shx;
3450 bx0 = gridi->c0[XX] + (ci_x )*gridi->sx + shx;
3451 bx1 = gridi->c0[XX] + (ci_x+1)*gridi->sx + shx;
3454 get_cell_range(bx0, bx1,
3455 gridj->ncx, gridj->c0[XX], gridj->sx, gridj->inv_sx,
3466 new_ci_entry(nbl, cell0_i+ci, shift, flags_i[ci]);
3470 new_sci_entry(nbl, cell0_i+ci, shift);
3473 if ((!c_pbcShiftBackward || (shift == CENTRAL &&
3477 /* Leave the pairs with i > j.
3478 * x is the major index, so skip half of it.
3485 set_icell_bb_simple(bb_i, ci, shx, shy, shz,
3491 set_icell_bbxxxx_supersub(pbb_i, ci, shx, shy, shz,
3494 set_icell_bb_supersub(bb_i, ci, shx, shy, shz,
3499 nbs->icell_set_x(cell0_i+ci, shx, shy, shz,
3500 nbat->xstride, nbat->x,
3503 for (int cx = cxf; cx <= cxl; cx++)
3506 if (gridj->c0[XX] + cx*gridj->sx > bx1)
3508 d2zx += gmx::square(gridj->c0[XX] + cx*gridj->sx - bx1);
3510 else if (gridj->c0[XX] + (cx+1)*gridj->sx < bx0)
3512 d2zx += gmx::square(gridj->c0[XX] + (cx+1)*gridj->sx - bx0);
3515 if (gridi == gridj &&
3517 (!c_pbcShiftBackward || shift == CENTRAL) &&
3520 /* Leave the pairs with i > j.
3521 * Skip half of y when i and j have the same x.
3530 for (int cy = cyf_x; cy <= cyl; cy++)
3532 const int columnStart = gridj->cxy_ind[cx*gridj->ncy + cy];
3533 const int columnEnd = gridj->cxy_ind[cx*gridj->ncy + cy + 1];
3536 if (gridj->c0[YY] + cy*gridj->sy > by1)
3538 d2zxy += gmx::square(gridj->c0[YY] + cy*gridj->sy - by1);
3540 else if (gridj->c0[YY] + (cy+1)*gridj->sy < by0)
3542 d2zxy += gmx::square(gridj->c0[YY] + (cy+1)*gridj->sy - by0);
3544 if (columnStart < columnEnd && d2zxy < rlist2)
3546 /* To improve efficiency in the common case
3547 * of a homogeneous particle distribution,
3548 * we estimate the index of the middle cell
3549 * in range (midCell). We search down and up
3550 * starting from this index.
3552 * Note that the bbcz_j array contains bounds
3553 * for i-clusters, thus for clusters of 4 atoms.
3554 * For the common case where the j-cluster size
3555 * is 8, we could step with a stride of 2,
3556 * but we do not do this because it would
3557 * complicate this code even more.
3559 int midCell = columnStart + static_cast<int>(bz1_frac*(columnEnd - columnStart));
3560 if (midCell >= columnEnd)
3562 midCell = columnEnd - 1;
3567 /* Find the lowest cell that can possibly
3569 * Check if we hit the bottom of the grid,
3570 * if the j-cell is below the i-cell and if so,
3571 * if it is within range.
3573 int firstCell = midCell;
3574 while (firstCell > columnStart &&
3575 (bbcz_j[firstCell*NNBSBB_D + 1] >= bz0 ||
3576 d2xy + gmx::square(bbcz_j[firstCell*NNBSBB_D + 1] - bz0) < rlist2))
3581 /* Find the highest cell that can possibly
3583 * Check if we hit the top of the grid,
3584 * if the j-cell is above the i-cell and if so,
3585 * if it is within range.
3587 int lastCell = midCell;
3588 while (lastCell < columnEnd - 1 &&
3589 (bbcz_j[lastCell*NNBSBB_D] <= bz1 ||
3590 d2xy + gmx::square(bbcz_j[lastCell*NNBSBB_D] - bz1) < rlist2))
3595 #define NBNXN_REFCODE 0
3598 /* Simple reference code, for debugging,
3599 * overrides the more complex code above.
3601 firstCell = columnEnd;
3603 for (int k = columnStart; k < columnEnd; k++)
3605 if (d2xy + gmx::square(bbcz_j[k*NNBSBB_D + 1] - bz0) < rlist2 &&
3610 if (d2xy + gmx::square(bbcz_j[k*NNBSBB_D] - bz1) < rlist2 &&
3621 /* We want each atom/cell pair only once,
3622 * only use cj >= ci.
3624 if (!c_pbcShiftBackward || shift == CENTRAL)
3626 firstCell = std::max(firstCell, ci);
3630 if (firstCell <= lastCell)
3632 GMX_ASSERT(firstCell >= columnStart && lastCell < columnEnd, "The range should reside within the current grid column");
3634 /* For f buffer flags with simple lists */
3635 ncj_old_j = nbl->ncj;
3639 /* We have a maximum of 2 j-clusters
3640 * per i-cluster sized cell.
3642 check_cell_list_space_simple(nbl, 2*(lastCell - firstCell + 1));
3646 check_cell_list_space_supersub(nbl, lastCell - firstCell + 1);
3649 switch (nb_kernel_type)
3651 case nbnxnk4x4_PlainC:
3652 makeClusterListSimple(gridj,
3653 nbl, ci, firstCell, lastCell,
3654 (gridi == gridj && shift == CENTRAL),
3657 &numDistanceChecks);
3659 #ifdef GMX_NBNXN_SIMD_4XN
3660 case nbnxnk4xN_SIMD_4xN:
3661 makeClusterListSimd4xn(gridj,
3662 nbl, ci, firstCell, lastCell,
3663 (gridi == gridj && shift == CENTRAL),
3666 &numDistanceChecks);
3669 #ifdef GMX_NBNXN_SIMD_2XNN
3670 case nbnxnk4xN_SIMD_2xNN:
3671 makeClusterListSimd2xnn(gridj,
3672 nbl, ci, firstCell, lastCell,
3673 (gridi == gridj && shift == CENTRAL),
3676 &numDistanceChecks);
3679 case nbnxnk8x8x8_PlainC:
3680 case nbnxnk8x8x8_GPU:
3681 for (cj = firstCell; cj <= lastCell; cj++)
3683 make_cluster_list_supersub(gridi, gridj,
3685 (gridi == gridj && shift == CENTRAL && ci == cj),
3686 nbat->xstride, nbat->x,
3688 &numDistanceChecks);
3693 if (bFBufferFlag && nbl->ncj > ncj_old_j)
3695 int cbf = nbl->cj[ncj_old_j].cj >> gridj_flag_shift;
3696 int cbl = nbl->cj[nbl->ncj-1].cj >> gridj_flag_shift;
3697 for (int cb = cbf; cb <= cbl; cb++)
3699 bitmask_init_bit(&gridj_flag[cb], th);
3703 nbl->ncjInUse += nbl->ncj - ncj_old_j;
3709 /* Set the exclusions for this ci list */
3712 set_ci_top_excls(nbs,
3714 shift == CENTRAL && gridi == gridj,
3717 &(nbl->ci[nbl->nci]),
3722 make_fep_list(nbs, nbat, nbl,
3723 shift == CENTRAL && gridi == gridj,
3724 &(nbl->ci[nbl->nci]),
3725 gridi, gridj, nbl_fep);
3730 set_sci_top_excls(nbs,
3732 shift == CENTRAL && gridi == gridj,
3734 &(nbl->sci[nbl->nsci]),
3739 make_fep_list_supersub(nbs, nbat, nbl,
3740 shift == CENTRAL && gridi == gridj,
3741 &(nbl->sci[nbl->nsci]),
3744 gridi, gridj, nbl_fep);
3748 /* Close this ci list */
3751 close_ci_entry_simple(nbl);
3755 close_ci_entry_supersub(nbl,
3757 progBal, nsubpair_tot_est,
3764 if (bFBufferFlag && nbl->ncj > ncj_old_i)
3766 bitmask_init_bit(&(work->buffer_flags.flag[(gridi->cell0+ci)>>gridi_flag_shift]), th);
3770 work->ndistc = numDistanceChecks;
3772 nbs_cycle_stop(&work->cc[enbsCCsearch]);
3774 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");
3778 fprintf(debug, "number of distance checks %d\n", numDistanceChecks);
3782 print_nblist_statistics_simple(debug, nbl, nbs, rlist);
3786 print_nblist_statistics_supersub(debug, nbl, nbs, rlist);
3791 fprintf(debug, "nbl FEP list pairs: %d\n", nbl_fep->nrj);
3796 static void reduce_buffer_flags(const nbnxn_search_t nbs,
3798 const nbnxn_buffer_flags_t *dest)
3800 for (int s = 0; s < nsrc; s++)
3802 gmx_bitmask_t * flag = nbs->work[s].buffer_flags.flag;
3804 for (int b = 0; b < dest->nflag; b++)
3806 bitmask_union(&(dest->flag[b]), flag[b]);
3811 static void print_reduction_cost(const nbnxn_buffer_flags_t *flags, int nout)
3813 int nelem, nkeep, ncopy, nred, out;
3814 gmx_bitmask_t mask_0;
3820 bitmask_init_bit(&mask_0, 0);
3821 for (int b = 0; b < flags->nflag; b++)
3823 if (bitmask_is_equal(flags->flag[b], mask_0))
3825 /* Only flag 0 is set, no copy of reduction required */
3829 else if (!bitmask_is_zero(flags->flag[b]))
3832 for (out = 0; out < nout; out++)
3834 if (bitmask_is_set(flags->flag[b], out))
3851 fprintf(debug, "nbnxn reduction: #flag %d #list %d elem %4.2f, keep %4.2f copy %4.2f red %4.2f\n",
3853 nelem/(double)(flags->nflag),
3854 nkeep/(double)(flags->nflag),
3855 ncopy/(double)(flags->nflag),
3856 nred/(double)(flags->nflag));
3859 /* Copies the list entries from src to dest when cjStart <= *cjGlobal < cjEnd.
3860 * *cjGlobal is updated with the cj count in src.
3861 * When setFlags==true, flag bit t is set in flag for all i and j clusters.
3863 template<bool setFlags>
3864 static void copySelectedListRange(const nbnxn_ci_t * gmx_restrict srcCi,
3865 const nbnxn_pairlist_t * gmx_restrict src,
3866 nbnxn_pairlist_t * gmx_restrict dest,
3867 gmx_bitmask_t *flag,
3868 int iFlagShift, int jFlagShift, int t)
3870 int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3872 if (dest->nci + 1 >= dest->ci_nalloc)
3874 nb_realloc_ci(dest, dest->nci + 1);
3876 check_cell_list_space_simple(dest, ncj);
3878 dest->ci[dest->nci] = *srcCi;
3879 dest->ci[dest->nci].cj_ind_start = dest->ncj;
3880 dest->ci[dest->nci].cj_ind_end = dest->ncj + ncj;
3884 bitmask_init_bit(&flag[srcCi->ci >> iFlagShift], t);
3887 for (int j = srcCi->cj_ind_start; j < srcCi->cj_ind_end; j++)
3889 dest->cj[dest->ncj++] = src->cj[j];
3893 /* NOTE: This is relatively expensive, since this
3894 * operation is done for all elements in the list,
3895 * whereas at list generation this is done only
3896 * once for each flag entry.
3898 bitmask_init_bit(&flag[src->cj[j].cj >> jFlagShift], t);
3905 /* This routine re-balances the pairlists such that all are nearly equally
3906 * sized. Only whole i-entries are moved between lists. These are moved
3907 * between the ends of the lists, such that the buffer reduction cost should
3908 * not change significantly.
3909 * Note that all original reduction flags are currently kept. This can lead
3910 * to reduction of parts of the force buffer that could be avoided. But since
3911 * the original lists are quite balanced, this will only give minor overhead.
3913 static void rebalanceSimpleLists(int numLists,
3914 nbnxn_pairlist_t * const * const srcSet,
3915 nbnxn_pairlist_t **destSet,
3916 nbnxn_search_work_t *searchWork)
3919 for (int s = 0; s < numLists; s++)
3921 ncjTotal += srcSet[s]->ncjInUse;
3923 int ncjTarget = (ncjTotal + numLists - 1)/numLists;
3925 #pragma omp parallel num_threads(numLists)
3927 int t = gmx_omp_get_thread_num();
3929 int cjStart = ncjTarget* t;
3930 int cjEnd = ncjTarget*(t + 1);
3932 /* The destination pair-list for task/thread t */
3933 nbnxn_pairlist_t *dest = destSet[t];
3935 clear_pairlist(dest);
3936 dest->bSimple = srcSet[0]->bSimple;
3937 dest->na_ci = srcSet[0]->na_ci;
3938 dest->na_cj = srcSet[0]->na_cj;
3940 /* Note that the flags in the work struct (still) contain flags
3941 * for all entries that are present in srcSet->nbl[t].
3943 gmx_bitmask_t *flag = searchWork[t].buffer_flags.flag;
3945 int iFlagShift = getBufferFlagShift(dest->na_ci);
3946 int jFlagShift = getBufferFlagShift(dest->na_cj);
3949 for (int s = 0; s < numLists && cjGlobal < cjEnd; s++)
3951 const nbnxn_pairlist_t *src = srcSet[s];
3953 if (cjGlobal + src->ncjInUse > cjStart)
3955 for (int i = 0; i < src->nci && cjGlobal < cjEnd; i++)
3957 const nbnxn_ci_t *srcCi = &src->ci[i];
3958 int ncj = srcCi->cj_ind_end - srcCi->cj_ind_start;
3959 if (cjGlobal >= cjStart)
3961 /* If the source list is not our own, we need to set
3962 * extra flags (the template bool parameter).
3966 copySelectedListRange
3969 flag, iFlagShift, jFlagShift, t);
3973 copySelectedListRange
3976 dest, flag, iFlagShift, jFlagShift, t);
3984 cjGlobal += src->ncjInUse;
3988 dest->ncjInUse = dest->ncj;
3992 int ncjTotalNew = 0;
3993 for (int s = 0; s < numLists; s++)
3995 ncjTotalNew += destSet[s]->ncjInUse;
3997 GMX_RELEASE_ASSERT(ncjTotalNew == ncjTotal, "The total size of the lists before and after rebalancing should match");
4001 /* Returns if the pairlists are so imbalanced that it is worth rebalancing. */
4002 static bool checkRebalanceSimpleLists(const nbnxn_pairlist_set_t *listSet)
4004 int numLists = listSet->nnbl;
4007 for (int s = 0; s < numLists; s++)
4009 ncjMax = std::max(ncjMax, listSet->nbl[s]->ncjInUse);
4010 ncjTotal += listSet->nbl[s]->ncjInUse;
4014 fprintf(debug, "Pair-list ncjMax %d ncjTotal %d\n", ncjMax, ncjTotal);
4016 /* The rebalancing adds 3% extra time to the search. Heuristically we
4017 * determined that under common conditions the non-bonded kernel balance
4018 * improvement will outweigh this when the imbalance is more than 3%.
4019 * But this will, obviously, depend on search vs kernel time and nstlist.
4021 const real rebalanceTolerance = 1.03;
4023 return numLists*ncjMax > ncjTotal*rebalanceTolerance;
4026 /* Perform a count (linear) sort to sort the smaller lists to the end.
4027 * This avoids load imbalance on the GPU, as large lists will be
4028 * scheduled and executed first and the smaller lists later.
4029 * Load balancing between multi-processors only happens at the end
4030 * and there smaller lists lead to more effective load balancing.
4031 * The sorting is done on the cj4 count, not on the actual pair counts.
4032 * Not only does this make the sort faster, but it also results in
4033 * better load balancing than using a list sorted on exact load.
4034 * This function swaps the pointer in the pair list to avoid a copy operation.
4036 static void sort_sci(nbnxn_pairlist_t *nbl)
4038 nbnxn_list_work_t *work;
4040 nbnxn_sci_t *sci_sort;
4042 if (nbl->ncj4 <= nbl->nsci)
4044 /* nsci = 0 or all sci have size 1, sorting won't change the order */
4050 /* We will distinguish differences up to double the average */
4051 m = (2*nbl->ncj4)/nbl->nsci;
4053 if (m + 1 > work->sort_nalloc)
4055 work->sort_nalloc = over_alloc_large(m + 1);
4056 srenew(work->sort, work->sort_nalloc);
4059 if (work->sci_sort_nalloc != nbl->sci_nalloc)
4061 work->sci_sort_nalloc = nbl->sci_nalloc;
4062 nbnxn_realloc_void((void **)&work->sci_sort,
4064 work->sci_sort_nalloc*sizeof(*work->sci_sort),
4065 nbl->alloc, nbl->free);
4068 /* Count the entries of each size */
4069 for (int i = 0; i <= m; i++)
4073 for (int s = 0; s < nbl->nsci; s++)
4075 int i = std::min(m, nbl->sci[s].cj4_ind_end - nbl->sci[s].cj4_ind_start);
4078 /* Calculate the offset for each count */
4081 for (int i = m - 1; i >= 0; i--)
4084 work->sort[i] = work->sort[i + 1] + s0;
4088 /* Sort entries directly into place */
4089 sci_sort = work->sci_sort;
4090 for (int s = 0; s < nbl->nsci; s++)
4092 int i = std::min(m, nbl->sci[s].cj4_ind_end - nbl->sci[s].cj4_ind_start);
4093 sci_sort[work->sort[i]++] = nbl->sci[s];
4096 /* Swap the sci pointers so we use the new, sorted list */
4097 work->sci_sort = nbl->sci;
4098 nbl->sci = sci_sort;
4101 /* Make a local or non-local pair-list, depending on iloc */
4102 void nbnxn_make_pairlist(const nbnxn_search_t nbs,
4103 nbnxn_atomdata_t *nbat,
4104 const t_blocka *excl,
4106 int min_ci_balanced,
4107 nbnxn_pairlist_set_t *nbl_list,
4112 nbnxn_grid_t *gridi, *gridj;
4115 int nsubpair_target;
4116 float nsubpair_tot_est;
4118 nbnxn_pairlist_t **nbl;
4120 gmx_bool CombineNBLists;
4122 int np_tot, np_noq, np_hlj, nap;
4124 /* Check if we are running hybrid GPU + CPU nbnxn mode */
4125 bGPUCPU = (!nbs->grid[0].bSimple && nbl_list->bSimple);
4127 nnbl = nbl_list->nnbl;
4128 nbl = nbl_list->nbl;
4129 CombineNBLists = nbl_list->bCombined;
4133 fprintf(debug, "ns making %d nblists\n", nnbl);
4136 nbat->bUseBufferFlags = (nbat->nout > 1);
4137 /* We should re-init the flags before making the first list */
4138 if (nbat->bUseBufferFlags && (LOCAL_I(iloc) || bGPUCPU))
4140 init_buffer_flags(&nbat->buffer_flags, nbat->natoms);
4143 if (nbl_list->bSimple)
4146 switch (nb_kernel_type)
4148 #ifdef GMX_NBNXN_SIMD_4XN
4149 case nbnxnk4xN_SIMD_4xN:
4150 nbs->icell_set_x = icell_set_x_simd_4xn;
4153 #ifdef GMX_NBNXN_SIMD_2XNN
4154 case nbnxnk4xN_SIMD_2xNN:
4155 nbs->icell_set_x = icell_set_x_simd_2xnn;
4159 nbs->icell_set_x = icell_set_x_simple;
4163 /* MSVC 2013 complains about switch statements without case */
4164 nbs->icell_set_x = icell_set_x_simple;
4169 nbs->icell_set_x = icell_set_x_supersub;
4174 /* Only zone (grid) 0 vs 0 */
4181 nzi = nbs->zones->nizone;
4184 if (!nbl_list->bSimple && min_ci_balanced > 0)
4186 get_nsubpair_target(nbs, iloc, rlist, min_ci_balanced,
4187 &nsubpair_target, &nsubpair_tot_est);
4191 nsubpair_target = 0;
4192 nsubpair_tot_est = 0;
4195 /* Clear all pair-lists */
4196 for (int th = 0; th < nnbl; th++)
4198 clear_pairlist(nbl[th]);
4202 clear_pairlist_fep(nbl_list->nbl_fep[th]);
4206 for (int zi = 0; zi < nzi; zi++)
4208 gridi = &nbs->grid[zi];
4210 if (NONLOCAL_I(iloc))
4212 zj0 = nbs->zones->izone[zi].j0;
4213 zj1 = nbs->zones->izone[zi].j1;
4219 for (int zj = zj0; zj < zj1; zj++)
4221 gridj = &nbs->grid[zj];
4225 fprintf(debug, "ns search grid %d vs %d\n", zi, zj);
4228 nbs_cycle_start(&nbs->cc[enbsCCsearch]);
4230 if (nbl[0]->bSimple && !gridi->bSimple)
4232 /* Hybrid list, determine blocking later */
4237 ci_block = get_ci_block_size(gridi, nbs->DomDec, nnbl);
4240 /* With GPU: generate progressively smaller lists for
4241 * load balancing for local only or non-local with 2 zones.
4243 progBal = (LOCAL_I(iloc) || nbs->zones->n <= 2);
4245 #pragma omp parallel for num_threads(nnbl) schedule(static)
4246 for (int th = 0; th < nnbl; th++)
4250 /* Re-init the thread-local work flag data before making
4251 * the first list (not an elegant conditional).
4253 if (nbat->bUseBufferFlags && ((zi == 0 && zj == 0) ||
4254 (bGPUCPU && zi == 0 && zj == 1)))
4256 init_buffer_flags(&nbs->work[th].buffer_flags, nbat->natoms);
4259 if (CombineNBLists && th > 0)
4261 clear_pairlist(nbl[th]);
4264 /* Divide the i super cell equally over the nblists */
4265 nbnxn_make_pairlist_part(nbs, gridi, gridj,
4266 &nbs->work[th], nbat, excl,
4270 nbat->bUseBufferFlags,
4272 progBal, nsubpair_tot_est,
4275 nbl_list->nbl_fep[th]);
4277 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4279 nbs_cycle_stop(&nbs->cc[enbsCCsearch]);
4284 for (int th = 0; th < nnbl; th++)
4286 inc_nrnb(nrnb, eNR_NBNXN_DIST2, nbs->work[th].ndistc);
4288 if (nbl_list->bSimple)
4290 np_tot += nbl[th]->ncj;
4291 np_noq += nbl[th]->work->ncj_noq;
4292 np_hlj += nbl[th]->work->ncj_hlj;
4296 /* This count ignores potential subsequent pair pruning */
4297 np_tot += nbl[th]->nci_tot;
4300 nap = nbl[0]->na_ci*nbl[0]->na_cj;
4301 nbl_list->natpair_ljq = (np_tot - np_noq)*nap - np_hlj*nap/2;
4302 nbl_list->natpair_lj = np_noq*nap;
4303 nbl_list->natpair_q = np_hlj*nap/2;
4305 if (CombineNBLists && nnbl > 1)
4307 nbs_cycle_start(&nbs->cc[enbsCCcombine]);
4309 combine_nblists(nnbl-1, nbl+1, nbl[0]);
4311 nbs_cycle_stop(&nbs->cc[enbsCCcombine]);
4316 if (nbl_list->bSimple)
4318 if (nnbl > 1 && checkRebalanceSimpleLists(nbl_list))
4320 rebalanceSimpleLists(nbl_list->nnbl, nbl_list->nbl, nbl_list->nbl_work, nbs->work);
4322 /* Swap the pointer of the sets of pair lists */
4323 nbnxn_pairlist_t **tmp = nbl_list->nbl;
4324 nbl_list->nbl = nbl_list->nbl_work;
4325 nbl_list->nbl_work = tmp;
4330 /* Sort the entries on size, large ones first */
4331 if (CombineNBLists || nnbl == 1)
4337 #pragma omp parallel for num_threads(nnbl) schedule(static)
4338 for (int th = 0; th < nnbl; th++)
4344 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4349 if (nbat->bUseBufferFlags)
4351 reduce_buffer_flags(nbs, nbl_list->nnbl, &nbat->buffer_flags);
4356 /* Balance the free-energy lists over all the threads */
4357 balance_fep_lists(nbs, nbl_list);
4360 /* Special performance logging stuff (env.var. GMX_NBNXN_CYCLE) */
4363 nbs->search_count++;
4365 if (nbs->print_cycles &&
4366 (!nbs->DomDec || !LOCAL_I(iloc)) &&
4367 nbs->search_count % 100 == 0)
4369 nbs_cycle_print(stderr, nbs);
4372 /* If we have more than one list, they either got rebalancing (CPU)
4373 * or combined (GPU), so we should dump the final result to debug.
4375 if (debug && nbl_list->nnbl > 1)
4377 if (nbl_list->bSimple)
4379 for (int t = 0; t < nbl_list->nnbl; t++)
4381 print_nblist_statistics_simple(debug, nbl_list->nbl[t], nbs, rlist);
4386 print_nblist_statistics_supersub(debug, nbl_list->nbl[0], nbs, rlist);
4394 if (nbl_list->bSimple)
4396 for (int t = 0; t < nbl_list->nnbl; t++)
4398 print_nblist_ci_cj(debug, nbl_list->nbl[t]);
4403 print_nblist_sci_cj(debug, nbl_list->nbl[0]);
4407 if (nbat->bUseBufferFlags)
4409 print_reduction_cost(&nbat->buffer_flags, nbl_list->nnbl);