1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2012, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
45 #include "nbnxn_consts.h"
46 /* nbnxn_internal.h included gmx_simd_macros.h */
47 #include "nbnxn_internal.h"
49 #include "gmx_simd_vec.h"
51 #include "nbnxn_atomdata.h"
52 #include "nbnxn_search.h"
53 #include "gmx_omp_nthreads.h"
56 #include "gromacs/fileio/gmxfio.h"
58 #ifdef NBNXN_SEARCH_BB_SIMD4
59 /* We use 4-wide SIMD for bounding box calculations */
62 /* Single precision BBs + coordinates, we can also load coordinates with SIMD */
63 #define NBNXN_SEARCH_SIMD4_FLOAT_X_BB
66 #if defined NBNXN_SEARCH_SIMD4_FLOAT_X_BB && (GPU_NSUBCELL == 4 || GPU_NSUBCELL == 8)
67 /* Store bounding boxes with x, y and z coordinates in packs of 4 */
68 #define NBNXN_PBB_SIMD4
71 /* The packed bounding box coordinate stride is always set to 4.
72 * With AVX we could use 8, but that turns out not to be faster.
75 #define STRIDE_PBB_2LOG 2
77 #endif /* NBNXN_SEARCH_BB_SIMD4 */
81 /* The functions below are macros as they are performance sensitive */
83 /* 4x4 list, pack=4: no complex conversion required */
84 /* i-cluster to j-cluster conversion */
85 #define CI_TO_CJ_J4(ci) (ci)
86 /* cluster index to coordinate array index conversion */
87 #define X_IND_CI_J4(ci) ((ci)*STRIDE_P4)
88 #define X_IND_CJ_J4(cj) ((cj)*STRIDE_P4)
90 /* 4x2 list, pack=4: j-cluster size is half the packing width */
91 /* i-cluster to j-cluster conversion */
92 #define CI_TO_CJ_J2(ci) ((ci)<<1)
93 /* cluster index to coordinate array index conversion */
94 #define X_IND_CI_J2(ci) ((ci)*STRIDE_P4)
95 #define X_IND_CJ_J2(cj) (((cj)>>1)*STRIDE_P4 + ((cj) & 1)*(PACK_X4>>1))
97 /* 4x8 list, pack=8: i-cluster size is half the packing width */
98 /* i-cluster to j-cluster conversion */
99 #define CI_TO_CJ_J8(ci) ((ci)>>1)
100 /* cluster index to coordinate array index conversion */
101 #define X_IND_CI_J8(ci) (((ci)>>1)*STRIDE_P8 + ((ci) & 1)*(PACK_X8>>1))
102 #define X_IND_CJ_J8(cj) ((cj)*STRIDE_P8)
104 /* The j-cluster size is matched to the SIMD width */
105 #if GMX_SIMD_WIDTH_HERE == 2
106 #define CI_TO_CJ_SIMD_4XN(ci) CI_TO_CJ_J2(ci)
107 #define X_IND_CI_SIMD_4XN(ci) X_IND_CI_J2(ci)
108 #define X_IND_CJ_SIMD_4XN(cj) X_IND_CJ_J2(cj)
110 #if GMX_SIMD_WIDTH_HERE == 4
111 #define CI_TO_CJ_SIMD_4XN(ci) CI_TO_CJ_J4(ci)
112 #define X_IND_CI_SIMD_4XN(ci) X_IND_CI_J4(ci)
113 #define X_IND_CJ_SIMD_4XN(cj) X_IND_CJ_J4(cj)
115 #if GMX_SIMD_WIDTH_HERE == 8
116 #define CI_TO_CJ_SIMD_4XN(ci) CI_TO_CJ_J8(ci)
117 #define X_IND_CI_SIMD_4XN(ci) X_IND_CI_J8(ci)
118 #define X_IND_CJ_SIMD_4XN(cj) X_IND_CJ_J8(cj)
119 /* Half SIMD with j-cluster size */
120 #define CI_TO_CJ_SIMD_2XNN(ci) CI_TO_CJ_J4(ci)
121 #define X_IND_CI_SIMD_2XNN(ci) X_IND_CI_J4(ci)
122 #define X_IND_CJ_SIMD_2XNN(cj) X_IND_CJ_J4(cj)
124 #if GMX_SIMD_WIDTH_HERE == 16
125 #define CI_TO_CJ_SIMD_2XNN(ci) CI_TO_CJ_J8(ci)
126 #define X_IND_CI_SIMD_2XNN(ci) X_IND_CI_J8(ci)
127 #define X_IND_CJ_SIMD_2XNN(cj) X_IND_CJ_J8(cj)
129 #error "unsupported GMX_NBNXN_SIMD_WIDTH"
135 #endif /* GMX_NBNXN_SIMD */
138 #ifdef NBNXN_SEARCH_BB_SIMD4
139 /* Store bounding boxes corners as quadruplets: xxxxyyyyzzzz */
141 /* Size of bounding box corners quadruplet */
142 #define NNBSBB_XXXX (NNBSBB_D*DIM*STRIDE_PBB)
145 /* We shift the i-particles backward for PBC.
146 * This leads to more conditionals than shifting forward.
147 * We do this to get more balanced pair lists.
149 #define NBNXN_SHIFT_BACKWARD
152 /* This define is a lazy way to avoid interdependence of the grid
153 * and searching data structures.
155 #define NBNXN_NA_SC_MAX (GPU_NSUBCELL*NBNXN_GPU_CLUSTER_SIZE)
158 static void nbs_cycle_clear(nbnxn_cycle_t *cc)
162 for (i = 0; i < enbsCCnr; i++)
169 static double Mcyc_av(const nbnxn_cycle_t *cc)
171 return (double)cc->c*1e-6/cc->count;
174 static void nbs_cycle_print(FILE *fp, const nbnxn_search_t nbs)
180 fprintf(fp, "ns %4d grid %4.1f search %4.1f red.f %5.3f",
181 nbs->cc[enbsCCgrid].count,
182 Mcyc_av(&nbs->cc[enbsCCgrid]),
183 Mcyc_av(&nbs->cc[enbsCCsearch]),
184 Mcyc_av(&nbs->cc[enbsCCreducef]));
186 if (nbs->nthread_max > 1)
188 if (nbs->cc[enbsCCcombine].count > 0)
190 fprintf(fp, " comb %5.2f",
191 Mcyc_av(&nbs->cc[enbsCCcombine]));
193 fprintf(fp, " s. th");
194 for (t = 0; t < nbs->nthread_max; t++)
196 fprintf(fp, " %4.1f",
197 Mcyc_av(&nbs->work[t].cc[enbsCCsearch]));
203 static void nbnxn_grid_init(nbnxn_grid_t * grid)
206 grid->cxy_ind = NULL;
207 grid->cxy_nalloc = 0;
213 static int get_2log(int n)
218 while ((1<<log2) < n)
224 gmx_fatal(FARGS, "nbnxn na_c (%d) is not a power of 2", n);
230 static int nbnxn_kernel_to_ci_size(int nb_kernel_type)
232 switch (nb_kernel_type)
234 case nbnxnk4x4_PlainC:
235 case nbnxnk4xN_SIMD_4xN:
236 case nbnxnk4xN_SIMD_2xNN:
237 return NBNXN_CPU_CLUSTER_I_SIZE;
238 case nbnxnk8x8x8_CUDA:
239 case nbnxnk8x8x8_PlainC:
240 /* The cluster size for super/sub lists is only set here.
241 * Any value should work for the pair-search and atomdata code.
242 * The kernels, of course, might require a particular value.
244 return NBNXN_GPU_CLUSTER_SIZE;
246 gmx_incons("unknown kernel type");
252 int nbnxn_kernel_to_cj_size(int nb_kernel_type)
254 int nbnxn_simd_width = 0;
257 #ifdef GMX_NBNXN_SIMD
258 nbnxn_simd_width = GMX_SIMD_WIDTH_HERE;
261 switch (nb_kernel_type)
263 case nbnxnk4x4_PlainC:
264 cj_size = NBNXN_CPU_CLUSTER_I_SIZE;
266 case nbnxnk4xN_SIMD_4xN:
267 cj_size = nbnxn_simd_width;
269 case nbnxnk4xN_SIMD_2xNN:
270 cj_size = nbnxn_simd_width/2;
272 case nbnxnk8x8x8_CUDA:
273 case nbnxnk8x8x8_PlainC:
274 cj_size = nbnxn_kernel_to_ci_size(nb_kernel_type);
277 gmx_incons("unknown kernel type");
283 static int ci_to_cj(int na_cj_2log, int ci)
287 case 2: return ci; break;
288 case 1: return (ci<<1); break;
289 case 3: return (ci>>1); break;
295 gmx_bool nbnxn_kernel_pairlist_simple(int nb_kernel_type)
297 if (nb_kernel_type == nbnxnkNotSet)
299 gmx_fatal(FARGS, "Non-bonded kernel type not set for Verlet-style pair-list.");
302 switch (nb_kernel_type)
304 case nbnxnk8x8x8_CUDA:
305 case nbnxnk8x8x8_PlainC:
308 case nbnxnk4x4_PlainC:
309 case nbnxnk4xN_SIMD_4xN:
310 case nbnxnk4xN_SIMD_2xNN:
314 gmx_incons("Invalid nonbonded kernel type passed!");
319 void nbnxn_init_search(nbnxn_search_t * nbs_ptr,
321 gmx_domdec_zones_t *zones,
330 nbs->DomDec = (n_dd_cells != NULL);
332 clear_ivec(nbs->dd_dim);
338 for (d = 0; d < DIM; d++)
340 if ((*n_dd_cells)[d] > 1)
343 /* Each grid matches a DD zone */
349 snew(nbs->grid, nbs->ngrid);
350 for (g = 0; g < nbs->ngrid; g++)
352 nbnxn_grid_init(&nbs->grid[g]);
355 nbs->cell_nalloc = 0;
359 nbs->nthread_max = nthread_max;
361 /* Initialize the work data structures for each thread */
362 snew(nbs->work, nbs->nthread_max);
363 for (t = 0; t < nbs->nthread_max; t++)
365 nbs->work[t].cxy_na = NULL;
366 nbs->work[t].cxy_na_nalloc = 0;
367 nbs->work[t].sort_work = NULL;
368 nbs->work[t].sort_work_nalloc = 0;
371 /* Initialize detailed nbsearch cycle counting */
372 nbs->print_cycles = (getenv("GMX_NBNXN_CYCLE") != 0);
373 nbs->search_count = 0;
374 nbs_cycle_clear(nbs->cc);
375 for (t = 0; t < nbs->nthread_max; t++)
377 nbs_cycle_clear(nbs->work[t].cc);
381 static real grid_atom_density(int n, rvec corner0, rvec corner1)
385 rvec_sub(corner1, corner0, size);
387 return n/(size[XX]*size[YY]*size[ZZ]);
390 static int set_grid_size_xy(const nbnxn_search_t nbs,
393 int n, rvec corner0, rvec corner1,
398 real adens, tlen, tlen_x, tlen_y, nc_max;
401 rvec_sub(corner1, corner0, size);
405 /* target cell length */
408 /* To minimize the zero interactions, we should make
409 * the largest of the i/j cell cubic.
411 na_c = max(grid->na_c, grid->na_cj);
413 /* Approximately cubic cells */
414 tlen = pow(na_c/atom_density, 1.0/3.0);
420 /* Approximately cubic sub cells */
421 tlen = pow(grid->na_c/atom_density, 1.0/3.0);
422 tlen_x = tlen*GPU_NSUBCELL_X;
423 tlen_y = tlen*GPU_NSUBCELL_Y;
425 /* We round ncx and ncy down, because we get less cell pairs
426 * in the nbsist when the fixed cell dimensions (x,y) are
427 * larger than the variable one (z) than the other way around.
429 grid->ncx = max(1, (int)(size[XX]/tlen_x));
430 grid->ncy = max(1, (int)(size[YY]/tlen_y));
438 grid->sx = size[XX]/grid->ncx;
439 grid->sy = size[YY]/grid->ncy;
440 grid->inv_sx = 1/grid->sx;
441 grid->inv_sy = 1/grid->sy;
445 /* This is a non-home zone, add an extra row of cells
446 * for particles communicated for bonded interactions.
447 * These can be beyond the cut-off. It doesn't matter where
448 * they end up on the grid, but for performance it's better
449 * if they don't end up in cells that can be within cut-off range.
455 /* We need one additional cell entry for particles moved by DD */
456 if (grid->ncx*grid->ncy+1 > grid->cxy_nalloc)
458 grid->cxy_nalloc = over_alloc_large(grid->ncx*grid->ncy+1);
459 srenew(grid->cxy_na, grid->cxy_nalloc);
460 srenew(grid->cxy_ind, grid->cxy_nalloc+1);
462 for (t = 0; t < nbs->nthread_max; t++)
464 if (grid->ncx*grid->ncy+1 > nbs->work[t].cxy_na_nalloc)
466 nbs->work[t].cxy_na_nalloc = over_alloc_large(grid->ncx*grid->ncy+1);
467 srenew(nbs->work[t].cxy_na, nbs->work[t].cxy_na_nalloc);
471 /* Worst case scenario of 1 atom in each last cell */
472 if (grid->na_cj <= grid->na_c)
474 nc_max = n/grid->na_sc + grid->ncx*grid->ncy;
478 nc_max = n/grid->na_sc + grid->ncx*grid->ncy*grid->na_cj/grid->na_c;
481 if (nc_max > grid->nc_nalloc)
483 grid->nc_nalloc = over_alloc_large(nc_max);
484 srenew(grid->nsubc, grid->nc_nalloc);
485 srenew(grid->bbcz, grid->nc_nalloc*NNBSBB_D);
487 sfree_aligned(grid->bb);
488 /* This snew also zeros the contents, this avoid possible
489 * floating exceptions in SIMD with the unused bb elements.
493 snew_aligned(grid->bb, grid->nc_nalloc, 16);
500 pbb_nalloc = grid->nc_nalloc*GPU_NSUBCELL/STRIDE_PBB*NNBSBB_XXXX;
501 snew_aligned(grid->pbb, pbb_nalloc, 16);
503 snew_aligned(grid->bb, grid->nc_nalloc*GPU_NSUBCELL, 16);
509 if (grid->na_cj == grid->na_c)
511 grid->bbj = grid->bb;
515 sfree_aligned(grid->bbj);
516 snew_aligned(grid->bbj, grid->nc_nalloc*grid->na_c/grid->na_cj, 16);
520 srenew(grid->flags, grid->nc_nalloc);
523 copy_rvec(corner0, grid->c0);
524 copy_rvec(corner1, grid->c1);
529 /* We need to sort paricles in grid columns on z-coordinate.
530 * As particle are very often distributed homogeneously, we a sorting
531 * algorithm similar to pigeonhole sort. We multiply the z-coordinate
532 * by a factor, cast to an int and try to store in that hole. If the hole
533 * is full, we move this or another particle. A second pass is needed to make
534 * contiguous elements. SORT_GRID_OVERSIZE is the ratio of holes to particles.
535 * 4 is the optimal value for homogeneous particle distribution and allows
536 * for an O(#particles) sort up till distributions were all particles are
537 * concentrated in 1/4 of the space. No NlogN fallback is implemented,
538 * as it can be expensive to detect imhomogeneous particle distributions.
539 * SGSF is the maximum ratio of holes used, in the worst case all particles
540 * end up in the last hole and we need #particles extra holes at the end.
542 #define SORT_GRID_OVERSIZE 4
543 #define SGSF (SORT_GRID_OVERSIZE + 1)
545 /* Sort particle index a on coordinates x along dim.
546 * Backwards tells if we want decreasing iso increasing coordinates.
547 * h0 is the minimum of the coordinate range.
548 * invh is the 1/length of the sorting range.
549 * n_per_h (>=n) is the expected average number of particles per 1/invh
550 * sort is the sorting work array.
551 * sort should have a size of at least n_per_h*SORT_GRID_OVERSIZE + n,
552 * or easier, allocate at least n*SGSF elements.
554 static void sort_atoms(int dim, gmx_bool Backwards,
555 int *a, int n, rvec *x,
556 real h0, real invh, int n_per_h,
560 int zi, zim, zi_min, zi_max;
572 gmx_incons("n > n_per_h");
576 /* Transform the inverse range height into the inverse hole height */
577 invh *= n_per_h*SORT_GRID_OVERSIZE;
579 /* Set nsort to the maximum possible number of holes used.
580 * In worst case all n elements end up in the last bin.
582 nsort = n_per_h*SORT_GRID_OVERSIZE + n;
584 /* Determine the index range used, so we can limit it for the second pass */
588 /* Sort the particles using a simple index sort */
589 for (i = 0; i < n; i++)
591 /* The cast takes care of float-point rounding effects below zero.
592 * This code assumes particles are less than 1/SORT_GRID_OVERSIZE
593 * times the box height out of the box.
595 zi = (int)((x[a[i]][dim] - h0)*invh);
598 /* As we can have rounding effect, we use > iso >= here */
599 if (zi < 0 || zi > n_per_h*SORT_GRID_OVERSIZE)
601 gmx_fatal(FARGS, "(int)((x[%d][%c]=%f - %f)*%f) = %d, not in 0 - %d*%d\n",
602 a[i], 'x'+dim, x[a[i]][dim], h0, invh, zi,
603 n_per_h, SORT_GRID_OVERSIZE);
607 /* Ideally this particle should go in sort cell zi,
608 * but that might already be in use,
609 * in that case find the first empty cell higher up
614 zi_min = min(zi_min, zi);
615 zi_max = max(zi_max, zi);
619 /* We have multiple atoms in the same sorting slot.
620 * Sort on real z for minimal bounding box size.
621 * There is an extra check for identical z to ensure
622 * well-defined output order, independent of input order
623 * to ensure binary reproducibility after restarts.
625 while (sort[zi] >= 0 && ( x[a[i]][dim] > x[sort[zi]][dim] ||
626 (x[a[i]][dim] == x[sort[zi]][dim] &&
634 /* Shift all elements by one slot until we find an empty slot */
637 while (sort[zim] >= 0)
645 zi_max = max(zi_max, zim);
648 zi_max = max(zi_max, zi);
655 for (zi = 0; zi < nsort; zi++)
666 for (zi = zi_max; zi >= zi_min; zi--)
677 gmx_incons("Lost particles while sorting");
682 #define R2F_D(x) ((float)((x) >= 0 ? ((1-GMX_FLOAT_EPS)*(x)) : ((1+GMX_FLOAT_EPS)*(x))))
683 #define R2F_U(x) ((float)((x) >= 0 ? ((1+GMX_FLOAT_EPS)*(x)) : ((1-GMX_FLOAT_EPS)*(x))))
689 /* Coordinate order x,y,z, bb order xyz0 */
690 static void calc_bounding_box(int na, int stride, const real *x, nbnxn_bb_t *bb)
693 real xl, xh, yl, yh, zl, zh;
703 for (j = 1; j < na; j++)
705 xl = min(xl, x[i+XX]);
706 xh = max(xh, x[i+XX]);
707 yl = min(yl, x[i+YY]);
708 yh = max(yh, x[i+YY]);
709 zl = min(zl, x[i+ZZ]);
710 zh = max(zh, x[i+ZZ]);
713 /* Note: possible double to float conversion here */
714 bb->lower[BB_X] = R2F_D(xl);
715 bb->lower[BB_Y] = R2F_D(yl);
716 bb->lower[BB_Z] = R2F_D(zl);
717 bb->upper[BB_X] = R2F_U(xh);
718 bb->upper[BB_Y] = R2F_U(yh);
719 bb->upper[BB_Z] = R2F_U(zh);
722 /* Packed coordinates, bb order xyz0 */
723 static void calc_bounding_box_x_x4(int na, const real *x, nbnxn_bb_t *bb)
726 real xl, xh, yl, yh, zl, zh;
734 for (j = 1; j < na; j++)
736 xl = min(xl, x[j+XX*PACK_X4]);
737 xh = max(xh, x[j+XX*PACK_X4]);
738 yl = min(yl, x[j+YY*PACK_X4]);
739 yh = max(yh, x[j+YY*PACK_X4]);
740 zl = min(zl, x[j+ZZ*PACK_X4]);
741 zh = max(zh, x[j+ZZ*PACK_X4]);
743 /* Note: possible double to float conversion here */
744 bb->lower[BB_X] = R2F_D(xl);
745 bb->lower[BB_Y] = R2F_D(yl);
746 bb->lower[BB_Z] = R2F_D(zl);
747 bb->upper[BB_X] = R2F_U(xh);
748 bb->upper[BB_Y] = R2F_U(yh);
749 bb->upper[BB_Z] = R2F_U(zh);
752 /* Packed coordinates, bb order xyz0 */
753 static void calc_bounding_box_x_x8(int na, const real *x, nbnxn_bb_t *bb)
756 real xl, xh, yl, yh, zl, zh;
764 for (j = 1; j < na; j++)
766 xl = min(xl, x[j+XX*PACK_X8]);
767 xh = max(xh, x[j+XX*PACK_X8]);
768 yl = min(yl, x[j+YY*PACK_X8]);
769 yh = max(yh, x[j+YY*PACK_X8]);
770 zl = min(zl, x[j+ZZ*PACK_X8]);
771 zh = max(zh, x[j+ZZ*PACK_X8]);
773 /* Note: possible double to float conversion here */
774 bb->lower[BB_X] = R2F_D(xl);
775 bb->lower[BB_Y] = R2F_D(yl);
776 bb->lower[BB_Z] = R2F_D(zl);
777 bb->upper[BB_X] = R2F_U(xh);
778 bb->upper[BB_Y] = R2F_U(yh);
779 bb->upper[BB_Z] = R2F_U(zh);
782 /* Packed coordinates, bb order xyz0 */
783 static void calc_bounding_box_x_x4_halves(int na, const real *x,
784 nbnxn_bb_t *bb, nbnxn_bb_t *bbj)
786 calc_bounding_box_x_x4(min(na, 2), x, bbj);
790 calc_bounding_box_x_x4(min(na-2, 2), x+(PACK_X4>>1), bbj+1);
794 /* Set the "empty" bounding box to the same as the first one,
795 * so we don't need to treat special cases in the rest of the code.
797 #ifdef NBNXN_SEARCH_BB_SIMD4
798 gmx_simd4_store_pr(&bbj[1].lower[0], gmx_simd4_load_bb_pr(&bbj[0].lower[0]));
799 gmx_simd4_store_pr(&bbj[1].upper[0], gmx_simd4_load_bb_pr(&bbj[0].upper[0]));
805 #ifdef NBNXN_SEARCH_BB_SIMD4
806 gmx_simd4_store_pr(&bb->lower[0],
807 gmx_simd4_min_pr(gmx_simd4_load_bb_pr(&bbj[0].lower[0]),
808 gmx_simd4_load_bb_pr(&bbj[1].lower[0])));
809 gmx_simd4_store_pr(&bb->upper[0],
810 gmx_simd4_max_pr(gmx_simd4_load_bb_pr(&bbj[0].upper[0]),
811 gmx_simd4_load_bb_pr(&bbj[1].upper[0])));
816 for (i = 0; i < NNBSBB_C; i++)
818 bb->lower[i] = min(bbj[0].lower[i], bbj[1].lower[i]);
819 bb->upper[i] = max(bbj[0].upper[i], bbj[1].upper[i]);
825 #ifdef NBNXN_SEARCH_BB_SIMD4
827 /* Coordinate order xyz, bb order xxxxyyyyzzzz */
828 static void calc_bounding_box_xxxx(int na, int stride, const real *x, float *bb)
831 real xl, xh, yl, yh, zl, zh;
841 for (j = 1; j < na; j++)
843 xl = min(xl, x[i+XX]);
844 xh = max(xh, x[i+XX]);
845 yl = min(yl, x[i+YY]);
846 yh = max(yh, x[i+YY]);
847 zl = min(zl, x[i+ZZ]);
848 zh = max(zh, x[i+ZZ]);
851 /* Note: possible double to float conversion here */
852 bb[0*STRIDE_PBB] = R2F_D(xl);
853 bb[1*STRIDE_PBB] = R2F_D(yl);
854 bb[2*STRIDE_PBB] = R2F_D(zl);
855 bb[3*STRIDE_PBB] = R2F_U(xh);
856 bb[4*STRIDE_PBB] = R2F_U(yh);
857 bb[5*STRIDE_PBB] = R2F_U(zh);
860 #endif /* NBNXN_SEARCH_BB_SIMD4 */
862 #ifdef NBNXN_SEARCH_SIMD4_FLOAT_X_BB
864 /* Coordinate order xyz?, bb order xyz0 */
865 static void calc_bounding_box_simd4(int na, const float *x, nbnxn_bb_t *bb)
867 gmx_simd4_pr bb_0_S, bb_1_S;
872 bb_0_S = gmx_simd4_load_bb_pr(x);
875 for (i = 1; i < na; i++)
877 x_S = gmx_simd4_load_bb_pr(x+i*NNBSBB_C);
878 bb_0_S = gmx_simd4_min_pr(bb_0_S, x_S);
879 bb_1_S = gmx_simd4_max_pr(bb_1_S, x_S);
882 gmx_simd4_store_pr(&bb->lower[0], bb_0_S);
883 gmx_simd4_store_pr(&bb->upper[0], bb_1_S);
886 /* Coordinate order xyz?, bb order xxxxyyyyzzzz */
887 static void calc_bounding_box_xxxx_simd4(int na, const float *x,
888 nbnxn_bb_t *bb_work_aligned,
891 calc_bounding_box_simd4(na, x, bb_work_aligned);
893 bb[0*STRIDE_PBB] = bb_work_aligned->lower[BB_X];
894 bb[1*STRIDE_PBB] = bb_work_aligned->lower[BB_Y];
895 bb[2*STRIDE_PBB] = bb_work_aligned->lower[BB_Z];
896 bb[3*STRIDE_PBB] = bb_work_aligned->upper[BB_X];
897 bb[4*STRIDE_PBB] = bb_work_aligned->upper[BB_Y];
898 bb[5*STRIDE_PBB] = bb_work_aligned->upper[BB_Z];
901 #endif /* NBNXN_SEARCH_SIMD4_FLOAT_X_BB */
904 /* Combines pairs of consecutive bounding boxes */
905 static void combine_bounding_box_pairs(nbnxn_grid_t *grid, const nbnxn_bb_t *bb)
907 int i, j, sc2, nc2, c2;
909 for (i = 0; i < grid->ncx*grid->ncy; i++)
911 /* Starting bb in a column is expected to be 2-aligned */
912 sc2 = grid->cxy_ind[i]>>1;
913 /* For odd numbers skip the last bb here */
914 nc2 = (grid->cxy_na[i]+3)>>(2+1);
915 for (c2 = sc2; c2 < sc2+nc2; c2++)
917 #ifdef NBNXN_SEARCH_BB_SIMD4
918 gmx_simd4_pr min_S, max_S;
920 min_S = gmx_simd4_min_pr(gmx_simd4_load_bb_pr(&bb[c2*2+0].lower[0]),
921 gmx_simd4_load_bb_pr(&bb[c2*2+1].lower[0]));
922 max_S = gmx_simd4_max_pr(gmx_simd4_load_bb_pr(&bb[c2*2+0].upper[0]),
923 gmx_simd4_load_bb_pr(&bb[c2*2+1].upper[0]));
924 gmx_simd4_store_pr(&grid->bbj[c2].lower[0], min_S);
925 gmx_simd4_store_pr(&grid->bbj[c2].upper[0], max_S);
927 for (j = 0; j < NNBSBB_C; j++)
929 grid->bbj[c2].lower[j] = min(bb[c2*2+0].lower[j],
930 bb[c2*2+1].lower[j]);
931 grid->bbj[c2].upper[j] = max(bb[c2*2+0].upper[j],
932 bb[c2*2+1].upper[j]);
936 if (((grid->cxy_na[i]+3)>>2) & 1)
938 /* The bb count in this column is odd: duplicate the last bb */
939 for (j = 0; j < NNBSBB_C; j++)
941 grid->bbj[c2].lower[j] = bb[c2*2].lower[j];
942 grid->bbj[c2].upper[j] = bb[c2*2].upper[j];
949 /* Prints the average bb size, used for debug output */
950 static void print_bbsizes_simple(FILE *fp,
951 const nbnxn_search_t nbs,
952 const nbnxn_grid_t *grid)
958 for (c = 0; c < grid->nc; c++)
960 for (d = 0; d < DIM; d++)
962 ba[d] += grid->bb[c].upper[d] - grid->bb[c].lower[d];
965 dsvmul(1.0/grid->nc, ba, ba);
967 fprintf(fp, "ns bb: %4.2f %4.2f %4.2f %4.2f %4.2f %4.2f rel %4.2f %4.2f %4.2f\n",
968 nbs->box[XX][XX]/grid->ncx,
969 nbs->box[YY][YY]/grid->ncy,
970 nbs->box[ZZ][ZZ]*grid->ncx*grid->ncy/grid->nc,
971 ba[XX], ba[YY], ba[ZZ],
972 ba[XX]*grid->ncx/nbs->box[XX][XX],
973 ba[YY]*grid->ncy/nbs->box[YY][YY],
974 ba[ZZ]*grid->nc/(grid->ncx*grid->ncy*nbs->box[ZZ][ZZ]));
977 /* Prints the average bb size, used for debug output */
978 static void print_bbsizes_supersub(FILE *fp,
979 const nbnxn_search_t nbs,
980 const nbnxn_grid_t *grid)
987 for (c = 0; c < grid->nc; c++)
990 for (s = 0; s < grid->nsubc[c]; s += STRIDE_PBB)
994 cs_w = (c*GPU_NSUBCELL + s)/STRIDE_PBB;
995 for (i = 0; i < STRIDE_PBB; i++)
997 for (d = 0; d < DIM; d++)
1000 grid->pbb[cs_w*NNBSBB_XXXX+(DIM+d)*STRIDE_PBB+i] -
1001 grid->pbb[cs_w*NNBSBB_XXXX+ d *STRIDE_PBB+i];
1006 for (s = 0; s < grid->nsubc[c]; s++)
1010 cs = c*GPU_NSUBCELL + s;
1011 for (d = 0; d < DIM; d++)
1013 ba[d] += grid->bb[cs].upper[d] - grid->bb[cs].lower[d];
1017 ns += grid->nsubc[c];
1019 dsvmul(1.0/ns, ba, ba);
1021 fprintf(fp, "ns bb: %4.2f %4.2f %4.2f %4.2f %4.2f %4.2f rel %4.2f %4.2f %4.2f\n",
1022 nbs->box[XX][XX]/(grid->ncx*GPU_NSUBCELL_X),
1023 nbs->box[YY][YY]/(grid->ncy*GPU_NSUBCELL_Y),
1024 nbs->box[ZZ][ZZ]*grid->ncx*grid->ncy/(grid->nc*GPU_NSUBCELL_Z),
1025 ba[XX], ba[YY], ba[ZZ],
1026 ba[XX]*grid->ncx*GPU_NSUBCELL_X/nbs->box[XX][XX],
1027 ba[YY]*grid->ncy*GPU_NSUBCELL_Y/nbs->box[YY][YY],
1028 ba[ZZ]*grid->nc*GPU_NSUBCELL_Z/(grid->ncx*grid->ncy*nbs->box[ZZ][ZZ]));
1031 /* Potentially sorts atoms on LJ coefficients !=0 and ==0.
1032 * Also sets interaction flags.
1034 void sort_on_lj(int na_c,
1035 int a0, int a1, const int *atinfo,
1039 int subc, s, a, n1, n2, a_lj_max, i, j;
1040 int sort1[NBNXN_NA_SC_MAX/GPU_NSUBCELL];
1041 int sort2[NBNXN_NA_SC_MAX/GPU_NSUBCELL];
1047 for (s = a0; s < a1; s += na_c)
1049 /* Make lists for this (sub-)cell on atoms with and without LJ */
1054 for (a = s; a < min(s+na_c, a1); a++)
1056 haveQ = haveQ || GET_CGINFO_HAS_Q(atinfo[order[a]]);
1058 if (GET_CGINFO_HAS_VDW(atinfo[order[a]]))
1060 sort1[n1++] = order[a];
1065 sort2[n2++] = order[a];
1069 /* If we don't have atom with LJ, there's nothing to sort */
1072 *flags |= NBNXN_CI_DO_LJ(subc);
1076 /* Only sort when strictly necessary. Ordering particles
1077 * Ordering particles can lead to less accurate summation
1078 * due to rounding, both for LJ and Coulomb interactions.
1080 if (2*(a_lj_max - s) >= na_c)
1082 for (i = 0; i < n1; i++)
1084 order[a0+i] = sort1[i];
1086 for (j = 0; j < n2; j++)
1088 order[a0+n1+j] = sort2[j];
1092 *flags |= NBNXN_CI_HALF_LJ(subc);
1097 *flags |= NBNXN_CI_DO_COUL(subc);
1103 /* Fill a pair search cell with atoms.
1104 * Potentially sorts atoms and sets the interaction flags.
1106 void fill_cell(const nbnxn_search_t nbs,
1108 nbnxn_atomdata_t *nbat,
1112 int sx, int sy, int sz,
1113 nbnxn_bb_t gmx_unused *bb_work_aligned)
1126 sort_on_lj(grid->na_c, a0, a1, atinfo, nbs->a,
1127 grid->flags+(a0>>grid->na_c_2log)-grid->cell0);
1130 /* Now we have sorted the atoms, set the cell indices */
1131 for (a = a0; a < a1; a++)
1133 nbs->cell[nbs->a[a]] = a;
1136 copy_rvec_to_nbat_real(nbs->a+a0, a1-a0, grid->na_c, x,
1137 nbat->XFormat, nbat->x, a0,
1140 if (nbat->XFormat == nbatX4)
1142 /* Store the bounding boxes as xyz.xyz. */
1143 offset = (a0 - grid->cell0*grid->na_sc) >> grid->na_c_2log;
1144 bb_ptr = grid->bb + offset;
1146 #if defined GMX_NBNXN_SIMD && GMX_SIMD_WIDTH_HERE == 2
1147 if (2*grid->na_cj == grid->na_c)
1149 calc_bounding_box_x_x4_halves(na, nbat->x+X4_IND_A(a0), bb_ptr,
1150 grid->bbj+offset*2);
1155 calc_bounding_box_x_x4(na, nbat->x+X4_IND_A(a0), bb_ptr);
1158 else if (nbat->XFormat == nbatX8)
1160 /* Store the bounding boxes as xyz.xyz. */
1161 offset = (a0 - grid->cell0*grid->na_sc) >> grid->na_c_2log;
1162 bb_ptr = grid->bb + offset;
1164 calc_bounding_box_x_x8(na, nbat->x+X8_IND_A(a0), bb_ptr);
1167 else if (!grid->bSimple)
1169 /* Store the bounding boxes in a format convenient
1170 * for SIMD4 calculations: xxxxyyyyzzzz...
1174 ((a0-grid->cell0*grid->na_sc)>>(grid->na_c_2log+STRIDE_PBB_2LOG))*NNBSBB_XXXX +
1175 (((a0-grid->cell0*grid->na_sc)>>grid->na_c_2log) & (STRIDE_PBB-1));
1177 #ifdef NBNXN_SEARCH_SIMD4_FLOAT_X_BB
1178 if (nbat->XFormat == nbatXYZQ)
1180 calc_bounding_box_xxxx_simd4(na, nbat->x+a0*nbat->xstride,
1181 bb_work_aligned, pbb_ptr);
1186 calc_bounding_box_xxxx(na, nbat->xstride, nbat->x+a0*nbat->xstride,
1191 fprintf(debug, "%2d %2d %2d bb %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f\n",
1193 pbb_ptr[0*STRIDE_PBB], pbb_ptr[3*STRIDE_PBB],
1194 pbb_ptr[1*STRIDE_PBB], pbb_ptr[4*STRIDE_PBB],
1195 pbb_ptr[2*STRIDE_PBB], pbb_ptr[5*STRIDE_PBB]);
1201 /* Store the bounding boxes as xyz.xyz. */
1202 bb_ptr = grid->bb+((a0-grid->cell0*grid->na_sc)>>grid->na_c_2log);
1204 calc_bounding_box(na, nbat->xstride, nbat->x+a0*nbat->xstride,
1210 bbo = (a0 - grid->cell0*grid->na_sc)/grid->na_c;
1211 fprintf(debug, "%2d %2d %2d bb %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f\n",
1213 grid->bb[bbo].lower[BB_X],
1214 grid->bb[bbo].lower[BB_Y],
1215 grid->bb[bbo].lower[BB_Z],
1216 grid->bb[bbo].upper[BB_X],
1217 grid->bb[bbo].upper[BB_Y],
1218 grid->bb[bbo].upper[BB_Z]);
1223 /* Spatially sort the atoms within one grid column */
1224 static void sort_columns_simple(const nbnxn_search_t nbs,
1230 nbnxn_atomdata_t *nbat,
1231 int cxy_start, int cxy_end,
1235 int cx, cy, cz, ncz, cfilled, c;
1236 int na, ash, ind, a;
1241 fprintf(debug, "cell0 %d sorting columns %d - %d, atoms %d - %d\n",
1242 grid->cell0, cxy_start, cxy_end, a0, a1);
1245 /* Sort the atoms within each x,y column in 3 dimensions */
1246 for (cxy = cxy_start; cxy < cxy_end; cxy++)
1249 cy = cxy - cx*grid->ncy;
1251 na = grid->cxy_na[cxy];
1252 ncz = grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy];
1253 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1255 /* Sort the atoms within each x,y column on z coordinate */
1256 sort_atoms(ZZ, FALSE,
1259 1.0/nbs->box[ZZ][ZZ], ncz*grid->na_sc,
1262 /* Fill the ncz cells in this column */
1263 cfilled = grid->cxy_ind[cxy];
1264 for (cz = 0; cz < ncz; cz++)
1266 c = grid->cxy_ind[cxy] + cz;
1268 ash_c = ash + cz*grid->na_sc;
1269 na_c = min(grid->na_sc, na-(ash_c-ash));
1271 fill_cell(nbs, grid, nbat,
1272 ash_c, ash_c+na_c, atinfo, x,
1273 grid->na_sc*cx + (dd_zone >> 2),
1274 grid->na_sc*cy + (dd_zone & 3),
1278 /* This copy to bbcz is not really necessary.
1279 * But it allows to use the same grid search code
1280 * for the simple and supersub cell setups.
1286 grid->bbcz[c*NNBSBB_D ] = grid->bb[cfilled].lower[BB_Z];
1287 grid->bbcz[c*NNBSBB_D+1] = grid->bb[cfilled].upper[BB_Z];
1290 /* Set the unused atom indices to -1 */
1291 for (ind = na; ind < ncz*grid->na_sc; ind++)
1293 nbs->a[ash+ind] = -1;
1298 /* Spatially sort the atoms within one grid column */
1299 static void sort_columns_supersub(const nbnxn_search_t nbs,
1305 nbnxn_atomdata_t *nbat,
1306 int cxy_start, int cxy_end,
1310 int cx, cy, cz = -1, c = -1, ncz;
1311 int na, ash, na_c, ind, a;
1312 int subdiv_z, sub_z, na_z, ash_z;
1313 int subdiv_y, sub_y, na_y, ash_y;
1314 int subdiv_x, sub_x, na_x, ash_x;
1316 /* cppcheck-suppress unassignedVariable */
1317 nbnxn_bb_t bb_work_array[2], *bb_work_aligned;
1319 bb_work_aligned = (nbnxn_bb_t *)(((size_t)(bb_work_array+1)) & (~((size_t)15)));
1323 fprintf(debug, "cell0 %d sorting columns %d - %d, atoms %d - %d\n",
1324 grid->cell0, cxy_start, cxy_end, a0, a1);
1327 subdiv_x = grid->na_c;
1328 subdiv_y = GPU_NSUBCELL_X*subdiv_x;
1329 subdiv_z = GPU_NSUBCELL_Y*subdiv_y;
1331 /* Sort the atoms within each x,y column in 3 dimensions */
1332 for (cxy = cxy_start; cxy < cxy_end; cxy++)
1335 cy = cxy - cx*grid->ncy;
1337 na = grid->cxy_na[cxy];
1338 ncz = grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy];
1339 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1341 /* Sort the atoms within each x,y column on z coordinate */
1342 sort_atoms(ZZ, FALSE,
1345 1.0/nbs->box[ZZ][ZZ], ncz*grid->na_sc,
1348 /* This loop goes over the supercells and subcells along z at once */
1349 for (sub_z = 0; sub_z < ncz*GPU_NSUBCELL_Z; sub_z++)
1351 ash_z = ash + sub_z*subdiv_z;
1352 na_z = min(subdiv_z, na-(ash_z-ash));
1354 /* We have already sorted on z */
1356 if (sub_z % GPU_NSUBCELL_Z == 0)
1358 cz = sub_z/GPU_NSUBCELL_Z;
1359 c = grid->cxy_ind[cxy] + cz;
1361 /* The number of atoms in this supercell */
1362 na_c = min(grid->na_sc, na-(ash_z-ash));
1364 grid->nsubc[c] = min(GPU_NSUBCELL, (na_c+grid->na_c-1)/grid->na_c);
1366 /* Store the z-boundaries of the super cell */
1367 grid->bbcz[c*NNBSBB_D ] = x[nbs->a[ash_z]][ZZ];
1368 grid->bbcz[c*NNBSBB_D+1] = x[nbs->a[ash_z+na_c-1]][ZZ];
1371 #if GPU_NSUBCELL_Y > 1
1372 /* Sort the atoms along y */
1373 sort_atoms(YY, (sub_z & 1),
1374 nbs->a+ash_z, na_z, x,
1375 grid->c0[YY]+cy*grid->sy,
1376 grid->inv_sy, subdiv_z,
1380 for (sub_y = 0; sub_y < GPU_NSUBCELL_Y; sub_y++)
1382 ash_y = ash_z + sub_y*subdiv_y;
1383 na_y = min(subdiv_y, na-(ash_y-ash));
1385 #if GPU_NSUBCELL_X > 1
1386 /* Sort the atoms along x */
1387 sort_atoms(XX, ((cz*GPU_NSUBCELL_Y + sub_y) & 1),
1388 nbs->a+ash_y, na_y, x,
1389 grid->c0[XX]+cx*grid->sx,
1390 grid->inv_sx, subdiv_y,
1394 for (sub_x = 0; sub_x < GPU_NSUBCELL_X; sub_x++)
1396 ash_x = ash_y + sub_x*subdiv_x;
1397 na_x = min(subdiv_x, na-(ash_x-ash));
1399 fill_cell(nbs, grid, nbat,
1400 ash_x, ash_x+na_x, atinfo, x,
1401 grid->na_c*(cx*GPU_NSUBCELL_X+sub_x) + (dd_zone >> 2),
1402 grid->na_c*(cy*GPU_NSUBCELL_Y+sub_y) + (dd_zone & 3),
1409 /* Set the unused atom indices to -1 */
1410 for (ind = na; ind < ncz*grid->na_sc; ind++)
1412 nbs->a[ash+ind] = -1;
1417 /* Determine in which grid column atoms should go */
1418 static void calc_column_indices(nbnxn_grid_t *grid,
1421 int dd_zone, const int *move,
1422 int thread, int nthread,
1429 /* We add one extra cell for particles which moved during DD */
1430 for (i = 0; i < grid->ncx*grid->ncy+1; i++)
1435 n0 = a0 + (int)((thread+0)*(a1 - a0))/nthread;
1436 n1 = a0 + (int)((thread+1)*(a1 - a0))/nthread;
1440 for (i = n0; i < n1; i++)
1442 if (move == NULL || move[i] >= 0)
1444 /* We need to be careful with rounding,
1445 * particles might be a few bits outside the local zone.
1446 * The int cast takes care of the lower bound,
1447 * we will explicitly take care of the upper bound.
1449 cx = (int)((x[i][XX] - grid->c0[XX])*grid->inv_sx);
1450 cy = (int)((x[i][YY] - grid->c0[YY])*grid->inv_sy);
1453 if (cx < 0 || cx > grid->ncx ||
1454 cy < 0 || cy > grid->ncy)
1457 "grid cell cx %d cy %d out of range (max %d %d)\n"
1458 "atom %f %f %f, grid->c0 %f %f",
1459 cx, cy, grid->ncx, grid->ncy,
1460 x[i][XX], x[i][YY], x[i][ZZ], grid->c0[XX], grid->c0[YY]);
1463 /* Take care of potential rouding issues */
1464 cx = min(cx, grid->ncx - 1);
1465 cy = min(cy, grid->ncy - 1);
1467 /* For the moment cell will contain only the, grid local,
1468 * x and y indices, not z.
1470 cell[i] = cx*grid->ncy + cy;
1474 /* Put this moved particle after the end of the grid,
1475 * so we can process it later without using conditionals.
1477 cell[i] = grid->ncx*grid->ncy;
1486 for (i = n0; i < n1; i++)
1488 cx = (int)((x[i][XX] - grid->c0[XX])*grid->inv_sx);
1489 cy = (int)((x[i][YY] - grid->c0[YY])*grid->inv_sy);
1491 /* For non-home zones there could be particles outside
1492 * the non-bonded cut-off range, which have been communicated
1493 * for bonded interactions only. For the result it doesn't
1494 * matter where these end up on the grid. For performance
1495 * we put them in an extra row at the border.
1498 cx = min(cx, grid->ncx - 1);
1500 cy = min(cy, grid->ncy - 1);
1502 /* For the moment cell will contain only the, grid local,
1503 * x and y indices, not z.
1505 cell[i] = cx*grid->ncy + cy;
1512 /* Determine in which grid cells the atoms should go */
1513 static void calc_cell_indices(const nbnxn_search_t nbs,
1520 nbnxn_atomdata_t *nbat)
1523 int cx, cy, cxy, ncz_max, ncz;
1524 int nthread, thread;
1525 int *cxy_na, cxy_na_i;
1527 nthread = gmx_omp_nthreads_get(emntPairsearch);
1529 #pragma omp parallel for num_threads(nthread) schedule(static)
1530 for (thread = 0; thread < nthread; thread++)
1532 calc_column_indices(grid, a0, a1, x, dd_zone, move, thread, nthread,
1533 nbs->cell, nbs->work[thread].cxy_na);
1536 /* Make the cell index as a function of x and y */
1539 grid->cxy_ind[0] = 0;
1540 for (i = 0; i < grid->ncx*grid->ncy+1; i++)
1542 /* We set ncz_max at the beginning of the loop iso at the end
1543 * to skip i=grid->ncx*grid->ncy which are moved particles
1544 * that do not need to be ordered on the grid.
1550 cxy_na_i = nbs->work[0].cxy_na[i];
1551 for (thread = 1; thread < nthread; thread++)
1553 cxy_na_i += nbs->work[thread].cxy_na[i];
1555 ncz = (cxy_na_i + grid->na_sc - 1)/grid->na_sc;
1556 if (nbat->XFormat == nbatX8)
1558 /* Make the number of cell a multiple of 2 */
1559 ncz = (ncz + 1) & ~1;
1561 grid->cxy_ind[i+1] = grid->cxy_ind[i] + ncz;
1562 /* Clear cxy_na, so we can reuse the array below */
1563 grid->cxy_na[i] = 0;
1565 grid->nc = grid->cxy_ind[grid->ncx*grid->ncy] - grid->cxy_ind[0];
1567 nbat->natoms = (grid->cell0 + grid->nc)*grid->na_sc;
1571 fprintf(debug, "ns na_sc %d na_c %d super-cells: %d x %d y %d z %.1f maxz %d\n",
1572 grid->na_sc, grid->na_c, grid->nc,
1573 grid->ncx, grid->ncy, grid->nc/((double)(grid->ncx*grid->ncy)),
1578 for (cy = 0; cy < grid->ncy; cy++)
1580 for (cx = 0; cx < grid->ncx; cx++)
1582 fprintf(debug, " %2d", grid->cxy_ind[i+1]-grid->cxy_ind[i]);
1585 fprintf(debug, "\n");
1590 /* Make sure the work array for sorting is large enough */
1591 if (ncz_max*grid->na_sc*SGSF > nbs->work[0].sort_work_nalloc)
1593 for (thread = 0; thread < nbs->nthread_max; thread++)
1595 nbs->work[thread].sort_work_nalloc =
1596 over_alloc_large(ncz_max*grid->na_sc*SGSF);
1597 srenew(nbs->work[thread].sort_work,
1598 nbs->work[thread].sort_work_nalloc);
1599 /* When not in use, all elements should be -1 */
1600 for (i = 0; i < nbs->work[thread].sort_work_nalloc; i++)
1602 nbs->work[thread].sort_work[i] = -1;
1607 /* Now we know the dimensions we can fill the grid.
1608 * This is the first, unsorted fill. We sort the columns after this.
1610 for (i = a0; i < a1; i++)
1612 /* At this point nbs->cell contains the local grid x,y indices */
1614 nbs->a[(grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc + grid->cxy_na[cxy]++] = i;
1619 /* Set the cell indices for the moved particles */
1620 n0 = grid->nc*grid->na_sc;
1621 n1 = grid->nc*grid->na_sc+grid->cxy_na[grid->ncx*grid->ncy];
1624 for (i = n0; i < n1; i++)
1626 nbs->cell[nbs->a[i]] = i;
1631 /* Sort the super-cell columns along z into the sub-cells. */
1632 #pragma omp parallel for num_threads(nbs->nthread_max) schedule(static)
1633 for (thread = 0; thread < nbs->nthread_max; thread++)
1637 sort_columns_simple(nbs, dd_zone, grid, a0, a1, atinfo, x, nbat,
1638 ((thread+0)*grid->ncx*grid->ncy)/nthread,
1639 ((thread+1)*grid->ncx*grid->ncy)/nthread,
1640 nbs->work[thread].sort_work);
1644 sort_columns_supersub(nbs, dd_zone, grid, a0, a1, atinfo, x, nbat,
1645 ((thread+0)*grid->ncx*grid->ncy)/nthread,
1646 ((thread+1)*grid->ncx*grid->ncy)/nthread,
1647 nbs->work[thread].sort_work);
1651 if (grid->bSimple && nbat->XFormat == nbatX8)
1653 combine_bounding_box_pairs(grid, grid->bb);
1658 grid->nsubc_tot = 0;
1659 for (i = 0; i < grid->nc; i++)
1661 grid->nsubc_tot += grid->nsubc[i];
1669 print_bbsizes_simple(debug, nbs, grid);
1673 fprintf(debug, "ns non-zero sub-cells: %d average atoms %.2f\n",
1674 grid->nsubc_tot, (a1-a0)/(double)grid->nsubc_tot);
1676 print_bbsizes_supersub(debug, nbs, grid);
1681 static void init_buffer_flags(nbnxn_buffer_flags_t *flags,
1686 flags->nflag = (natoms + NBNXN_BUFFERFLAG_SIZE - 1)/NBNXN_BUFFERFLAG_SIZE;
1687 if (flags->nflag > flags->flag_nalloc)
1689 flags->flag_nalloc = over_alloc_large(flags->nflag);
1690 srenew(flags->flag, flags->flag_nalloc);
1692 for (b = 0; b < flags->nflag; b++)
1698 /* Sets up a grid and puts the atoms on the grid.
1699 * This function only operates on one domain of the domain decompostion.
1700 * Note that without domain decomposition there is only one domain.
1702 void nbnxn_put_on_grid(nbnxn_search_t nbs,
1703 int ePBC, matrix box,
1705 rvec corner0, rvec corner1,
1710 int nmoved, int *move,
1712 nbnxn_atomdata_t *nbat)
1716 int nc_max_grid, nc_max;
1718 grid = &nbs->grid[dd_zone];
1720 nbs_cycle_start(&nbs->cc[enbsCCgrid]);
1722 grid->bSimple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
1724 grid->na_c = nbnxn_kernel_to_ci_size(nb_kernel_type);
1725 grid->na_cj = nbnxn_kernel_to_cj_size(nb_kernel_type);
1726 grid->na_sc = (grid->bSimple ? 1 : GPU_NSUBCELL)*grid->na_c;
1727 grid->na_c_2log = get_2log(grid->na_c);
1729 nbat->na_c = grid->na_c;
1738 (nbs->grid[dd_zone-1].cell0 + nbs->grid[dd_zone-1].nc)*
1739 nbs->grid[dd_zone-1].na_sc/grid->na_sc;
1747 copy_mat(box, nbs->box);
1749 if (atom_density >= 0)
1751 grid->atom_density = atom_density;
1755 grid->atom_density = grid_atom_density(n-nmoved, corner0, corner1);
1760 nbs->natoms_local = a1 - nmoved;
1761 /* We assume that nbnxn_put_on_grid is called first
1762 * for the local atoms (dd_zone=0).
1764 nbs->natoms_nonlocal = a1 - nmoved;
1768 nbs->natoms_nonlocal = max(nbs->natoms_nonlocal, a1);
1771 nc_max_grid = set_grid_size_xy(nbs, grid,
1772 dd_zone, n-nmoved, corner0, corner1,
1773 nbs->grid[0].atom_density);
1775 nc_max = grid->cell0 + nc_max_grid;
1777 if (a1 > nbs->cell_nalloc)
1779 nbs->cell_nalloc = over_alloc_large(a1);
1780 srenew(nbs->cell, nbs->cell_nalloc);
1783 /* To avoid conditionals we store the moved particles at the end of a,
1784 * make sure we have enough space.
1786 if (nc_max*grid->na_sc + nmoved > nbs->a_nalloc)
1788 nbs->a_nalloc = over_alloc_large(nc_max*grid->na_sc + nmoved);
1789 srenew(nbs->a, nbs->a_nalloc);
1792 /* We need padding up to a multiple of the buffer flag size: simply add */
1793 if (nc_max*grid->na_sc + NBNXN_BUFFERFLAG_SIZE > nbat->nalloc)
1795 nbnxn_atomdata_realloc(nbat, nc_max*grid->na_sc+NBNXN_BUFFERFLAG_SIZE);
1798 calc_cell_indices(nbs, dd_zone, grid, a0, a1, atinfo, x, move, nbat);
1802 nbat->natoms_local = nbat->natoms;
1805 nbs_cycle_stop(&nbs->cc[enbsCCgrid]);
1808 /* Calls nbnxn_put_on_grid for all non-local domains */
1809 void nbnxn_put_on_grid_nonlocal(nbnxn_search_t nbs,
1810 const gmx_domdec_zones_t *zones,
1814 nbnxn_atomdata_t *nbat)
1819 for (zone = 1; zone < zones->n; zone++)
1821 for (d = 0; d < DIM; d++)
1823 c0[d] = zones->size[zone].bb_x0[d];
1824 c1[d] = zones->size[zone].bb_x1[d];
1827 nbnxn_put_on_grid(nbs, nbs->ePBC, NULL,
1829 zones->cg_range[zone],
1830 zones->cg_range[zone+1],
1840 /* Add simple grid type information to the local super/sub grid */
1841 void nbnxn_grid_add_simple(nbnxn_search_t nbs,
1842 nbnxn_atomdata_t *nbat)
1849 grid = &nbs->grid[0];
1853 gmx_incons("nbnxn_grid_simple called with a simple grid");
1856 ncd = grid->na_sc/NBNXN_CPU_CLUSTER_I_SIZE;
1858 if (grid->nc*ncd > grid->nc_nalloc_simple)
1860 grid->nc_nalloc_simple = over_alloc_large(grid->nc*ncd);
1861 srenew(grid->bbcz_simple, grid->nc_nalloc_simple*NNBSBB_D);
1862 srenew(grid->bb_simple, grid->nc_nalloc_simple);
1863 srenew(grid->flags_simple, grid->nc_nalloc_simple);
1866 sfree_aligned(grid->bbj);
1867 snew_aligned(grid->bbj, grid->nc_nalloc_simple/2, 16);
1871 bbcz = grid->bbcz_simple;
1872 bb = grid->bb_simple;
1874 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntPairsearch)) schedule(static)
1875 for (sc = 0; sc < grid->nc; sc++)
1879 for (c = 0; c < ncd; c++)
1883 na = NBNXN_CPU_CLUSTER_I_SIZE;
1885 nbat->type[tx*NBNXN_CPU_CLUSTER_I_SIZE+na-1] == nbat->ntype-1)
1892 switch (nbat->XFormat)
1895 /* PACK_X4==NBNXN_CPU_CLUSTER_I_SIZE, so this is simple */
1896 calc_bounding_box_x_x4(na, nbat->x+tx*STRIDE_P4,
1900 /* PACK_X8>NBNXN_CPU_CLUSTER_I_SIZE, more complicated */
1901 calc_bounding_box_x_x8(na, nbat->x+X8_IND_A(tx*NBNXN_CPU_CLUSTER_I_SIZE),
1905 calc_bounding_box(na, nbat->xstride,
1906 nbat->x+tx*NBNXN_CPU_CLUSTER_I_SIZE*nbat->xstride,
1910 bbcz[tx*NNBSBB_D+0] = bb[tx].lower[BB_Z];
1911 bbcz[tx*NNBSBB_D+1] = bb[tx].upper[BB_Z];
1913 /* No interaction optimization yet here */
1914 grid->flags_simple[tx] = NBNXN_CI_DO_LJ(0) | NBNXN_CI_DO_COUL(0);
1918 grid->flags_simple[tx] = 0;
1923 if (grid->bSimple && nbat->XFormat == nbatX8)
1925 combine_bounding_box_pairs(grid, grid->bb_simple);
1929 void nbnxn_get_ncells(nbnxn_search_t nbs, int *ncx, int *ncy)
1931 *ncx = nbs->grid[0].ncx;
1932 *ncy = nbs->grid[0].ncy;
1935 void nbnxn_get_atomorder(nbnxn_search_t nbs, int **a, int *n)
1937 const nbnxn_grid_t *grid;
1939 grid = &nbs->grid[0];
1941 /* Return the atom order for the home cell (index 0) */
1944 *n = grid->cxy_ind[grid->ncx*grid->ncy]*grid->na_sc;
1947 void nbnxn_set_atomorder(nbnxn_search_t nbs)
1950 int ao, cx, cy, cxy, cz, j;
1952 /* Set the atom order for the home cell (index 0) */
1953 grid = &nbs->grid[0];
1956 for (cx = 0; cx < grid->ncx; cx++)
1958 for (cy = 0; cy < grid->ncy; cy++)
1960 cxy = cx*grid->ncy + cy;
1961 j = grid->cxy_ind[cxy]*grid->na_sc;
1962 for (cz = 0; cz < grid->cxy_na[cxy]; cz++)
1973 /* Determines the cell range along one dimension that
1974 * the bounding box b0 - b1 sees.
1976 static void get_cell_range(real b0, real b1,
1977 int nc, real c0, real s, real invs,
1978 real d2, real r2, int *cf, int *cl)
1980 *cf = max((int)((b0 - c0)*invs), 0);
1982 while (*cf > 0 && d2 + sqr((b0 - c0) - (*cf-1+1)*s) < r2)
1987 *cl = min((int)((b1 - c0)*invs), nc-1);
1988 while (*cl < nc-1 && d2 + sqr((*cl+1)*s - (b1 - c0)) < r2)
1994 /* Reference code calculating the distance^2 between two bounding boxes */
1995 static float box_dist2(float bx0, float bx1, float by0,
1996 float by1, float bz0, float bz1,
1997 const nbnxn_bb_t *bb)
2000 float dl, dh, dm, dm0;
2004 dl = bx0 - bb->upper[BB_X];
2005 dh = bb->lower[BB_X] - bx1;
2010 dl = by0 - bb->upper[BB_Y];
2011 dh = bb->lower[BB_Y] - by1;
2016 dl = bz0 - bb->upper[BB_Z];
2017 dh = bb->lower[BB_Z] - bz1;
2025 /* Plain C code calculating the distance^2 between two bounding boxes */
2026 static float subc_bb_dist2(int si, const nbnxn_bb_t *bb_i_ci,
2027 int csj, const nbnxn_bb_t *bb_j_all)
2029 const nbnxn_bb_t *bb_i, *bb_j;
2031 float dl, dh, dm, dm0;
2033 bb_i = bb_i_ci + si;
2034 bb_j = bb_j_all + csj;
2038 dl = bb_i->lower[BB_X] - bb_j->upper[BB_X];
2039 dh = bb_j->lower[BB_X] - bb_i->upper[BB_X];
2044 dl = bb_i->lower[BB_Y] - bb_j->upper[BB_Y];
2045 dh = bb_j->lower[BB_Y] - bb_i->upper[BB_Y];
2050 dl = bb_i->lower[BB_Z] - bb_j->upper[BB_Z];
2051 dh = bb_j->lower[BB_Z] - bb_i->upper[BB_Z];
2059 #ifdef NBNXN_SEARCH_BB_SIMD4
2061 /* 4-wide SIMD code for bb distance for bb format xyz0 */
2062 static float subc_bb_dist2_simd4(int si, const nbnxn_bb_t *bb_i_ci,
2063 int csj, const nbnxn_bb_t *bb_j_all)
2065 gmx_simd4_pr bb_i_S0, bb_i_S1;
2066 gmx_simd4_pr bb_j_S0, bb_j_S1;
2072 bb_i_S0 = gmx_simd4_load_bb_pr(&bb_i_ci[si].lower[0]);
2073 bb_i_S1 = gmx_simd4_load_bb_pr(&bb_i_ci[si].upper[0]);
2074 bb_j_S0 = gmx_simd4_load_bb_pr(&bb_j_all[csj].lower[0]);
2075 bb_j_S1 = gmx_simd4_load_bb_pr(&bb_j_all[csj].upper[0]);
2077 dl_S = gmx_simd4_sub_pr(bb_i_S0, bb_j_S1);
2078 dh_S = gmx_simd4_sub_pr(bb_j_S0, bb_i_S1);
2080 dm_S = gmx_simd4_max_pr(dl_S, dh_S);
2081 dm0_S = gmx_simd4_max_pr(dm_S, gmx_simd4_setzero_pr());
2083 return gmx_simd4_dotproduct3(dm0_S, dm0_S);
2086 /* Calculate bb bounding distances of bb_i[si,...,si+3] and store them in d2 */
2087 #define SUBC_BB_DIST2_SIMD4_XXXX_INNER(si, bb_i, d2) \
2091 gmx_simd4_pr dx_0, dy_0, dz_0; \
2092 gmx_simd4_pr dx_1, dy_1, dz_1; \
2094 gmx_simd4_pr mx, my, mz; \
2095 gmx_simd4_pr m0x, m0y, m0z; \
2097 gmx_simd4_pr d2x, d2y, d2z; \
2098 gmx_simd4_pr d2s, d2t; \
2100 shi = si*NNBSBB_D*DIM; \
2102 xi_l = gmx_simd4_load_bb_pr(bb_i+shi+0*STRIDE_PBB); \
2103 yi_l = gmx_simd4_load_bb_pr(bb_i+shi+1*STRIDE_PBB); \
2104 zi_l = gmx_simd4_load_bb_pr(bb_i+shi+2*STRIDE_PBB); \
2105 xi_h = gmx_simd4_load_bb_pr(bb_i+shi+3*STRIDE_PBB); \
2106 yi_h = gmx_simd4_load_bb_pr(bb_i+shi+4*STRIDE_PBB); \
2107 zi_h = gmx_simd4_load_bb_pr(bb_i+shi+5*STRIDE_PBB); \
2109 dx_0 = gmx_simd4_sub_pr(xi_l, xj_h); \
2110 dy_0 = gmx_simd4_sub_pr(yi_l, yj_h); \
2111 dz_0 = gmx_simd4_sub_pr(zi_l, zj_h); \
2113 dx_1 = gmx_simd4_sub_pr(xj_l, xi_h); \
2114 dy_1 = gmx_simd4_sub_pr(yj_l, yi_h); \
2115 dz_1 = gmx_simd4_sub_pr(zj_l, zi_h); \
2117 mx = gmx_simd4_max_pr(dx_0, dx_1); \
2118 my = gmx_simd4_max_pr(dy_0, dy_1); \
2119 mz = gmx_simd4_max_pr(dz_0, dz_1); \
2121 m0x = gmx_simd4_max_pr(mx, zero); \
2122 m0y = gmx_simd4_max_pr(my, zero); \
2123 m0z = gmx_simd4_max_pr(mz, zero); \
2125 d2x = gmx_simd4_mul_pr(m0x, m0x); \
2126 d2y = gmx_simd4_mul_pr(m0y, m0y); \
2127 d2z = gmx_simd4_mul_pr(m0z, m0z); \
2129 d2s = gmx_simd4_add_pr(d2x, d2y); \
2130 d2t = gmx_simd4_add_pr(d2s, d2z); \
2132 gmx_simd4_store_pr(d2+si, d2t); \
2135 /* 4-wide SIMD code for nsi bb distances for bb format xxxxyyyyzzzz */
2136 static void subc_bb_dist2_simd4_xxxx(const float *bb_j,
2137 int nsi, const float *bb_i,
2140 gmx_simd4_pr xj_l, yj_l, zj_l;
2141 gmx_simd4_pr xj_h, yj_h, zj_h;
2142 gmx_simd4_pr xi_l, yi_l, zi_l;
2143 gmx_simd4_pr xi_h, yi_h, zi_h;
2147 zero = gmx_simd4_setzero_pr();
2149 xj_l = gmx_simd4_set1_pr(bb_j[0*STRIDE_PBB]);
2150 yj_l = gmx_simd4_set1_pr(bb_j[1*STRIDE_PBB]);
2151 zj_l = gmx_simd4_set1_pr(bb_j[2*STRIDE_PBB]);
2152 xj_h = gmx_simd4_set1_pr(bb_j[3*STRIDE_PBB]);
2153 yj_h = gmx_simd4_set1_pr(bb_j[4*STRIDE_PBB]);
2154 zj_h = gmx_simd4_set1_pr(bb_j[5*STRIDE_PBB]);
2156 /* Here we "loop" over si (0,STRIDE_PBB) from 0 to nsi with step STRIDE_PBB.
2157 * But as we know the number of iterations is 1 or 2, we unroll manually.
2159 SUBC_BB_DIST2_SIMD4_XXXX_INNER(0, bb_i, d2);
2160 if (STRIDE_PBB < nsi)
2162 SUBC_BB_DIST2_SIMD4_XXXX_INNER(STRIDE_PBB, bb_i, d2);
2166 #endif /* NBNXN_SEARCH_BB_SIMD4 */
2168 /* Plain C function which determines if any atom pair between two cells
2169 * is within distance sqrt(rl2).
2171 static gmx_bool subc_in_range_x(int na_c,
2172 int si, const real *x_i,
2173 int csj, int stride, const real *x_j,
2179 for (i = 0; i < na_c; i++)
2181 i0 = (si*na_c + i)*DIM;
2182 for (j = 0; j < na_c; j++)
2184 j0 = (csj*na_c + j)*stride;
2186 d2 = sqr(x_i[i0 ] - x_j[j0 ]) +
2187 sqr(x_i[i0+1] - x_j[j0+1]) +
2188 sqr(x_i[i0+2] - x_j[j0+2]);
2200 #ifdef NBNXN_SEARCH_SIMD4_FLOAT_X_BB
2201 /* When we make seperate single/double precision SIMD vector operation
2202 * include files, this function should be moved there (also using FMA).
2204 static inline gmx_simd4_pr
2205 gmx_simd4_calc_rsq_pr(gmx_simd4_pr x, gmx_simd4_pr y, gmx_simd4_pr z)
2207 return gmx_simd4_add_pr( gmx_simd4_add_pr( gmx_simd4_mul_pr(x, x), gmx_simd4_mul_pr(y, y) ), gmx_simd4_mul_pr(z, z) );
2210 /* 4-wide SIMD function which determines if any atom pair between two cells,
2211 * both with 8 atoms, is within distance sqrt(rl2).
2212 * Using 8-wide AVX is not faster on Intel Sandy Bridge.
2214 static gmx_bool subc_in_range_simd4(int na_c,
2215 int si, const real *x_i,
2216 int csj, int stride, const real *x_j,
2219 gmx_simd4_pr ix_S0, iy_S0, iz_S0;
2220 gmx_simd4_pr ix_S1, iy_S1, iz_S1;
2227 rc2_S = gmx_simd4_set1_pr(rl2);
2229 dim_stride = NBNXN_GPU_CLUSTER_SIZE/STRIDE_PBB*DIM;
2230 ix_S0 = gmx_simd4_load_bb_pr(x_i+(si*dim_stride+0)*STRIDE_PBB);
2231 iy_S0 = gmx_simd4_load_bb_pr(x_i+(si*dim_stride+1)*STRIDE_PBB);
2232 iz_S0 = gmx_simd4_load_bb_pr(x_i+(si*dim_stride+2)*STRIDE_PBB);
2233 ix_S1 = gmx_simd4_load_bb_pr(x_i+(si*dim_stride+3)*STRIDE_PBB);
2234 iy_S1 = gmx_simd4_load_bb_pr(x_i+(si*dim_stride+4)*STRIDE_PBB);
2235 iz_S1 = gmx_simd4_load_bb_pr(x_i+(si*dim_stride+5)*STRIDE_PBB);
2237 /* We loop from the outer to the inner particles to maximize
2238 * the chance that we find a pair in range quickly and return.
2244 gmx_simd4_pr jx0_S, jy0_S, jz0_S;
2245 gmx_simd4_pr jx1_S, jy1_S, jz1_S;
2247 gmx_simd4_pr dx_S0, dy_S0, dz_S0;
2248 gmx_simd4_pr dx_S1, dy_S1, dz_S1;
2249 gmx_simd4_pr dx_S2, dy_S2, dz_S2;
2250 gmx_simd4_pr dx_S3, dy_S3, dz_S3;
2252 gmx_simd4_pr rsq_S0;
2253 gmx_simd4_pr rsq_S1;
2254 gmx_simd4_pr rsq_S2;
2255 gmx_simd4_pr rsq_S3;
2257 gmx_simd4_pb wco_S0;
2258 gmx_simd4_pb wco_S1;
2259 gmx_simd4_pb wco_S2;
2260 gmx_simd4_pb wco_S3;
2261 gmx_simd4_pb wco_any_S01, wco_any_S23, wco_any_S;
2263 jx0_S = gmx_simd4_set1_pr(x_j[j0*stride+0]);
2264 jy0_S = gmx_simd4_set1_pr(x_j[j0*stride+1]);
2265 jz0_S = gmx_simd4_set1_pr(x_j[j0*stride+2]);
2267 jx1_S = gmx_simd4_set1_pr(x_j[j1*stride+0]);
2268 jy1_S = gmx_simd4_set1_pr(x_j[j1*stride+1]);
2269 jz1_S = gmx_simd4_set1_pr(x_j[j1*stride+2]);
2271 /* Calculate distance */
2272 dx_S0 = gmx_simd4_sub_pr(ix_S0, jx0_S);
2273 dy_S0 = gmx_simd4_sub_pr(iy_S0, jy0_S);
2274 dz_S0 = gmx_simd4_sub_pr(iz_S0, jz0_S);
2275 dx_S1 = gmx_simd4_sub_pr(ix_S1, jx0_S);
2276 dy_S1 = gmx_simd4_sub_pr(iy_S1, jy0_S);
2277 dz_S1 = gmx_simd4_sub_pr(iz_S1, jz0_S);
2278 dx_S2 = gmx_simd4_sub_pr(ix_S0, jx1_S);
2279 dy_S2 = gmx_simd4_sub_pr(iy_S0, jy1_S);
2280 dz_S2 = gmx_simd4_sub_pr(iz_S0, jz1_S);
2281 dx_S3 = gmx_simd4_sub_pr(ix_S1, jx1_S);
2282 dy_S3 = gmx_simd4_sub_pr(iy_S1, jy1_S);
2283 dz_S3 = gmx_simd4_sub_pr(iz_S1, jz1_S);
2285 /* rsq = dx*dx+dy*dy+dz*dz */
2286 rsq_S0 = gmx_simd4_calc_rsq_pr(dx_S0, dy_S0, dz_S0);
2287 rsq_S1 = gmx_simd4_calc_rsq_pr(dx_S1, dy_S1, dz_S1);
2288 rsq_S2 = gmx_simd4_calc_rsq_pr(dx_S2, dy_S2, dz_S2);
2289 rsq_S3 = gmx_simd4_calc_rsq_pr(dx_S3, dy_S3, dz_S3);
2291 wco_S0 = gmx_simd4_cmplt_pr(rsq_S0, rc2_S);
2292 wco_S1 = gmx_simd4_cmplt_pr(rsq_S1, rc2_S);
2293 wco_S2 = gmx_simd4_cmplt_pr(rsq_S2, rc2_S);
2294 wco_S3 = gmx_simd4_cmplt_pr(rsq_S3, rc2_S);
2296 wco_any_S01 = gmx_simd4_or_pb(wco_S0, wco_S1);
2297 wco_any_S23 = gmx_simd4_or_pb(wco_S2, wco_S3);
2298 wco_any_S = gmx_simd4_or_pb(wco_any_S01, wco_any_S23);
2300 if (gmx_simd4_anytrue_pb(wco_any_S))
2314 /* Returns the j sub-cell for index cj_ind */
2315 static int nbl_cj(const nbnxn_pairlist_t *nbl, int cj_ind)
2317 return nbl->cj4[cj_ind >> NBNXN_GPU_JGROUP_SIZE_2LOG].cj[cj_ind & (NBNXN_GPU_JGROUP_SIZE - 1)];
2320 /* Returns the i-interaction mask of the j sub-cell for index cj_ind */
2321 static unsigned nbl_imask0(const nbnxn_pairlist_t *nbl, int cj_ind)
2323 return nbl->cj4[cj_ind >> NBNXN_GPU_JGROUP_SIZE_2LOG].imei[0].imask;
2326 /* Ensures there is enough space for extra extra exclusion masks */
2327 static void check_excl_space(nbnxn_pairlist_t *nbl, int extra)
2329 if (nbl->nexcl+extra > nbl->excl_nalloc)
2331 nbl->excl_nalloc = over_alloc_small(nbl->nexcl+extra);
2332 nbnxn_realloc_void((void **)&nbl->excl,
2333 nbl->nexcl*sizeof(*nbl->excl),
2334 nbl->excl_nalloc*sizeof(*nbl->excl),
2335 nbl->alloc, nbl->free);
2339 /* Ensures there is enough space for ncell extra j-cells in the list */
2340 static void check_subcell_list_space_simple(nbnxn_pairlist_t *nbl,
2345 cj_max = nbl->ncj + ncell;
2347 if (cj_max > nbl->cj_nalloc)
2349 nbl->cj_nalloc = over_alloc_small(cj_max);
2350 nbnxn_realloc_void((void **)&nbl->cj,
2351 nbl->ncj*sizeof(*nbl->cj),
2352 nbl->cj_nalloc*sizeof(*nbl->cj),
2353 nbl->alloc, nbl->free);
2357 /* Ensures there is enough space for ncell extra j-subcells in the list */
2358 static void check_subcell_list_space_supersub(nbnxn_pairlist_t *nbl,
2361 int ncj4_max, j4, j, w, t;
2364 #define WARP_SIZE 32
2366 /* We can have maximally nsupercell*GPU_NSUBCELL sj lists */
2367 /* We can store 4 j-subcell - i-supercell pairs in one struct.
2368 * since we round down, we need one extra entry.
2370 ncj4_max = ((nbl->work->cj_ind + nsupercell*GPU_NSUBCELL + NBNXN_GPU_JGROUP_SIZE - 1) >> NBNXN_GPU_JGROUP_SIZE_2LOG);
2372 if (ncj4_max > nbl->cj4_nalloc)
2374 nbl->cj4_nalloc = over_alloc_small(ncj4_max);
2375 nbnxn_realloc_void((void **)&nbl->cj4,
2376 nbl->work->cj4_init*sizeof(*nbl->cj4),
2377 nbl->cj4_nalloc*sizeof(*nbl->cj4),
2378 nbl->alloc, nbl->free);
2381 if (ncj4_max > nbl->work->cj4_init)
2383 for (j4 = nbl->work->cj4_init; j4 < ncj4_max; j4++)
2385 /* No i-subcells and no excl's in the list initially */
2386 for (w = 0; w < NWARP; w++)
2388 nbl->cj4[j4].imei[w].imask = 0U;
2389 nbl->cj4[j4].imei[w].excl_ind = 0;
2393 nbl->work->cj4_init = ncj4_max;
2397 /* Set all excl masks for one GPU warp no exclusions */
2398 static void set_no_excls(nbnxn_excl_t *excl)
2402 for (t = 0; t < WARP_SIZE; t++)
2404 /* Turn all interaction bits on */
2405 excl->pair[t] = NBNXN_INTERACTION_MASK_ALL;
2409 /* Initializes a single nbnxn_pairlist_t data structure */
2410 static void nbnxn_init_pairlist(nbnxn_pairlist_t *nbl,
2412 nbnxn_alloc_t *alloc,
2417 nbl->alloc = nbnxn_alloc_aligned;
2425 nbl->free = nbnxn_free_aligned;
2432 nbl->bSimple = bSimple;
2443 /* We need one element extra in sj, so alloc initially with 1 */
2444 nbl->cj4_nalloc = 0;
2451 nbl->excl_nalloc = 0;
2453 check_excl_space(nbl, 1);
2455 set_no_excls(&nbl->excl[0]);
2461 snew_aligned(nbl->work->bb_ci, 1, NBNXN_SEARCH_BB_MEM_ALIGN);
2466 snew_aligned(nbl->work->pbb_ci, GPU_NSUBCELL/STRIDE_PBB*NNBSBB_XXXX, NBNXN_SEARCH_BB_MEM_ALIGN);
2468 snew_aligned(nbl->work->bb_ci, GPU_NSUBCELL, NBNXN_SEARCH_BB_MEM_ALIGN);
2471 snew_aligned(nbl->work->x_ci, NBNXN_NA_SC_MAX*DIM, NBNXN_SEARCH_BB_MEM_ALIGN);
2472 #ifdef GMX_NBNXN_SIMD
2473 snew_aligned(nbl->work->x_ci_simd_4xn, 1, NBNXN_MEM_ALIGN);
2474 snew_aligned(nbl->work->x_ci_simd_2xnn, 1, NBNXN_MEM_ALIGN);
2476 snew_aligned(nbl->work->d2, GPU_NSUBCELL, NBNXN_SEARCH_BB_MEM_ALIGN);
2478 nbl->work->sort = NULL;
2479 nbl->work->sort_nalloc = 0;
2480 nbl->work->sci_sort = NULL;
2481 nbl->work->sci_sort_nalloc = 0;
2484 void nbnxn_init_pairlist_set(nbnxn_pairlist_set_t *nbl_list,
2485 gmx_bool bSimple, gmx_bool bCombined,
2486 nbnxn_alloc_t *alloc,
2491 nbl_list->bSimple = bSimple;
2492 nbl_list->bCombined = bCombined;
2494 nbl_list->nnbl = gmx_omp_nthreads_get(emntNonbonded);
2496 if (!nbl_list->bCombined &&
2497 nbl_list->nnbl > NBNXN_BUFFERFLAG_MAX_THREADS)
2499 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.",
2500 nbl_list->nnbl, NBNXN_BUFFERFLAG_MAX_THREADS, NBNXN_BUFFERFLAG_MAX_THREADS);
2503 snew(nbl_list->nbl, nbl_list->nnbl);
2504 /* Execute in order to avoid memory interleaving between threads */
2505 #pragma omp parallel for num_threads(nbl_list->nnbl) schedule(static)
2506 for (i = 0; i < nbl_list->nnbl; i++)
2508 /* Allocate the nblist data structure locally on each thread
2509 * to optimize memory access for NUMA architectures.
2511 snew(nbl_list->nbl[i], 1);
2513 /* Only list 0 is used on the GPU, use normal allocation for i>0 */
2516 nbnxn_init_pairlist(nbl_list->nbl[i], nbl_list->bSimple, alloc, free);
2520 nbnxn_init_pairlist(nbl_list->nbl[i], nbl_list->bSimple, NULL, NULL);
2525 /* Print statistics of a pair list, used for debug output */
2526 static void print_nblist_statistics_simple(FILE *fp, const nbnxn_pairlist_t *nbl,
2527 const nbnxn_search_t nbs, real rl)
2529 const nbnxn_grid_t *grid;
2534 /* This code only produces correct statistics with domain decomposition */
2535 grid = &nbs->grid[0];
2537 fprintf(fp, "nbl nci %d ncj %d\n",
2538 nbl->nci, nbl->ncj);
2539 fprintf(fp, "nbl na_sc %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
2540 nbl->na_sc, rl, nbl->ncj, nbl->ncj/(double)grid->nc,
2541 nbl->ncj/(double)grid->nc*grid->na_sc,
2542 nbl->ncj/(double)grid->nc*grid->na_sc/(0.5*4.0/3.0*M_PI*rl*rl*rl*grid->nc*grid->na_sc/det(nbs->box)));
2544 fprintf(fp, "nbl average j cell list length %.1f\n",
2545 0.25*nbl->ncj/(double)nbl->nci);
2547 for (s = 0; s < SHIFTS; s++)
2552 for (i = 0; i < nbl->nci; i++)
2554 cs[nbl->ci[i].shift & NBNXN_CI_SHIFT] +=
2555 nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start;
2557 j = nbl->ci[i].cj_ind_start;
2558 while (j < nbl->ci[i].cj_ind_end &&
2559 nbl->cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
2565 fprintf(fp, "nbl cell pairs, total: %d excl: %d %.1f%%\n",
2566 nbl->ncj, npexcl, 100*npexcl/(double)nbl->ncj);
2567 for (s = 0; s < SHIFTS; s++)
2571 fprintf(fp, "nbl shift %2d ncj %3d\n", s, cs[s]);
2576 /* Print statistics of a pair lists, used for debug output */
2577 static void print_nblist_statistics_supersub(FILE *fp, const nbnxn_pairlist_t *nbl,
2578 const nbnxn_search_t nbs, real rl)
2580 const nbnxn_grid_t *grid;
2581 int i, j4, j, si, b;
2582 int c[GPU_NSUBCELL+1];
2584 /* This code only produces correct statistics with domain decomposition */
2585 grid = &nbs->grid[0];
2587 fprintf(fp, "nbl nsci %d ncj4 %d nsi %d excl4 %d\n",
2588 nbl->nsci, nbl->ncj4, nbl->nci_tot, nbl->nexcl);
2589 fprintf(fp, "nbl na_c %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
2590 nbl->na_ci, rl, nbl->nci_tot, nbl->nci_tot/(double)grid->nsubc_tot,
2591 nbl->nci_tot/(double)grid->nsubc_tot*grid->na_c,
2592 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/det(nbs->box)));
2594 fprintf(fp, "nbl average j super cell list length %.1f\n",
2595 0.25*nbl->ncj4/(double)nbl->nsci);
2596 fprintf(fp, "nbl average i sub cell list length %.1f\n",
2597 nbl->nci_tot/((double)nbl->ncj4));
2599 for (si = 0; si <= GPU_NSUBCELL; si++)
2603 for (i = 0; i < nbl->nsci; i++)
2605 for (j4 = nbl->sci[i].cj4_ind_start; j4 < nbl->sci[i].cj4_ind_end; j4++)
2607 for (j = 0; j < NBNXN_GPU_JGROUP_SIZE; j++)
2610 for (si = 0; si < GPU_NSUBCELL; si++)
2612 if (nbl->cj4[j4].imei[0].imask & (1U << (j*GPU_NSUBCELL + si)))
2621 for (b = 0; b <= GPU_NSUBCELL; b++)
2623 fprintf(fp, "nbl j-list #i-subcell %d %7d %4.1f\n",
2624 b, c[b], 100.0*c[b]/(double)(nbl->ncj4*NBNXN_GPU_JGROUP_SIZE));
2628 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp */
2629 static void low_get_nbl_exclusions(nbnxn_pairlist_t *nbl, int cj4,
2630 int warp, nbnxn_excl_t **excl)
2632 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
2634 /* No exclusions set, make a new list entry */
2635 nbl->cj4[cj4].imei[warp].excl_ind = nbl->nexcl;
2637 *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
2638 set_no_excls(*excl);
2642 /* We already have some exclusions, new ones can be added to the list */
2643 *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
2647 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp,
2648 * allocates extra memory, if necessary.
2650 static void get_nbl_exclusions_1(nbnxn_pairlist_t *nbl, int cj4,
2651 int warp, nbnxn_excl_t **excl)
2653 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
2655 /* We need to make a new list entry, check if we have space */
2656 check_excl_space(nbl, 1);
2658 low_get_nbl_exclusions(nbl, cj4, warp, excl);
2661 /* Returns pointers to the exclusion mask for cj4-unit cj4 for both warps,
2662 * allocates extra memory, if necessary.
2664 static void get_nbl_exclusions_2(nbnxn_pairlist_t *nbl, int cj4,
2665 nbnxn_excl_t **excl_w0,
2666 nbnxn_excl_t **excl_w1)
2668 /* Check for space we might need */
2669 check_excl_space(nbl, 2);
2671 low_get_nbl_exclusions(nbl, cj4, 0, excl_w0);
2672 low_get_nbl_exclusions(nbl, cj4, 1, excl_w1);
2675 /* Sets the self exclusions i=j and pair exclusions i>j */
2676 static void set_self_and_newton_excls_supersub(nbnxn_pairlist_t *nbl,
2677 int cj4_ind, int sj_offset,
2680 nbnxn_excl_t *excl[2];
2683 /* Here we only set the set self and double pair exclusions */
2685 get_nbl_exclusions_2(nbl, cj4_ind, &excl[0], &excl[1]);
2687 /* Only minor < major bits set */
2688 for (ej = 0; ej < nbl->na_ci; ej++)
2691 for (ei = ej; ei < nbl->na_ci; ei++)
2693 excl[w]->pair[(ej & (NBNXN_GPU_JGROUP_SIZE-1))*nbl->na_ci + ei] &=
2694 ~(1U << (sj_offset*GPU_NSUBCELL + si));
2699 /* Returns a diagonal or off-diagonal interaction mask for plain C lists */
2700 static unsigned int get_imask(gmx_bool rdiag, int ci, int cj)
2702 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
2705 /* Returns a diagonal or off-diagonal interaction mask for cj-size=2 */
2706 static unsigned int get_imask_simd_j2(gmx_bool rdiag, int ci, int cj)
2708 return (rdiag && ci*2 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_0 :
2709 (rdiag && ci*2+1 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_1 :
2710 NBNXN_INTERACTION_MASK_ALL));
2713 /* Returns a diagonal or off-diagonal interaction mask for cj-size=4 */
2714 static unsigned int get_imask_simd_j4(gmx_bool rdiag, int ci, int cj)
2716 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
2719 /* Returns a diagonal or off-diagonal interaction mask for cj-size=8 */
2720 static unsigned int get_imask_simd_j8(gmx_bool rdiag, int ci, int cj)
2722 return (rdiag && ci == cj*2 ? NBNXN_INTERACTION_MASK_DIAG_J8_0 :
2723 (rdiag && ci == cj*2+1 ? NBNXN_INTERACTION_MASK_DIAG_J8_1 :
2724 NBNXN_INTERACTION_MASK_ALL));
2727 #ifdef GMX_NBNXN_SIMD
2728 #if GMX_SIMD_WIDTH_HERE == 2
2729 #define get_imask_simd_4xn get_imask_simd_j2
2731 #if GMX_SIMD_WIDTH_HERE == 4
2732 #define get_imask_simd_4xn get_imask_simd_j4
2734 #if GMX_SIMD_WIDTH_HERE == 8
2735 #define get_imask_simd_4xn get_imask_simd_j8
2736 #define get_imask_simd_2xnn get_imask_simd_j4
2738 #if GMX_SIMD_WIDTH_HERE == 16
2739 #define get_imask_simd_2xnn get_imask_simd_j8
2743 /* Plain C code for making a pair list of cell ci vs cell cjf-cjl.
2744 * Checks bounding box distances and possibly atom pair distances.
2746 static void make_cluster_list_simple(const nbnxn_grid_t *gridj,
2747 nbnxn_pairlist_t *nbl,
2748 int ci, int cjf, int cjl,
2749 gmx_bool remove_sub_diag,
2751 real rl2, float rbb2,
2754 const nbnxn_list_work_t *work;
2756 const nbnxn_bb_t *bb_ci;
2761 int cjf_gl, cjl_gl, cj;
2765 bb_ci = nbl->work->bb_ci;
2766 x_ci = nbl->work->x_ci;
2769 while (!InRange && cjf <= cjl)
2771 d2 = subc_bb_dist2(0, bb_ci, cjf, gridj->bb);
2774 /* Check if the distance is within the distance where
2775 * we use only the bounding box distance rbb,
2776 * or within the cut-off and there is at least one atom pair
2777 * within the cut-off.
2787 cjf_gl = gridj->cell0 + cjf;
2788 for (i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
2790 for (j = 0; j < NBNXN_CPU_CLUSTER_I_SIZE; j++)
2792 InRange = InRange ||
2793 (sqr(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
2794 sqr(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
2795 sqr(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rl2);
2798 *ndistc += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
2811 while (!InRange && cjl > cjf)
2813 d2 = subc_bb_dist2(0, bb_ci, cjl, gridj->bb);
2816 /* Check if the distance is within the distance where
2817 * we use only the bounding box distance rbb,
2818 * or within the cut-off and there is at least one atom pair
2819 * within the cut-off.
2829 cjl_gl = gridj->cell0 + cjl;
2830 for (i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
2832 for (j = 0; j < NBNXN_CPU_CLUSTER_I_SIZE; j++)
2834 InRange = InRange ||
2835 (sqr(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
2836 sqr(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
2837 sqr(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rl2);
2840 *ndistc += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
2850 for (cj = cjf; cj <= cjl; cj++)
2852 /* Store cj and the interaction mask */
2853 nbl->cj[nbl->ncj].cj = gridj->cell0 + cj;
2854 nbl->cj[nbl->ncj].excl = get_imask(remove_sub_diag, ci, cj);
2857 /* Increase the closing index in i super-cell list */
2858 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
2862 #ifdef GMX_NBNXN_SIMD_4XN
2863 #include "nbnxn_search_simd_4xn.h"
2865 #ifdef GMX_NBNXN_SIMD_2XNN
2866 #include "nbnxn_search_simd_2xnn.h"
2869 /* Plain C or SIMD4 code for making a pair list of super-cell sci vs scj.
2870 * Checks bounding box distances and possibly atom pair distances.
2872 static void make_cluster_list_supersub(const nbnxn_grid_t *gridi,
2873 const nbnxn_grid_t *gridj,
2874 nbnxn_pairlist_t *nbl,
2876 gmx_bool sci_equals_scj,
2877 int stride, const real *x,
2878 real rl2, float rbb2,
2883 int cjo, ci1, ci, cj, cj_gl;
2884 int cj4_ind, cj_offset;
2888 const float *pbb_ci;
2890 const nbnxn_bb_t *bb_ci;
2895 #define PRUNE_LIST_CPU_ONE
2896 #ifdef PRUNE_LIST_CPU_ONE
2900 d2l = nbl->work->d2;
2903 pbb_ci = nbl->work->pbb_ci;
2905 bb_ci = nbl->work->bb_ci;
2907 x_ci = nbl->work->x_ci;
2911 for (cjo = 0; cjo < gridj->nsubc[scj]; cjo++)
2913 cj4_ind = (nbl->work->cj_ind >> NBNXN_GPU_JGROUP_SIZE_2LOG);
2914 cj_offset = nbl->work->cj_ind - cj4_ind*NBNXN_GPU_JGROUP_SIZE;
2915 cj4 = &nbl->cj4[cj4_ind];
2917 cj = scj*GPU_NSUBCELL + cjo;
2919 cj_gl = gridj->cell0*GPU_NSUBCELL + cj;
2921 /* Initialize this j-subcell i-subcell list */
2922 cj4->cj[cj_offset] = cj_gl;
2931 ci1 = gridi->nsubc[sci];
2935 /* Determine all ci1 bb distances in one call with SIMD4 */
2936 subc_bb_dist2_simd4_xxxx(gridj->pbb+(cj>>STRIDE_PBB_2LOG)*NNBSBB_XXXX+(cj & (STRIDE_PBB-1)),
2942 /* We use a fixed upper-bound instead of ci1 to help optimization */
2943 for (ci = 0; ci < GPU_NSUBCELL; ci++)
2950 #ifndef NBNXN_BBXXXX
2951 /* Determine the bb distance between ci and cj */
2952 d2l[ci] = subc_bb_dist2(ci, bb_ci, cj, gridj->bb);
2957 #ifdef PRUNE_LIST_CPU_ALL
2958 /* Check if the distance is within the distance where
2959 * we use only the bounding box distance rbb,
2960 * or within the cut-off and there is at least one atom pair
2961 * within the cut-off. This check is very costly.
2963 *ndistc += na_c*na_c;
2966 #ifdef NBNXN_PBB_SIMD4
2971 (na_c, ci, x_ci, cj_gl, stride, x, rl2)))
2973 /* Check if the distance between the two bounding boxes
2974 * in within the pair-list cut-off.
2979 /* Flag this i-subcell to be taken into account */
2980 imask |= (1U << (cj_offset*GPU_NSUBCELL+ci));
2982 #ifdef PRUNE_LIST_CPU_ONE
2990 #ifdef PRUNE_LIST_CPU_ONE
2991 /* If we only found 1 pair, check if any atoms are actually
2992 * within the cut-off, so we could get rid of it.
2994 if (npair == 1 && d2l[ci_last] >= rbb2)
2996 /* Avoid using function pointers here, as it's slower */
2998 #ifdef NBNXN_PBB_SIMD4
2999 !subc_in_range_simd4
3003 (na_c, ci_last, x_ci, cj_gl, stride, x, rl2))
3005 imask &= ~(1U << (cj_offset*GPU_NSUBCELL+ci_last));
3013 /* We have a useful sj entry, close it now */
3015 /* Set the exclucions for the ci== sj entry.
3016 * Here we don't bother to check if this entry is actually flagged,
3017 * as it will nearly always be in the list.
3021 set_self_and_newton_excls_supersub(nbl, cj4_ind, cj_offset, cjo);
3024 /* Copy the cluster interaction mask to the list */
3025 for (w = 0; w < NWARP; w++)
3027 cj4->imei[w].imask |= imask;
3030 nbl->work->cj_ind++;
3032 /* Keep the count */
3033 nbl->nci_tot += npair;
3035 /* Increase the closing index in i super-cell list */
3036 nbl->sci[nbl->nsci].cj4_ind_end =
3037 ((nbl->work->cj_ind+NBNXN_GPU_JGROUP_SIZE-1) >> NBNXN_GPU_JGROUP_SIZE_2LOG);
3042 /* Set all atom-pair exclusions from the topology stored in excl
3043 * as masks in the pair-list for simple list i-entry nbl_ci
3045 static void set_ci_top_excls(const nbnxn_search_t nbs,
3046 nbnxn_pairlist_t *nbl,
3047 gmx_bool diagRemoved,
3050 const nbnxn_ci_t *nbl_ci,
3051 const t_blocka *excl)
3055 int cj_ind_first, cj_ind_last;
3056 int cj_first, cj_last;
3058 int i, ai, aj, si, eind, ge, se;
3059 int found, cj_ind_0, cj_ind_1, cj_ind_m;
3063 nbnxn_excl_t *nbl_excl;
3064 int inner_i, inner_e;
3068 if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
3076 cj_ind_first = nbl_ci->cj_ind_start;
3077 cj_ind_last = nbl->ncj - 1;
3079 cj_first = nbl->cj[cj_ind_first].cj;
3080 cj_last = nbl->cj[cj_ind_last].cj;
3082 /* Determine how many contiguous j-cells we have starting
3083 * from the first i-cell. This number can be used to directly
3084 * calculate j-cell indices for excluded atoms.
3087 if (na_ci_2log == na_cj_2log)
3089 while (cj_ind_first + ndirect <= cj_ind_last &&
3090 nbl->cj[cj_ind_first+ndirect].cj == ci + ndirect)
3095 #ifdef NBNXN_SEARCH_BB_SIMD4
3098 while (cj_ind_first + ndirect <= cj_ind_last &&
3099 nbl->cj[cj_ind_first+ndirect].cj == ci_to_cj(na_cj_2log, ci) + ndirect)
3106 /* Loop over the atoms in the i super-cell */
3107 for (i = 0; i < nbl->na_sc; i++)
3109 ai = nbs->a[ci*nbl->na_sc+i];
3112 si = (i>>na_ci_2log);
3114 /* Loop over the topology-based exclusions for this i-atom */
3115 for (eind = excl->index[ai]; eind < excl->index[ai+1]; eind++)
3121 /* The self exclusion are already set, save some time */
3127 /* Without shifts we only calculate interactions j>i
3128 * for one-way pair-lists.
3130 if (diagRemoved && ge <= ci*nbl->na_sc + i)
3135 se = (ge >> na_cj_2log);
3137 /* Could the cluster se be in our list? */
3138 if (se >= cj_first && se <= cj_last)
3140 if (se < cj_first + ndirect)
3142 /* We can calculate cj_ind directly from se */
3143 found = cj_ind_first + se - cj_first;
3147 /* Search for se using bisection */
3149 cj_ind_0 = cj_ind_first + ndirect;
3150 cj_ind_1 = cj_ind_last + 1;
3151 while (found == -1 && cj_ind_0 < cj_ind_1)
3153 cj_ind_m = (cj_ind_0 + cj_ind_1)>>1;
3155 cj_m = nbl->cj[cj_ind_m].cj;
3163 cj_ind_1 = cj_ind_m;
3167 cj_ind_0 = cj_ind_m + 1;
3174 inner_i = i - (si << na_ci_2log);
3175 inner_e = ge - (se << na_cj_2log);
3177 nbl->cj[found].excl &= ~(1U<<((inner_i<<na_cj_2log) + inner_e));
3178 /* The next code line is usually not needed. We do not want to version
3179 * away the above line, because there is logic that relies on being
3180 * able to detect easily whether any exclusions exist. */
3181 #if (defined GMX_CPU_ACCELERATION_IBM_QPX)
3182 nbl->cj[found].interaction_mask_indices[inner_i] &= ~(1U << inner_e);
3191 /* Set all atom-pair exclusions from the topology stored in excl
3192 * as masks in the pair-list for i-super-cell entry nbl_sci
3194 static void set_sci_top_excls(const nbnxn_search_t nbs,
3195 nbnxn_pairlist_t *nbl,
3196 gmx_bool diagRemoved,
3198 const nbnxn_sci_t *nbl_sci,
3199 const t_blocka *excl)
3204 int cj_ind_first, cj_ind_last;
3205 int cj_first, cj_last;
3207 int i, ai, aj, si, eind, ge, se;
3208 int found, cj_ind_0, cj_ind_1, cj_ind_m;
3212 nbnxn_excl_t *nbl_excl;
3213 int inner_i, inner_e, w;
3219 if (nbl_sci->cj4_ind_end == nbl_sci->cj4_ind_start)
3227 cj_ind_first = nbl_sci->cj4_ind_start*NBNXN_GPU_JGROUP_SIZE;
3228 cj_ind_last = nbl->work->cj_ind - 1;
3230 cj_first = nbl->cj4[nbl_sci->cj4_ind_start].cj[0];
3231 cj_last = nbl_cj(nbl, cj_ind_last);
3233 /* Determine how many contiguous j-clusters we have starting
3234 * from the first i-cluster. This number can be used to directly
3235 * calculate j-cluster indices for excluded atoms.
3238 while (cj_ind_first + ndirect <= cj_ind_last &&
3239 nbl_cj(nbl, cj_ind_first+ndirect) == sci*GPU_NSUBCELL + ndirect)
3244 /* Loop over the atoms in the i super-cell */
3245 for (i = 0; i < nbl->na_sc; i++)
3247 ai = nbs->a[sci*nbl->na_sc+i];
3250 si = (i>>na_c_2log);
3252 /* Loop over the topology-based exclusions for this i-atom */
3253 for (eind = excl->index[ai]; eind < excl->index[ai+1]; eind++)
3259 /* The self exclusion are already set, save some time */
3265 /* Without shifts we only calculate interactions j>i
3266 * for one-way pair-lists.
3268 if (diagRemoved && ge <= sci*nbl->na_sc + i)
3274 /* Could the cluster se be in our list? */
3275 if (se >= cj_first && se <= cj_last)
3277 if (se < cj_first + ndirect)
3279 /* We can calculate cj_ind directly from se */
3280 found = cj_ind_first + se - cj_first;
3284 /* Search for se using bisection */
3286 cj_ind_0 = cj_ind_first + ndirect;
3287 cj_ind_1 = cj_ind_last + 1;
3288 while (found == -1 && cj_ind_0 < cj_ind_1)
3290 cj_ind_m = (cj_ind_0 + cj_ind_1)>>1;
3292 cj_m = nbl_cj(nbl, cj_ind_m);
3300 cj_ind_1 = cj_ind_m;
3304 cj_ind_0 = cj_ind_m + 1;
3311 inner_i = i - si*na_c;
3312 inner_e = ge - se*na_c;
3314 /* Macro for getting the index of atom a within a cluster */
3315 #define AMODCJ4(a) ((a) & (NBNXN_GPU_JGROUP_SIZE - 1))
3316 /* Macro for converting an atom number to a cluster number */
3317 #define A2CJ4(a) ((a) >> NBNXN_GPU_JGROUP_SIZE_2LOG)
3318 /* Macro for getting the index of an i-atom within a warp */
3319 #define AMODWI(a) ((a) & (NBNXN_GPU_CLUSTER_SIZE/2 - 1))
3321 if (nbl_imask0(nbl, found) & (1U << (AMODCJ4(found)*GPU_NSUBCELL + si)))
3325 get_nbl_exclusions_1(nbl, A2CJ4(found), w, &nbl_excl);
3327 nbl_excl->pair[AMODWI(inner_e)*nbl->na_ci+inner_i] &=
3328 ~(1U << (AMODCJ4(found)*GPU_NSUBCELL + si));
3341 /* Reallocate the simple ci list for at least n entries */
3342 static void nb_realloc_ci(nbnxn_pairlist_t *nbl, int n)
3344 nbl->ci_nalloc = over_alloc_small(n);
3345 nbnxn_realloc_void((void **)&nbl->ci,
3346 nbl->nci*sizeof(*nbl->ci),
3347 nbl->ci_nalloc*sizeof(*nbl->ci),
3348 nbl->alloc, nbl->free);
3351 /* Reallocate the super-cell sci list for at least n entries */
3352 static void nb_realloc_sci(nbnxn_pairlist_t *nbl, int n)
3354 nbl->sci_nalloc = over_alloc_small(n);
3355 nbnxn_realloc_void((void **)&nbl->sci,
3356 nbl->nsci*sizeof(*nbl->sci),
3357 nbl->sci_nalloc*sizeof(*nbl->sci),
3358 nbl->alloc, nbl->free);
3361 /* Make a new ci entry at index nbl->nci */
3362 static void new_ci_entry(nbnxn_pairlist_t *nbl, int ci, int shift, int flags)
3364 if (nbl->nci + 1 > nbl->ci_nalloc)
3366 nb_realloc_ci(nbl, nbl->nci+1);
3368 nbl->ci[nbl->nci].ci = ci;
3369 nbl->ci[nbl->nci].shift = shift;
3370 /* Store the interaction flags along with the shift */
3371 nbl->ci[nbl->nci].shift |= flags;
3372 nbl->ci[nbl->nci].cj_ind_start = nbl->ncj;
3373 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
3376 /* Make a new sci entry at index nbl->nsci */
3377 static void new_sci_entry(nbnxn_pairlist_t *nbl, int sci, int shift)
3379 if (nbl->nsci + 1 > nbl->sci_nalloc)
3381 nb_realloc_sci(nbl, nbl->nsci+1);
3383 nbl->sci[nbl->nsci].sci = sci;
3384 nbl->sci[nbl->nsci].shift = shift;
3385 nbl->sci[nbl->nsci].cj4_ind_start = nbl->ncj4;
3386 nbl->sci[nbl->nsci].cj4_ind_end = nbl->ncj4;
3389 /* Sort the simple j-list cj on exclusions.
3390 * Entries with exclusions will all be sorted to the beginning of the list.
3392 static void sort_cj_excl(nbnxn_cj_t *cj, int ncj,
3393 nbnxn_list_work_t *work)
3397 if (ncj > work->cj_nalloc)
3399 work->cj_nalloc = over_alloc_large(ncj);
3400 srenew(work->cj, work->cj_nalloc);
3403 /* Make a list of the j-cells involving exclusions */
3405 for (j = 0; j < ncj; j++)
3407 if (cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
3409 work->cj[jnew++] = cj[j];
3412 /* Check if there are exclusions at all or not just the first entry */
3413 if (!((jnew == 0) ||
3414 (jnew == 1 && cj[0].excl != NBNXN_INTERACTION_MASK_ALL)))
3416 for (j = 0; j < ncj; j++)
3418 if (cj[j].excl == NBNXN_INTERACTION_MASK_ALL)
3420 work->cj[jnew++] = cj[j];
3423 for (j = 0; j < ncj; j++)
3425 cj[j] = work->cj[j];
3430 /* Close this simple list i entry */
3431 static void close_ci_entry_simple(nbnxn_pairlist_t *nbl)
3435 /* All content of the new ci entry have already been filled correctly,
3436 * we only need to increase the count here (for non empty lists).
3438 jlen = nbl->ci[nbl->nci].cj_ind_end - nbl->ci[nbl->nci].cj_ind_start;
3441 sort_cj_excl(nbl->cj+nbl->ci[nbl->nci].cj_ind_start, jlen, nbl->work);
3443 /* The counts below are used for non-bonded pair/flop counts
3444 * and should therefore match the available kernel setups.
3446 if (!(nbl->ci[nbl->nci].shift & NBNXN_CI_DO_COUL(0)))
3448 nbl->work->ncj_noq += jlen;
3450 else if ((nbl->ci[nbl->nci].shift & NBNXN_CI_HALF_LJ(0)) ||
3451 !(nbl->ci[nbl->nci].shift & NBNXN_CI_DO_LJ(0)))
3453 nbl->work->ncj_hlj += jlen;
3460 /* Split sci entry for load balancing on the GPU.
3461 * Splitting ensures we have enough lists to fully utilize the whole GPU.
3462 * With progBal we generate progressively smaller lists, which improves
3463 * load balancing. As we only know the current count on our own thread,
3464 * we will need to estimate the current total amount of i-entries.
3465 * As the lists get concatenated later, this estimate depends
3466 * both on nthread and our own thread index.
3468 static void split_sci_entry(nbnxn_pairlist_t *nbl,
3469 int nsp_max_av, gmx_bool progBal, int nc_bal,
3470 int thread, int nthread)
3474 int cj4_start, cj4_end, j4len, cj4;
3476 int nsp, nsp_sci, nsp_cj4, nsp_cj4_e, nsp_cj4_p;
3481 /* Estimate the total numbers of ci's of the nblist combined
3482 * over all threads using the target number of ci's.
3484 nsci_est = nc_bal*thread/nthread + nbl->nsci;
3486 /* The first ci blocks should be larger, to avoid overhead.
3487 * The last ci blocks should be smaller, to improve load balancing.
3490 nsp_max_av*nc_bal*3/(2*(nsci_est - 1 + nc_bal)));
3494 nsp_max = nsp_max_av;
3497 cj4_start = nbl->sci[nbl->nsci-1].cj4_ind_start;
3498 cj4_end = nbl->sci[nbl->nsci-1].cj4_ind_end;
3499 j4len = cj4_end - cj4_start;
3501 if (j4len > 1 && j4len*GPU_NSUBCELL*NBNXN_GPU_JGROUP_SIZE > nsp_max)
3503 /* Remove the last ci entry and process the cj4's again */
3511 for (cj4 = cj4_start; cj4 < cj4_end; cj4++)
3513 nsp_cj4_p = nsp_cj4;
3514 /* Count the number of cluster pairs in this cj4 group */
3516 for (p = 0; p < GPU_NSUBCELL*NBNXN_GPU_JGROUP_SIZE; p++)
3518 nsp_cj4 += (nbl->cj4[cj4].imei[0].imask >> p) & 1;
3521 if (nsp_cj4 > 0 && nsp + nsp_cj4 > nsp_max)
3523 /* Split the list at cj4 */
3524 nbl->sci[sci].cj4_ind_end = cj4;
3525 /* Create a new sci entry */
3528 if (nbl->nsci+1 > nbl->sci_nalloc)
3530 nb_realloc_sci(nbl, nbl->nsci+1);
3532 nbl->sci[sci].sci = nbl->sci[nbl->nsci-1].sci;
3533 nbl->sci[sci].shift = nbl->sci[nbl->nsci-1].shift;
3534 nbl->sci[sci].cj4_ind_start = cj4;
3536 nsp_cj4_e = nsp_cj4_p;
3542 /* Put the remaining cj4's in the last sci entry */
3543 nbl->sci[sci].cj4_ind_end = cj4_end;
3545 /* Possibly balance out the last two sci's
3546 * by moving the last cj4 of the second last sci.
3548 if (nsp_sci - nsp_cj4_e >= nsp + nsp_cj4_e)
3550 nbl->sci[sci-1].cj4_ind_end--;
3551 nbl->sci[sci].cj4_ind_start--;
3558 /* Clost this super/sub list i entry */
3559 static void close_ci_entry_supersub(nbnxn_pairlist_t *nbl,
3561 gmx_bool progBal, int nc_bal,
3562 int thread, int nthread)
3567 /* All content of the new ci entry have already been filled correctly,
3568 * we only need to increase the count here (for non empty lists).
3570 j4len = nbl->sci[nbl->nsci].cj4_ind_end - nbl->sci[nbl->nsci].cj4_ind_start;
3573 /* We can only have complete blocks of 4 j-entries in a list,
3574 * so round the count up before closing.
3576 nbl->ncj4 = ((nbl->work->cj_ind + NBNXN_GPU_JGROUP_SIZE - 1) >> NBNXN_GPU_JGROUP_SIZE_2LOG);
3577 nbl->work->cj_ind = nbl->ncj4*NBNXN_GPU_JGROUP_SIZE;
3583 /* Measure the size of the new entry and potentially split it */
3584 split_sci_entry(nbl, nsp_max_av, progBal, nc_bal, thread, nthread);
3589 /* Syncs the working array before adding another grid pair to the list */
3590 static void sync_work(nbnxn_pairlist_t *nbl)
3594 nbl->work->cj_ind = nbl->ncj4*NBNXN_GPU_JGROUP_SIZE;
3595 nbl->work->cj4_init = nbl->ncj4;
3599 /* Clears an nbnxn_pairlist_t data structure */
3600 static void clear_pairlist(nbnxn_pairlist_t *nbl)
3609 nbl->work->ncj_noq = 0;
3610 nbl->work->ncj_hlj = 0;
3613 /* Sets a simple list i-cell bounding box, including PBC shift */
3614 static gmx_inline void set_icell_bb_simple(const nbnxn_bb_t *bb, int ci,
3615 real shx, real shy, real shz,
3618 bb_ci->lower[BB_X] = bb[ci].lower[BB_X] + shx;
3619 bb_ci->lower[BB_Y] = bb[ci].lower[BB_Y] + shy;
3620 bb_ci->lower[BB_Z] = bb[ci].lower[BB_Z] + shz;
3621 bb_ci->upper[BB_X] = bb[ci].upper[BB_X] + shx;
3622 bb_ci->upper[BB_Y] = bb[ci].upper[BB_Y] + shy;
3623 bb_ci->upper[BB_Z] = bb[ci].upper[BB_Z] + shz;
3627 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
3628 static void set_icell_bbxxxx_supersub(const float *bb, int ci,
3629 real shx, real shy, real shz,
3634 ia = ci*(GPU_NSUBCELL>>STRIDE_PBB_2LOG)*NNBSBB_XXXX;
3635 for (m = 0; m < (GPU_NSUBCELL>>STRIDE_PBB_2LOG)*NNBSBB_XXXX; m += NNBSBB_XXXX)
3637 for (i = 0; i < STRIDE_PBB; i++)
3639 bb_ci[m+0*STRIDE_PBB+i] = bb[ia+m+0*STRIDE_PBB+i] + shx;
3640 bb_ci[m+1*STRIDE_PBB+i] = bb[ia+m+1*STRIDE_PBB+i] + shy;
3641 bb_ci[m+2*STRIDE_PBB+i] = bb[ia+m+2*STRIDE_PBB+i] + shz;
3642 bb_ci[m+3*STRIDE_PBB+i] = bb[ia+m+3*STRIDE_PBB+i] + shx;
3643 bb_ci[m+4*STRIDE_PBB+i] = bb[ia+m+4*STRIDE_PBB+i] + shy;
3644 bb_ci[m+5*STRIDE_PBB+i] = bb[ia+m+5*STRIDE_PBB+i] + shz;
3650 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
3651 static void set_icell_bb_supersub(const nbnxn_bb_t *bb, int ci,
3652 real shx, real shy, real shz,
3657 for (i = 0; i < GPU_NSUBCELL; i++)
3659 set_icell_bb_simple(bb, ci*GPU_NSUBCELL+i,
3665 /* Copies PBC shifted i-cell atom coordinates x,y,z to working array */
3666 static void icell_set_x_simple(int ci,
3667 real shx, real shy, real shz,
3668 int gmx_unused na_c,
3669 int stride, const real *x,
3670 nbnxn_list_work_t *work)
3674 ia = ci*NBNXN_CPU_CLUSTER_I_SIZE;
3676 for (i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE; i++)
3678 work->x_ci[i*STRIDE_XYZ+XX] = x[(ia+i)*stride+XX] + shx;
3679 work->x_ci[i*STRIDE_XYZ+YY] = x[(ia+i)*stride+YY] + shy;
3680 work->x_ci[i*STRIDE_XYZ+ZZ] = x[(ia+i)*stride+ZZ] + shz;
3684 /* Copies PBC shifted super-cell atom coordinates x,y,z to working array */
3685 static void icell_set_x_supersub(int ci,
3686 real shx, real shy, real shz,
3688 int stride, const real *x,
3689 nbnxn_list_work_t *work)
3696 ia = ci*GPU_NSUBCELL*na_c;
3697 for (i = 0; i < GPU_NSUBCELL*na_c; i++)
3699 x_ci[i*DIM + XX] = x[(ia+i)*stride + XX] + shx;
3700 x_ci[i*DIM + YY] = x[(ia+i)*stride + YY] + shy;
3701 x_ci[i*DIM + ZZ] = x[(ia+i)*stride + ZZ] + shz;
3705 #ifdef NBNXN_SEARCH_BB_SIMD4
3706 /* Copies PBC shifted super-cell packed atom coordinates to working array */
3707 static void icell_set_x_supersub_simd4(int ci,
3708 real shx, real shy, real shz,
3710 int stride, const real *x,
3711 nbnxn_list_work_t *work)
3713 int si, io, ia, i, j;
3718 for (si = 0; si < GPU_NSUBCELL; si++)
3720 for (i = 0; i < na_c; i += STRIDE_PBB)
3723 ia = ci*GPU_NSUBCELL*na_c + io;
3724 for (j = 0; j < STRIDE_PBB; j++)
3726 x_ci[io*DIM + j + XX*STRIDE_PBB] = x[(ia+j)*stride+XX] + shx;
3727 x_ci[io*DIM + j + YY*STRIDE_PBB] = x[(ia+j)*stride+YY] + shy;
3728 x_ci[io*DIM + j + ZZ*STRIDE_PBB] = x[(ia+j)*stride+ZZ] + shz;
3735 static real nbnxn_rlist_inc_nonloc_fac = 0.6;
3737 /* Due to the cluster size the effective pair-list is longer than
3738 * that of a simple atom pair-list. This function gives the extra distance.
3740 real nbnxn_get_rlist_effective_inc(int cluster_size, real atom_density)
3742 return ((0.5 + nbnxn_rlist_inc_nonloc_fac)*sqr(((cluster_size) - 1.0)/(cluster_size))*pow((cluster_size)/(atom_density), 1.0/3.0));
3745 /* Estimates the interaction volume^2 for non-local interactions */
3746 static real nonlocal_vol2(const gmx_domdec_zones_t *zones, rvec ls, real r)
3755 /* Here we simply add up the volumes of 1, 2 or 3 1D decomposition
3756 * not home interaction volume^2. As these volumes are not additive,
3757 * this is an overestimate, but it would only be significant in the limit
3758 * of small cells, where we anyhow need to split the lists into
3759 * as small parts as possible.
3762 for (z = 0; z < zones->n; z++)
3764 if (zones->shift[z][XX] + zones->shift[z][YY] + zones->shift[z][ZZ] == 1)
3769 for (d = 0; d < DIM; d++)
3771 if (zones->shift[z][d] == 0)
3775 za *= zones->size[z].x1[d] - zones->size[z].x0[d];
3779 /* 4 octants of a sphere */
3780 vold_est = 0.25*M_PI*r*r*r*r;
3781 /* 4 quarter pie slices on the edges */
3782 vold_est += 4*cl*M_PI/6.0*r*r*r;
3783 /* One rectangular volume on a face */
3784 vold_est += ca*0.5*r*r;
3786 vol2_est_tot += vold_est*za;
3790 return vol2_est_tot;
3793 /* Estimates the average size of a full j-list for super/sub setup */
3794 static int get_nsubpair_max(const nbnxn_search_t nbs,
3797 int min_ci_balanced)
3799 const nbnxn_grid_t *grid;
3801 real xy_diag2, r_eff_sup, vol_est, nsp_est, nsp_est_nl;
3804 grid = &nbs->grid[0];
3806 ls[XX] = (grid->c1[XX] - grid->c0[XX])/(grid->ncx*GPU_NSUBCELL_X);
3807 ls[YY] = (grid->c1[YY] - grid->c0[YY])/(grid->ncy*GPU_NSUBCELL_Y);
3808 ls[ZZ] = (grid->c1[ZZ] - grid->c0[ZZ])*grid->ncx*grid->ncy/(grid->nc*GPU_NSUBCELL_Z);
3810 /* The average squared length of the diagonal of a sub cell */
3811 xy_diag2 = ls[XX]*ls[XX] + ls[YY]*ls[YY] + ls[ZZ]*ls[ZZ];
3813 /* The formulas below are a heuristic estimate of the average nsj per si*/
3814 r_eff_sup = rlist + nbnxn_rlist_inc_nonloc_fac*sqr((grid->na_c - 1.0)/grid->na_c)*sqrt(xy_diag2/3);
3816 if (!nbs->DomDec || nbs->zones->n == 1)
3823 sqr(grid->atom_density/grid->na_c)*
3824 nonlocal_vol2(nbs->zones, ls, r_eff_sup);
3829 /* Sub-cell interacts with itself */
3830 vol_est = ls[XX]*ls[YY]*ls[ZZ];
3831 /* 6/2 rectangular volume on the faces */
3832 vol_est += (ls[XX]*ls[YY] + ls[XX]*ls[ZZ] + ls[YY]*ls[ZZ])*r_eff_sup;
3833 /* 12/2 quarter pie slices on the edges */
3834 vol_est += 2*(ls[XX] + ls[YY] + ls[ZZ])*0.25*M_PI*sqr(r_eff_sup);
3835 /* 4 octants of a sphere */
3836 vol_est += 0.5*4.0/3.0*M_PI*pow(r_eff_sup, 3);
3838 nsp_est = grid->nsubc_tot*vol_est*grid->atom_density/grid->na_c;
3840 /* Subtract the non-local pair count */
3841 nsp_est -= nsp_est_nl;
3845 fprintf(debug, "nsp_est local %5.1f non-local %5.1f\n",
3846 nsp_est, nsp_est_nl);
3851 nsp_est = nsp_est_nl;
3854 if (min_ci_balanced <= 0 || grid->nc >= min_ci_balanced || grid->nc == 0)
3856 /* We don't need to worry */
3861 /* Thus the (average) maximum j-list size should be as follows */
3862 nsubpair_max = max(1, (int)(nsp_est/min_ci_balanced+0.5));
3864 /* Since the target value is a maximum (this avoids high outliers,
3865 * which lead to load imbalance), not average, we add half the
3866 * number of pairs in a cj4 block to get the average about right.
3868 nsubpair_max += GPU_NSUBCELL*NBNXN_GPU_JGROUP_SIZE/2;
3873 fprintf(debug, "nbl nsp estimate %.1f, nsubpair_max %d\n",
3874 nsp_est, nsubpair_max);
3877 return nsubpair_max;
3880 /* Debug list print function */
3881 static void print_nblist_ci_cj(FILE *fp, const nbnxn_pairlist_t *nbl)
3885 for (i = 0; i < nbl->nci; i++)
3887 fprintf(fp, "ci %4d shift %2d ncj %3d\n",
3888 nbl->ci[i].ci, nbl->ci[i].shift,
3889 nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start);
3891 for (j = nbl->ci[i].cj_ind_start; j < nbl->ci[i].cj_ind_end; j++)
3893 fprintf(fp, " cj %5d imask %x\n",
3900 /* Debug list print function */
3901 static void print_nblist_sci_cj(FILE *fp, const nbnxn_pairlist_t *nbl)
3903 int i, j4, j, ncp, si;
3905 for (i = 0; i < nbl->nsci; i++)
3907 fprintf(fp, "ci %4d shift %2d ncj4 %2d\n",
3908 nbl->sci[i].sci, nbl->sci[i].shift,
3909 nbl->sci[i].cj4_ind_end - nbl->sci[i].cj4_ind_start);
3912 for (j4 = nbl->sci[i].cj4_ind_start; j4 < nbl->sci[i].cj4_ind_end; j4++)
3914 for (j = 0; j < NBNXN_GPU_JGROUP_SIZE; j++)
3916 fprintf(fp, " sj %5d imask %x\n",
3918 nbl->cj4[j4].imei[0].imask);
3919 for (si = 0; si < GPU_NSUBCELL; si++)
3921 if (nbl->cj4[j4].imei[0].imask & (1U << (j*GPU_NSUBCELL + si)))
3928 fprintf(fp, "ci %4d shift %2d ncj4 %2d ncp %3d\n",
3929 nbl->sci[i].sci, nbl->sci[i].shift,
3930 nbl->sci[i].cj4_ind_end - nbl->sci[i].cj4_ind_start,
3935 /* Combine pair lists *nbl generated on multiple threads nblc */
3936 static void combine_nblists(int nnbl, nbnxn_pairlist_t **nbl,
3937 nbnxn_pairlist_t *nblc)
3939 int nsci, ncj4, nexcl;
3944 gmx_incons("combine_nblists does not support simple lists");
3949 nexcl = nblc->nexcl;
3950 for (i = 0; i < nnbl; i++)
3952 nsci += nbl[i]->nsci;
3953 ncj4 += nbl[i]->ncj4;
3954 nexcl += nbl[i]->nexcl;
3957 if (nsci > nblc->sci_nalloc)
3959 nb_realloc_sci(nblc, nsci);
3961 if (ncj4 > nblc->cj4_nalloc)
3963 nblc->cj4_nalloc = over_alloc_small(ncj4);
3964 nbnxn_realloc_void((void **)&nblc->cj4,
3965 nblc->ncj4*sizeof(*nblc->cj4),
3966 nblc->cj4_nalloc*sizeof(*nblc->cj4),
3967 nblc->alloc, nblc->free);
3969 if (nexcl > nblc->excl_nalloc)
3971 nblc->excl_nalloc = over_alloc_small(nexcl);
3972 nbnxn_realloc_void((void **)&nblc->excl,
3973 nblc->nexcl*sizeof(*nblc->excl),
3974 nblc->excl_nalloc*sizeof(*nblc->excl),
3975 nblc->alloc, nblc->free);
3978 /* Each thread should copy its own data to the combined arrays,
3979 * as otherwise data will go back and forth between different caches.
3981 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntPairsearch)) schedule(static)
3982 for (n = 0; n < nnbl; n++)
3989 const nbnxn_pairlist_t *nbli;
3991 /* Determine the offset in the combined data for our thread */
3992 sci_offset = nblc->nsci;
3993 cj4_offset = nblc->ncj4;
3994 ci_offset = nblc->nci_tot;
3995 excl_offset = nblc->nexcl;
3997 for (i = 0; i < n; i++)
3999 sci_offset += nbl[i]->nsci;
4000 cj4_offset += nbl[i]->ncj4;
4001 ci_offset += nbl[i]->nci_tot;
4002 excl_offset += nbl[i]->nexcl;
4007 for (i = 0; i < nbli->nsci; i++)
4009 nblc->sci[sci_offset+i] = nbli->sci[i];
4010 nblc->sci[sci_offset+i].cj4_ind_start += cj4_offset;
4011 nblc->sci[sci_offset+i].cj4_ind_end += cj4_offset;
4014 for (j4 = 0; j4 < nbli->ncj4; j4++)
4016 nblc->cj4[cj4_offset+j4] = nbli->cj4[j4];
4017 nblc->cj4[cj4_offset+j4].imei[0].excl_ind += excl_offset;
4018 nblc->cj4[cj4_offset+j4].imei[1].excl_ind += excl_offset;
4021 for (j4 = 0; j4 < nbli->nexcl; j4++)
4023 nblc->excl[excl_offset+j4] = nbli->excl[j4];
4027 for (n = 0; n < nnbl; n++)
4029 nblc->nsci += nbl[n]->nsci;
4030 nblc->ncj4 += nbl[n]->ncj4;
4031 nblc->nci_tot += nbl[n]->nci_tot;
4032 nblc->nexcl += nbl[n]->nexcl;
4036 /* Returns the next ci to be processes by our thread */
4037 static gmx_bool next_ci(const nbnxn_grid_t *grid,
4039 int nth, int ci_block,
4040 int *ci_x, int *ci_y,
4046 if (*ci_b == ci_block)
4048 /* Jump to the next block assigned to this task */
4049 *ci += (nth - 1)*ci_block;
4053 if (*ci >= grid->nc*conv)
4058 while (*ci >= grid->cxy_ind[*ci_x*grid->ncy + *ci_y + 1]*conv)
4061 if (*ci_y == grid->ncy)
4071 /* Returns the distance^2 for which we put cell pairs in the list
4072 * without checking atom pair distances. This is usually < rlist^2.
4074 static float boundingbox_only_distance2(const nbnxn_grid_t *gridi,
4075 const nbnxn_grid_t *gridj,
4079 /* If the distance between two sub-cell bounding boxes is less
4080 * than this distance, do not check the distance between
4081 * all particle pairs in the sub-cell, since then it is likely
4082 * that the box pair has atom pairs within the cut-off.
4083 * We use the nblist cut-off minus 0.5 times the average x/y diagonal
4084 * spacing of the sub-cells. Around 40% of the checked pairs are pruned.
4085 * Using more than 0.5 gains at most 0.5%.
4086 * If forces are calculated more than twice, the performance gain
4087 * in the force calculation outweighs the cost of checking.
4088 * Note that with subcell lists, the atom-pair distance check
4089 * is only performed when only 1 out of 8 sub-cells in within range,
4090 * this is because the GPU is much faster than the cpu.
4095 bbx = 0.5*(gridi->sx + gridj->sx);
4096 bby = 0.5*(gridi->sy + gridj->sy);
4099 bbx /= GPU_NSUBCELL_X;
4100 bby /= GPU_NSUBCELL_Y;
4103 rbb2 = sqr(max(0, rlist - 0.5*sqrt(bbx*bbx + bby*bby)));
4108 return (float)((1+GMX_FLOAT_EPS)*rbb2);
4112 static int get_ci_block_size(const nbnxn_grid_t *gridi,
4113 gmx_bool bDomDec, int nth)
4115 const int ci_block_enum = 5;
4116 const int ci_block_denom = 11;
4117 const int ci_block_min_atoms = 16;
4120 /* Here we decide how to distribute the blocks over the threads.
4121 * We use prime numbers to try to avoid that the grid size becomes
4122 * a multiple of the number of threads, which would lead to some
4123 * threads getting "inner" pairs and others getting boundary pairs,
4124 * which in turns will lead to load imbalance between threads.
4125 * Set the block size as 5/11/ntask times the average number of cells
4126 * in a y,z slab. This should ensure a quite uniform distribution
4127 * of the grid parts of the different thread along all three grid
4128 * zone boundaries with 3D domain decomposition. At the same time
4129 * the blocks will not become too small.
4131 ci_block = (gridi->nc*ci_block_enum)/(ci_block_denom*gridi->ncx*nth);
4133 /* Ensure the blocks are not too small: avoids cache invalidation */
4134 if (ci_block*gridi->na_sc < ci_block_min_atoms)
4136 ci_block = (ci_block_min_atoms + gridi->na_sc - 1)/gridi->na_sc;
4139 /* Without domain decomposition
4140 * or with less than 3 blocks per task, divide in nth blocks.
4142 if (!bDomDec || ci_block*3*nth > gridi->nc)
4144 ci_block = (gridi->nc + nth - 1)/nth;
4150 /* Generates the part of pair-list nbl assigned to our thread */
4151 static void nbnxn_make_pairlist_part(const nbnxn_search_t nbs,
4152 const nbnxn_grid_t *gridi,
4153 const nbnxn_grid_t *gridj,
4154 nbnxn_search_work_t *work,
4155 const nbnxn_atomdata_t *nbat,
4156 const t_blocka *excl,
4160 gmx_bool bFBufferFlag,
4163 int min_ci_balanced,
4165 nbnxn_pairlist_t *nbl)
4172 int ci_b, ci, ci_x, ci_y, ci_xy, cj;
4178 int conv_i, cell0_i;
4179 const nbnxn_bb_t *bb_i = NULL;
4181 const float *pbb_i = NULL;
4183 const float *bbcz_i, *bbcz_j;
4185 real bx0, bx1, by0, by1, bz0, bz1;
4187 real d2cx, d2z, d2z_cx, d2z_cy, d2zx, d2zxy, d2xy;
4188 int cxf, cxl, cyf, cyf_x, cyl;
4190 int c0, c1, cs, cf, cl;
4193 int gridi_flag_shift = 0, gridj_flag_shift = 0;
4194 unsigned *gridj_flag = NULL;
4195 int ncj_old_i, ncj_old_j;
4197 nbs_cycle_start(&work->cc[enbsCCsearch]);
4199 if (gridj->bSimple != nbl->bSimple)
4201 gmx_incons("Grid incompatible with pair-list");
4205 nbl->na_sc = gridj->na_sc;
4206 nbl->na_ci = gridj->na_c;
4207 nbl->na_cj = nbnxn_kernel_to_cj_size(nb_kernel_type);
4208 na_cj_2log = get_2log(nbl->na_cj);
4214 /* Determine conversion of clusters to flag blocks */
4215 gridi_flag_shift = 0;
4216 while ((nbl->na_ci<<gridi_flag_shift) < NBNXN_BUFFERFLAG_SIZE)
4220 gridj_flag_shift = 0;
4221 while ((nbl->na_cj<<gridj_flag_shift) < NBNXN_BUFFERFLAG_SIZE)
4226 gridj_flag = work->buffer_flags.flag;
4229 copy_mat(nbs->box, box);
4231 rl2 = nbl->rlist*nbl->rlist;
4233 rbb2 = boundingbox_only_distance2(gridi, gridj, nbl->rlist, nbl->bSimple);
4237 fprintf(debug, "nbl bounding box only distance %f\n", sqrt(rbb2));
4240 /* Set the shift range */
4241 for (d = 0; d < DIM; d++)
4243 /* Check if we need periodicity shifts.
4244 * Without PBC or with domain decomposition we don't need them.
4246 if (d >= ePBC2npbcdim(nbs->ePBC) || nbs->dd_dim[d])
4253 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
4264 if (nbl->bSimple && !gridi->bSimple)
4266 conv_i = gridi->na_sc/gridj->na_sc;
4267 bb_i = gridi->bb_simple;
4268 bbcz_i = gridi->bbcz_simple;
4269 flags_i = gridi->flags_simple;
4284 /* We use the normal bounding box format for both grid types */
4287 bbcz_i = gridi->bbcz;
4288 flags_i = gridi->flags;
4290 cell0_i = gridi->cell0*conv_i;
4292 bbcz_j = gridj->bbcz;
4296 /* Blocks of the conversion factor - 1 give a large repeat count
4297 * combined with a small block size. This should result in good
4298 * load balancing for both small and large domains.
4300 ci_block = conv_i - 1;
4304 fprintf(debug, "nbl nc_i %d col.av. %.1f ci_block %d\n",
4305 gridi->nc, gridi->nc/(double)(gridi->ncx*gridi->ncy), ci_block);
4311 /* Initially ci_b and ci to 1 before where we want them to start,
4312 * as they will both be incremented in next_ci.
4315 ci = th*ci_block - 1;
4318 while (next_ci(gridi, conv_i, nth, ci_block, &ci_x, &ci_y, &ci_b, &ci))
4320 if (nbl->bSimple && flags_i[ci] == 0)
4325 ncj_old_i = nbl->ncj;
4328 if (gridj != gridi && shp[XX] == 0)
4332 bx1 = bb_i[ci].upper[BB_X];
4336 bx1 = gridi->c0[XX] + (ci_x+1)*gridi->sx;
4338 if (bx1 < gridj->c0[XX])
4340 d2cx = sqr(gridj->c0[XX] - bx1);
4349 ci_xy = ci_x*gridi->ncy + ci_y;
4351 /* Loop over shift vectors in three dimensions */
4352 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
4354 shz = tz*box[ZZ][ZZ];
4356 bz0 = bbcz_i[ci*NNBSBB_D ] + shz;
4357 bz1 = bbcz_i[ci*NNBSBB_D+1] + shz;
4369 d2z = sqr(bz0 - box[ZZ][ZZ]);
4372 d2z_cx = d2z + d2cx;
4380 bz1/((real)(gridi->cxy_ind[ci_xy+1] - gridi->cxy_ind[ci_xy]));
4385 /* The check with bz1_frac close to or larger than 1 comes later */
4387 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
4389 shy = ty*box[YY][YY] + tz*box[ZZ][YY];
4393 by0 = bb_i[ci].lower[BB_Y] + shy;
4394 by1 = bb_i[ci].upper[BB_Y] + shy;
4398 by0 = gridi->c0[YY] + (ci_y )*gridi->sy + shy;
4399 by1 = gridi->c0[YY] + (ci_y+1)*gridi->sy + shy;
4402 get_cell_range(by0, by1,
4403 gridj->ncy, gridj->c0[YY], gridj->sy, gridj->inv_sy,
4413 if (by1 < gridj->c0[YY])
4415 d2z_cy += sqr(gridj->c0[YY] - by1);
4417 else if (by0 > gridj->c1[YY])
4419 d2z_cy += sqr(by0 - gridj->c1[YY]);
4422 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
4424 shift = XYZ2IS(tx, ty, tz);
4426 #ifdef NBNXN_SHIFT_BACKWARD
4427 if (gridi == gridj && shift > CENTRAL)
4433 shx = tx*box[XX][XX] + ty*box[YY][XX] + tz*box[ZZ][XX];
4437 bx0 = bb_i[ci].lower[BB_X] + shx;
4438 bx1 = bb_i[ci].upper[BB_X] + shx;
4442 bx0 = gridi->c0[XX] + (ci_x )*gridi->sx + shx;
4443 bx1 = gridi->c0[XX] + (ci_x+1)*gridi->sx + shx;
4446 get_cell_range(bx0, bx1,
4447 gridj->ncx, gridj->c0[XX], gridj->sx, gridj->inv_sx,
4458 new_ci_entry(nbl, cell0_i+ci, shift, flags_i[ci]);
4462 new_sci_entry(nbl, cell0_i+ci, shift);
4465 #ifndef NBNXN_SHIFT_BACKWARD
4468 if (shift == CENTRAL && gridi == gridj &&
4472 /* Leave the pairs with i > j.
4473 * x is the major index, so skip half of it.
4480 set_icell_bb_simple(bb_i, ci, shx, shy, shz,
4486 set_icell_bbxxxx_supersub(pbb_i, ci, shx, shy, shz,
4489 set_icell_bb_supersub(bb_i, ci, shx, shy, shz,
4494 nbs->icell_set_x(cell0_i+ci, shx, shy, shz,
4495 gridi->na_c, nbat->xstride, nbat->x,
4498 for (cx = cxf; cx <= cxl; cx++)
4501 if (gridj->c0[XX] + cx*gridj->sx > bx1)
4503 d2zx += sqr(gridj->c0[XX] + cx*gridj->sx - bx1);
4505 else if (gridj->c0[XX] + (cx+1)*gridj->sx < bx0)
4507 d2zx += sqr(gridj->c0[XX] + (cx+1)*gridj->sx - bx0);
4510 #ifndef NBNXN_SHIFT_BACKWARD
4511 if (gridi == gridj &&
4512 cx == 0 && cyf < ci_y)
4514 if (gridi == gridj &&
4515 cx == 0 && shift == CENTRAL && cyf < ci_y)
4518 /* Leave the pairs with i > j.
4519 * Skip half of y when i and j have the same x.
4528 for (cy = cyf_x; cy <= cyl; cy++)
4530 c0 = gridj->cxy_ind[cx*gridj->ncy+cy];
4531 c1 = gridj->cxy_ind[cx*gridj->ncy+cy+1];
4532 #ifdef NBNXN_SHIFT_BACKWARD
4533 if (gridi == gridj &&
4534 shift == CENTRAL && c0 < ci)
4541 if (gridj->c0[YY] + cy*gridj->sy > by1)
4543 d2zxy += sqr(gridj->c0[YY] + cy*gridj->sy - by1);
4545 else if (gridj->c0[YY] + (cy+1)*gridj->sy < by0)
4547 d2zxy += sqr(gridj->c0[YY] + (cy+1)*gridj->sy - by0);
4549 if (c1 > c0 && d2zxy < rl2)
4551 cs = c0 + (int)(bz1_frac*(c1 - c0));
4559 /* Find the lowest cell that can possibly
4564 (bbcz_j[cf*NNBSBB_D+1] >= bz0 ||
4565 d2xy + sqr(bbcz_j[cf*NNBSBB_D+1] - bz0) < rl2))
4570 /* Find the highest cell that can possibly
4575 (bbcz_j[cl*NNBSBB_D] <= bz1 ||
4576 d2xy + sqr(bbcz_j[cl*NNBSBB_D] - bz1) < rl2))
4581 #ifdef NBNXN_REFCODE
4583 /* Simple reference code, for debugging,
4584 * overrides the more complex code above.
4589 for (k = c0; k < c1; k++)
4591 if (box_dist2(bx0, bx1, by0, by1, bz0, bz1, bb+k) < rl2 &&
4596 if (box_dist2(bx0, bx1, by0, by1, bz0, bz1, bb+k) < rl2 &&
4607 /* We want each atom/cell pair only once,
4608 * only use cj >= ci.
4610 #ifndef NBNXN_SHIFT_BACKWARD
4613 if (shift == CENTRAL)
4622 /* For f buffer flags with simple lists */
4623 ncj_old_j = nbl->ncj;
4625 switch (nb_kernel_type)
4627 case nbnxnk4x4_PlainC:
4628 check_subcell_list_space_simple(nbl, cl-cf+1);
4630 make_cluster_list_simple(gridj,
4632 (gridi == gridj && shift == CENTRAL),
4637 #ifdef GMX_NBNXN_SIMD_4XN
4638 case nbnxnk4xN_SIMD_4xN:
4639 check_subcell_list_space_simple(nbl, ci_to_cj(na_cj_2log, cl-cf)+2);
4640 make_cluster_list_simd_4xn(gridj,
4642 (gridi == gridj && shift == CENTRAL),
4648 #ifdef GMX_NBNXN_SIMD_2XNN
4649 case nbnxnk4xN_SIMD_2xNN:
4650 check_subcell_list_space_simple(nbl, ci_to_cj(na_cj_2log, cl-cf)+2);
4651 make_cluster_list_simd_2xnn(gridj,
4653 (gridi == gridj && shift == CENTRAL),
4659 case nbnxnk8x8x8_PlainC:
4660 case nbnxnk8x8x8_CUDA:
4661 check_subcell_list_space_supersub(nbl, cl-cf+1);
4662 for (cj = cf; cj <= cl; cj++)
4664 make_cluster_list_supersub(gridi, gridj,
4666 (gridi == gridj && shift == CENTRAL && ci == cj),
4667 nbat->xstride, nbat->x,
4673 ncpcheck += cl - cf + 1;
4675 if (bFBufferFlag && nbl->ncj > ncj_old_j)
4679 cbf = nbl->cj[ncj_old_j].cj >> gridj_flag_shift;
4680 cbl = nbl->cj[nbl->ncj-1].cj >> gridj_flag_shift;
4681 for (cb = cbf; cb <= cbl; cb++)
4683 gridj_flag[cb] = 1U<<th;
4691 /* Set the exclusions for this ci list */
4694 set_ci_top_excls(nbs,
4696 shift == CENTRAL && gridi == gridj,
4699 &(nbl->ci[nbl->nci]),
4704 set_sci_top_excls(nbs,
4706 shift == CENTRAL && gridi == gridj,
4708 &(nbl->sci[nbl->nsci]),
4712 /* Close this ci list */
4715 close_ci_entry_simple(nbl);
4719 close_ci_entry_supersub(nbl,
4721 progBal, min_ci_balanced,
4728 if (bFBufferFlag && nbl->ncj > ncj_old_i)
4730 work->buffer_flags.flag[(gridi->cell0+ci)>>gridi_flag_shift] = 1U<<th;
4734 work->ndistc = ndistc;
4736 nbs_cycle_stop(&work->cc[enbsCCsearch]);
4740 fprintf(debug, "number of distance checks %d\n", ndistc);
4741 fprintf(debug, "ncpcheck %s %d\n", gridi == gridj ? "local" : "non-local",
4746 print_nblist_statistics_simple(debug, nbl, nbs, rlist);
4750 print_nblist_statistics_supersub(debug, nbl, nbs, rlist);
4756 static void reduce_buffer_flags(const nbnxn_search_t nbs,
4758 const nbnxn_buffer_flags_t *dest)
4761 const unsigned *flag;
4763 for (s = 0; s < nsrc; s++)
4765 flag = nbs->work[s].buffer_flags.flag;
4767 for (b = 0; b < dest->nflag; b++)
4769 dest->flag[b] |= flag[b];
4774 static void print_reduction_cost(const nbnxn_buffer_flags_t *flags, int nout)
4776 int nelem, nkeep, ncopy, nred, b, c, out;
4782 for (b = 0; b < flags->nflag; b++)
4784 if (flags->flag[b] == 1)
4786 /* Only flag 0 is set, no copy of reduction required */
4790 else if (flags->flag[b] > 0)
4793 for (out = 0; out < nout; out++)
4795 if (flags->flag[b] & (1U<<out))
4812 fprintf(debug, "nbnxn reduction: #flag %d #list %d elem %4.2f, keep %4.2f copy %4.2f red %4.2f\n",
4814 nelem/(double)(flags->nflag),
4815 nkeep/(double)(flags->nflag),
4816 ncopy/(double)(flags->nflag),
4817 nred/(double)(flags->nflag));
4820 /* Perform a count (linear) sort to sort the smaller lists to the end.
4821 * This avoids load imbalance on the GPU, as large lists will be
4822 * scheduled and executed first and the smaller lists later.
4823 * Load balancing between multi-processors only happens at the end
4824 * and there smaller lists lead to more effective load balancing.
4825 * The sorting is done on the cj4 count, not on the actual pair counts.
4826 * Not only does this make the sort faster, but it also results in
4827 * better load balancing than using a list sorted on exact load.
4828 * This function swaps the pointer in the pair list to avoid a copy operation.
4830 static void sort_sci(nbnxn_pairlist_t *nbl)
4832 nbnxn_list_work_t *work;
4833 int m, i, s, s0, s1;
4834 nbnxn_sci_t *sci_sort;
4836 if (nbl->ncj4 <= nbl->nsci)
4838 /* nsci = 0 or all sci have size 1, sorting won't change the order */
4844 /* We will distinguish differences up to double the average */
4845 m = (2*nbl->ncj4)/nbl->nsci;
4847 if (m + 1 > work->sort_nalloc)
4849 work->sort_nalloc = over_alloc_large(m + 1);
4850 srenew(work->sort, work->sort_nalloc);
4853 if (work->sci_sort_nalloc != nbl->sci_nalloc)
4855 work->sci_sort_nalloc = nbl->sci_nalloc;
4856 nbnxn_realloc_void((void **)&work->sci_sort,
4858 work->sci_sort_nalloc*sizeof(*work->sci_sort),
4859 nbl->alloc, nbl->free);
4862 /* Count the entries of each size */
4863 for (i = 0; i <= m; i++)
4867 for (s = 0; s < nbl->nsci; s++)
4869 i = min(m, nbl->sci[s].cj4_ind_end - nbl->sci[s].cj4_ind_start);
4872 /* Calculate the offset for each count */
4875 for (i = m - 1; i >= 0; i--)
4878 work->sort[i] = work->sort[i + 1] + s0;
4882 /* Sort entries directly into place */
4883 sci_sort = work->sci_sort;
4884 for (s = 0; s < nbl->nsci; s++)
4886 i = min(m, nbl->sci[s].cj4_ind_end - nbl->sci[s].cj4_ind_start);
4887 sci_sort[work->sort[i]++] = nbl->sci[s];
4890 /* Swap the sci pointers so we use the new, sorted list */
4891 work->sci_sort = nbl->sci;
4892 nbl->sci = sci_sort;
4895 /* Make a local or non-local pair-list, depending on iloc */
4896 void nbnxn_make_pairlist(const nbnxn_search_t nbs,
4897 nbnxn_atomdata_t *nbat,
4898 const t_blocka *excl,
4900 int min_ci_balanced,
4901 nbnxn_pairlist_set_t *nbl_list,
4906 nbnxn_grid_t *gridi, *gridj;
4908 int nzi, zi, zj0, zj1, zj;
4912 nbnxn_pairlist_t **nbl;
4914 gmx_bool CombineNBLists;
4916 int np_tot, np_noq, np_hlj, nap;
4918 /* Check if we are running hybrid GPU + CPU nbnxn mode */
4919 bGPUCPU = (!nbs->grid[0].bSimple && nbl_list->bSimple);
4921 nnbl = nbl_list->nnbl;
4922 nbl = nbl_list->nbl;
4923 CombineNBLists = nbl_list->bCombined;
4927 fprintf(debug, "ns making %d nblists\n", nnbl);
4930 nbat->bUseBufferFlags = (nbat->nout > 1);
4931 /* We should re-init the flags before making the first list */
4932 if (nbat->bUseBufferFlags && (LOCAL_I(iloc) || bGPUCPU))
4934 init_buffer_flags(&nbat->buffer_flags, nbat->natoms);
4937 if (nbl_list->bSimple)
4939 switch (nb_kernel_type)
4941 #ifdef GMX_NBNXN_SIMD_4XN
4942 case nbnxnk4xN_SIMD_4xN:
4943 nbs->icell_set_x = icell_set_x_simd_4xn;
4946 #ifdef GMX_NBNXN_SIMD_2XNN
4947 case nbnxnk4xN_SIMD_2xNN:
4948 nbs->icell_set_x = icell_set_x_simd_2xnn;
4952 nbs->icell_set_x = icell_set_x_simple;
4958 #ifdef NBNXN_SEARCH_BB_SIMD4
4959 nbs->icell_set_x = icell_set_x_supersub_simd4;
4961 nbs->icell_set_x = icell_set_x_supersub;
4967 /* Only zone (grid) 0 vs 0 */
4974 nzi = nbs->zones->nizone;
4977 if (!nbl_list->bSimple && min_ci_balanced > 0)
4979 nsubpair_max = get_nsubpair_max(nbs, iloc, rlist, min_ci_balanced);
4986 /* Clear all pair-lists */
4987 for (th = 0; th < nnbl; th++)
4989 clear_pairlist(nbl[th]);
4992 for (zi = 0; zi < nzi; zi++)
4994 gridi = &nbs->grid[zi];
4996 if (NONLOCAL_I(iloc))
4998 zj0 = nbs->zones->izone[zi].j0;
4999 zj1 = nbs->zones->izone[zi].j1;
5005 for (zj = zj0; zj < zj1; zj++)
5007 gridj = &nbs->grid[zj];
5011 fprintf(debug, "ns search grid %d vs %d\n", zi, zj);
5014 nbs_cycle_start(&nbs->cc[enbsCCsearch]);
5016 if (nbl[0]->bSimple && !gridi->bSimple)
5018 /* Hybrid list, determine blocking later */
5023 ci_block = get_ci_block_size(gridi, nbs->DomDec, nnbl);
5026 #pragma omp parallel for num_threads(nnbl) schedule(static)
5027 for (th = 0; th < nnbl; th++)
5029 /* Re-init the thread-local work flag data before making
5030 * the first list (not an elegant conditional).
5032 if (nbat->bUseBufferFlags && ((zi == 0 && zj == 0) ||
5033 (bGPUCPU && zi == 0 && zj == 1)))
5035 init_buffer_flags(&nbs->work[th].buffer_flags, nbat->natoms);
5038 if (CombineNBLists && th > 0)
5040 clear_pairlist(nbl[th]);
5043 /* With GPU: generate progressively smaller lists for
5044 * load balancing for local only or non-local with 2 zones.
5046 progBal = (LOCAL_I(iloc) || nbs->zones->n <= 2);
5048 /* Divide the i super cell equally over the nblists */
5049 nbnxn_make_pairlist_part(nbs, gridi, gridj,
5050 &nbs->work[th], nbat, excl,
5054 nbat->bUseBufferFlags,
5056 progBal, min_ci_balanced,
5060 nbs_cycle_stop(&nbs->cc[enbsCCsearch]);
5065 for (th = 0; th < nnbl; th++)
5067 inc_nrnb(nrnb, eNR_NBNXN_DIST2, nbs->work[th].ndistc);
5069 if (nbl_list->bSimple)
5071 np_tot += nbl[th]->ncj;
5072 np_noq += nbl[th]->work->ncj_noq;
5073 np_hlj += nbl[th]->work->ncj_hlj;
5077 /* This count ignores potential subsequent pair pruning */
5078 np_tot += nbl[th]->nci_tot;
5081 nap = nbl[0]->na_ci*nbl[0]->na_cj;
5082 nbl_list->natpair_ljq = (np_tot - np_noq)*nap - np_hlj*nap/2;
5083 nbl_list->natpair_lj = np_noq*nap;
5084 nbl_list->natpair_q = np_hlj*nap/2;
5086 if (CombineNBLists && nnbl > 1)
5088 nbs_cycle_start(&nbs->cc[enbsCCcombine]);
5090 combine_nblists(nnbl-1, nbl+1, nbl[0]);
5092 nbs_cycle_stop(&nbs->cc[enbsCCcombine]);
5097 if (!nbl_list->bSimple)
5099 /* Sort the entries on size, large ones first */
5100 if (CombineNBLists || nnbl == 1)
5106 #pragma omp parallel for num_threads(nnbl) schedule(static)
5107 for (th = 0; th < nnbl; th++)
5114 if (nbat->bUseBufferFlags)
5116 reduce_buffer_flags(nbs, nnbl, &nbat->buffer_flags);
5119 /* Special performance logging stuff (env.var. GMX_NBNXN_CYCLE) */
5122 nbs->search_count++;
5124 if (nbs->print_cycles &&
5125 (!nbs->DomDec || (nbs->DomDec && !LOCAL_I(iloc))) &&
5126 nbs->search_count % 100 == 0)
5128 nbs_cycle_print(stderr, nbs);
5131 if (debug && (CombineNBLists && nnbl > 1))
5133 if (nbl[0]->bSimple)
5135 print_nblist_statistics_simple(debug, nbl[0], nbs, rlist);
5139 print_nblist_statistics_supersub(debug, nbl[0], nbs, rlist);
5147 if (nbl[0]->bSimple)
5149 print_nblist_ci_cj(debug, nbl[0]);
5153 print_nblist_sci_cj(debug, nbl[0]);
5157 if (nbat->bUseBufferFlags)
5159 print_reduction_cost(&nbat->buffer_flags, nnbl);