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_cyclecounter.h"
54 #include "gromacs/fileio/gmxfio.h"
55 #include "gmx_omp_nthreads.h"
59 #ifdef NBNXN_SEARCH_BB_SIMD4
60 /* We use 4-wide SIMD for bounding box calculations */
63 /* Single precision BBs + coordinates, we can also load coordinates with SIMD */
64 #define NBNXN_SEARCH_SIMD4_FLOAT_X_BB
67 #if defined NBNXN_SEARCH_SIMD4_FLOAT_X_BB && (GPU_NSUBCELL == 4 || GPU_NSUBCELL == 8)
68 /* Store bounding boxes with x, y and z coordinates in packs of 4 */
69 #define NBNXN_PBB_SIMD4
72 /* The packed bounding box coordinate stride is always set to 4.
73 * With AVX we could use 8, but that turns out not to be faster.
76 #define STRIDE_PBB_2LOG 2
78 #endif /* NBNXN_SEARCH_BB_SIMD4 */
82 /* The functions below are macros as they are performance sensitive */
84 /* 4x4 list, pack=4: no complex conversion required */
85 /* i-cluster to j-cluster conversion */
86 #define CI_TO_CJ_J4(ci) (ci)
87 /* cluster index to coordinate array index conversion */
88 #define X_IND_CI_J4(ci) ((ci)*STRIDE_P4)
89 #define X_IND_CJ_J4(cj) ((cj)*STRIDE_P4)
91 /* 4x2 list, pack=4: j-cluster size is half the packing width */
92 /* i-cluster to j-cluster conversion */
93 #define CI_TO_CJ_J2(ci) ((ci)<<1)
94 /* cluster index to coordinate array index conversion */
95 #define X_IND_CI_J2(ci) ((ci)*STRIDE_P4)
96 #define X_IND_CJ_J2(cj) (((cj)>>1)*STRIDE_P4 + ((cj) & 1)*(PACK_X4>>1))
98 /* 4x8 list, pack=8: i-cluster size is half the packing width */
99 /* i-cluster to j-cluster conversion */
100 #define CI_TO_CJ_J8(ci) ((ci)>>1)
101 /* cluster index to coordinate array index conversion */
102 #define X_IND_CI_J8(ci) (((ci)>>1)*STRIDE_P8 + ((ci) & 1)*(PACK_X8>>1))
103 #define X_IND_CJ_J8(cj) ((cj)*STRIDE_P8)
105 /* The j-cluster size is matched to the SIMD width */
106 #if GMX_SIMD_WIDTH_HERE == 2
107 #define CI_TO_CJ_SIMD_4XN(ci) CI_TO_CJ_J2(ci)
108 #define X_IND_CI_SIMD_4XN(ci) X_IND_CI_J2(ci)
109 #define X_IND_CJ_SIMD_4XN(cj) X_IND_CJ_J2(cj)
111 #if GMX_SIMD_WIDTH_HERE == 4
112 #define CI_TO_CJ_SIMD_4XN(ci) CI_TO_CJ_J4(ci)
113 #define X_IND_CI_SIMD_4XN(ci) X_IND_CI_J4(ci)
114 #define X_IND_CJ_SIMD_4XN(cj) X_IND_CJ_J4(cj)
116 #if GMX_SIMD_WIDTH_HERE == 8
117 #define CI_TO_CJ_SIMD_4XN(ci) CI_TO_CJ_J8(ci)
118 #define X_IND_CI_SIMD_4XN(ci) X_IND_CI_J8(ci)
119 #define X_IND_CJ_SIMD_4XN(cj) X_IND_CJ_J8(cj)
120 /* Half SIMD with j-cluster size */
121 #define CI_TO_CJ_SIMD_2XNN(ci) CI_TO_CJ_J4(ci)
122 #define X_IND_CI_SIMD_2XNN(ci) X_IND_CI_J4(ci)
123 #define X_IND_CJ_SIMD_2XNN(cj) X_IND_CJ_J4(cj)
125 #if GMX_SIMD_WIDTH_HERE == 16
126 #define CI_TO_CJ_SIMD_2XNN(ci) CI_TO_CJ_J8(ci)
127 #define X_IND_CI_SIMD_2XNN(ci) X_IND_CI_J8(ci)
128 #define X_IND_CJ_SIMD_2XNN(cj) X_IND_CJ_J8(cj)
130 #error "unsupported GMX_NBNXN_SIMD_WIDTH"
136 #endif /* GMX_NBNXN_SIMD */
139 #ifdef NBNXN_SEARCH_BB_SIMD4
140 /* Store bounding boxes corners as quadruplets: xxxxyyyyzzzz */
142 /* Size of bounding box corners quadruplet */
143 #define NNBSBB_XXXX (NNBSBB_D*DIM*STRIDE_PBB)
146 /* We shift the i-particles backward for PBC.
147 * This leads to more conditionals than shifting forward.
148 * We do this to get more balanced pair lists.
150 #define NBNXN_SHIFT_BACKWARD
153 /* This define is a lazy way to avoid interdependence of the grid
154 * and searching data structures.
156 #define NBNXN_NA_SC_MAX (GPU_NSUBCELL*NBNXN_GPU_CLUSTER_SIZE)
159 static void nbs_cycle_clear(nbnxn_cycle_t *cc)
163 for (i = 0; i < enbsCCnr; i++)
170 static double Mcyc_av(const nbnxn_cycle_t *cc)
172 return (double)cc->c*1e-6/cc->count;
175 static void nbs_cycle_print(FILE *fp, const nbnxn_search_t nbs)
181 fprintf(fp, "ns %4d grid %4.1f search %4.1f red.f %5.3f",
182 nbs->cc[enbsCCgrid].count,
183 Mcyc_av(&nbs->cc[enbsCCgrid]),
184 Mcyc_av(&nbs->cc[enbsCCsearch]),
185 Mcyc_av(&nbs->cc[enbsCCreducef]));
187 if (nbs->nthread_max > 1)
189 if (nbs->cc[enbsCCcombine].count > 0)
191 fprintf(fp, " comb %5.2f",
192 Mcyc_av(&nbs->cc[enbsCCcombine]));
194 fprintf(fp, " s. th");
195 for (t = 0; t < nbs->nthread_max; t++)
197 fprintf(fp, " %4.1f",
198 Mcyc_av(&nbs->work[t].cc[enbsCCsearch]));
204 static void nbnxn_grid_init(nbnxn_grid_t * grid)
207 grid->cxy_ind = NULL;
208 grid->cxy_nalloc = 0;
214 static int get_2log(int n)
219 while ((1<<log2) < n)
225 gmx_fatal(FARGS, "nbnxn na_c (%d) is not a power of 2", n);
231 static int nbnxn_kernel_to_ci_size(int nb_kernel_type)
233 switch (nb_kernel_type)
235 case nbnxnk4x4_PlainC:
236 case nbnxnk4xN_SIMD_4xN:
237 case nbnxnk4xN_SIMD_2xNN:
238 return NBNXN_CPU_CLUSTER_I_SIZE;
239 case nbnxnk8x8x8_CUDA:
240 case nbnxnk8x8x8_PlainC:
241 /* The cluster size for super/sub lists is only set here.
242 * Any value should work for the pair-search and atomdata code.
243 * The kernels, of course, might require a particular value.
245 return NBNXN_GPU_CLUSTER_SIZE;
247 gmx_incons("unknown kernel type");
253 int nbnxn_kernel_to_cj_size(int nb_kernel_type)
255 int nbnxn_simd_width = 0;
258 #ifdef GMX_NBNXN_SIMD
259 nbnxn_simd_width = GMX_SIMD_WIDTH_HERE;
262 switch (nb_kernel_type)
264 case nbnxnk4x4_PlainC:
265 cj_size = NBNXN_CPU_CLUSTER_I_SIZE;
267 case nbnxnk4xN_SIMD_4xN:
268 cj_size = nbnxn_simd_width;
270 case nbnxnk4xN_SIMD_2xNN:
271 cj_size = nbnxn_simd_width/2;
273 case nbnxnk8x8x8_CUDA:
274 case nbnxnk8x8x8_PlainC:
275 cj_size = nbnxn_kernel_to_ci_size(nb_kernel_type);
278 gmx_incons("unknown kernel type");
284 static int ci_to_cj(int na_cj_2log, int ci)
288 case 2: return ci; break;
289 case 1: return (ci<<1); break;
290 case 3: return (ci>>1); break;
296 gmx_bool nbnxn_kernel_pairlist_simple(int nb_kernel_type)
298 if (nb_kernel_type == nbnxnkNotSet)
300 gmx_fatal(FARGS, "Non-bonded kernel type not set for Verlet-style pair-list.");
303 switch (nb_kernel_type)
305 case nbnxnk8x8x8_CUDA:
306 case nbnxnk8x8x8_PlainC:
309 case nbnxnk4x4_PlainC:
310 case nbnxnk4xN_SIMD_4xN:
311 case nbnxnk4xN_SIMD_2xNN:
315 gmx_incons("Invalid nonbonded kernel type passed!");
320 void nbnxn_init_search(nbnxn_search_t * nbs_ptr,
322 gmx_domdec_zones_t *zones,
331 nbs->DomDec = (n_dd_cells != NULL);
333 clear_ivec(nbs->dd_dim);
339 for (d = 0; d < DIM; d++)
341 if ((*n_dd_cells)[d] > 1)
344 /* Each grid matches a DD zone */
350 snew(nbs->grid, nbs->ngrid);
351 for (g = 0; g < nbs->ngrid; g++)
353 nbnxn_grid_init(&nbs->grid[g]);
356 nbs->cell_nalloc = 0;
360 nbs->nthread_max = nthread_max;
362 /* Initialize the work data structures for each thread */
363 snew(nbs->work, nbs->nthread_max);
364 for (t = 0; t < nbs->nthread_max; t++)
366 nbs->work[t].cxy_na = NULL;
367 nbs->work[t].cxy_na_nalloc = 0;
368 nbs->work[t].sort_work = NULL;
369 nbs->work[t].sort_work_nalloc = 0;
372 /* Initialize detailed nbsearch cycle counting */
373 nbs->print_cycles = (getenv("GMX_NBNXN_CYCLE") != 0);
374 nbs->search_count = 0;
375 nbs_cycle_clear(nbs->cc);
376 for (t = 0; t < nbs->nthread_max; t++)
378 nbs_cycle_clear(nbs->work[t].cc);
382 static real grid_atom_density(int n, rvec corner0, rvec corner1)
386 rvec_sub(corner1, corner0, size);
388 return n/(size[XX]*size[YY]*size[ZZ]);
391 static int set_grid_size_xy(const nbnxn_search_t nbs,
394 int n, rvec corner0, rvec corner1,
399 real adens, tlen, tlen_x, tlen_y, nc_max;
402 rvec_sub(corner1, corner0, size);
406 /* target cell length */
409 /* To minimize the zero interactions, we should make
410 * the largest of the i/j cell cubic.
412 na_c = max(grid->na_c, grid->na_cj);
414 /* Approximately cubic cells */
415 tlen = pow(na_c/atom_density, 1.0/3.0);
421 /* Approximately cubic sub cells */
422 tlen = pow(grid->na_c/atom_density, 1.0/3.0);
423 tlen_x = tlen*GPU_NSUBCELL_X;
424 tlen_y = tlen*GPU_NSUBCELL_Y;
426 /* We round ncx and ncy down, because we get less cell pairs
427 * in the nbsist when the fixed cell dimensions (x,y) are
428 * larger than the variable one (z) than the other way around.
430 grid->ncx = max(1, (int)(size[XX]/tlen_x));
431 grid->ncy = max(1, (int)(size[YY]/tlen_y));
439 grid->sx = size[XX]/grid->ncx;
440 grid->sy = size[YY]/grid->ncy;
441 grid->inv_sx = 1/grid->sx;
442 grid->inv_sy = 1/grid->sy;
446 /* This is a non-home zone, add an extra row of cells
447 * for particles communicated for bonded interactions.
448 * These can be beyond the cut-off. It doesn't matter where
449 * they end up on the grid, but for performance it's better
450 * if they don't end up in cells that can be within cut-off range.
456 /* We need one additional cell entry for particles moved by DD */
457 if (grid->ncx*grid->ncy+1 > grid->cxy_nalloc)
459 grid->cxy_nalloc = over_alloc_large(grid->ncx*grid->ncy+1);
460 srenew(grid->cxy_na, grid->cxy_nalloc);
461 srenew(grid->cxy_ind, grid->cxy_nalloc+1);
463 for (t = 0; t < nbs->nthread_max; t++)
465 if (grid->ncx*grid->ncy+1 > nbs->work[t].cxy_na_nalloc)
467 nbs->work[t].cxy_na_nalloc = over_alloc_large(grid->ncx*grid->ncy+1);
468 srenew(nbs->work[t].cxy_na, nbs->work[t].cxy_na_nalloc);
472 /* Worst case scenario of 1 atom in each last cell */
473 if (grid->na_cj <= grid->na_c)
475 nc_max = n/grid->na_sc + grid->ncx*grid->ncy;
479 nc_max = n/grid->na_sc + grid->ncx*grid->ncy*grid->na_cj/grid->na_c;
482 if (nc_max > grid->nc_nalloc)
484 grid->nc_nalloc = over_alloc_large(nc_max);
485 srenew(grid->nsubc, grid->nc_nalloc);
486 srenew(grid->bbcz, grid->nc_nalloc*NNBSBB_D);
488 sfree_aligned(grid->bb);
489 /* This snew also zeros the contents, this avoid possible
490 * floating exceptions in SIMD with the unused bb elements.
494 snew_aligned(grid->bb, grid->nc_nalloc, 16);
501 pbb_nalloc = grid->nc_nalloc*GPU_NSUBCELL/STRIDE_PBB*NNBSBB_XXXX;
502 snew_aligned(grid->pbb, pbb_nalloc, 16);
504 snew_aligned(grid->bb, grid->nc_nalloc*GPU_NSUBCELL, 16);
510 if (grid->na_cj == grid->na_c)
512 grid->bbj = grid->bb;
516 sfree_aligned(grid->bbj);
517 snew_aligned(grid->bbj, grid->nc_nalloc*grid->na_c/grid->na_cj, 16);
521 srenew(grid->flags, grid->nc_nalloc);
524 copy_rvec(corner0, grid->c0);
525 copy_rvec(corner1, grid->c1);
530 /* We need to sort paricles in grid columns on z-coordinate.
531 * As particle are very often distributed homogeneously, we a sorting
532 * algorithm similar to pigeonhole sort. We multiply the z-coordinate
533 * by a factor, cast to an int and try to store in that hole. If the hole
534 * is full, we move this or another particle. A second pass is needed to make
535 * contiguous elements. SORT_GRID_OVERSIZE is the ratio of holes to particles.
536 * 4 is the optimal value for homogeneous particle distribution and allows
537 * for an O(#particles) sort up till distributions were all particles are
538 * concentrated in 1/4 of the space. No NlogN fallback is implemented,
539 * as it can be expensive to detect imhomogeneous particle distributions.
540 * SGSF is the maximum ratio of holes used, in the worst case all particles
541 * end up in the last hole and we need #particles extra holes at the end.
543 #define SORT_GRID_OVERSIZE 4
544 #define SGSF (SORT_GRID_OVERSIZE + 1)
546 /* Sort particle index a on coordinates x along dim.
547 * Backwards tells if we want decreasing iso increasing coordinates.
548 * h0 is the minimum of the coordinate range.
549 * invh is the 1/length of the sorting range.
550 * n_per_h (>=n) is the expected average number of particles per 1/invh
551 * sort is the sorting work array.
552 * sort should have a size of at least n_per_h*SORT_GRID_OVERSIZE + n,
553 * or easier, allocate at least n*SGSF elements.
555 static void sort_atoms(int dim, gmx_bool Backwards,
556 int *a, int n, rvec *x,
557 real h0, real invh, int n_per_h,
561 int zi, zim, zi_min, zi_max;
573 gmx_incons("n > n_per_h");
577 /* Transform the inverse range height into the inverse hole height */
578 invh *= n_per_h*SORT_GRID_OVERSIZE;
580 /* Set nsort to the maximum possible number of holes used.
581 * In worst case all n elements end up in the last bin.
583 nsort = n_per_h*SORT_GRID_OVERSIZE + n;
585 /* Determine the index range used, so we can limit it for the second pass */
589 /* Sort the particles using a simple index sort */
590 for (i = 0; i < n; i++)
592 /* The cast takes care of float-point rounding effects below zero.
593 * This code assumes particles are less than 1/SORT_GRID_OVERSIZE
594 * times the box height out of the box.
596 zi = (int)((x[a[i]][dim] - h0)*invh);
599 /* As we can have rounding effect, we use > iso >= here */
600 if (zi < 0 || zi > n_per_h*SORT_GRID_OVERSIZE)
602 gmx_fatal(FARGS, "(int)((x[%d][%c]=%f - %f)*%f) = %d, not in 0 - %d*%d\n",
603 a[i], 'x'+dim, x[a[i]][dim], h0, invh, zi,
604 n_per_h, SORT_GRID_OVERSIZE);
608 /* Ideally this particle should go in sort cell zi,
609 * but that might already be in use,
610 * in that case find the first empty cell higher up
615 zi_min = min(zi_min, zi);
616 zi_max = max(zi_max, zi);
620 /* We have multiple atoms in the same sorting slot.
621 * Sort on real z for minimal bounding box size.
622 * There is an extra check for identical z to ensure
623 * well-defined output order, independent of input order
624 * to ensure binary reproducibility after restarts.
626 while (sort[zi] >= 0 && ( x[a[i]][dim] > x[sort[zi]][dim] ||
627 (x[a[i]][dim] == x[sort[zi]][dim] &&
635 /* Shift all elements by one slot until we find an empty slot */
638 while (sort[zim] >= 0)
646 zi_max = max(zi_max, zim);
649 zi_max = max(zi_max, zi);
656 for (zi = 0; zi < nsort; zi++)
667 for (zi = zi_max; zi >= zi_min; zi--)
678 gmx_incons("Lost particles while sorting");
683 #define R2F_D(x) ((float)((x) >= 0 ? ((1-GMX_FLOAT_EPS)*(x)) : ((1+GMX_FLOAT_EPS)*(x))))
684 #define R2F_U(x) ((float)((x) >= 0 ? ((1+GMX_FLOAT_EPS)*(x)) : ((1-GMX_FLOAT_EPS)*(x))))
690 /* Coordinate order x,y,z, bb order xyz0 */
691 static void calc_bounding_box(int na, int stride, const real *x, nbnxn_bb_t *bb)
694 real xl, xh, yl, yh, zl, zh;
704 for (j = 1; j < na; j++)
706 xl = min(xl, x[i+XX]);
707 xh = max(xh, x[i+XX]);
708 yl = min(yl, x[i+YY]);
709 yh = max(yh, x[i+YY]);
710 zl = min(zl, x[i+ZZ]);
711 zh = max(zh, x[i+ZZ]);
714 /* Note: possible double to float conversion here */
715 bb->lower[BB_X] = R2F_D(xl);
716 bb->lower[BB_Y] = R2F_D(yl);
717 bb->lower[BB_Z] = R2F_D(zl);
718 bb->upper[BB_X] = R2F_U(xh);
719 bb->upper[BB_Y] = R2F_U(yh);
720 bb->upper[BB_Z] = R2F_U(zh);
723 /* Packed coordinates, bb order xyz0 */
724 static void calc_bounding_box_x_x4(int na, const real *x, nbnxn_bb_t *bb)
727 real xl, xh, yl, yh, zl, zh;
735 for (j = 1; j < na; j++)
737 xl = min(xl, x[j+XX*PACK_X4]);
738 xh = max(xh, x[j+XX*PACK_X4]);
739 yl = min(yl, x[j+YY*PACK_X4]);
740 yh = max(yh, x[j+YY*PACK_X4]);
741 zl = min(zl, x[j+ZZ*PACK_X4]);
742 zh = max(zh, x[j+ZZ*PACK_X4]);
744 /* Note: possible double to float conversion here */
745 bb->lower[BB_X] = R2F_D(xl);
746 bb->lower[BB_Y] = R2F_D(yl);
747 bb->lower[BB_Z] = R2F_D(zl);
748 bb->upper[BB_X] = R2F_U(xh);
749 bb->upper[BB_Y] = R2F_U(yh);
750 bb->upper[BB_Z] = R2F_U(zh);
753 /* Packed coordinates, bb order xyz0 */
754 static void calc_bounding_box_x_x8(int na, const real *x, nbnxn_bb_t *bb)
757 real xl, xh, yl, yh, zl, zh;
765 for (j = 1; j < na; j++)
767 xl = min(xl, x[j+XX*PACK_X8]);
768 xh = max(xh, x[j+XX*PACK_X8]);
769 yl = min(yl, x[j+YY*PACK_X8]);
770 yh = max(yh, x[j+YY*PACK_X8]);
771 zl = min(zl, x[j+ZZ*PACK_X8]);
772 zh = max(zh, x[j+ZZ*PACK_X8]);
774 /* Note: possible double to float conversion here */
775 bb->lower[BB_X] = R2F_D(xl);
776 bb->lower[BB_Y] = R2F_D(yl);
777 bb->lower[BB_Z] = R2F_D(zl);
778 bb->upper[BB_X] = R2F_U(xh);
779 bb->upper[BB_Y] = R2F_U(yh);
780 bb->upper[BB_Z] = R2F_U(zh);
783 /* Packed coordinates, bb order xyz0 */
784 static void calc_bounding_box_x_x4_halves(int na, const real *x,
785 nbnxn_bb_t *bb, nbnxn_bb_t *bbj)
787 calc_bounding_box_x_x4(min(na, 2), x, bbj);
791 calc_bounding_box_x_x4(min(na-2, 2), x+(PACK_X4>>1), bbj+1);
795 /* Set the "empty" bounding box to the same as the first one,
796 * so we don't need to treat special cases in the rest of the code.
798 #ifdef NBNXN_SEARCH_BB_SIMD4
799 gmx_simd4_store_pr(&bbj[1].lower[0], gmx_simd4_load_bb_pr(&bbj[0].lower[0]));
800 gmx_simd4_store_pr(&bbj[1].upper[0], gmx_simd4_load_bb_pr(&bbj[0].upper[0]));
806 #ifdef NBNXN_SEARCH_BB_SIMD4
807 gmx_simd4_store_pr(&bb->lower[0],
808 gmx_simd4_min_pr(gmx_simd4_load_bb_pr(&bbj[0].lower[0]),
809 gmx_simd4_load_bb_pr(&bbj[1].lower[0])));
810 gmx_simd4_store_pr(&bb->upper[0],
811 gmx_simd4_max_pr(gmx_simd4_load_bb_pr(&bbj[0].upper[0]),
812 gmx_simd4_load_bb_pr(&bbj[1].upper[0])));
817 for (i = 0; i < NNBSBB_C; i++)
819 bb->lower[i] = min(bbj[0].lower[i], bbj[1].lower[i]);
820 bb->upper[i] = max(bbj[0].upper[i], bbj[1].upper[i]);
826 #ifdef NBNXN_SEARCH_BB_SIMD4
828 /* Coordinate order xyz, bb order xxxxyyyyzzzz */
829 static void calc_bounding_box_xxxx(int na, int stride, const real *x, float *bb)
832 real xl, xh, yl, yh, zl, zh;
842 for (j = 1; j < na; j++)
844 xl = min(xl, x[i+XX]);
845 xh = max(xh, x[i+XX]);
846 yl = min(yl, x[i+YY]);
847 yh = max(yh, x[i+YY]);
848 zl = min(zl, x[i+ZZ]);
849 zh = max(zh, x[i+ZZ]);
852 /* Note: possible double to float conversion here */
853 bb[0*STRIDE_PBB] = R2F_D(xl);
854 bb[1*STRIDE_PBB] = R2F_D(yl);
855 bb[2*STRIDE_PBB] = R2F_D(zl);
856 bb[3*STRIDE_PBB] = R2F_U(xh);
857 bb[4*STRIDE_PBB] = R2F_U(yh);
858 bb[5*STRIDE_PBB] = R2F_U(zh);
861 #endif /* NBNXN_SEARCH_BB_SIMD4 */
863 #ifdef NBNXN_SEARCH_SIMD4_FLOAT_X_BB
865 /* Coordinate order xyz?, bb order xyz0 */
866 static void calc_bounding_box_simd4(int na, const float *x, nbnxn_bb_t *bb)
868 gmx_simd4_pr bb_0_S, bb_1_S;
873 bb_0_S = gmx_simd4_load_bb_pr(x);
876 for (i = 1; i < na; i++)
878 x_S = gmx_simd4_load_bb_pr(x+i*NNBSBB_C);
879 bb_0_S = gmx_simd4_min_pr(bb_0_S, x_S);
880 bb_1_S = gmx_simd4_max_pr(bb_1_S, x_S);
883 gmx_simd4_store_pr(&bb->lower[0], bb_0_S);
884 gmx_simd4_store_pr(&bb->upper[0], bb_1_S);
887 /* Coordinate order xyz?, bb order xxxxyyyyzzzz */
888 static void calc_bounding_box_xxxx_simd4(int na, const float *x,
889 nbnxn_bb_t *bb_work_aligned,
892 calc_bounding_box_simd4(na, x, bb_work_aligned);
894 bb[0*STRIDE_PBB] = bb_work_aligned->lower[BB_X];
895 bb[1*STRIDE_PBB] = bb_work_aligned->lower[BB_Y];
896 bb[2*STRIDE_PBB] = bb_work_aligned->lower[BB_Z];
897 bb[3*STRIDE_PBB] = bb_work_aligned->upper[BB_X];
898 bb[4*STRIDE_PBB] = bb_work_aligned->upper[BB_Y];
899 bb[5*STRIDE_PBB] = bb_work_aligned->upper[BB_Z];
902 #endif /* NBNXN_SEARCH_SIMD4_FLOAT_X_BB */
905 /* Combines pairs of consecutive bounding boxes */
906 static void combine_bounding_box_pairs(nbnxn_grid_t *grid, const nbnxn_bb_t *bb)
908 int i, j, sc2, nc2, c2;
910 for (i = 0; i < grid->ncx*grid->ncy; i++)
912 /* Starting bb in a column is expected to be 2-aligned */
913 sc2 = grid->cxy_ind[i]>>1;
914 /* For odd numbers skip the last bb here */
915 nc2 = (grid->cxy_na[i]+3)>>(2+1);
916 for (c2 = sc2; c2 < sc2+nc2; c2++)
918 #ifdef NBNXN_SEARCH_BB_SIMD4
919 gmx_simd4_pr min_S, max_S;
921 min_S = gmx_simd4_min_pr(gmx_simd4_load_bb_pr(&bb[c2*2+0].lower[0]),
922 gmx_simd4_load_bb_pr(&bb[c2*2+1].lower[0]));
923 max_S = gmx_simd4_max_pr(gmx_simd4_load_bb_pr(&bb[c2*2+0].upper[0]),
924 gmx_simd4_load_bb_pr(&bb[c2*2+1].upper[0]));
925 gmx_simd4_store_pr(&grid->bbj[c2].lower[0], min_S);
926 gmx_simd4_store_pr(&grid->bbj[c2].upper[0], max_S);
928 for (j = 0; j < NNBSBB_C; j++)
930 grid->bbj[c2].lower[j] = min(bb[c2*2+0].lower[j],
931 bb[c2*2+1].lower[j]);
932 grid->bbj[c2].upper[j] = max(bb[c2*2+0].upper[j],
933 bb[c2*2+1].upper[j]);
937 if (((grid->cxy_na[i]+3)>>2) & 1)
939 /* The bb count in this column is odd: duplicate the last bb */
940 for (j = 0; j < NNBSBB_C; j++)
942 grid->bbj[c2].lower[j] = bb[c2*2].lower[j];
943 grid->bbj[c2].upper[j] = bb[c2*2].upper[j];
950 /* Prints the average bb size, used for debug output */
951 static void print_bbsizes_simple(FILE *fp,
952 const nbnxn_search_t nbs,
953 const nbnxn_grid_t *grid)
959 for (c = 0; c < grid->nc; c++)
961 for (d = 0; d < DIM; d++)
963 ba[d] += grid->bb[c].upper[d] - grid->bb[c].lower[d];
966 dsvmul(1.0/grid->nc, ba, ba);
968 fprintf(fp, "ns bb: %4.2f %4.2f %4.2f %4.2f %4.2f %4.2f rel %4.2f %4.2f %4.2f\n",
969 nbs->box[XX][XX]/grid->ncx,
970 nbs->box[YY][YY]/grid->ncy,
971 nbs->box[ZZ][ZZ]*grid->ncx*grid->ncy/grid->nc,
972 ba[XX], ba[YY], ba[ZZ],
973 ba[XX]*grid->ncx/nbs->box[XX][XX],
974 ba[YY]*grid->ncy/nbs->box[YY][YY],
975 ba[ZZ]*grid->nc/(grid->ncx*grid->ncy*nbs->box[ZZ][ZZ]));
978 /* Prints the average bb size, used for debug output */
979 static void print_bbsizes_supersub(FILE *fp,
980 const nbnxn_search_t nbs,
981 const nbnxn_grid_t *grid)
988 for (c = 0; c < grid->nc; c++)
991 for (s = 0; s < grid->nsubc[c]; s += STRIDE_PBB)
995 cs_w = (c*GPU_NSUBCELL + s)/STRIDE_PBB;
996 for (i = 0; i < STRIDE_PBB; i++)
998 for (d = 0; d < DIM; d++)
1001 grid->pbb[cs_w*NNBSBB_XXXX+(DIM+d)*STRIDE_PBB+i] -
1002 grid->pbb[cs_w*NNBSBB_XXXX+ d *STRIDE_PBB+i];
1007 for (s = 0; s < grid->nsubc[c]; s++)
1011 cs = c*GPU_NSUBCELL + s;
1012 for (d = 0; d < DIM; d++)
1014 ba[d] += grid->bb[cs].upper[d] - grid->bb[cs].lower[d];
1018 ns += grid->nsubc[c];
1020 dsvmul(1.0/ns, ba, ba);
1022 fprintf(fp, "ns bb: %4.2f %4.2f %4.2f %4.2f %4.2f %4.2f rel %4.2f %4.2f %4.2f\n",
1023 nbs->box[XX][XX]/(grid->ncx*GPU_NSUBCELL_X),
1024 nbs->box[YY][YY]/(grid->ncy*GPU_NSUBCELL_Y),
1025 nbs->box[ZZ][ZZ]*grid->ncx*grid->ncy/(grid->nc*GPU_NSUBCELL_Z),
1026 ba[XX], ba[YY], ba[ZZ],
1027 ba[XX]*grid->ncx*GPU_NSUBCELL_X/nbs->box[XX][XX],
1028 ba[YY]*grid->ncy*GPU_NSUBCELL_Y/nbs->box[YY][YY],
1029 ba[ZZ]*grid->nc*GPU_NSUBCELL_Z/(grid->ncx*grid->ncy*nbs->box[ZZ][ZZ]));
1032 /* Potentially sorts atoms on LJ coefficients !=0 and ==0.
1033 * Also sets interaction flags.
1035 void sort_on_lj(int na_c,
1036 int a0, int a1, const int *atinfo,
1040 int subc, s, a, n1, n2, a_lj_max, i, j;
1041 int sort1[NBNXN_NA_SC_MAX/GPU_NSUBCELL];
1042 int sort2[NBNXN_NA_SC_MAX/GPU_NSUBCELL];
1048 for (s = a0; s < a1; s += na_c)
1050 /* Make lists for this (sub-)cell on atoms with and without LJ */
1055 for (a = s; a < min(s+na_c, a1); a++)
1057 haveQ = haveQ || GET_CGINFO_HAS_Q(atinfo[order[a]]);
1059 if (GET_CGINFO_HAS_VDW(atinfo[order[a]]))
1061 sort1[n1++] = order[a];
1066 sort2[n2++] = order[a];
1070 /* If we don't have atom with LJ, there's nothing to sort */
1073 *flags |= NBNXN_CI_DO_LJ(subc);
1077 /* Only sort when strictly necessary. Ordering particles
1078 * Ordering particles can lead to less accurate summation
1079 * due to rounding, both for LJ and Coulomb interactions.
1081 if (2*(a_lj_max - s) >= na_c)
1083 for (i = 0; i < n1; i++)
1085 order[a0+i] = sort1[i];
1087 for (j = 0; j < n2; j++)
1089 order[a0+n1+j] = sort2[j];
1093 *flags |= NBNXN_CI_HALF_LJ(subc);
1098 *flags |= NBNXN_CI_DO_COUL(subc);
1104 /* Fill a pair search cell with atoms.
1105 * Potentially sorts atoms and sets the interaction flags.
1107 void fill_cell(const nbnxn_search_t nbs,
1109 nbnxn_atomdata_t *nbat,
1113 int sx, int sy, int sz,
1114 nbnxn_bb_t gmx_unused *bb_work_aligned)
1127 sort_on_lj(grid->na_c, a0, a1, atinfo, nbs->a,
1128 grid->flags+(a0>>grid->na_c_2log)-grid->cell0);
1131 /* Now we have sorted the atoms, set the cell indices */
1132 for (a = a0; a < a1; a++)
1134 nbs->cell[nbs->a[a]] = a;
1137 copy_rvec_to_nbat_real(nbs->a+a0, a1-a0, grid->na_c, x,
1138 nbat->XFormat, nbat->x, a0,
1141 if (nbat->XFormat == nbatX4)
1143 /* Store the bounding boxes as xyz.xyz. */
1144 offset = (a0 - grid->cell0*grid->na_sc) >> grid->na_c_2log;
1145 bb_ptr = grid->bb + offset;
1147 #if defined GMX_NBNXN_SIMD && GMX_SIMD_WIDTH_HERE == 2
1148 if (2*grid->na_cj == grid->na_c)
1150 calc_bounding_box_x_x4_halves(na, nbat->x+X4_IND_A(a0), bb_ptr,
1151 grid->bbj+offset*2);
1156 calc_bounding_box_x_x4(na, nbat->x+X4_IND_A(a0), bb_ptr);
1159 else if (nbat->XFormat == nbatX8)
1161 /* Store the bounding boxes as xyz.xyz. */
1162 offset = (a0 - grid->cell0*grid->na_sc) >> grid->na_c_2log;
1163 bb_ptr = grid->bb + offset;
1165 calc_bounding_box_x_x8(na, nbat->x+X8_IND_A(a0), bb_ptr);
1168 else if (!grid->bSimple)
1170 /* Store the bounding boxes in a format convenient
1171 * for SIMD4 calculations: xxxxyyyyzzzz...
1175 ((a0-grid->cell0*grid->na_sc)>>(grid->na_c_2log+STRIDE_PBB_2LOG))*NNBSBB_XXXX +
1176 (((a0-grid->cell0*grid->na_sc)>>grid->na_c_2log) & (STRIDE_PBB-1));
1178 #ifdef NBNXN_SEARCH_SIMD4_FLOAT_X_BB
1179 if (nbat->XFormat == nbatXYZQ)
1181 calc_bounding_box_xxxx_simd4(na, nbat->x+a0*nbat->xstride,
1182 bb_work_aligned, pbb_ptr);
1187 calc_bounding_box_xxxx(na, nbat->xstride, nbat->x+a0*nbat->xstride,
1192 fprintf(debug, "%2d %2d %2d bb %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f\n",
1194 pbb_ptr[0*STRIDE_PBB], pbb_ptr[3*STRIDE_PBB],
1195 pbb_ptr[1*STRIDE_PBB], pbb_ptr[4*STRIDE_PBB],
1196 pbb_ptr[2*STRIDE_PBB], pbb_ptr[5*STRIDE_PBB]);
1202 /* Store the bounding boxes as xyz.xyz. */
1203 bb_ptr = grid->bb+((a0-grid->cell0*grid->na_sc)>>grid->na_c_2log);
1205 calc_bounding_box(na, nbat->xstride, nbat->x+a0*nbat->xstride,
1211 bbo = (a0 - grid->cell0*grid->na_sc)/grid->na_c;
1212 fprintf(debug, "%2d %2d %2d bb %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f\n",
1214 grid->bb[bbo].lower[BB_X],
1215 grid->bb[bbo].lower[BB_Y],
1216 grid->bb[bbo].lower[BB_Z],
1217 grid->bb[bbo].upper[BB_X],
1218 grid->bb[bbo].upper[BB_Y],
1219 grid->bb[bbo].upper[BB_Z]);
1224 /* Spatially sort the atoms within one grid column */
1225 static void sort_columns_simple(const nbnxn_search_t nbs,
1231 nbnxn_atomdata_t *nbat,
1232 int cxy_start, int cxy_end,
1236 int cx, cy, cz, ncz, cfilled, c;
1237 int na, ash, ind, a;
1242 fprintf(debug, "cell0 %d sorting columns %d - %d, atoms %d - %d\n",
1243 grid->cell0, cxy_start, cxy_end, a0, a1);
1246 /* Sort the atoms within each x,y column in 3 dimensions */
1247 for (cxy = cxy_start; cxy < cxy_end; cxy++)
1250 cy = cxy - cx*grid->ncy;
1252 na = grid->cxy_na[cxy];
1253 ncz = grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy];
1254 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1256 /* Sort the atoms within each x,y column on z coordinate */
1257 sort_atoms(ZZ, FALSE,
1260 1.0/nbs->box[ZZ][ZZ], ncz*grid->na_sc,
1263 /* Fill the ncz cells in this column */
1264 cfilled = grid->cxy_ind[cxy];
1265 for (cz = 0; cz < ncz; cz++)
1267 c = grid->cxy_ind[cxy] + cz;
1269 ash_c = ash + cz*grid->na_sc;
1270 na_c = min(grid->na_sc, na-(ash_c-ash));
1272 fill_cell(nbs, grid, nbat,
1273 ash_c, ash_c+na_c, atinfo, x,
1274 grid->na_sc*cx + (dd_zone >> 2),
1275 grid->na_sc*cy + (dd_zone & 3),
1279 /* This copy to bbcz is not really necessary.
1280 * But it allows to use the same grid search code
1281 * for the simple and supersub cell setups.
1287 grid->bbcz[c*NNBSBB_D ] = grid->bb[cfilled].lower[BB_Z];
1288 grid->bbcz[c*NNBSBB_D+1] = grid->bb[cfilled].upper[BB_Z];
1291 /* Set the unused atom indices to -1 */
1292 for (ind = na; ind < ncz*grid->na_sc; ind++)
1294 nbs->a[ash+ind] = -1;
1299 /* Spatially sort the atoms within one grid column */
1300 static void sort_columns_supersub(const nbnxn_search_t nbs,
1306 nbnxn_atomdata_t *nbat,
1307 int cxy_start, int cxy_end,
1311 int cx, cy, cz = -1, c = -1, ncz;
1312 int na, ash, na_c, ind, a;
1313 int subdiv_z, sub_z, na_z, ash_z;
1314 int subdiv_y, sub_y, na_y, ash_y;
1315 int subdiv_x, sub_x, na_x, ash_x;
1317 /* cppcheck-suppress unassignedVariable */
1318 nbnxn_bb_t bb_work_array[2], *bb_work_aligned;
1320 bb_work_aligned = (nbnxn_bb_t *)(((size_t)(bb_work_array+1)) & (~((size_t)15)));
1324 fprintf(debug, "cell0 %d sorting columns %d - %d, atoms %d - %d\n",
1325 grid->cell0, cxy_start, cxy_end, a0, a1);
1328 subdiv_x = grid->na_c;
1329 subdiv_y = GPU_NSUBCELL_X*subdiv_x;
1330 subdiv_z = GPU_NSUBCELL_Y*subdiv_y;
1332 /* Sort the atoms within each x,y column in 3 dimensions */
1333 for (cxy = cxy_start; cxy < cxy_end; cxy++)
1336 cy = cxy - cx*grid->ncy;
1338 na = grid->cxy_na[cxy];
1339 ncz = grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy];
1340 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1342 /* Sort the atoms within each x,y column on z coordinate */
1343 sort_atoms(ZZ, FALSE,
1346 1.0/nbs->box[ZZ][ZZ], ncz*grid->na_sc,
1349 /* This loop goes over the supercells and subcells along z at once */
1350 for (sub_z = 0; sub_z < ncz*GPU_NSUBCELL_Z; sub_z++)
1352 ash_z = ash + sub_z*subdiv_z;
1353 na_z = min(subdiv_z, na-(ash_z-ash));
1355 /* We have already sorted on z */
1357 if (sub_z % GPU_NSUBCELL_Z == 0)
1359 cz = sub_z/GPU_NSUBCELL_Z;
1360 c = grid->cxy_ind[cxy] + cz;
1362 /* The number of atoms in this supercell */
1363 na_c = min(grid->na_sc, na-(ash_z-ash));
1365 grid->nsubc[c] = min(GPU_NSUBCELL, (na_c+grid->na_c-1)/grid->na_c);
1367 /* Store the z-boundaries of the super cell */
1368 grid->bbcz[c*NNBSBB_D ] = x[nbs->a[ash_z]][ZZ];
1369 grid->bbcz[c*NNBSBB_D+1] = x[nbs->a[ash_z+na_c-1]][ZZ];
1372 #if GPU_NSUBCELL_Y > 1
1373 /* Sort the atoms along y */
1374 sort_atoms(YY, (sub_z & 1),
1375 nbs->a+ash_z, na_z, x,
1376 grid->c0[YY]+cy*grid->sy,
1377 grid->inv_sy, subdiv_z,
1381 for (sub_y = 0; sub_y < GPU_NSUBCELL_Y; sub_y++)
1383 ash_y = ash_z + sub_y*subdiv_y;
1384 na_y = min(subdiv_y, na-(ash_y-ash));
1386 #if GPU_NSUBCELL_X > 1
1387 /* Sort the atoms along x */
1388 sort_atoms(XX, ((cz*GPU_NSUBCELL_Y + sub_y) & 1),
1389 nbs->a+ash_y, na_y, x,
1390 grid->c0[XX]+cx*grid->sx,
1391 grid->inv_sx, subdiv_y,
1395 for (sub_x = 0; sub_x < GPU_NSUBCELL_X; sub_x++)
1397 ash_x = ash_y + sub_x*subdiv_x;
1398 na_x = min(subdiv_x, na-(ash_x-ash));
1400 fill_cell(nbs, grid, nbat,
1401 ash_x, ash_x+na_x, atinfo, x,
1402 grid->na_c*(cx*GPU_NSUBCELL_X+sub_x) + (dd_zone >> 2),
1403 grid->na_c*(cy*GPU_NSUBCELL_Y+sub_y) + (dd_zone & 3),
1410 /* Set the unused atom indices to -1 */
1411 for (ind = na; ind < ncz*grid->na_sc; ind++)
1413 nbs->a[ash+ind] = -1;
1418 /* Determine in which grid column atoms should go */
1419 static void calc_column_indices(nbnxn_grid_t *grid,
1422 int dd_zone, const int *move,
1423 int thread, int nthread,
1430 /* We add one extra cell for particles which moved during DD */
1431 for (i = 0; i < grid->ncx*grid->ncy+1; i++)
1436 n0 = a0 + (int)((thread+0)*(a1 - a0))/nthread;
1437 n1 = a0 + (int)((thread+1)*(a1 - a0))/nthread;
1441 for (i = n0; i < n1; i++)
1443 if (move == NULL || move[i] >= 0)
1445 /* We need to be careful with rounding,
1446 * particles might be a few bits outside the local zone.
1447 * The int cast takes care of the lower bound,
1448 * we will explicitly take care of the upper bound.
1450 cx = (int)((x[i][XX] - grid->c0[XX])*grid->inv_sx);
1451 cy = (int)((x[i][YY] - grid->c0[YY])*grid->inv_sy);
1454 if (cx < 0 || cx > grid->ncx ||
1455 cy < 0 || cy > grid->ncy)
1458 "grid cell cx %d cy %d out of range (max %d %d)\n"
1459 "atom %f %f %f, grid->c0 %f %f",
1460 cx, cy, grid->ncx, grid->ncy,
1461 x[i][XX], x[i][YY], x[i][ZZ], grid->c0[XX], grid->c0[YY]);
1464 /* Take care of potential rouding issues */
1465 cx = min(cx, grid->ncx - 1);
1466 cy = min(cy, grid->ncy - 1);
1468 /* For the moment cell will contain only the, grid local,
1469 * x and y indices, not z.
1471 cell[i] = cx*grid->ncy + cy;
1475 /* Put this moved particle after the end of the grid,
1476 * so we can process it later without using conditionals.
1478 cell[i] = grid->ncx*grid->ncy;
1487 for (i = n0; i < n1; i++)
1489 cx = (int)((x[i][XX] - grid->c0[XX])*grid->inv_sx);
1490 cy = (int)((x[i][YY] - grid->c0[YY])*grid->inv_sy);
1492 /* For non-home zones there could be particles outside
1493 * the non-bonded cut-off range, which have been communicated
1494 * for bonded interactions only. For the result it doesn't
1495 * matter where these end up on the grid. For performance
1496 * we put them in an extra row at the border.
1499 cx = min(cx, grid->ncx - 1);
1501 cy = min(cy, grid->ncy - 1);
1503 /* For the moment cell will contain only the, grid local,
1504 * x and y indices, not z.
1506 cell[i] = cx*grid->ncy + cy;
1513 /* Determine in which grid cells the atoms should go */
1514 static void calc_cell_indices(const nbnxn_search_t nbs,
1521 nbnxn_atomdata_t *nbat)
1524 int cx, cy, cxy, ncz_max, ncz;
1525 int nthread, thread;
1526 int *cxy_na, cxy_na_i;
1528 nthread = gmx_omp_nthreads_get(emntPairsearch);
1530 #pragma omp parallel for num_threads(nthread) schedule(static)
1531 for (thread = 0; thread < nthread; thread++)
1533 calc_column_indices(grid, a0, a1, x, dd_zone, move, thread, nthread,
1534 nbs->cell, nbs->work[thread].cxy_na);
1537 /* Make the cell index as a function of x and y */
1540 grid->cxy_ind[0] = 0;
1541 for (i = 0; i < grid->ncx*grid->ncy+1; i++)
1543 /* We set ncz_max at the beginning of the loop iso at the end
1544 * to skip i=grid->ncx*grid->ncy which are moved particles
1545 * that do not need to be ordered on the grid.
1551 cxy_na_i = nbs->work[0].cxy_na[i];
1552 for (thread = 1; thread < nthread; thread++)
1554 cxy_na_i += nbs->work[thread].cxy_na[i];
1556 ncz = (cxy_na_i + grid->na_sc - 1)/grid->na_sc;
1557 if (nbat->XFormat == nbatX8)
1559 /* Make the number of cell a multiple of 2 */
1560 ncz = (ncz + 1) & ~1;
1562 grid->cxy_ind[i+1] = grid->cxy_ind[i] + ncz;
1563 /* Clear cxy_na, so we can reuse the array below */
1564 grid->cxy_na[i] = 0;
1566 grid->nc = grid->cxy_ind[grid->ncx*grid->ncy] - grid->cxy_ind[0];
1568 nbat->natoms = (grid->cell0 + grid->nc)*grid->na_sc;
1572 fprintf(debug, "ns na_sc %d na_c %d super-cells: %d x %d y %d z %.1f maxz %d\n",
1573 grid->na_sc, grid->na_c, grid->nc,
1574 grid->ncx, grid->ncy, grid->nc/((double)(grid->ncx*grid->ncy)),
1579 for (cy = 0; cy < grid->ncy; cy++)
1581 for (cx = 0; cx < grid->ncx; cx++)
1583 fprintf(debug, " %2d", grid->cxy_ind[i+1]-grid->cxy_ind[i]);
1586 fprintf(debug, "\n");
1591 /* Make sure the work array for sorting is large enough */
1592 if (ncz_max*grid->na_sc*SGSF > nbs->work[0].sort_work_nalloc)
1594 for (thread = 0; thread < nbs->nthread_max; thread++)
1596 nbs->work[thread].sort_work_nalloc =
1597 over_alloc_large(ncz_max*grid->na_sc*SGSF);
1598 srenew(nbs->work[thread].sort_work,
1599 nbs->work[thread].sort_work_nalloc);
1600 /* When not in use, all elements should be -1 */
1601 for (i = 0; i < nbs->work[thread].sort_work_nalloc; i++)
1603 nbs->work[thread].sort_work[i] = -1;
1608 /* Now we know the dimensions we can fill the grid.
1609 * This is the first, unsorted fill. We sort the columns after this.
1611 for (i = a0; i < a1; i++)
1613 /* At this point nbs->cell contains the local grid x,y indices */
1615 nbs->a[(grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc + grid->cxy_na[cxy]++] = i;
1620 /* Set the cell indices for the moved particles */
1621 n0 = grid->nc*grid->na_sc;
1622 n1 = grid->nc*grid->na_sc+grid->cxy_na[grid->ncx*grid->ncy];
1625 for (i = n0; i < n1; i++)
1627 nbs->cell[nbs->a[i]] = i;
1632 /* Sort the super-cell columns along z into the sub-cells. */
1633 #pragma omp parallel for num_threads(nbs->nthread_max) schedule(static)
1634 for (thread = 0; thread < nbs->nthread_max; thread++)
1638 sort_columns_simple(nbs, dd_zone, grid, a0, a1, atinfo, x, nbat,
1639 ((thread+0)*grid->ncx*grid->ncy)/nthread,
1640 ((thread+1)*grid->ncx*grid->ncy)/nthread,
1641 nbs->work[thread].sort_work);
1645 sort_columns_supersub(nbs, dd_zone, grid, a0, a1, atinfo, x, nbat,
1646 ((thread+0)*grid->ncx*grid->ncy)/nthread,
1647 ((thread+1)*grid->ncx*grid->ncy)/nthread,
1648 nbs->work[thread].sort_work);
1652 if (grid->bSimple && nbat->XFormat == nbatX8)
1654 combine_bounding_box_pairs(grid, grid->bb);
1659 grid->nsubc_tot = 0;
1660 for (i = 0; i < grid->nc; i++)
1662 grid->nsubc_tot += grid->nsubc[i];
1670 print_bbsizes_simple(debug, nbs, grid);
1674 fprintf(debug, "ns non-zero sub-cells: %d average atoms %.2f\n",
1675 grid->nsubc_tot, (a1-a0)/(double)grid->nsubc_tot);
1677 print_bbsizes_supersub(debug, nbs, grid);
1682 static void init_buffer_flags(nbnxn_buffer_flags_t *flags,
1687 flags->nflag = (natoms + NBNXN_BUFFERFLAG_SIZE - 1)/NBNXN_BUFFERFLAG_SIZE;
1688 if (flags->nflag > flags->flag_nalloc)
1690 flags->flag_nalloc = over_alloc_large(flags->nflag);
1691 srenew(flags->flag, flags->flag_nalloc);
1693 for (b = 0; b < flags->nflag; b++)
1699 /* Sets up a grid and puts the atoms on the grid.
1700 * This function only operates on one domain of the domain decompostion.
1701 * Note that without domain decomposition there is only one domain.
1703 void nbnxn_put_on_grid(nbnxn_search_t nbs,
1704 int ePBC, matrix box,
1706 rvec corner0, rvec corner1,
1711 int nmoved, int *move,
1713 nbnxn_atomdata_t *nbat)
1717 int nc_max_grid, nc_max;
1719 grid = &nbs->grid[dd_zone];
1721 nbs_cycle_start(&nbs->cc[enbsCCgrid]);
1723 grid->bSimple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
1725 grid->na_c = nbnxn_kernel_to_ci_size(nb_kernel_type);
1726 grid->na_cj = nbnxn_kernel_to_cj_size(nb_kernel_type);
1727 grid->na_sc = (grid->bSimple ? 1 : GPU_NSUBCELL)*grid->na_c;
1728 grid->na_c_2log = get_2log(grid->na_c);
1730 nbat->na_c = grid->na_c;
1739 (nbs->grid[dd_zone-1].cell0 + nbs->grid[dd_zone-1].nc)*
1740 nbs->grid[dd_zone-1].na_sc/grid->na_sc;
1748 copy_mat(box, nbs->box);
1750 if (atom_density >= 0)
1752 grid->atom_density = atom_density;
1756 grid->atom_density = grid_atom_density(n-nmoved, corner0, corner1);
1761 nbs->natoms_local = a1 - nmoved;
1762 /* We assume that nbnxn_put_on_grid is called first
1763 * for the local atoms (dd_zone=0).
1765 nbs->natoms_nonlocal = a1 - nmoved;
1769 nbs->natoms_nonlocal = max(nbs->natoms_nonlocal, a1);
1772 nc_max_grid = set_grid_size_xy(nbs, grid,
1773 dd_zone, n-nmoved, corner0, corner1,
1774 nbs->grid[0].atom_density);
1776 nc_max = grid->cell0 + nc_max_grid;
1778 if (a1 > nbs->cell_nalloc)
1780 nbs->cell_nalloc = over_alloc_large(a1);
1781 srenew(nbs->cell, nbs->cell_nalloc);
1784 /* To avoid conditionals we store the moved particles at the end of a,
1785 * make sure we have enough space.
1787 if (nc_max*grid->na_sc + nmoved > nbs->a_nalloc)
1789 nbs->a_nalloc = over_alloc_large(nc_max*grid->na_sc + nmoved);
1790 srenew(nbs->a, nbs->a_nalloc);
1793 /* We need padding up to a multiple of the buffer flag size: simply add */
1794 if (nc_max*grid->na_sc + NBNXN_BUFFERFLAG_SIZE > nbat->nalloc)
1796 nbnxn_atomdata_realloc(nbat, nc_max*grid->na_sc+NBNXN_BUFFERFLAG_SIZE);
1799 calc_cell_indices(nbs, dd_zone, grid, a0, a1, atinfo, x, move, nbat);
1803 nbat->natoms_local = nbat->natoms;
1806 nbs_cycle_stop(&nbs->cc[enbsCCgrid]);
1809 /* Calls nbnxn_put_on_grid for all non-local domains */
1810 void nbnxn_put_on_grid_nonlocal(nbnxn_search_t nbs,
1811 const gmx_domdec_zones_t *zones,
1815 nbnxn_atomdata_t *nbat)
1820 for (zone = 1; zone < zones->n; zone++)
1822 for (d = 0; d < DIM; d++)
1824 c0[d] = zones->size[zone].bb_x0[d];
1825 c1[d] = zones->size[zone].bb_x1[d];
1828 nbnxn_put_on_grid(nbs, nbs->ePBC, NULL,
1830 zones->cg_range[zone],
1831 zones->cg_range[zone+1],
1841 /* Add simple grid type information to the local super/sub grid */
1842 void nbnxn_grid_add_simple(nbnxn_search_t nbs,
1843 nbnxn_atomdata_t *nbat)
1850 grid = &nbs->grid[0];
1854 gmx_incons("nbnxn_grid_simple called with a simple grid");
1857 ncd = grid->na_sc/NBNXN_CPU_CLUSTER_I_SIZE;
1859 if (grid->nc*ncd > grid->nc_nalloc_simple)
1861 grid->nc_nalloc_simple = over_alloc_large(grid->nc*ncd);
1862 srenew(grid->bbcz_simple, grid->nc_nalloc_simple*NNBSBB_D);
1863 srenew(grid->bb_simple, grid->nc_nalloc_simple);
1864 srenew(grid->flags_simple, grid->nc_nalloc_simple);
1867 sfree_aligned(grid->bbj);
1868 snew_aligned(grid->bbj, grid->nc_nalloc_simple/2, 16);
1872 bbcz = grid->bbcz_simple;
1873 bb = grid->bb_simple;
1875 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntPairsearch)) schedule(static)
1876 for (sc = 0; sc < grid->nc; sc++)
1880 for (c = 0; c < ncd; c++)
1884 na = NBNXN_CPU_CLUSTER_I_SIZE;
1886 nbat->type[tx*NBNXN_CPU_CLUSTER_I_SIZE+na-1] == nbat->ntype-1)
1893 switch (nbat->XFormat)
1896 /* PACK_X4==NBNXN_CPU_CLUSTER_I_SIZE, so this is simple */
1897 calc_bounding_box_x_x4(na, nbat->x+tx*STRIDE_P4,
1901 /* PACK_X8>NBNXN_CPU_CLUSTER_I_SIZE, more complicated */
1902 calc_bounding_box_x_x8(na, nbat->x+X8_IND_A(tx*NBNXN_CPU_CLUSTER_I_SIZE),
1906 calc_bounding_box(na, nbat->xstride,
1907 nbat->x+tx*NBNXN_CPU_CLUSTER_I_SIZE*nbat->xstride,
1911 bbcz[tx*NNBSBB_D+0] = bb[tx].lower[BB_Z];
1912 bbcz[tx*NNBSBB_D+1] = bb[tx].upper[BB_Z];
1914 /* No interaction optimization yet here */
1915 grid->flags_simple[tx] = NBNXN_CI_DO_LJ(0) | NBNXN_CI_DO_COUL(0);
1919 grid->flags_simple[tx] = 0;
1924 if (grid->bSimple && nbat->XFormat == nbatX8)
1926 combine_bounding_box_pairs(grid, grid->bb_simple);
1930 void nbnxn_get_ncells(nbnxn_search_t nbs, int *ncx, int *ncy)
1932 *ncx = nbs->grid[0].ncx;
1933 *ncy = nbs->grid[0].ncy;
1936 void nbnxn_get_atomorder(nbnxn_search_t nbs, int **a, int *n)
1938 const nbnxn_grid_t *grid;
1940 grid = &nbs->grid[0];
1942 /* Return the atom order for the home cell (index 0) */
1945 *n = grid->cxy_ind[grid->ncx*grid->ncy]*grid->na_sc;
1948 void nbnxn_set_atomorder(nbnxn_search_t nbs)
1951 int ao, cx, cy, cxy, cz, j;
1953 /* Set the atom order for the home cell (index 0) */
1954 grid = &nbs->grid[0];
1957 for (cx = 0; cx < grid->ncx; cx++)
1959 for (cy = 0; cy < grid->ncy; cy++)
1961 cxy = cx*grid->ncy + cy;
1962 j = grid->cxy_ind[cxy]*grid->na_sc;
1963 for (cz = 0; cz < grid->cxy_na[cxy]; cz++)
1974 /* Determines the cell range along one dimension that
1975 * the bounding box b0 - b1 sees.
1977 static void get_cell_range(real b0, real b1,
1978 int nc, real c0, real s, real invs,
1979 real d2, real r2, int *cf, int *cl)
1981 *cf = max((int)((b0 - c0)*invs), 0);
1983 while (*cf > 0 && d2 + sqr((b0 - c0) - (*cf-1+1)*s) < r2)
1988 *cl = min((int)((b1 - c0)*invs), nc-1);
1989 while (*cl < nc-1 && d2 + sqr((*cl+1)*s - (b1 - c0)) < r2)
1995 /* Reference code calculating the distance^2 between two bounding boxes */
1996 static float box_dist2(float bx0, float bx1, float by0,
1997 float by1, float bz0, float bz1,
1998 const nbnxn_bb_t *bb)
2001 float dl, dh, dm, dm0;
2005 dl = bx0 - bb->upper[BB_X];
2006 dh = bb->lower[BB_X] - bx1;
2011 dl = by0 - bb->upper[BB_Y];
2012 dh = bb->lower[BB_Y] - by1;
2017 dl = bz0 - bb->upper[BB_Z];
2018 dh = bb->lower[BB_Z] - bz1;
2026 /* Plain C code calculating the distance^2 between two bounding boxes */
2027 static float subc_bb_dist2(int si, const nbnxn_bb_t *bb_i_ci,
2028 int csj, const nbnxn_bb_t *bb_j_all)
2030 const nbnxn_bb_t *bb_i, *bb_j;
2032 float dl, dh, dm, dm0;
2034 bb_i = bb_i_ci + si;
2035 bb_j = bb_j_all + csj;
2039 dl = bb_i->lower[BB_X] - bb_j->upper[BB_X];
2040 dh = bb_j->lower[BB_X] - bb_i->upper[BB_X];
2045 dl = bb_i->lower[BB_Y] - bb_j->upper[BB_Y];
2046 dh = bb_j->lower[BB_Y] - bb_i->upper[BB_Y];
2051 dl = bb_i->lower[BB_Z] - bb_j->upper[BB_Z];
2052 dh = bb_j->lower[BB_Z] - bb_i->upper[BB_Z];
2060 #ifdef NBNXN_SEARCH_BB_SIMD4
2062 /* 4-wide SIMD code for bb distance for bb format xyz0 */
2063 static float subc_bb_dist2_simd4(int si, const nbnxn_bb_t *bb_i_ci,
2064 int csj, const nbnxn_bb_t *bb_j_all)
2066 gmx_simd4_pr bb_i_S0, bb_i_S1;
2067 gmx_simd4_pr bb_j_S0, bb_j_S1;
2073 bb_i_S0 = gmx_simd4_load_bb_pr(&bb_i_ci[si].lower[0]);
2074 bb_i_S1 = gmx_simd4_load_bb_pr(&bb_i_ci[si].upper[0]);
2075 bb_j_S0 = gmx_simd4_load_bb_pr(&bb_j_all[csj].lower[0]);
2076 bb_j_S1 = gmx_simd4_load_bb_pr(&bb_j_all[csj].upper[0]);
2078 dl_S = gmx_simd4_sub_pr(bb_i_S0, bb_j_S1);
2079 dh_S = gmx_simd4_sub_pr(bb_j_S0, bb_i_S1);
2081 dm_S = gmx_simd4_max_pr(dl_S, dh_S);
2082 dm0_S = gmx_simd4_max_pr(dm_S, gmx_simd4_setzero_pr());
2084 return gmx_simd4_dotproduct3(dm0_S, dm0_S);
2087 /* Calculate bb bounding distances of bb_i[si,...,si+3] and store them in d2 */
2088 #define SUBC_BB_DIST2_SIMD4_XXXX_INNER(si, bb_i, d2) \
2092 gmx_simd4_pr dx_0, dy_0, dz_0; \
2093 gmx_simd4_pr dx_1, dy_1, dz_1; \
2095 gmx_simd4_pr mx, my, mz; \
2096 gmx_simd4_pr m0x, m0y, m0z; \
2098 gmx_simd4_pr d2x, d2y, d2z; \
2099 gmx_simd4_pr d2s, d2t; \
2101 shi = si*NNBSBB_D*DIM; \
2103 xi_l = gmx_simd4_load_bb_pr(bb_i+shi+0*STRIDE_PBB); \
2104 yi_l = gmx_simd4_load_bb_pr(bb_i+shi+1*STRIDE_PBB); \
2105 zi_l = gmx_simd4_load_bb_pr(bb_i+shi+2*STRIDE_PBB); \
2106 xi_h = gmx_simd4_load_bb_pr(bb_i+shi+3*STRIDE_PBB); \
2107 yi_h = gmx_simd4_load_bb_pr(bb_i+shi+4*STRIDE_PBB); \
2108 zi_h = gmx_simd4_load_bb_pr(bb_i+shi+5*STRIDE_PBB); \
2110 dx_0 = gmx_simd4_sub_pr(xi_l, xj_h); \
2111 dy_0 = gmx_simd4_sub_pr(yi_l, yj_h); \
2112 dz_0 = gmx_simd4_sub_pr(zi_l, zj_h); \
2114 dx_1 = gmx_simd4_sub_pr(xj_l, xi_h); \
2115 dy_1 = gmx_simd4_sub_pr(yj_l, yi_h); \
2116 dz_1 = gmx_simd4_sub_pr(zj_l, zi_h); \
2118 mx = gmx_simd4_max_pr(dx_0, dx_1); \
2119 my = gmx_simd4_max_pr(dy_0, dy_1); \
2120 mz = gmx_simd4_max_pr(dz_0, dz_1); \
2122 m0x = gmx_simd4_max_pr(mx, zero); \
2123 m0y = gmx_simd4_max_pr(my, zero); \
2124 m0z = gmx_simd4_max_pr(mz, zero); \
2126 d2x = gmx_simd4_mul_pr(m0x, m0x); \
2127 d2y = gmx_simd4_mul_pr(m0y, m0y); \
2128 d2z = gmx_simd4_mul_pr(m0z, m0z); \
2130 d2s = gmx_simd4_add_pr(d2x, d2y); \
2131 d2t = gmx_simd4_add_pr(d2s, d2z); \
2133 gmx_simd4_store_pr(d2+si, d2t); \
2136 /* 4-wide SIMD code for nsi bb distances for bb format xxxxyyyyzzzz */
2137 static void subc_bb_dist2_simd4_xxxx(const float *bb_j,
2138 int nsi, const float *bb_i,
2141 gmx_simd4_pr xj_l, yj_l, zj_l;
2142 gmx_simd4_pr xj_h, yj_h, zj_h;
2143 gmx_simd4_pr xi_l, yi_l, zi_l;
2144 gmx_simd4_pr xi_h, yi_h, zi_h;
2148 zero = gmx_simd4_setzero_pr();
2150 xj_l = gmx_simd4_set1_pr(bb_j[0*STRIDE_PBB]);
2151 yj_l = gmx_simd4_set1_pr(bb_j[1*STRIDE_PBB]);
2152 zj_l = gmx_simd4_set1_pr(bb_j[2*STRIDE_PBB]);
2153 xj_h = gmx_simd4_set1_pr(bb_j[3*STRIDE_PBB]);
2154 yj_h = gmx_simd4_set1_pr(bb_j[4*STRIDE_PBB]);
2155 zj_h = gmx_simd4_set1_pr(bb_j[5*STRIDE_PBB]);
2157 /* Here we "loop" over si (0,STRIDE_PBB) from 0 to nsi with step STRIDE_PBB.
2158 * But as we know the number of iterations is 1 or 2, we unroll manually.
2160 SUBC_BB_DIST2_SIMD4_XXXX_INNER(0, bb_i, d2);
2161 if (STRIDE_PBB < nsi)
2163 SUBC_BB_DIST2_SIMD4_XXXX_INNER(STRIDE_PBB, bb_i, d2);
2167 #endif /* NBNXN_SEARCH_BB_SIMD4 */
2169 /* Plain C function which determines if any atom pair between two cells
2170 * is within distance sqrt(rl2).
2172 static gmx_bool subc_in_range_x(int na_c,
2173 int si, const real *x_i,
2174 int csj, int stride, const real *x_j,
2180 for (i = 0; i < na_c; i++)
2182 i0 = (si*na_c + i)*DIM;
2183 for (j = 0; j < na_c; j++)
2185 j0 = (csj*na_c + j)*stride;
2187 d2 = sqr(x_i[i0 ] - x_j[j0 ]) +
2188 sqr(x_i[i0+1] - x_j[j0+1]) +
2189 sqr(x_i[i0+2] - x_j[j0+2]);
2201 #ifdef NBNXN_SEARCH_SIMD4_FLOAT_X_BB
2202 /* When we make seperate single/double precision SIMD vector operation
2203 * include files, this function should be moved there (also using FMA).
2205 static inline gmx_simd4_pr
2206 gmx_simd4_calc_rsq_pr(gmx_simd4_pr x, gmx_simd4_pr y, gmx_simd4_pr z)
2208 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) );
2211 /* 4-wide SIMD function which determines if any atom pair between two cells,
2212 * both with 8 atoms, is within distance sqrt(rl2).
2213 * Using 8-wide AVX is not faster on Intel Sandy Bridge.
2215 static gmx_bool subc_in_range_simd4(int na_c,
2216 int si, const real *x_i,
2217 int csj, int stride, const real *x_j,
2220 gmx_simd4_pr ix_S0, iy_S0, iz_S0;
2221 gmx_simd4_pr ix_S1, iy_S1, iz_S1;
2228 rc2_S = gmx_simd4_set1_pr(rl2);
2230 dim_stride = NBNXN_GPU_CLUSTER_SIZE/STRIDE_PBB*DIM;
2231 ix_S0 = gmx_simd4_load_bb_pr(x_i+(si*dim_stride+0)*STRIDE_PBB);
2232 iy_S0 = gmx_simd4_load_bb_pr(x_i+(si*dim_stride+1)*STRIDE_PBB);
2233 iz_S0 = gmx_simd4_load_bb_pr(x_i+(si*dim_stride+2)*STRIDE_PBB);
2234 ix_S1 = gmx_simd4_load_bb_pr(x_i+(si*dim_stride+3)*STRIDE_PBB);
2235 iy_S1 = gmx_simd4_load_bb_pr(x_i+(si*dim_stride+4)*STRIDE_PBB);
2236 iz_S1 = gmx_simd4_load_bb_pr(x_i+(si*dim_stride+5)*STRIDE_PBB);
2238 /* We loop from the outer to the inner particles to maximize
2239 * the chance that we find a pair in range quickly and return.
2245 gmx_simd4_pr jx0_S, jy0_S, jz0_S;
2246 gmx_simd4_pr jx1_S, jy1_S, jz1_S;
2248 gmx_simd4_pr dx_S0, dy_S0, dz_S0;
2249 gmx_simd4_pr dx_S1, dy_S1, dz_S1;
2250 gmx_simd4_pr dx_S2, dy_S2, dz_S2;
2251 gmx_simd4_pr dx_S3, dy_S3, dz_S3;
2253 gmx_simd4_pr rsq_S0;
2254 gmx_simd4_pr rsq_S1;
2255 gmx_simd4_pr rsq_S2;
2256 gmx_simd4_pr rsq_S3;
2258 gmx_simd4_pb wco_S0;
2259 gmx_simd4_pb wco_S1;
2260 gmx_simd4_pb wco_S2;
2261 gmx_simd4_pb wco_S3;
2262 gmx_simd4_pb wco_any_S01, wco_any_S23, wco_any_S;
2264 jx0_S = gmx_simd4_set1_pr(x_j[j0*stride+0]);
2265 jy0_S = gmx_simd4_set1_pr(x_j[j0*stride+1]);
2266 jz0_S = gmx_simd4_set1_pr(x_j[j0*stride+2]);
2268 jx1_S = gmx_simd4_set1_pr(x_j[j1*stride+0]);
2269 jy1_S = gmx_simd4_set1_pr(x_j[j1*stride+1]);
2270 jz1_S = gmx_simd4_set1_pr(x_j[j1*stride+2]);
2272 /* Calculate distance */
2273 dx_S0 = gmx_simd4_sub_pr(ix_S0, jx0_S);
2274 dy_S0 = gmx_simd4_sub_pr(iy_S0, jy0_S);
2275 dz_S0 = gmx_simd4_sub_pr(iz_S0, jz0_S);
2276 dx_S1 = gmx_simd4_sub_pr(ix_S1, jx0_S);
2277 dy_S1 = gmx_simd4_sub_pr(iy_S1, jy0_S);
2278 dz_S1 = gmx_simd4_sub_pr(iz_S1, jz0_S);
2279 dx_S2 = gmx_simd4_sub_pr(ix_S0, jx1_S);
2280 dy_S2 = gmx_simd4_sub_pr(iy_S0, jy1_S);
2281 dz_S2 = gmx_simd4_sub_pr(iz_S0, jz1_S);
2282 dx_S3 = gmx_simd4_sub_pr(ix_S1, jx1_S);
2283 dy_S3 = gmx_simd4_sub_pr(iy_S1, jy1_S);
2284 dz_S3 = gmx_simd4_sub_pr(iz_S1, jz1_S);
2286 /* rsq = dx*dx+dy*dy+dz*dz */
2287 rsq_S0 = gmx_simd4_calc_rsq_pr(dx_S0, dy_S0, dz_S0);
2288 rsq_S1 = gmx_simd4_calc_rsq_pr(dx_S1, dy_S1, dz_S1);
2289 rsq_S2 = gmx_simd4_calc_rsq_pr(dx_S2, dy_S2, dz_S2);
2290 rsq_S3 = gmx_simd4_calc_rsq_pr(dx_S3, dy_S3, dz_S3);
2292 wco_S0 = gmx_simd4_cmplt_pr(rsq_S0, rc2_S);
2293 wco_S1 = gmx_simd4_cmplt_pr(rsq_S1, rc2_S);
2294 wco_S2 = gmx_simd4_cmplt_pr(rsq_S2, rc2_S);
2295 wco_S3 = gmx_simd4_cmplt_pr(rsq_S3, rc2_S);
2297 wco_any_S01 = gmx_simd4_or_pb(wco_S0, wco_S1);
2298 wco_any_S23 = gmx_simd4_or_pb(wco_S2, wco_S3);
2299 wco_any_S = gmx_simd4_or_pb(wco_any_S01, wco_any_S23);
2301 if (gmx_simd4_anytrue_pb(wco_any_S))
2315 /* Returns the j sub-cell for index cj_ind */
2316 static int nbl_cj(const nbnxn_pairlist_t *nbl, int cj_ind)
2318 return nbl->cj4[cj_ind >> NBNXN_GPU_JGROUP_SIZE_2LOG].cj[cj_ind & (NBNXN_GPU_JGROUP_SIZE - 1)];
2321 /* Returns the i-interaction mask of the j sub-cell for index cj_ind */
2322 static unsigned nbl_imask0(const nbnxn_pairlist_t *nbl, int cj_ind)
2324 return nbl->cj4[cj_ind >> NBNXN_GPU_JGROUP_SIZE_2LOG].imei[0].imask;
2327 /* Ensures there is enough space for extra extra exclusion masks */
2328 static void check_excl_space(nbnxn_pairlist_t *nbl, int extra)
2330 if (nbl->nexcl+extra > nbl->excl_nalloc)
2332 nbl->excl_nalloc = over_alloc_small(nbl->nexcl+extra);
2333 nbnxn_realloc_void((void **)&nbl->excl,
2334 nbl->nexcl*sizeof(*nbl->excl),
2335 nbl->excl_nalloc*sizeof(*nbl->excl),
2336 nbl->alloc, nbl->free);
2340 /* Ensures there is enough space for ncell extra j-cells in the list */
2341 static void check_subcell_list_space_simple(nbnxn_pairlist_t *nbl,
2346 cj_max = nbl->ncj + ncell;
2348 if (cj_max > nbl->cj_nalloc)
2350 nbl->cj_nalloc = over_alloc_small(cj_max);
2351 nbnxn_realloc_void((void **)&nbl->cj,
2352 nbl->ncj*sizeof(*nbl->cj),
2353 nbl->cj_nalloc*sizeof(*nbl->cj),
2354 nbl->alloc, nbl->free);
2358 /* Ensures there is enough space for ncell extra j-subcells in the list */
2359 static void check_subcell_list_space_supersub(nbnxn_pairlist_t *nbl,
2362 int ncj4_max, j4, j, w, t;
2365 #define WARP_SIZE 32
2367 /* We can have maximally nsupercell*GPU_NSUBCELL sj lists */
2368 /* We can store 4 j-subcell - i-supercell pairs in one struct.
2369 * since we round down, we need one extra entry.
2371 ncj4_max = ((nbl->work->cj_ind + nsupercell*GPU_NSUBCELL + NBNXN_GPU_JGROUP_SIZE - 1) >> NBNXN_GPU_JGROUP_SIZE_2LOG);
2373 if (ncj4_max > nbl->cj4_nalloc)
2375 nbl->cj4_nalloc = over_alloc_small(ncj4_max);
2376 nbnxn_realloc_void((void **)&nbl->cj4,
2377 nbl->work->cj4_init*sizeof(*nbl->cj4),
2378 nbl->cj4_nalloc*sizeof(*nbl->cj4),
2379 nbl->alloc, nbl->free);
2382 if (ncj4_max > nbl->work->cj4_init)
2384 for (j4 = nbl->work->cj4_init; j4 < ncj4_max; j4++)
2386 /* No i-subcells and no excl's in the list initially */
2387 for (w = 0; w < NWARP; w++)
2389 nbl->cj4[j4].imei[w].imask = 0U;
2390 nbl->cj4[j4].imei[w].excl_ind = 0;
2394 nbl->work->cj4_init = ncj4_max;
2398 /* Set all excl masks for one GPU warp no exclusions */
2399 static void set_no_excls(nbnxn_excl_t *excl)
2403 for (t = 0; t < WARP_SIZE; t++)
2405 /* Turn all interaction bits on */
2406 excl->pair[t] = NBNXN_INTERACTION_MASK_ALL;
2410 /* Initializes a single nbnxn_pairlist_t data structure */
2411 static void nbnxn_init_pairlist(nbnxn_pairlist_t *nbl,
2413 nbnxn_alloc_t *alloc,
2418 nbl->alloc = nbnxn_alloc_aligned;
2426 nbl->free = nbnxn_free_aligned;
2433 nbl->bSimple = bSimple;
2444 /* We need one element extra in sj, so alloc initially with 1 */
2445 nbl->cj4_nalloc = 0;
2452 nbl->excl_nalloc = 0;
2454 check_excl_space(nbl, 1);
2456 set_no_excls(&nbl->excl[0]);
2462 snew_aligned(nbl->work->bb_ci, 1, NBNXN_SEARCH_BB_MEM_ALIGN);
2467 snew_aligned(nbl->work->pbb_ci, GPU_NSUBCELL/STRIDE_PBB*NNBSBB_XXXX, NBNXN_SEARCH_BB_MEM_ALIGN);
2469 snew_aligned(nbl->work->bb_ci, GPU_NSUBCELL, NBNXN_SEARCH_BB_MEM_ALIGN);
2472 snew_aligned(nbl->work->x_ci, NBNXN_NA_SC_MAX*DIM, NBNXN_SEARCH_BB_MEM_ALIGN);
2473 #ifdef GMX_NBNXN_SIMD
2474 snew_aligned(nbl->work->x_ci_simd_4xn, 1, NBNXN_MEM_ALIGN);
2475 snew_aligned(nbl->work->x_ci_simd_2xnn, 1, NBNXN_MEM_ALIGN);
2477 snew_aligned(nbl->work->d2, GPU_NSUBCELL, NBNXN_SEARCH_BB_MEM_ALIGN);
2479 nbl->work->sort = NULL;
2480 nbl->work->sort_nalloc = 0;
2481 nbl->work->sci_sort = NULL;
2482 nbl->work->sci_sort_nalloc = 0;
2485 void nbnxn_init_pairlist_set(nbnxn_pairlist_set_t *nbl_list,
2486 gmx_bool bSimple, gmx_bool bCombined,
2487 nbnxn_alloc_t *alloc,
2492 nbl_list->bSimple = bSimple;
2493 nbl_list->bCombined = bCombined;
2495 nbl_list->nnbl = gmx_omp_nthreads_get(emntNonbonded);
2497 if (!nbl_list->bCombined &&
2498 nbl_list->nnbl > NBNXN_BUFFERFLAG_MAX_THREADS)
2500 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.",
2501 nbl_list->nnbl, NBNXN_BUFFERFLAG_MAX_THREADS, NBNXN_BUFFERFLAG_MAX_THREADS);
2504 snew(nbl_list->nbl, nbl_list->nnbl);
2505 /* Execute in order to avoid memory interleaving between threads */
2506 #pragma omp parallel for num_threads(nbl_list->nnbl) schedule(static)
2507 for (i = 0; i < nbl_list->nnbl; i++)
2509 /* Allocate the nblist data structure locally on each thread
2510 * to optimize memory access for NUMA architectures.
2512 snew(nbl_list->nbl[i], 1);
2514 /* Only list 0 is used on the GPU, use normal allocation for i>0 */
2517 nbnxn_init_pairlist(nbl_list->nbl[i], nbl_list->bSimple, alloc, free);
2521 nbnxn_init_pairlist(nbl_list->nbl[i], nbl_list->bSimple, NULL, NULL);
2526 /* Print statistics of a pair list, used for debug output */
2527 static void print_nblist_statistics_simple(FILE *fp, const nbnxn_pairlist_t *nbl,
2528 const nbnxn_search_t nbs, real rl)
2530 const nbnxn_grid_t *grid;
2535 /* This code only produces correct statistics with domain decomposition */
2536 grid = &nbs->grid[0];
2538 fprintf(fp, "nbl nci %d ncj %d\n",
2539 nbl->nci, nbl->ncj);
2540 fprintf(fp, "nbl na_sc %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
2541 nbl->na_sc, rl, nbl->ncj, nbl->ncj/(double)grid->nc,
2542 nbl->ncj/(double)grid->nc*grid->na_sc,
2543 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)));
2545 fprintf(fp, "nbl average j cell list length %.1f\n",
2546 0.25*nbl->ncj/(double)nbl->nci);
2548 for (s = 0; s < SHIFTS; s++)
2553 for (i = 0; i < nbl->nci; i++)
2555 cs[nbl->ci[i].shift & NBNXN_CI_SHIFT] +=
2556 nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start;
2558 j = nbl->ci[i].cj_ind_start;
2559 while (j < nbl->ci[i].cj_ind_end &&
2560 nbl->cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
2566 fprintf(fp, "nbl cell pairs, total: %d excl: %d %.1f%%\n",
2567 nbl->ncj, npexcl, 100*npexcl/(double)nbl->ncj);
2568 for (s = 0; s < SHIFTS; s++)
2572 fprintf(fp, "nbl shift %2d ncj %3d\n", s, cs[s]);
2577 /* Print statistics of a pair lists, used for debug output */
2578 static void print_nblist_statistics_supersub(FILE *fp, const nbnxn_pairlist_t *nbl,
2579 const nbnxn_search_t nbs, real rl)
2581 const nbnxn_grid_t *grid;
2582 int i, j4, j, si, b;
2583 int c[GPU_NSUBCELL+1];
2585 /* This code only produces correct statistics with domain decomposition */
2586 grid = &nbs->grid[0];
2588 fprintf(fp, "nbl nsci %d ncj4 %d nsi %d excl4 %d\n",
2589 nbl->nsci, nbl->ncj4, nbl->nci_tot, nbl->nexcl);
2590 fprintf(fp, "nbl na_c %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
2591 nbl->na_ci, rl, nbl->nci_tot, nbl->nci_tot/(double)grid->nsubc_tot,
2592 nbl->nci_tot/(double)grid->nsubc_tot*grid->na_c,
2593 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)));
2595 fprintf(fp, "nbl average j super cell list length %.1f\n",
2596 0.25*nbl->ncj4/(double)nbl->nsci);
2597 fprintf(fp, "nbl average i sub cell list length %.1f\n",
2598 nbl->nci_tot/((double)nbl->ncj4));
2600 for (si = 0; si <= GPU_NSUBCELL; si++)
2604 for (i = 0; i < nbl->nsci; i++)
2606 for (j4 = nbl->sci[i].cj4_ind_start; j4 < nbl->sci[i].cj4_ind_end; j4++)
2608 for (j = 0; j < NBNXN_GPU_JGROUP_SIZE; j++)
2611 for (si = 0; si < GPU_NSUBCELL; si++)
2613 if (nbl->cj4[j4].imei[0].imask & (1U << (j*GPU_NSUBCELL + si)))
2622 for (b = 0; b <= GPU_NSUBCELL; b++)
2624 fprintf(fp, "nbl j-list #i-subcell %d %7d %4.1f\n",
2625 b, c[b], 100.0*c[b]/(double)(nbl->ncj4*NBNXN_GPU_JGROUP_SIZE));
2629 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp */
2630 static void low_get_nbl_exclusions(nbnxn_pairlist_t *nbl, int cj4,
2631 int warp, nbnxn_excl_t **excl)
2633 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
2635 /* No exclusions set, make a new list entry */
2636 nbl->cj4[cj4].imei[warp].excl_ind = nbl->nexcl;
2638 *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
2639 set_no_excls(*excl);
2643 /* We already have some exclusions, new ones can be added to the list */
2644 *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
2648 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp,
2649 * allocates extra memory, if necessary.
2651 static void get_nbl_exclusions_1(nbnxn_pairlist_t *nbl, int cj4,
2652 int warp, nbnxn_excl_t **excl)
2654 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
2656 /* We need to make a new list entry, check if we have space */
2657 check_excl_space(nbl, 1);
2659 low_get_nbl_exclusions(nbl, cj4, warp, excl);
2662 /* Returns pointers to the exclusion mask for cj4-unit cj4 for both warps,
2663 * allocates extra memory, if necessary.
2665 static void get_nbl_exclusions_2(nbnxn_pairlist_t *nbl, int cj4,
2666 nbnxn_excl_t **excl_w0,
2667 nbnxn_excl_t **excl_w1)
2669 /* Check for space we might need */
2670 check_excl_space(nbl, 2);
2672 low_get_nbl_exclusions(nbl, cj4, 0, excl_w0);
2673 low_get_nbl_exclusions(nbl, cj4, 1, excl_w1);
2676 /* Sets the self exclusions i=j and pair exclusions i>j */
2677 static void set_self_and_newton_excls_supersub(nbnxn_pairlist_t *nbl,
2678 int cj4_ind, int sj_offset,
2681 nbnxn_excl_t *excl[2];
2684 /* Here we only set the set self and double pair exclusions */
2686 get_nbl_exclusions_2(nbl, cj4_ind, &excl[0], &excl[1]);
2688 /* Only minor < major bits set */
2689 for (ej = 0; ej < nbl->na_ci; ej++)
2692 for (ei = ej; ei < nbl->na_ci; ei++)
2694 excl[w]->pair[(ej & (NBNXN_GPU_JGROUP_SIZE-1))*nbl->na_ci + ei] &=
2695 ~(1U << (sj_offset*GPU_NSUBCELL + si));
2700 /* Returns a diagonal or off-diagonal interaction mask for plain C lists */
2701 static unsigned int get_imask(gmx_bool rdiag, int ci, int cj)
2703 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
2706 /* Returns a diagonal or off-diagonal interaction mask for cj-size=2 */
2707 static unsigned int get_imask_simd_j2(gmx_bool rdiag, int ci, int cj)
2709 return (rdiag && ci*2 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_0 :
2710 (rdiag && ci*2+1 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_1 :
2711 NBNXN_INTERACTION_MASK_ALL));
2714 /* Returns a diagonal or off-diagonal interaction mask for cj-size=4 */
2715 static unsigned int get_imask_simd_j4(gmx_bool rdiag, int ci, int cj)
2717 return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
2720 /* Returns a diagonal or off-diagonal interaction mask for cj-size=8 */
2721 static unsigned int get_imask_simd_j8(gmx_bool rdiag, int ci, int cj)
2723 return (rdiag && ci == cj*2 ? NBNXN_INTERACTION_MASK_DIAG_J8_0 :
2724 (rdiag && ci == cj*2+1 ? NBNXN_INTERACTION_MASK_DIAG_J8_1 :
2725 NBNXN_INTERACTION_MASK_ALL));
2728 #ifdef GMX_NBNXN_SIMD
2729 #if GMX_SIMD_WIDTH_HERE == 2
2730 #define get_imask_simd_4xn get_imask_simd_j2
2732 #if GMX_SIMD_WIDTH_HERE == 4
2733 #define get_imask_simd_4xn get_imask_simd_j4
2735 #if GMX_SIMD_WIDTH_HERE == 8
2736 #define get_imask_simd_4xn get_imask_simd_j8
2737 #define get_imask_simd_2xnn get_imask_simd_j4
2739 #if GMX_SIMD_WIDTH_HERE == 16
2740 #define get_imask_simd_2xnn get_imask_simd_j8
2744 /* Plain C code for making a pair list of cell ci vs cell cjf-cjl.
2745 * Checks bounding box distances and possibly atom pair distances.
2747 static void make_cluster_list_simple(const nbnxn_grid_t *gridj,
2748 nbnxn_pairlist_t *nbl,
2749 int ci, int cjf, int cjl,
2750 gmx_bool remove_sub_diag,
2752 real rl2, float rbb2,
2755 const nbnxn_list_work_t *work;
2757 const nbnxn_bb_t *bb_ci;
2762 int cjf_gl, cjl_gl, cj;
2766 bb_ci = nbl->work->bb_ci;
2767 x_ci = nbl->work->x_ci;
2770 while (!InRange && cjf <= cjl)
2772 d2 = subc_bb_dist2(0, bb_ci, cjf, gridj->bb);
2775 /* Check if the distance is within the distance where
2776 * we use only the bounding box distance rbb,
2777 * or within the cut-off and there is at least one atom pair
2778 * within the cut-off.
2788 cjf_gl = gridj->cell0 + cjf;
2789 for (i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
2791 for (j = 0; j < NBNXN_CPU_CLUSTER_I_SIZE; j++)
2793 InRange = InRange ||
2794 (sqr(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
2795 sqr(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
2796 sqr(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rl2);
2799 *ndistc += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
2812 while (!InRange && cjl > cjf)
2814 d2 = subc_bb_dist2(0, bb_ci, cjl, gridj->bb);
2817 /* Check if the distance is within the distance where
2818 * we use only the bounding box distance rbb,
2819 * or within the cut-off and there is at least one atom pair
2820 * within the cut-off.
2830 cjl_gl = gridj->cell0 + cjl;
2831 for (i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
2833 for (j = 0; j < NBNXN_CPU_CLUSTER_I_SIZE; j++)
2835 InRange = InRange ||
2836 (sqr(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
2837 sqr(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
2838 sqr(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rl2);
2841 *ndistc += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
2851 for (cj = cjf; cj <= cjl; cj++)
2853 /* Store cj and the interaction mask */
2854 nbl->cj[nbl->ncj].cj = gridj->cell0 + cj;
2855 nbl->cj[nbl->ncj].excl = get_imask(remove_sub_diag, ci, cj);
2858 /* Increase the closing index in i super-cell list */
2859 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
2863 #ifdef GMX_NBNXN_SIMD_4XN
2864 #include "nbnxn_search_simd_4xn.h"
2866 #ifdef GMX_NBNXN_SIMD_2XNN
2867 #include "nbnxn_search_simd_2xnn.h"
2870 /* Plain C or SIMD4 code for making a pair list of super-cell sci vs scj.
2871 * Checks bounding box distances and possibly atom pair distances.
2873 static void make_cluster_list_supersub(const nbnxn_grid_t *gridi,
2874 const nbnxn_grid_t *gridj,
2875 nbnxn_pairlist_t *nbl,
2877 gmx_bool sci_equals_scj,
2878 int stride, const real *x,
2879 real rl2, float rbb2,
2884 int cjo, ci1, ci, cj, cj_gl;
2885 int cj4_ind, cj_offset;
2889 const float *pbb_ci;
2891 const nbnxn_bb_t *bb_ci;
2896 #define PRUNE_LIST_CPU_ONE
2897 #ifdef PRUNE_LIST_CPU_ONE
2901 d2l = nbl->work->d2;
2904 pbb_ci = nbl->work->pbb_ci;
2906 bb_ci = nbl->work->bb_ci;
2908 x_ci = nbl->work->x_ci;
2912 for (cjo = 0; cjo < gridj->nsubc[scj]; cjo++)
2914 cj4_ind = (nbl->work->cj_ind >> NBNXN_GPU_JGROUP_SIZE_2LOG);
2915 cj_offset = nbl->work->cj_ind - cj4_ind*NBNXN_GPU_JGROUP_SIZE;
2916 cj4 = &nbl->cj4[cj4_ind];
2918 cj = scj*GPU_NSUBCELL + cjo;
2920 cj_gl = gridj->cell0*GPU_NSUBCELL + cj;
2922 /* Initialize this j-subcell i-subcell list */
2923 cj4->cj[cj_offset] = cj_gl;
2932 ci1 = gridi->nsubc[sci];
2936 /* Determine all ci1 bb distances in one call with SIMD4 */
2937 subc_bb_dist2_simd4_xxxx(gridj->pbb+(cj>>STRIDE_PBB_2LOG)*NNBSBB_XXXX+(cj & (STRIDE_PBB-1)),
2943 /* We use a fixed upper-bound instead of ci1 to help optimization */
2944 for (ci = 0; ci < GPU_NSUBCELL; ci++)
2951 #ifndef NBNXN_BBXXXX
2952 /* Determine the bb distance between ci and cj */
2953 d2l[ci] = subc_bb_dist2(ci, bb_ci, cj, gridj->bb);
2958 #ifdef PRUNE_LIST_CPU_ALL
2959 /* Check if the distance is within the distance where
2960 * we use only the bounding box distance rbb,
2961 * or within the cut-off and there is at least one atom pair
2962 * within the cut-off. This check is very costly.
2964 *ndistc += na_c*na_c;
2967 #ifdef NBNXN_PBB_SIMD4
2972 (na_c, ci, x_ci, cj_gl, stride, x, rl2)))
2974 /* Check if the distance between the two bounding boxes
2975 * in within the pair-list cut-off.
2980 /* Flag this i-subcell to be taken into account */
2981 imask |= (1U << (cj_offset*GPU_NSUBCELL+ci));
2983 #ifdef PRUNE_LIST_CPU_ONE
2991 #ifdef PRUNE_LIST_CPU_ONE
2992 /* If we only found 1 pair, check if any atoms are actually
2993 * within the cut-off, so we could get rid of it.
2995 if (npair == 1 && d2l[ci_last] >= rbb2)
2997 /* Avoid using function pointers here, as it's slower */
2999 #ifdef NBNXN_PBB_SIMD4
3000 !subc_in_range_simd4
3004 (na_c, ci_last, x_ci, cj_gl, stride, x, rl2))
3006 imask &= ~(1U << (cj_offset*GPU_NSUBCELL+ci_last));
3014 /* We have a useful sj entry, close it now */
3016 /* Set the exclucions for the ci== sj entry.
3017 * Here we don't bother to check if this entry is actually flagged,
3018 * as it will nearly always be in the list.
3022 set_self_and_newton_excls_supersub(nbl, cj4_ind, cj_offset, cjo);
3025 /* Copy the cluster interaction mask to the list */
3026 for (w = 0; w < NWARP; w++)
3028 cj4->imei[w].imask |= imask;
3031 nbl->work->cj_ind++;
3033 /* Keep the count */
3034 nbl->nci_tot += npair;
3036 /* Increase the closing index in i super-cell list */
3037 nbl->sci[nbl->nsci].cj4_ind_end =
3038 ((nbl->work->cj_ind+NBNXN_GPU_JGROUP_SIZE-1) >> NBNXN_GPU_JGROUP_SIZE_2LOG);
3043 /* Set all atom-pair exclusions from the topology stored in excl
3044 * as masks in the pair-list for simple list i-entry nbl_ci
3046 static void set_ci_top_excls(const nbnxn_search_t nbs,
3047 nbnxn_pairlist_t *nbl,
3048 gmx_bool diagRemoved,
3051 const nbnxn_ci_t *nbl_ci,
3052 const t_blocka *excl)
3056 int cj_ind_first, cj_ind_last;
3057 int cj_first, cj_last;
3059 int i, ai, aj, si, eind, ge, se;
3060 int found, cj_ind_0, cj_ind_1, cj_ind_m;
3064 nbnxn_excl_t *nbl_excl;
3065 int inner_i, inner_e;
3069 if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
3077 cj_ind_first = nbl_ci->cj_ind_start;
3078 cj_ind_last = nbl->ncj - 1;
3080 cj_first = nbl->cj[cj_ind_first].cj;
3081 cj_last = nbl->cj[cj_ind_last].cj;
3083 /* Determine how many contiguous j-cells we have starting
3084 * from the first i-cell. This number can be used to directly
3085 * calculate j-cell indices for excluded atoms.
3088 if (na_ci_2log == na_cj_2log)
3090 while (cj_ind_first + ndirect <= cj_ind_last &&
3091 nbl->cj[cj_ind_first+ndirect].cj == ci + ndirect)
3096 #ifdef NBNXN_SEARCH_BB_SIMD4
3099 while (cj_ind_first + ndirect <= cj_ind_last &&
3100 nbl->cj[cj_ind_first+ndirect].cj == ci_to_cj(na_cj_2log, ci) + ndirect)
3107 /* Loop over the atoms in the i super-cell */
3108 for (i = 0; i < nbl->na_sc; i++)
3110 ai = nbs->a[ci*nbl->na_sc+i];
3113 si = (i>>na_ci_2log);
3115 /* Loop over the topology-based exclusions for this i-atom */
3116 for (eind = excl->index[ai]; eind < excl->index[ai+1]; eind++)
3122 /* The self exclusion are already set, save some time */
3128 /* Without shifts we only calculate interactions j>i
3129 * for one-way pair-lists.
3131 if (diagRemoved && ge <= ci*nbl->na_sc + i)
3136 se = (ge >> na_cj_2log);
3138 /* Could the cluster se be in our list? */
3139 if (se >= cj_first && se <= cj_last)
3141 if (se < cj_first + ndirect)
3143 /* We can calculate cj_ind directly from se */
3144 found = cj_ind_first + se - cj_first;
3148 /* Search for se using bisection */
3150 cj_ind_0 = cj_ind_first + ndirect;
3151 cj_ind_1 = cj_ind_last + 1;
3152 while (found == -1 && cj_ind_0 < cj_ind_1)
3154 cj_ind_m = (cj_ind_0 + cj_ind_1)>>1;
3156 cj_m = nbl->cj[cj_ind_m].cj;
3164 cj_ind_1 = cj_ind_m;
3168 cj_ind_0 = cj_ind_m + 1;
3175 inner_i = i - (si << na_ci_2log);
3176 inner_e = ge - (se << na_cj_2log);
3178 nbl->cj[found].excl &= ~(1U<<((inner_i<<na_cj_2log) + inner_e));
3179 /* The next code line is usually not needed. We do not want to version
3180 * away the above line, because there is logic that relies on being
3181 * able to detect easily whether any exclusions exist. */
3182 #if (defined GMX_CPU_ACCELERATION_IBM_QPX)
3183 nbl->cj[found].interaction_mask_indices[inner_i] &= ~(1U << inner_e);
3192 /* Set all atom-pair exclusions from the topology stored in excl
3193 * as masks in the pair-list for i-super-cell entry nbl_sci
3195 static void set_sci_top_excls(const nbnxn_search_t nbs,
3196 nbnxn_pairlist_t *nbl,
3197 gmx_bool diagRemoved,
3199 const nbnxn_sci_t *nbl_sci,
3200 const t_blocka *excl)
3205 int cj_ind_first, cj_ind_last;
3206 int cj_first, cj_last;
3208 int i, ai, aj, si, eind, ge, se;
3209 int found, cj_ind_0, cj_ind_1, cj_ind_m;
3213 nbnxn_excl_t *nbl_excl;
3214 int inner_i, inner_e, w;
3220 if (nbl_sci->cj4_ind_end == nbl_sci->cj4_ind_start)
3228 cj_ind_first = nbl_sci->cj4_ind_start*NBNXN_GPU_JGROUP_SIZE;
3229 cj_ind_last = nbl->work->cj_ind - 1;
3231 cj_first = nbl->cj4[nbl_sci->cj4_ind_start].cj[0];
3232 cj_last = nbl_cj(nbl, cj_ind_last);
3234 /* Determine how many contiguous j-clusters we have starting
3235 * from the first i-cluster. This number can be used to directly
3236 * calculate j-cluster indices for excluded atoms.
3239 while (cj_ind_first + ndirect <= cj_ind_last &&
3240 nbl_cj(nbl, cj_ind_first+ndirect) == sci*GPU_NSUBCELL + ndirect)
3245 /* Loop over the atoms in the i super-cell */
3246 for (i = 0; i < nbl->na_sc; i++)
3248 ai = nbs->a[sci*nbl->na_sc+i];
3251 si = (i>>na_c_2log);
3253 /* Loop over the topology-based exclusions for this i-atom */
3254 for (eind = excl->index[ai]; eind < excl->index[ai+1]; eind++)
3260 /* The self exclusion are already set, save some time */
3266 /* Without shifts we only calculate interactions j>i
3267 * for one-way pair-lists.
3269 if (diagRemoved && ge <= sci*nbl->na_sc + i)
3275 /* Could the cluster se be in our list? */
3276 if (se >= cj_first && se <= cj_last)
3278 if (se < cj_first + ndirect)
3280 /* We can calculate cj_ind directly from se */
3281 found = cj_ind_first + se - cj_first;
3285 /* Search for se using bisection */
3287 cj_ind_0 = cj_ind_first + ndirect;
3288 cj_ind_1 = cj_ind_last + 1;
3289 while (found == -1 && cj_ind_0 < cj_ind_1)
3291 cj_ind_m = (cj_ind_0 + cj_ind_1)>>1;
3293 cj_m = nbl_cj(nbl, cj_ind_m);
3301 cj_ind_1 = cj_ind_m;
3305 cj_ind_0 = cj_ind_m + 1;
3312 inner_i = i - si*na_c;
3313 inner_e = ge - se*na_c;
3315 /* Macro for getting the index of atom a within a cluster */
3316 #define AMODCJ4(a) ((a) & (NBNXN_GPU_JGROUP_SIZE - 1))
3317 /* Macro for converting an atom number to a cluster number */
3318 #define A2CJ4(a) ((a) >> NBNXN_GPU_JGROUP_SIZE_2LOG)
3319 /* Macro for getting the index of an i-atom within a warp */
3320 #define AMODWI(a) ((a) & (NBNXN_GPU_CLUSTER_SIZE/2 - 1))
3322 if (nbl_imask0(nbl, found) & (1U << (AMODCJ4(found)*GPU_NSUBCELL + si)))
3326 get_nbl_exclusions_1(nbl, A2CJ4(found), w, &nbl_excl);
3328 nbl_excl->pair[AMODWI(inner_e)*nbl->na_ci+inner_i] &=
3329 ~(1U << (AMODCJ4(found)*GPU_NSUBCELL + si));
3342 /* Reallocate the simple ci list for at least n entries */
3343 static void nb_realloc_ci(nbnxn_pairlist_t *nbl, int n)
3345 nbl->ci_nalloc = over_alloc_small(n);
3346 nbnxn_realloc_void((void **)&nbl->ci,
3347 nbl->nci*sizeof(*nbl->ci),
3348 nbl->ci_nalloc*sizeof(*nbl->ci),
3349 nbl->alloc, nbl->free);
3352 /* Reallocate the super-cell sci list for at least n entries */
3353 static void nb_realloc_sci(nbnxn_pairlist_t *nbl, int n)
3355 nbl->sci_nalloc = over_alloc_small(n);
3356 nbnxn_realloc_void((void **)&nbl->sci,
3357 nbl->nsci*sizeof(*nbl->sci),
3358 nbl->sci_nalloc*sizeof(*nbl->sci),
3359 nbl->alloc, nbl->free);
3362 /* Make a new ci entry at index nbl->nci */
3363 static void new_ci_entry(nbnxn_pairlist_t *nbl, int ci, int shift, int flags)
3365 if (nbl->nci + 1 > nbl->ci_nalloc)
3367 nb_realloc_ci(nbl, nbl->nci+1);
3369 nbl->ci[nbl->nci].ci = ci;
3370 nbl->ci[nbl->nci].shift = shift;
3371 /* Store the interaction flags along with the shift */
3372 nbl->ci[nbl->nci].shift |= flags;
3373 nbl->ci[nbl->nci].cj_ind_start = nbl->ncj;
3374 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
3377 /* Make a new sci entry at index nbl->nsci */
3378 static void new_sci_entry(nbnxn_pairlist_t *nbl, int sci, int shift)
3380 if (nbl->nsci + 1 > nbl->sci_nalloc)
3382 nb_realloc_sci(nbl, nbl->nsci+1);
3384 nbl->sci[nbl->nsci].sci = sci;
3385 nbl->sci[nbl->nsci].shift = shift;
3386 nbl->sci[nbl->nsci].cj4_ind_start = nbl->ncj4;
3387 nbl->sci[nbl->nsci].cj4_ind_end = nbl->ncj4;
3390 /* Sort the simple j-list cj on exclusions.
3391 * Entries with exclusions will all be sorted to the beginning of the list.
3393 static void sort_cj_excl(nbnxn_cj_t *cj, int ncj,
3394 nbnxn_list_work_t *work)
3398 if (ncj > work->cj_nalloc)
3400 work->cj_nalloc = over_alloc_large(ncj);
3401 srenew(work->cj, work->cj_nalloc);
3404 /* Make a list of the j-cells involving exclusions */
3406 for (j = 0; j < ncj; j++)
3408 if (cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
3410 work->cj[jnew++] = cj[j];
3413 /* Check if there are exclusions at all or not just the first entry */
3414 if (!((jnew == 0) ||
3415 (jnew == 1 && cj[0].excl != NBNXN_INTERACTION_MASK_ALL)))
3417 for (j = 0; j < ncj; j++)
3419 if (cj[j].excl == NBNXN_INTERACTION_MASK_ALL)
3421 work->cj[jnew++] = cj[j];
3424 for (j = 0; j < ncj; j++)
3426 cj[j] = work->cj[j];
3431 /* Close this simple list i entry */
3432 static void close_ci_entry_simple(nbnxn_pairlist_t *nbl)
3436 /* All content of the new ci entry have already been filled correctly,
3437 * we only need to increase the count here (for non empty lists).
3439 jlen = nbl->ci[nbl->nci].cj_ind_end - nbl->ci[nbl->nci].cj_ind_start;
3442 sort_cj_excl(nbl->cj+nbl->ci[nbl->nci].cj_ind_start, jlen, nbl->work);
3444 /* The counts below are used for non-bonded pair/flop counts
3445 * and should therefore match the available kernel setups.
3447 if (!(nbl->ci[nbl->nci].shift & NBNXN_CI_DO_COUL(0)))
3449 nbl->work->ncj_noq += jlen;
3451 else if ((nbl->ci[nbl->nci].shift & NBNXN_CI_HALF_LJ(0)) ||
3452 !(nbl->ci[nbl->nci].shift & NBNXN_CI_DO_LJ(0)))
3454 nbl->work->ncj_hlj += jlen;
3461 /* Split sci entry for load balancing on the GPU.
3462 * Splitting ensures we have enough lists to fully utilize the whole GPU.
3463 * With progBal we generate progressively smaller lists, which improves
3464 * load balancing. As we only know the current count on our own thread,
3465 * we will need to estimate the current total amount of i-entries.
3466 * As the lists get concatenated later, this estimate depends
3467 * both on nthread and our own thread index.
3469 static void split_sci_entry(nbnxn_pairlist_t *nbl,
3470 int nsp_max_av, gmx_bool progBal, int nc_bal,
3471 int thread, int nthread)
3475 int cj4_start, cj4_end, j4len, cj4;
3477 int nsp, nsp_sci, nsp_cj4, nsp_cj4_e, nsp_cj4_p;
3482 /* Estimate the total numbers of ci's of the nblist combined
3483 * over all threads using the target number of ci's.
3485 nsci_est = nc_bal*thread/nthread + nbl->nsci;
3487 /* The first ci blocks should be larger, to avoid overhead.
3488 * The last ci blocks should be smaller, to improve load balancing.
3491 nsp_max_av*nc_bal*3/(2*(nsci_est - 1 + nc_bal)));
3495 nsp_max = nsp_max_av;
3498 cj4_start = nbl->sci[nbl->nsci-1].cj4_ind_start;
3499 cj4_end = nbl->sci[nbl->nsci-1].cj4_ind_end;
3500 j4len = cj4_end - cj4_start;
3502 if (j4len > 1 && j4len*GPU_NSUBCELL*NBNXN_GPU_JGROUP_SIZE > nsp_max)
3504 /* Remove the last ci entry and process the cj4's again */
3512 for (cj4 = cj4_start; cj4 < cj4_end; cj4++)
3514 nsp_cj4_p = nsp_cj4;
3515 /* Count the number of cluster pairs in this cj4 group */
3517 for (p = 0; p < GPU_NSUBCELL*NBNXN_GPU_JGROUP_SIZE; p++)
3519 nsp_cj4 += (nbl->cj4[cj4].imei[0].imask >> p) & 1;
3522 if (nsp_cj4 > 0 && nsp + nsp_cj4 > nsp_max)
3524 /* Split the list at cj4 */
3525 nbl->sci[sci].cj4_ind_end = cj4;
3526 /* Create a new sci entry */
3529 if (nbl->nsci+1 > nbl->sci_nalloc)
3531 nb_realloc_sci(nbl, nbl->nsci+1);
3533 nbl->sci[sci].sci = nbl->sci[nbl->nsci-1].sci;
3534 nbl->sci[sci].shift = nbl->sci[nbl->nsci-1].shift;
3535 nbl->sci[sci].cj4_ind_start = cj4;
3537 nsp_cj4_e = nsp_cj4_p;
3543 /* Put the remaining cj4's in the last sci entry */
3544 nbl->sci[sci].cj4_ind_end = cj4_end;
3546 /* Possibly balance out the last two sci's
3547 * by moving the last cj4 of the second last sci.
3549 if (nsp_sci - nsp_cj4_e >= nsp + nsp_cj4_e)
3551 nbl->sci[sci-1].cj4_ind_end--;
3552 nbl->sci[sci].cj4_ind_start--;
3559 /* Clost this super/sub list i entry */
3560 static void close_ci_entry_supersub(nbnxn_pairlist_t *nbl,
3562 gmx_bool progBal, int nc_bal,
3563 int thread, int nthread)
3568 /* All content of the new ci entry have already been filled correctly,
3569 * we only need to increase the count here (for non empty lists).
3571 j4len = nbl->sci[nbl->nsci].cj4_ind_end - nbl->sci[nbl->nsci].cj4_ind_start;
3574 /* We can only have complete blocks of 4 j-entries in a list,
3575 * so round the count up before closing.
3577 nbl->ncj4 = ((nbl->work->cj_ind + NBNXN_GPU_JGROUP_SIZE - 1) >> NBNXN_GPU_JGROUP_SIZE_2LOG);
3578 nbl->work->cj_ind = nbl->ncj4*NBNXN_GPU_JGROUP_SIZE;
3584 /* Measure the size of the new entry and potentially split it */
3585 split_sci_entry(nbl, nsp_max_av, progBal, nc_bal, thread, nthread);
3590 /* Syncs the working array before adding another grid pair to the list */
3591 static void sync_work(nbnxn_pairlist_t *nbl)
3595 nbl->work->cj_ind = nbl->ncj4*NBNXN_GPU_JGROUP_SIZE;
3596 nbl->work->cj4_init = nbl->ncj4;
3600 /* Clears an nbnxn_pairlist_t data structure */
3601 static void clear_pairlist(nbnxn_pairlist_t *nbl)
3610 nbl->work->ncj_noq = 0;
3611 nbl->work->ncj_hlj = 0;
3614 /* Sets a simple list i-cell bounding box, including PBC shift */
3615 static gmx_inline void set_icell_bb_simple(const nbnxn_bb_t *bb, int ci,
3616 real shx, real shy, real shz,
3619 bb_ci->lower[BB_X] = bb[ci].lower[BB_X] + shx;
3620 bb_ci->lower[BB_Y] = bb[ci].lower[BB_Y] + shy;
3621 bb_ci->lower[BB_Z] = bb[ci].lower[BB_Z] + shz;
3622 bb_ci->upper[BB_X] = bb[ci].upper[BB_X] + shx;
3623 bb_ci->upper[BB_Y] = bb[ci].upper[BB_Y] + shy;
3624 bb_ci->upper[BB_Z] = bb[ci].upper[BB_Z] + shz;
3628 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
3629 static void set_icell_bbxxxx_supersub(const float *bb, int ci,
3630 real shx, real shy, real shz,
3635 ia = ci*(GPU_NSUBCELL>>STRIDE_PBB_2LOG)*NNBSBB_XXXX;
3636 for (m = 0; m < (GPU_NSUBCELL>>STRIDE_PBB_2LOG)*NNBSBB_XXXX; m += NNBSBB_XXXX)
3638 for (i = 0; i < STRIDE_PBB; i++)
3640 bb_ci[m+0*STRIDE_PBB+i] = bb[ia+m+0*STRIDE_PBB+i] + shx;
3641 bb_ci[m+1*STRIDE_PBB+i] = bb[ia+m+1*STRIDE_PBB+i] + shy;
3642 bb_ci[m+2*STRIDE_PBB+i] = bb[ia+m+2*STRIDE_PBB+i] + shz;
3643 bb_ci[m+3*STRIDE_PBB+i] = bb[ia+m+3*STRIDE_PBB+i] + shx;
3644 bb_ci[m+4*STRIDE_PBB+i] = bb[ia+m+4*STRIDE_PBB+i] + shy;
3645 bb_ci[m+5*STRIDE_PBB+i] = bb[ia+m+5*STRIDE_PBB+i] + shz;
3651 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
3652 static void set_icell_bb_supersub(const nbnxn_bb_t *bb, int ci,
3653 real shx, real shy, real shz,
3658 for (i = 0; i < GPU_NSUBCELL; i++)
3660 set_icell_bb_simple(bb, ci*GPU_NSUBCELL+i,
3666 /* Copies PBC shifted i-cell atom coordinates x,y,z to working array */
3667 static void icell_set_x_simple(int ci,
3668 real shx, real shy, real shz,
3669 int gmx_unused na_c,
3670 int stride, const real *x,
3671 nbnxn_list_work_t *work)
3675 ia = ci*NBNXN_CPU_CLUSTER_I_SIZE;
3677 for (i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE; i++)
3679 work->x_ci[i*STRIDE_XYZ+XX] = x[(ia+i)*stride+XX] + shx;
3680 work->x_ci[i*STRIDE_XYZ+YY] = x[(ia+i)*stride+YY] + shy;
3681 work->x_ci[i*STRIDE_XYZ+ZZ] = x[(ia+i)*stride+ZZ] + shz;
3685 /* Copies PBC shifted super-cell atom coordinates x,y,z to working array */
3686 static void icell_set_x_supersub(int ci,
3687 real shx, real shy, real shz,
3689 int stride, const real *x,
3690 nbnxn_list_work_t *work)
3697 ia = ci*GPU_NSUBCELL*na_c;
3698 for (i = 0; i < GPU_NSUBCELL*na_c; i++)
3700 x_ci[i*DIM + XX] = x[(ia+i)*stride + XX] + shx;
3701 x_ci[i*DIM + YY] = x[(ia+i)*stride + YY] + shy;
3702 x_ci[i*DIM + ZZ] = x[(ia+i)*stride + ZZ] + shz;
3706 #ifdef NBNXN_SEARCH_BB_SIMD4
3707 /* Copies PBC shifted super-cell packed atom coordinates to working array */
3708 static void icell_set_x_supersub_simd4(int ci,
3709 real shx, real shy, real shz,
3711 int stride, const real *x,
3712 nbnxn_list_work_t *work)
3714 int si, io, ia, i, j;
3719 for (si = 0; si < GPU_NSUBCELL; si++)
3721 for (i = 0; i < na_c; i += STRIDE_PBB)
3724 ia = ci*GPU_NSUBCELL*na_c + io;
3725 for (j = 0; j < STRIDE_PBB; j++)
3727 x_ci[io*DIM + j + XX*STRIDE_PBB] = x[(ia+j)*stride+XX] + shx;
3728 x_ci[io*DIM + j + YY*STRIDE_PBB] = x[(ia+j)*stride+YY] + shy;
3729 x_ci[io*DIM + j + ZZ*STRIDE_PBB] = x[(ia+j)*stride+ZZ] + shz;
3736 static real nbnxn_rlist_inc_nonloc_fac = 0.6;
3738 /* Due to the cluster size the effective pair-list is longer than
3739 * that of a simple atom pair-list. This function gives the extra distance.
3741 real nbnxn_get_rlist_effective_inc(int cluster_size, real atom_density)
3743 return ((0.5 + nbnxn_rlist_inc_nonloc_fac)*sqr(((cluster_size) - 1.0)/(cluster_size))*pow((cluster_size)/(atom_density), 1.0/3.0));
3746 /* Estimates the interaction volume^2 for non-local interactions */
3747 static real nonlocal_vol2(const gmx_domdec_zones_t *zones, rvec ls, real r)
3756 /* Here we simply add up the volumes of 1, 2 or 3 1D decomposition
3757 * not home interaction volume^2. As these volumes are not additive,
3758 * this is an overestimate, but it would only be significant in the limit
3759 * of small cells, where we anyhow need to split the lists into
3760 * as small parts as possible.
3763 for (z = 0; z < zones->n; z++)
3765 if (zones->shift[z][XX] + zones->shift[z][YY] + zones->shift[z][ZZ] == 1)
3770 for (d = 0; d < DIM; d++)
3772 if (zones->shift[z][d] == 0)
3776 za *= zones->size[z].x1[d] - zones->size[z].x0[d];
3780 /* 4 octants of a sphere */
3781 vold_est = 0.25*M_PI*r*r*r*r;
3782 /* 4 quarter pie slices on the edges */
3783 vold_est += 4*cl*M_PI/6.0*r*r*r;
3784 /* One rectangular volume on a face */
3785 vold_est += ca*0.5*r*r;
3787 vol2_est_tot += vold_est*za;
3791 return vol2_est_tot;
3794 /* Estimates the average size of a full j-list for super/sub setup */
3795 static int get_nsubpair_max(const nbnxn_search_t nbs,
3798 int min_ci_balanced)
3800 const nbnxn_grid_t *grid;
3802 real xy_diag2, r_eff_sup, vol_est, nsp_est, nsp_est_nl;
3805 grid = &nbs->grid[0];
3807 ls[XX] = (grid->c1[XX] - grid->c0[XX])/(grid->ncx*GPU_NSUBCELL_X);
3808 ls[YY] = (grid->c1[YY] - grid->c0[YY])/(grid->ncy*GPU_NSUBCELL_Y);
3809 ls[ZZ] = (grid->c1[ZZ] - grid->c0[ZZ])*grid->ncx*grid->ncy/(grid->nc*GPU_NSUBCELL_Z);
3811 /* The average squared length of the diagonal of a sub cell */
3812 xy_diag2 = ls[XX]*ls[XX] + ls[YY]*ls[YY] + ls[ZZ]*ls[ZZ];
3814 /* The formulas below are a heuristic estimate of the average nsj per si*/
3815 r_eff_sup = rlist + nbnxn_rlist_inc_nonloc_fac*sqr((grid->na_c - 1.0)/grid->na_c)*sqrt(xy_diag2/3);
3817 if (!nbs->DomDec || nbs->zones->n == 1)
3824 sqr(grid->atom_density/grid->na_c)*
3825 nonlocal_vol2(nbs->zones, ls, r_eff_sup);
3830 /* Sub-cell interacts with itself */
3831 vol_est = ls[XX]*ls[YY]*ls[ZZ];
3832 /* 6/2 rectangular volume on the faces */
3833 vol_est += (ls[XX]*ls[YY] + ls[XX]*ls[ZZ] + ls[YY]*ls[ZZ])*r_eff_sup;
3834 /* 12/2 quarter pie slices on the edges */
3835 vol_est += 2*(ls[XX] + ls[YY] + ls[ZZ])*0.25*M_PI*sqr(r_eff_sup);
3836 /* 4 octants of a sphere */
3837 vol_est += 0.5*4.0/3.0*M_PI*pow(r_eff_sup, 3);
3839 nsp_est = grid->nsubc_tot*vol_est*grid->atom_density/grid->na_c;
3841 /* Subtract the non-local pair count */
3842 nsp_est -= nsp_est_nl;
3846 fprintf(debug, "nsp_est local %5.1f non-local %5.1f\n",
3847 nsp_est, nsp_est_nl);
3852 nsp_est = nsp_est_nl;
3855 if (min_ci_balanced <= 0 || grid->nc >= min_ci_balanced || grid->nc == 0)
3857 /* We don't need to worry */
3862 /* Thus the (average) maximum j-list size should be as follows */
3863 nsubpair_max = max(1, (int)(nsp_est/min_ci_balanced+0.5));
3865 /* Since the target value is a maximum (this avoids high outliers,
3866 * which lead to load imbalance), not average, we add half the
3867 * number of pairs in a cj4 block to get the average about right.
3869 nsubpair_max += GPU_NSUBCELL*NBNXN_GPU_JGROUP_SIZE/2;
3874 fprintf(debug, "nbl nsp estimate %.1f, nsubpair_max %d\n",
3875 nsp_est, nsubpair_max);
3878 return nsubpair_max;
3881 /* Debug list print function */
3882 static void print_nblist_ci_cj(FILE *fp, const nbnxn_pairlist_t *nbl)
3886 for (i = 0; i < nbl->nci; i++)
3888 fprintf(fp, "ci %4d shift %2d ncj %3d\n",
3889 nbl->ci[i].ci, nbl->ci[i].shift,
3890 nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start);
3892 for (j = nbl->ci[i].cj_ind_start; j < nbl->ci[i].cj_ind_end; j++)
3894 fprintf(fp, " cj %5d imask %x\n",
3901 /* Debug list print function */
3902 static void print_nblist_sci_cj(FILE *fp, const nbnxn_pairlist_t *nbl)
3904 int i, j4, j, ncp, si;
3906 for (i = 0; i < nbl->nsci; i++)
3908 fprintf(fp, "ci %4d shift %2d ncj4 %2d\n",
3909 nbl->sci[i].sci, nbl->sci[i].shift,
3910 nbl->sci[i].cj4_ind_end - nbl->sci[i].cj4_ind_start);
3913 for (j4 = nbl->sci[i].cj4_ind_start; j4 < nbl->sci[i].cj4_ind_end; j4++)
3915 for (j = 0; j < NBNXN_GPU_JGROUP_SIZE; j++)
3917 fprintf(fp, " sj %5d imask %x\n",
3919 nbl->cj4[j4].imei[0].imask);
3920 for (si = 0; si < GPU_NSUBCELL; si++)
3922 if (nbl->cj4[j4].imei[0].imask & (1U << (j*GPU_NSUBCELL + si)))
3929 fprintf(fp, "ci %4d shift %2d ncj4 %2d ncp %3d\n",
3930 nbl->sci[i].sci, nbl->sci[i].shift,
3931 nbl->sci[i].cj4_ind_end - nbl->sci[i].cj4_ind_start,
3936 /* Combine pair lists *nbl generated on multiple threads nblc */
3937 static void combine_nblists(int nnbl, nbnxn_pairlist_t **nbl,
3938 nbnxn_pairlist_t *nblc)
3940 int nsci, ncj4, nexcl;
3945 gmx_incons("combine_nblists does not support simple lists");
3950 nexcl = nblc->nexcl;
3951 for (i = 0; i < nnbl; i++)
3953 nsci += nbl[i]->nsci;
3954 ncj4 += nbl[i]->ncj4;
3955 nexcl += nbl[i]->nexcl;
3958 if (nsci > nblc->sci_nalloc)
3960 nb_realloc_sci(nblc, nsci);
3962 if (ncj4 > nblc->cj4_nalloc)
3964 nblc->cj4_nalloc = over_alloc_small(ncj4);
3965 nbnxn_realloc_void((void **)&nblc->cj4,
3966 nblc->ncj4*sizeof(*nblc->cj4),
3967 nblc->cj4_nalloc*sizeof(*nblc->cj4),
3968 nblc->alloc, nblc->free);
3970 if (nexcl > nblc->excl_nalloc)
3972 nblc->excl_nalloc = over_alloc_small(nexcl);
3973 nbnxn_realloc_void((void **)&nblc->excl,
3974 nblc->nexcl*sizeof(*nblc->excl),
3975 nblc->excl_nalloc*sizeof(*nblc->excl),
3976 nblc->alloc, nblc->free);
3979 /* Each thread should copy its own data to the combined arrays,
3980 * as otherwise data will go back and forth between different caches.
3982 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntPairsearch)) schedule(static)
3983 for (n = 0; n < nnbl; n++)
3990 const nbnxn_pairlist_t *nbli;
3992 /* Determine the offset in the combined data for our thread */
3993 sci_offset = nblc->nsci;
3994 cj4_offset = nblc->ncj4;
3995 ci_offset = nblc->nci_tot;
3996 excl_offset = nblc->nexcl;
3998 for (i = 0; i < n; i++)
4000 sci_offset += nbl[i]->nsci;
4001 cj4_offset += nbl[i]->ncj4;
4002 ci_offset += nbl[i]->nci_tot;
4003 excl_offset += nbl[i]->nexcl;
4008 for (i = 0; i < nbli->nsci; i++)
4010 nblc->sci[sci_offset+i] = nbli->sci[i];
4011 nblc->sci[sci_offset+i].cj4_ind_start += cj4_offset;
4012 nblc->sci[sci_offset+i].cj4_ind_end += cj4_offset;
4015 for (j4 = 0; j4 < nbli->ncj4; j4++)
4017 nblc->cj4[cj4_offset+j4] = nbli->cj4[j4];
4018 nblc->cj4[cj4_offset+j4].imei[0].excl_ind += excl_offset;
4019 nblc->cj4[cj4_offset+j4].imei[1].excl_ind += excl_offset;
4022 for (j4 = 0; j4 < nbli->nexcl; j4++)
4024 nblc->excl[excl_offset+j4] = nbli->excl[j4];
4028 for (n = 0; n < nnbl; n++)
4030 nblc->nsci += nbl[n]->nsci;
4031 nblc->ncj4 += nbl[n]->ncj4;
4032 nblc->nci_tot += nbl[n]->nci_tot;
4033 nblc->nexcl += nbl[n]->nexcl;
4037 /* Returns the next ci to be processes by our thread */
4038 static gmx_bool next_ci(const nbnxn_grid_t *grid,
4040 int nth, int ci_block,
4041 int *ci_x, int *ci_y,
4047 if (*ci_b == ci_block)
4049 /* Jump to the next block assigned to this task */
4050 *ci += (nth - 1)*ci_block;
4054 if (*ci >= grid->nc*conv)
4059 while (*ci >= grid->cxy_ind[*ci_x*grid->ncy + *ci_y + 1]*conv)
4062 if (*ci_y == grid->ncy)
4072 /* Returns the distance^2 for which we put cell pairs in the list
4073 * without checking atom pair distances. This is usually < rlist^2.
4075 static float boundingbox_only_distance2(const nbnxn_grid_t *gridi,
4076 const nbnxn_grid_t *gridj,
4080 /* If the distance between two sub-cell bounding boxes is less
4081 * than this distance, do not check the distance between
4082 * all particle pairs in the sub-cell, since then it is likely
4083 * that the box pair has atom pairs within the cut-off.
4084 * We use the nblist cut-off minus 0.5 times the average x/y diagonal
4085 * spacing of the sub-cells. Around 40% of the checked pairs are pruned.
4086 * Using more than 0.5 gains at most 0.5%.
4087 * If forces are calculated more than twice, the performance gain
4088 * in the force calculation outweighs the cost of checking.
4089 * Note that with subcell lists, the atom-pair distance check
4090 * is only performed when only 1 out of 8 sub-cells in within range,
4091 * this is because the GPU is much faster than the cpu.
4096 bbx = 0.5*(gridi->sx + gridj->sx);
4097 bby = 0.5*(gridi->sy + gridj->sy);
4100 bbx /= GPU_NSUBCELL_X;
4101 bby /= GPU_NSUBCELL_Y;
4104 rbb2 = sqr(max(0, rlist - 0.5*sqrt(bbx*bbx + bby*bby)));
4109 return (float)((1+GMX_FLOAT_EPS)*rbb2);
4113 static int get_ci_block_size(const nbnxn_grid_t *gridi,
4114 gmx_bool bDomDec, int nth)
4116 const int ci_block_enum = 5;
4117 const int ci_block_denom = 11;
4118 const int ci_block_min_atoms = 16;
4121 /* Here we decide how to distribute the blocks over the threads.
4122 * We use prime numbers to try to avoid that the grid size becomes
4123 * a multiple of the number of threads, which would lead to some
4124 * threads getting "inner" pairs and others getting boundary pairs,
4125 * which in turns will lead to load imbalance between threads.
4126 * Set the block size as 5/11/ntask times the average number of cells
4127 * in a y,z slab. This should ensure a quite uniform distribution
4128 * of the grid parts of the different thread along all three grid
4129 * zone boundaries with 3D domain decomposition. At the same time
4130 * the blocks will not become too small.
4132 ci_block = (gridi->nc*ci_block_enum)/(ci_block_denom*gridi->ncx*nth);
4134 /* Ensure the blocks are not too small: avoids cache invalidation */
4135 if (ci_block*gridi->na_sc < ci_block_min_atoms)
4137 ci_block = (ci_block_min_atoms + gridi->na_sc - 1)/gridi->na_sc;
4140 /* Without domain decomposition
4141 * or with less than 3 blocks per task, divide in nth blocks.
4143 if (!bDomDec || ci_block*3*nth > gridi->nc)
4145 ci_block = (gridi->nc + nth - 1)/nth;
4151 /* Generates the part of pair-list nbl assigned to our thread */
4152 static void nbnxn_make_pairlist_part(const nbnxn_search_t nbs,
4153 const nbnxn_grid_t *gridi,
4154 const nbnxn_grid_t *gridj,
4155 nbnxn_search_work_t *work,
4156 const nbnxn_atomdata_t *nbat,
4157 const t_blocka *excl,
4161 gmx_bool bFBufferFlag,
4164 int min_ci_balanced,
4166 nbnxn_pairlist_t *nbl)
4173 int ci_b, ci, ci_x, ci_y, ci_xy, cj;
4179 int conv_i, cell0_i;
4180 const nbnxn_bb_t *bb_i = NULL;
4182 const float *pbb_i = NULL;
4184 const float *bbcz_i, *bbcz_j;
4186 real bx0, bx1, by0, by1, bz0, bz1;
4188 real d2cx, d2z, d2z_cx, d2z_cy, d2zx, d2zxy, d2xy;
4189 int cxf, cxl, cyf, cyf_x, cyl;
4191 int c0, c1, cs, cf, cl;
4194 int gridi_flag_shift = 0, gridj_flag_shift = 0;
4195 unsigned *gridj_flag = NULL;
4196 int ncj_old_i, ncj_old_j;
4198 nbs_cycle_start(&work->cc[enbsCCsearch]);
4200 if (gridj->bSimple != nbl->bSimple)
4202 gmx_incons("Grid incompatible with pair-list");
4206 nbl->na_sc = gridj->na_sc;
4207 nbl->na_ci = gridj->na_c;
4208 nbl->na_cj = nbnxn_kernel_to_cj_size(nb_kernel_type);
4209 na_cj_2log = get_2log(nbl->na_cj);
4215 /* Determine conversion of clusters to flag blocks */
4216 gridi_flag_shift = 0;
4217 while ((nbl->na_ci<<gridi_flag_shift) < NBNXN_BUFFERFLAG_SIZE)
4221 gridj_flag_shift = 0;
4222 while ((nbl->na_cj<<gridj_flag_shift) < NBNXN_BUFFERFLAG_SIZE)
4227 gridj_flag = work->buffer_flags.flag;
4230 copy_mat(nbs->box, box);
4232 rl2 = nbl->rlist*nbl->rlist;
4234 rbb2 = boundingbox_only_distance2(gridi, gridj, nbl->rlist, nbl->bSimple);
4238 fprintf(debug, "nbl bounding box only distance %f\n", sqrt(rbb2));
4241 /* Set the shift range */
4242 for (d = 0; d < DIM; d++)
4244 /* Check if we need periodicity shifts.
4245 * Without PBC or with domain decomposition we don't need them.
4247 if (d >= ePBC2npbcdim(nbs->ePBC) || nbs->dd_dim[d])
4254 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
4265 if (nbl->bSimple && !gridi->bSimple)
4267 conv_i = gridi->na_sc/gridj->na_sc;
4268 bb_i = gridi->bb_simple;
4269 bbcz_i = gridi->bbcz_simple;
4270 flags_i = gridi->flags_simple;
4285 /* We use the normal bounding box format for both grid types */
4288 bbcz_i = gridi->bbcz;
4289 flags_i = gridi->flags;
4291 cell0_i = gridi->cell0*conv_i;
4293 bbcz_j = gridj->bbcz;
4297 /* Blocks of the conversion factor - 1 give a large repeat count
4298 * combined with a small block size. This should result in good
4299 * load balancing for both small and large domains.
4301 ci_block = conv_i - 1;
4305 fprintf(debug, "nbl nc_i %d col.av. %.1f ci_block %d\n",
4306 gridi->nc, gridi->nc/(double)(gridi->ncx*gridi->ncy), ci_block);
4312 /* Initially ci_b and ci to 1 before where we want them to start,
4313 * as they will both be incremented in next_ci.
4316 ci = th*ci_block - 1;
4319 while (next_ci(gridi, conv_i, nth, ci_block, &ci_x, &ci_y, &ci_b, &ci))
4321 if (nbl->bSimple && flags_i[ci] == 0)
4326 ncj_old_i = nbl->ncj;
4329 if (gridj != gridi && shp[XX] == 0)
4333 bx1 = bb_i[ci].upper[BB_X];
4337 bx1 = gridi->c0[XX] + (ci_x+1)*gridi->sx;
4339 if (bx1 < gridj->c0[XX])
4341 d2cx = sqr(gridj->c0[XX] - bx1);
4350 ci_xy = ci_x*gridi->ncy + ci_y;
4352 /* Loop over shift vectors in three dimensions */
4353 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
4355 shz = tz*box[ZZ][ZZ];
4357 bz0 = bbcz_i[ci*NNBSBB_D ] + shz;
4358 bz1 = bbcz_i[ci*NNBSBB_D+1] + shz;
4370 d2z = sqr(bz0 - box[ZZ][ZZ]);
4373 d2z_cx = d2z + d2cx;
4381 bz1/((real)(gridi->cxy_ind[ci_xy+1] - gridi->cxy_ind[ci_xy]));
4386 /* The check with bz1_frac close to or larger than 1 comes later */
4388 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
4390 shy = ty*box[YY][YY] + tz*box[ZZ][YY];
4394 by0 = bb_i[ci].lower[BB_Y] + shy;
4395 by1 = bb_i[ci].upper[BB_Y] + shy;
4399 by0 = gridi->c0[YY] + (ci_y )*gridi->sy + shy;
4400 by1 = gridi->c0[YY] + (ci_y+1)*gridi->sy + shy;
4403 get_cell_range(by0, by1,
4404 gridj->ncy, gridj->c0[YY], gridj->sy, gridj->inv_sy,
4414 if (by1 < gridj->c0[YY])
4416 d2z_cy += sqr(gridj->c0[YY] - by1);
4418 else if (by0 > gridj->c1[YY])
4420 d2z_cy += sqr(by0 - gridj->c1[YY]);
4423 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
4425 shift = XYZ2IS(tx, ty, tz);
4427 #ifdef NBNXN_SHIFT_BACKWARD
4428 if (gridi == gridj && shift > CENTRAL)
4434 shx = tx*box[XX][XX] + ty*box[YY][XX] + tz*box[ZZ][XX];
4438 bx0 = bb_i[ci].lower[BB_X] + shx;
4439 bx1 = bb_i[ci].upper[BB_X] + shx;
4443 bx0 = gridi->c0[XX] + (ci_x )*gridi->sx + shx;
4444 bx1 = gridi->c0[XX] + (ci_x+1)*gridi->sx + shx;
4447 get_cell_range(bx0, bx1,
4448 gridj->ncx, gridj->c0[XX], gridj->sx, gridj->inv_sx,
4459 new_ci_entry(nbl, cell0_i+ci, shift, flags_i[ci]);
4463 new_sci_entry(nbl, cell0_i+ci, shift);
4466 #ifndef NBNXN_SHIFT_BACKWARD
4469 if (shift == CENTRAL && gridi == gridj &&
4473 /* Leave the pairs with i > j.
4474 * x is the major index, so skip half of it.
4481 set_icell_bb_simple(bb_i, ci, shx, shy, shz,
4487 set_icell_bbxxxx_supersub(pbb_i, ci, shx, shy, shz,
4490 set_icell_bb_supersub(bb_i, ci, shx, shy, shz,
4495 nbs->icell_set_x(cell0_i+ci, shx, shy, shz,
4496 gridi->na_c, nbat->xstride, nbat->x,
4499 for (cx = cxf; cx <= cxl; cx++)
4502 if (gridj->c0[XX] + cx*gridj->sx > bx1)
4504 d2zx += sqr(gridj->c0[XX] + cx*gridj->sx - bx1);
4506 else if (gridj->c0[XX] + (cx+1)*gridj->sx < bx0)
4508 d2zx += sqr(gridj->c0[XX] + (cx+1)*gridj->sx - bx0);
4511 #ifndef NBNXN_SHIFT_BACKWARD
4512 if (gridi == gridj &&
4513 cx == 0 && cyf < ci_y)
4515 if (gridi == gridj &&
4516 cx == 0 && shift == CENTRAL && cyf < ci_y)
4519 /* Leave the pairs with i > j.
4520 * Skip half of y when i and j have the same x.
4529 for (cy = cyf_x; cy <= cyl; cy++)
4531 c0 = gridj->cxy_ind[cx*gridj->ncy+cy];
4532 c1 = gridj->cxy_ind[cx*gridj->ncy+cy+1];
4533 #ifdef NBNXN_SHIFT_BACKWARD
4534 if (gridi == gridj &&
4535 shift == CENTRAL && c0 < ci)
4542 if (gridj->c0[YY] + cy*gridj->sy > by1)
4544 d2zxy += sqr(gridj->c0[YY] + cy*gridj->sy - by1);
4546 else if (gridj->c0[YY] + (cy+1)*gridj->sy < by0)
4548 d2zxy += sqr(gridj->c0[YY] + (cy+1)*gridj->sy - by0);
4550 if (c1 > c0 && d2zxy < rl2)
4552 cs = c0 + (int)(bz1_frac*(c1 - c0));
4560 /* Find the lowest cell that can possibly
4565 (bbcz_j[cf*NNBSBB_D+1] >= bz0 ||
4566 d2xy + sqr(bbcz_j[cf*NNBSBB_D+1] - bz0) < rl2))
4571 /* Find the highest cell that can possibly
4576 (bbcz_j[cl*NNBSBB_D] <= bz1 ||
4577 d2xy + sqr(bbcz_j[cl*NNBSBB_D] - bz1) < rl2))
4582 #ifdef NBNXN_REFCODE
4584 /* Simple reference code, for debugging,
4585 * overrides the more complex code above.
4590 for (k = c0; k < c1; k++)
4592 if (box_dist2(bx0, bx1, by0, by1, bz0, bz1, bb+k) < rl2 &&
4597 if (box_dist2(bx0, bx1, by0, by1, bz0, bz1, bb+k) < rl2 &&
4608 /* We want each atom/cell pair only once,
4609 * only use cj >= ci.
4611 #ifndef NBNXN_SHIFT_BACKWARD
4614 if (shift == CENTRAL)
4623 /* For f buffer flags with simple lists */
4624 ncj_old_j = nbl->ncj;
4626 switch (nb_kernel_type)
4628 case nbnxnk4x4_PlainC:
4629 check_subcell_list_space_simple(nbl, cl-cf+1);
4631 make_cluster_list_simple(gridj,
4633 (gridi == gridj && shift == CENTRAL),
4638 #ifdef GMX_NBNXN_SIMD_4XN
4639 case nbnxnk4xN_SIMD_4xN:
4640 check_subcell_list_space_simple(nbl, ci_to_cj(na_cj_2log, cl-cf)+2);
4641 make_cluster_list_simd_4xn(gridj,
4643 (gridi == gridj && shift == CENTRAL),
4649 #ifdef GMX_NBNXN_SIMD_2XNN
4650 case nbnxnk4xN_SIMD_2xNN:
4651 check_subcell_list_space_simple(nbl, ci_to_cj(na_cj_2log, cl-cf)+2);
4652 make_cluster_list_simd_2xnn(gridj,
4654 (gridi == gridj && shift == CENTRAL),
4660 case nbnxnk8x8x8_PlainC:
4661 case nbnxnk8x8x8_CUDA:
4662 check_subcell_list_space_supersub(nbl, cl-cf+1);
4663 for (cj = cf; cj <= cl; cj++)
4665 make_cluster_list_supersub(gridi, gridj,
4667 (gridi == gridj && shift == CENTRAL && ci == cj),
4668 nbat->xstride, nbat->x,
4674 ncpcheck += cl - cf + 1;
4676 if (bFBufferFlag && nbl->ncj > ncj_old_j)
4680 cbf = nbl->cj[ncj_old_j].cj >> gridj_flag_shift;
4681 cbl = nbl->cj[nbl->ncj-1].cj >> gridj_flag_shift;
4682 for (cb = cbf; cb <= cbl; cb++)
4684 gridj_flag[cb] = 1U<<th;
4692 /* Set the exclusions for this ci list */
4695 set_ci_top_excls(nbs,
4697 shift == CENTRAL && gridi == gridj,
4700 &(nbl->ci[nbl->nci]),
4705 set_sci_top_excls(nbs,
4707 shift == CENTRAL && gridi == gridj,
4709 &(nbl->sci[nbl->nsci]),
4713 /* Close this ci list */
4716 close_ci_entry_simple(nbl);
4720 close_ci_entry_supersub(nbl,
4722 progBal, min_ci_balanced,
4729 if (bFBufferFlag && nbl->ncj > ncj_old_i)
4731 work->buffer_flags.flag[(gridi->cell0+ci)>>gridi_flag_shift] = 1U<<th;
4735 work->ndistc = ndistc;
4737 nbs_cycle_stop(&work->cc[enbsCCsearch]);
4741 fprintf(debug, "number of distance checks %d\n", ndistc);
4742 fprintf(debug, "ncpcheck %s %d\n", gridi == gridj ? "local" : "non-local",
4747 print_nblist_statistics_simple(debug, nbl, nbs, rlist);
4751 print_nblist_statistics_supersub(debug, nbl, nbs, rlist);
4757 static void reduce_buffer_flags(const nbnxn_search_t nbs,
4759 const nbnxn_buffer_flags_t *dest)
4762 const unsigned *flag;
4764 for (s = 0; s < nsrc; s++)
4766 flag = nbs->work[s].buffer_flags.flag;
4768 for (b = 0; b < dest->nflag; b++)
4770 dest->flag[b] |= flag[b];
4775 static void print_reduction_cost(const nbnxn_buffer_flags_t *flags, int nout)
4777 int nelem, nkeep, ncopy, nred, b, c, out;
4783 for (b = 0; b < flags->nflag; b++)
4785 if (flags->flag[b] == 1)
4787 /* Only flag 0 is set, no copy of reduction required */
4791 else if (flags->flag[b] > 0)
4794 for (out = 0; out < nout; out++)
4796 if (flags->flag[b] & (1U<<out))
4813 fprintf(debug, "nbnxn reduction: #flag %d #list %d elem %4.2f, keep %4.2f copy %4.2f red %4.2f\n",
4815 nelem/(double)(flags->nflag),
4816 nkeep/(double)(flags->nflag),
4817 ncopy/(double)(flags->nflag),
4818 nred/(double)(flags->nflag));
4821 /* Perform a count (linear) sort to sort the smaller lists to the end.
4822 * This avoids load imbalance on the GPU, as large lists will be
4823 * scheduled and executed first and the smaller lists later.
4824 * Load balancing between multi-processors only happens at the end
4825 * and there smaller lists lead to more effective load balancing.
4826 * The sorting is done on the cj4 count, not on the actual pair counts.
4827 * Not only does this make the sort faster, but it also results in
4828 * better load balancing than using a list sorted on exact load.
4829 * This function swaps the pointer in the pair list to avoid a copy operation.
4831 static void sort_sci(nbnxn_pairlist_t *nbl)
4833 nbnxn_list_work_t *work;
4834 int m, i, s, s0, s1;
4835 nbnxn_sci_t *sci_sort;
4837 if (nbl->ncj4 <= nbl->nsci)
4839 /* nsci = 0 or all sci have size 1, sorting won't change the order */
4845 /* We will distinguish differences up to double the average */
4846 m = (2*nbl->ncj4)/nbl->nsci;
4848 if (m + 1 > work->sort_nalloc)
4850 work->sort_nalloc = over_alloc_large(m + 1);
4851 srenew(work->sort, work->sort_nalloc);
4854 if (work->sci_sort_nalloc != nbl->sci_nalloc)
4856 work->sci_sort_nalloc = nbl->sci_nalloc;
4857 nbnxn_realloc_void((void **)&work->sci_sort,
4859 work->sci_sort_nalloc*sizeof(*work->sci_sort),
4860 nbl->alloc, nbl->free);
4863 /* Count the entries of each size */
4864 for (i = 0; i <= m; i++)
4868 for (s = 0; s < nbl->nsci; s++)
4870 i = min(m, nbl->sci[s].cj4_ind_end - nbl->sci[s].cj4_ind_start);
4873 /* Calculate the offset for each count */
4876 for (i = m - 1; i >= 0; i--)
4879 work->sort[i] = work->sort[i + 1] + s0;
4883 /* Sort entries directly into place */
4884 sci_sort = work->sci_sort;
4885 for (s = 0; s < nbl->nsci; s++)
4887 i = min(m, nbl->sci[s].cj4_ind_end - nbl->sci[s].cj4_ind_start);
4888 sci_sort[work->sort[i]++] = nbl->sci[s];
4891 /* Swap the sci pointers so we use the new, sorted list */
4892 work->sci_sort = nbl->sci;
4893 nbl->sci = sci_sort;
4896 /* Make a local or non-local pair-list, depending on iloc */
4897 void nbnxn_make_pairlist(const nbnxn_search_t nbs,
4898 nbnxn_atomdata_t *nbat,
4899 const t_blocka *excl,
4901 int min_ci_balanced,
4902 nbnxn_pairlist_set_t *nbl_list,
4907 nbnxn_grid_t *gridi, *gridj;
4909 int nzi, zi, zj0, zj1, zj;
4913 nbnxn_pairlist_t **nbl;
4915 gmx_bool CombineNBLists;
4917 int np_tot, np_noq, np_hlj, nap;
4919 /* Check if we are running hybrid GPU + CPU nbnxn mode */
4920 bGPUCPU = (!nbs->grid[0].bSimple && nbl_list->bSimple);
4922 nnbl = nbl_list->nnbl;
4923 nbl = nbl_list->nbl;
4924 CombineNBLists = nbl_list->bCombined;
4928 fprintf(debug, "ns making %d nblists\n", nnbl);
4931 nbat->bUseBufferFlags = (nbat->nout > 1);
4932 /* We should re-init the flags before making the first list */
4933 if (nbat->bUseBufferFlags && (LOCAL_I(iloc) || bGPUCPU))
4935 init_buffer_flags(&nbat->buffer_flags, nbat->natoms);
4938 if (nbl_list->bSimple)
4940 switch (nb_kernel_type)
4942 #ifdef GMX_NBNXN_SIMD_4XN
4943 case nbnxnk4xN_SIMD_4xN:
4944 nbs->icell_set_x = icell_set_x_simd_4xn;
4947 #ifdef GMX_NBNXN_SIMD_2XNN
4948 case nbnxnk4xN_SIMD_2xNN:
4949 nbs->icell_set_x = icell_set_x_simd_2xnn;
4953 nbs->icell_set_x = icell_set_x_simple;
4959 #ifdef NBNXN_SEARCH_BB_SIMD4
4960 nbs->icell_set_x = icell_set_x_supersub_simd4;
4962 nbs->icell_set_x = icell_set_x_supersub;
4968 /* Only zone (grid) 0 vs 0 */
4975 nzi = nbs->zones->nizone;
4978 if (!nbl_list->bSimple && min_ci_balanced > 0)
4980 nsubpair_max = get_nsubpair_max(nbs, iloc, rlist, min_ci_balanced);
4987 /* Clear all pair-lists */
4988 for (th = 0; th < nnbl; th++)
4990 clear_pairlist(nbl[th]);
4993 for (zi = 0; zi < nzi; zi++)
4995 gridi = &nbs->grid[zi];
4997 if (NONLOCAL_I(iloc))
4999 zj0 = nbs->zones->izone[zi].j0;
5000 zj1 = nbs->zones->izone[zi].j1;
5006 for (zj = zj0; zj < zj1; zj++)
5008 gridj = &nbs->grid[zj];
5012 fprintf(debug, "ns search grid %d vs %d\n", zi, zj);
5015 nbs_cycle_start(&nbs->cc[enbsCCsearch]);
5017 if (nbl[0]->bSimple && !gridi->bSimple)
5019 /* Hybrid list, determine blocking later */
5024 ci_block = get_ci_block_size(gridi, nbs->DomDec, nnbl);
5027 #pragma omp parallel for num_threads(nnbl) schedule(static)
5028 for (th = 0; th < nnbl; th++)
5030 /* Re-init the thread-local work flag data before making
5031 * the first list (not an elegant conditional).
5033 if (nbat->bUseBufferFlags && ((zi == 0 && zj == 0) ||
5034 (bGPUCPU && zi == 0 && zj == 1)))
5036 init_buffer_flags(&nbs->work[th].buffer_flags, nbat->natoms);
5039 if (CombineNBLists && th > 0)
5041 clear_pairlist(nbl[th]);
5044 /* With GPU: generate progressively smaller lists for
5045 * load balancing for local only or non-local with 2 zones.
5047 progBal = (LOCAL_I(iloc) || nbs->zones->n <= 2);
5049 /* Divide the i super cell equally over the nblists */
5050 nbnxn_make_pairlist_part(nbs, gridi, gridj,
5051 &nbs->work[th], nbat, excl,
5055 nbat->bUseBufferFlags,
5057 progBal, min_ci_balanced,
5061 nbs_cycle_stop(&nbs->cc[enbsCCsearch]);
5066 for (th = 0; th < nnbl; th++)
5068 inc_nrnb(nrnb, eNR_NBNXN_DIST2, nbs->work[th].ndistc);
5070 if (nbl_list->bSimple)
5072 np_tot += nbl[th]->ncj;
5073 np_noq += nbl[th]->work->ncj_noq;
5074 np_hlj += nbl[th]->work->ncj_hlj;
5078 /* This count ignores potential subsequent pair pruning */
5079 np_tot += nbl[th]->nci_tot;
5082 nap = nbl[0]->na_ci*nbl[0]->na_cj;
5083 nbl_list->natpair_ljq = (np_tot - np_noq)*nap - np_hlj*nap/2;
5084 nbl_list->natpair_lj = np_noq*nap;
5085 nbl_list->natpair_q = np_hlj*nap/2;
5087 if (CombineNBLists && nnbl > 1)
5089 nbs_cycle_start(&nbs->cc[enbsCCcombine]);
5091 combine_nblists(nnbl-1, nbl+1, nbl[0]);
5093 nbs_cycle_stop(&nbs->cc[enbsCCcombine]);
5098 if (!nbl_list->bSimple)
5100 /* Sort the entries on size, large ones first */
5101 if (CombineNBLists || nnbl == 1)
5107 #pragma omp parallel for num_threads(nnbl) schedule(static)
5108 for (th = 0; th < nnbl; th++)
5115 if (nbat->bUseBufferFlags)
5117 reduce_buffer_flags(nbs, nnbl, &nbat->buffer_flags);
5120 /* Special performance logging stuff (env.var. GMX_NBNXN_CYCLE) */
5123 nbs->search_count++;
5125 if (nbs->print_cycles &&
5126 (!nbs->DomDec || (nbs->DomDec && !LOCAL_I(iloc))) &&
5127 nbs->search_count % 100 == 0)
5129 nbs_cycle_print(stderr, nbs);
5132 if (debug && (CombineNBLists && nnbl > 1))
5134 if (nbl[0]->bSimple)
5136 print_nblist_statistics_simple(debug, nbl[0], nbs, rlist);
5140 print_nblist_statistics_supersub(debug, nbl[0], nbs, rlist);
5148 if (nbl[0]->bSimple)
5150 print_nblist_ci_cj(debug, nbl[0]);
5154 print_nblist_sci_cj(debug, nbl[0]);
5158 if (nbat->bUseBufferFlags)
5160 print_reduction_cost(&nbat->buffer_flags, nnbl);