2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2012, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
51 #include "nbnxn_consts.h"
52 #include "nbnxn_internal.h"
53 #include "nbnxn_atomdata.h"
54 #include "nbnxn_search.h"
55 #include "gmx_cyclecounter.h"
57 #include "gmx_omp_nthreads.h"
61 /* Pair search box lower and upper corner in x,y,z.
62 * Store this in 4 iso 3 reals, which is useful with SSE.
63 * To avoid complicating the code we also use 4 without SSE.
66 #define NNBSBB_B (2*NNBSBB_C)
67 /* Pair search box lower and upper bound in z only. */
69 /* Pair search box lower and upper corner x,y,z indices */
78 #ifdef NBNXN_SEARCH_SSE
81 #define NBNXN_SEARCH_SSE_SINGLE
82 #include "gmx_x86_simd_single.h"
84 #include "gmx_x86_simd_double.h"
87 #if defined NBNXN_SEARCH_SSE_SINGLE && GPU_NSUBCELL == 8
91 /* The width of SSE/AVX128 with single precision for bounding boxes with GPU.
92 * Here AVX-256 turns out to be slightly slower than AVX-128.
95 #define STRIDE_8BB_2LOG 2
98 /* The functions below are macros as they are performance sensitive */
100 /* 4x4 list, pack=4: no complex conversion required */
101 /* i-cluster to j-cluster conversion */
102 #define CI_TO_CJ_J4(ci) (ci)
103 /* cluster index to coordinate array index conversion */
104 #define X_IND_CI_J4(ci) ((ci)*STRIDE_P4)
105 #define X_IND_CJ_J4(cj) ((cj)*STRIDE_P4)
107 /* 4x2 list, pack=4: j-cluster size is half the packing width */
108 /* i-cluster to j-cluster conversion */
109 #define CI_TO_CJ_J2(ci) ((ci)<<1)
110 /* cluster index to coordinate array index conversion */
111 #define X_IND_CI_J2(ci) ((ci)*STRIDE_P4)
112 #define X_IND_CJ_J2(cj) (((cj)>>1)*STRIDE_P4 + ((cj) & 1)*(PACK_X4>>1))
114 /* 4x8 list, pack=8: i-cluster size is half the packing width */
115 /* i-cluster to j-cluster conversion */
116 #define CI_TO_CJ_J8(ci) ((ci)>>1)
117 /* cluster index to coordinate array index conversion */
118 #define X_IND_CI_J8(ci) (((ci)>>1)*STRIDE_P8 + ((ci) & 1)*(PACK_X8>>1))
119 #define X_IND_CJ_J8(cj) ((cj)*STRIDE_P8)
121 /* The j-cluster size is matched to the SIMD width */
123 /* 128 bits can hold 4 floats */
124 #define CI_TO_CJ_S128(ci) CI_TO_CJ_J4(ci)
125 #define X_IND_CI_S128(ci) X_IND_CI_J4(ci)
126 #define X_IND_CJ_S128(cj) X_IND_CJ_J4(cj)
127 /* 256 bits can hold 8 floats */
128 #define CI_TO_CJ_S256(ci) CI_TO_CJ_J8(ci)
129 #define X_IND_CI_S256(ci) X_IND_CI_J8(ci)
130 #define X_IND_CJ_S256(cj) X_IND_CJ_J8(cj)
132 /* 128 bits can hold 2 doubles */
133 #define CI_TO_CJ_S128(ci) CI_TO_CJ_J2(ci)
134 #define X_IND_CI_S128(ci) X_IND_CI_J2(ci)
135 #define X_IND_CJ_S128(cj) X_IND_CJ_J2(cj)
136 /* 256 bits can hold 4 doubles */
137 #define CI_TO_CJ_S256(ci) CI_TO_CJ_J4(ci)
138 #define X_IND_CI_S256(ci) X_IND_CI_J4(ci)
139 #define X_IND_CJ_S256(cj) X_IND_CJ_J4(cj)
142 #endif /* NBNXN_SEARCH_SSE */
145 /* Interaction masks for 4xN atom interactions.
146 * Bit i*CJ_SIZE + j tells if atom i and j interact.
148 /* All interaction mask is the same for all kernels */
149 #define NBNXN_INT_MASK_ALL 0xffffffff
150 /* 4x4 kernel diagonal mask */
151 #define NBNXN_INT_MASK_DIAG 0x08ce
152 /* 4x2 kernel diagonal masks */
153 #define NBNXN_INT_MASK_DIAG_J2_0 0x0002
154 #define NBNXN_INT_MASK_DIAG_J2_1 0x002F
155 /* 4x8 kernel diagonal masks */
156 #define NBNXN_INT_MASK_DIAG_J8_0 0xf0f8fcfe
157 #define NBNXN_INT_MASK_DIAG_J8_1 0x0080c0e0
160 #ifdef NBNXN_SEARCH_SSE
161 /* Store bounding boxes corners as quadruplets: xxxxyyyyzzzz */
163 /* Size of bounding box corners quadruplet */
164 #define NNBSBB_XXXX (NNBSBB_D*DIM*STRIDE_8BB)
167 /* We shift the i-particles backward for PBC.
168 * This leads to more conditionals than shifting forward.
169 * We do this to get more balanced pair lists.
171 #define NBNXN_SHIFT_BACKWARD
174 /* This define is a lazy way to avoid interdependence of the grid
175 * and searching data structures.
177 #define NBNXN_NA_SC_MAX (GPU_NSUBCELL*NBNXN_GPU_CLUSTER_SIZE)
180 static void nbs_cycle_clear(nbnxn_cycle_t *cc)
184 for(i=0; i<enbsCCnr; i++)
191 static double Mcyc_av(const nbnxn_cycle_t *cc)
193 return (double)cc->c*1e-6/cc->count;
196 static void nbs_cycle_print(FILE *fp,const nbnxn_search_t nbs)
202 fprintf(fp,"ns %4d grid %4.1f search %4.1f red.f %5.3f",
203 nbs->cc[enbsCCgrid].count,
204 Mcyc_av(&nbs->cc[enbsCCgrid]),
205 Mcyc_av(&nbs->cc[enbsCCsearch]),
206 Mcyc_av(&nbs->cc[enbsCCreducef]));
208 if (nbs->nthread_max > 1)
210 if (nbs->cc[enbsCCcombine].count > 0)
212 fprintf(fp," comb %5.2f",
213 Mcyc_av(&nbs->cc[enbsCCcombine]));
215 fprintf(fp," s. th");
216 for(t=0; t<nbs->nthread_max; t++)
219 Mcyc_av(&nbs->work[t].cc[enbsCCsearch]));
225 static void nbnxn_grid_init(nbnxn_grid_t * grid)
228 grid->cxy_ind = NULL;
229 grid->cxy_nalloc = 0;
235 static int get_2log(int n)
240 while ((1<<log2) < n)
246 gmx_fatal(FARGS,"nbnxn na_c (%d) is not a power of 2",n);
252 static int nbnxn_kernel_to_ci_size(int nb_kernel_type)
254 switch (nb_kernel_type)
257 case nbk4xN_X86_SIMD128:
258 case nbk4xN_X86_SIMD256:
259 return NBNXN_CPU_CLUSTER_I_SIZE;
261 case nbk8x8x8_PlainC:
262 /* The cluster size for super/sub lists is only set here.
263 * Any value should work for the pair-search and atomdata code.
264 * The kernels, of course, might require a particular value.
266 return NBNXN_GPU_CLUSTER_SIZE;
268 gmx_incons("unknown kernel type");
274 int nbnxn_kernel_to_cj_size(int nb_kernel_type)
276 switch (nb_kernel_type)
279 return NBNXN_CPU_CLUSTER_I_SIZE;
280 case nbk4xN_X86_SIMD128:
281 /* Number of reals that fit in SIMD (128 bits = 16 bytes) */
282 return 16/sizeof(real);
283 case nbk4xN_X86_SIMD256:
284 /* Number of reals that fit in SIMD (256 bits = 32 bytes) */
285 return 32/sizeof(real);
287 case nbk8x8x8_PlainC:
288 return nbnxn_kernel_to_ci_size(nb_kernel_type);
290 gmx_incons("unknown kernel type");
296 static int ci_to_cj(int na_cj_2log,int ci)
300 case 2: return ci; break;
301 case 1: return (ci<<1); break;
302 case 3: return (ci>>1); break;
308 gmx_bool nbnxn_kernel_pairlist_simple(int nb_kernel_type)
310 if (nb_kernel_type == nbkNotSet)
312 gmx_fatal(FARGS, "Non-bonded kernel type not set for Verlet-style pair-list.");
315 switch (nb_kernel_type)
318 case nbk8x8x8_PlainC:
322 case nbk4xN_X86_SIMD128:
323 case nbk4xN_X86_SIMD256:
327 gmx_incons("Invalid nonbonded kernel type passed!");
332 void nbnxn_init_search(nbnxn_search_t * nbs_ptr,
334 gmx_domdec_zones_t *zones,
343 nbs->DomDec = (n_dd_cells != NULL);
345 clear_ivec(nbs->dd_dim);
353 if ((*n_dd_cells)[d] > 1)
356 /* Each grid matches a DD zone */
362 snew(nbs->grid,nbs->ngrid);
363 for(g=0; g<nbs->ngrid; g++)
365 nbnxn_grid_init(&nbs->grid[g]);
368 nbs->cell_nalloc = 0;
372 nbs->nthread_max = nthread_max;
374 /* Initialize the work data structures for each thread */
375 snew(nbs->work,nbs->nthread_max);
376 for(t=0; t<nbs->nthread_max; t++)
378 nbs->work[t].cxy_na = NULL;
379 nbs->work[t].cxy_na_nalloc = 0;
380 nbs->work[t].sort_work = NULL;
381 nbs->work[t].sort_work_nalloc = 0;
384 /* Initialize detailed nbsearch cycle counting */
385 nbs->print_cycles = (getenv("GMX_NBNXN_CYCLE") != 0);
386 nbs->search_count = 0;
387 nbs_cycle_clear(nbs->cc);
388 for(t=0; t<nbs->nthread_max; t++)
390 nbs_cycle_clear(nbs->work[t].cc);
394 static real grid_atom_density(int n,rvec corner0,rvec corner1)
398 rvec_sub(corner1,corner0,size);
400 return n/(size[XX]*size[YY]*size[ZZ]);
403 static int set_grid_size_xy(const nbnxn_search_t nbs,
405 int n,rvec corner0,rvec corner1,
411 real adens,tlen,tlen_x,tlen_y,nc_max;
414 rvec_sub(corner1,corner0,size);
418 /* target cell length */
421 /* To minimize the zero interactions, we should make
422 * the largest of the i/j cell cubic.
424 na_c = max(grid->na_c,grid->na_cj);
426 /* Approximately cubic cells */
427 tlen = pow(na_c/atom_density,1.0/3.0);
433 /* Approximately cubic sub cells */
434 tlen = pow(grid->na_c/atom_density,1.0/3.0);
435 tlen_x = tlen*GPU_NSUBCELL_X;
436 tlen_y = tlen*GPU_NSUBCELL_Y;
438 /* We round ncx and ncy down, because we get less cell pairs
439 * in the nbsist when the fixed cell dimensions (x,y) are
440 * larger than the variable one (z) than the other way around.
442 grid->ncx = max(1,(int)(size[XX]/tlen_x));
443 grid->ncy = max(1,(int)(size[YY]/tlen_y));
451 /* We need one additional cell entry for particles moved by DD */
452 if (grid->ncx*grid->ncy+1 > grid->cxy_nalloc)
454 grid->cxy_nalloc = over_alloc_large(grid->ncx*grid->ncy+1);
455 srenew(grid->cxy_na,grid->cxy_nalloc);
456 srenew(grid->cxy_ind,grid->cxy_nalloc+1);
458 for(t=0; t<nbs->nthread_max; t++)
460 if (grid->ncx*grid->ncy+1 > nbs->work[t].cxy_na_nalloc)
462 nbs->work[t].cxy_na_nalloc = over_alloc_large(grid->ncx*grid->ncy+1);
463 srenew(nbs->work[t].cxy_na,nbs->work[t].cxy_na_nalloc);
467 /* Worst case scenario of 1 atom in each last cell */
468 if (grid->na_cj <= grid->na_c)
470 nc_max = n/grid->na_sc + grid->ncx*grid->ncy;
474 nc_max = n/grid->na_sc + grid->ncx*grid->ncy*grid->na_cj/grid->na_c;
477 if (nc_max > grid->nc_nalloc)
481 grid->nc_nalloc = over_alloc_large(nc_max);
482 srenew(grid->nsubc,grid->nc_nalloc);
483 srenew(grid->bbcz,grid->nc_nalloc*NNBSBB_D);
485 bb_nalloc = grid->nc_nalloc*GPU_NSUBCELL/STRIDE_8BB*NNBSBB_XXXX;
487 bb_nalloc = grid->nc_nalloc*GPU_NSUBCELL*NNBSBB_B;
489 sfree_aligned(grid->bb);
490 /* This snew also zeros the contents, this avoid possible
491 * floating exceptions in SSE with the unused bb elements.
493 snew_aligned(grid->bb,bb_nalloc,16);
497 if (grid->na_cj == grid->na_c)
499 grid->bbj = grid->bb;
503 sfree_aligned(grid->bbj);
504 snew_aligned(grid->bbj,bb_nalloc*grid->na_c/grid->na_cj,16);
508 srenew(grid->flags,grid->nc_nalloc);
511 copy_rvec(corner0,grid->c0);
512 copy_rvec(corner1,grid->c1);
513 grid->sx = size[XX]/grid->ncx;
514 grid->sy = size[YY]/grid->ncy;
515 grid->inv_sx = 1/grid->sx;
516 grid->inv_sy = 1/grid->sy;
521 #define SORT_GRID_OVERSIZE 2
522 #define SGSF (SORT_GRID_OVERSIZE + 1)
524 static void sort_atoms(int dim,gmx_bool Backwards,
525 int *a,int n,rvec *x,
526 real h0,real invh,int nsort,int *sort)
538 /* For small oversize factors clearing the whole area is fastest.
539 * For large oversize we should clear the used elements after use.
541 for(i=0; i<nsort; i++)
545 /* Sort the particles using a simple index sort */
548 /* The cast takes care of float-point rounding effects below zero.
549 * This code assumes particles are less than 1/SORT_GRID_OVERSIZE
550 * times the box height out of the box.
552 zi = (int)((x[a[i]][dim] - h0)*invh);
554 #ifdef DEBUG_NBNXN_GRIDDING
555 if (zi < 0 || zi >= nsort)
557 gmx_fatal(FARGS,"(int)((x[%d][%c]=%f - %f)*%f) = %d, not in 0 - %d\n",
558 a[i],'x'+dim,x[a[i]][dim],h0,invh,zi,nsort);
562 /* Ideally this particle should go in sort cell zi,
563 * but that might already be in use,
564 * in that case find the first empty cell higher up
572 /* We have multiple atoms in the same sorting slot.
573 * Sort on real z for minimal bounding box size.
574 * There is an extra check for identical z to ensure
575 * well-defined output order, independent of input order
576 * to ensure binary reproducibility after restarts.
578 while(sort[zi] >= 0 && ( x[a[i]][dim] > x[sort[zi]][dim] ||
579 (x[a[i]][dim] == x[sort[zi]][dim] &&
587 /* Shift all elements by one slot until we find an empty slot */
590 while (sort[zim] >= 0)
606 for(zi=0; zi<nsort; zi++)
616 for(zi=nsort-1; zi>=0; zi--)
626 gmx_incons("Lost particles while sorting");
631 #define R2F_D(x) ((float)((x) >= 0 ? ((1-GMX_FLOAT_EPS)*(x)) : ((1+GMX_FLOAT_EPS)*(x))))
632 #define R2F_U(x) ((float)((x) >= 0 ? ((1+GMX_FLOAT_EPS)*(x)) : ((1-GMX_FLOAT_EPS)*(x))))
638 /* Coordinate order x,y,z, bb order xyz0 */
639 static void calc_bounding_box(int na,int stride,const real *x,float *bb)
642 real xl,xh,yl,yh,zl,zh;
654 xl = min(xl,x[i+XX]);
655 xh = max(xh,x[i+XX]);
656 yl = min(yl,x[i+YY]);
657 yh = max(yh,x[i+YY]);
658 zl = min(zl,x[i+ZZ]);
659 zh = max(zh,x[i+ZZ]);
662 /* Note: possible double to float conversion here */
663 bb[BBL_X] = R2F_D(xl);
664 bb[BBL_Y] = R2F_D(yl);
665 bb[BBL_Z] = R2F_D(zl);
666 bb[BBU_X] = R2F_U(xh);
667 bb[BBU_Y] = R2F_U(yh);
668 bb[BBU_Z] = R2F_U(zh);
671 /* Packed coordinates, bb order xyz0 */
672 static void calc_bounding_box_x_x4(int na,const real *x,float *bb)
675 real xl,xh,yl,yh,zl,zh;
685 xl = min(xl,x[j+XX*PACK_X4]);
686 xh = max(xh,x[j+XX*PACK_X4]);
687 yl = min(yl,x[j+YY*PACK_X4]);
688 yh = max(yh,x[j+YY*PACK_X4]);
689 zl = min(zl,x[j+ZZ*PACK_X4]);
690 zh = max(zh,x[j+ZZ*PACK_X4]);
692 /* Note: possible double to float conversion here */
693 bb[BBL_X] = R2F_D(xl);
694 bb[BBL_Y] = R2F_D(yl);
695 bb[BBL_Z] = R2F_D(zl);
696 bb[BBU_X] = R2F_U(xh);
697 bb[BBU_Y] = R2F_U(yh);
698 bb[BBU_Z] = R2F_U(zh);
701 /* Packed coordinates, bb order xyz0 */
702 static void calc_bounding_box_x_x8(int na,const real *x,float *bb)
705 real xl,xh,yl,yh,zl,zh;
715 xl = min(xl,x[j+XX*PACK_X8]);
716 xh = max(xh,x[j+XX*PACK_X8]);
717 yl = min(yl,x[j+YY*PACK_X8]);
718 yh = max(yh,x[j+YY*PACK_X8]);
719 zl = min(zl,x[j+ZZ*PACK_X8]);
720 zh = max(zh,x[j+ZZ*PACK_X8]);
722 /* Note: possible double to float conversion here */
723 bb[BBL_X] = R2F_D(xl);
724 bb[BBL_Y] = R2F_D(yl);
725 bb[BBL_Z] = R2F_D(zl);
726 bb[BBU_X] = R2F_U(xh);
727 bb[BBU_Y] = R2F_U(yh);
728 bb[BBU_Z] = R2F_U(zh);
731 #ifdef NBNXN_SEARCH_SSE
733 /* Packed coordinates, bb order xyz0 */
734 static void calc_bounding_box_x_x4_halves(int na,const real *x,
735 float *bb,float *bbj)
737 calc_bounding_box_x_x4(min(na,2),x,bbj);
741 calc_bounding_box_x_x4(min(na-2,2),x+(PACK_X4>>1),bbj+NNBSBB_B);
745 /* Set the "empty" bounding box to the same as the first one,
746 * so we don't need to treat special cases in the rest of the code.
748 _mm_store_ps(bbj+NNBSBB_B ,_mm_load_ps(bbj));
749 _mm_store_ps(bbj+NNBSBB_B+NNBSBB_C,_mm_load_ps(bbj+NNBSBB_C));
752 _mm_store_ps(bb ,_mm_min_ps(_mm_load_ps(bbj),
753 _mm_load_ps(bbj+NNBSBB_B)));
754 _mm_store_ps(bb+NNBSBB_C,_mm_max_ps(_mm_load_ps(bbj+NNBSBB_C),
755 _mm_load_ps(bbj+NNBSBB_B+NNBSBB_C)));
758 /* Coordinate order xyz, bb order xxxxyyyyzzzz */
759 static void calc_bounding_box_xxxx(int na,int stride,const real *x,float *bb)
762 real xl,xh,yl,yh,zl,zh;
774 xl = min(xl,x[i+XX]);
775 xh = max(xh,x[i+XX]);
776 yl = min(yl,x[i+YY]);
777 yh = max(yh,x[i+YY]);
778 zl = min(zl,x[i+ZZ]);
779 zh = max(zh,x[i+ZZ]);
782 /* Note: possible double to float conversion here */
783 bb[0*STRIDE_8BB] = R2F_D(xl);
784 bb[1*STRIDE_8BB] = R2F_D(yl);
785 bb[2*STRIDE_8BB] = R2F_D(zl);
786 bb[3*STRIDE_8BB] = R2F_U(xh);
787 bb[4*STRIDE_8BB] = R2F_U(yh);
788 bb[5*STRIDE_8BB] = R2F_U(zh);
791 #endif /* NBNXN_SEARCH_SSE */
793 #ifdef NBNXN_SEARCH_SSE_SINGLE
795 /* Coordinate order xyz?, bb order xyz0 */
796 static void calc_bounding_box_sse(int na,const float *x,float *bb)
798 __m128 bb_0_SSE,bb_1_SSE;
803 bb_0_SSE = _mm_load_ps(x);
808 x_SSE = _mm_load_ps(x+i*NNBSBB_C);
809 bb_0_SSE = _mm_min_ps(bb_0_SSE,x_SSE);
810 bb_1_SSE = _mm_max_ps(bb_1_SSE,x_SSE);
813 _mm_store_ps(bb ,bb_0_SSE);
814 _mm_store_ps(bb+4,bb_1_SSE);
817 /* Coordinate order xyz?, bb order xxxxyyyyzzzz */
818 static void calc_bounding_box_xxxx_sse(int na,const float *x,
822 calc_bounding_box_sse(na,x,bb_work);
824 bb[0*STRIDE_8BB] = bb_work[BBL_X];
825 bb[1*STRIDE_8BB] = bb_work[BBL_Y];
826 bb[2*STRIDE_8BB] = bb_work[BBL_Z];
827 bb[3*STRIDE_8BB] = bb_work[BBU_X];
828 bb[4*STRIDE_8BB] = bb_work[BBU_Y];
829 bb[5*STRIDE_8BB] = bb_work[BBU_Z];
832 #endif /* NBNXN_SEARCH_SSE_SINGLE */
834 #ifdef NBNXN_SEARCH_SSE
836 /* Combines pairs of consecutive bounding boxes */
837 static void combine_bounding_box_pairs(nbnxn_grid_t *grid,const float *bb)
840 __m128 min_SSE,max_SSE;
842 for(i=0; i<grid->ncx*grid->ncy; i++)
844 /* Starting bb in a column is expected to be 2-aligned */
845 sc2 = grid->cxy_ind[i]>>1;
846 /* For odd numbers skip the last bb here */
847 nc2 = (grid->cxy_na[i]+3)>>(2+1);
848 for(c2=sc2; c2<sc2+nc2; c2++)
850 min_SSE = _mm_min_ps(_mm_load_ps(bb+(c2*4+0)*NNBSBB_C),
851 _mm_load_ps(bb+(c2*4+2)*NNBSBB_C));
852 max_SSE = _mm_max_ps(_mm_load_ps(bb+(c2*4+1)*NNBSBB_C),
853 _mm_load_ps(bb+(c2*4+3)*NNBSBB_C));
854 _mm_store_ps(grid->bbj+(c2*2+0)*NNBSBB_C,min_SSE);
855 _mm_store_ps(grid->bbj+(c2*2+1)*NNBSBB_C,max_SSE);
857 if (((grid->cxy_na[i]+3)>>2) & 1)
859 /* Copy the last bb for odd bb count in this column */
860 for(j=0; j<NNBSBB_C; j++)
862 grid->bbj[(c2*2+0)*NNBSBB_C+j] = bb[(c2*4+0)*NNBSBB_C+j];
863 grid->bbj[(c2*2+1)*NNBSBB_C+j] = bb[(c2*4+1)*NNBSBB_C+j];
872 /* Prints the average bb size, used for debug output */
873 static void print_bbsizes_simple(FILE *fp,
874 const nbnxn_search_t nbs,
875 const nbnxn_grid_t *grid)
881 for(c=0; c<grid->nc; c++)
885 ba[d] += grid->bb[c*NNBSBB_B+NNBSBB_C+d] - grid->bb[c*NNBSBB_B+d];
888 dsvmul(1.0/grid->nc,ba,ba);
890 fprintf(fp,"ns bb: %4.2f %4.2f %4.2f %4.2f %4.2f %4.2f rel %4.2f %4.2f %4.2f\n",
891 nbs->box[XX][XX]/grid->ncx,
892 nbs->box[YY][YY]/grid->ncy,
893 nbs->box[ZZ][ZZ]*grid->ncx*grid->ncy/grid->nc,
894 ba[XX],ba[YY],ba[ZZ],
895 ba[XX]*grid->ncx/nbs->box[XX][XX],
896 ba[YY]*grid->ncy/nbs->box[YY][YY],
897 ba[ZZ]*grid->nc/(grid->ncx*grid->ncy*nbs->box[ZZ][ZZ]));
900 /* Prints the average bb size, used for debug output */
901 static void print_bbsizes_supersub(FILE *fp,
902 const nbnxn_search_t nbs,
903 const nbnxn_grid_t *grid)
910 for(c=0; c<grid->nc; c++)
913 for(s=0; s<grid->nsubc[c]; s+=STRIDE_8BB)
917 cs_w = (c*GPU_NSUBCELL + s)/STRIDE_8BB;
918 for(i=0; i<STRIDE_8BB; i++)
923 grid->bb[cs_w*NNBSBB_XXXX+(DIM+d)*STRIDE_8BB+i] -
924 grid->bb[cs_w*NNBSBB_XXXX+ d *STRIDE_8BB+i];
929 for(s=0; s<grid->nsubc[c]; s++)
933 cs = c*GPU_NSUBCELL + s;
937 grid->bb[cs*NNBSBB_B+NNBSBB_C+d] -
938 grid->bb[cs*NNBSBB_B +d];
942 ns += grid->nsubc[c];
944 dsvmul(1.0/ns,ba,ba);
946 fprintf(fp,"ns bb: %4.2f %4.2f %4.2f %4.2f %4.2f %4.2f rel %4.2f %4.2f %4.2f\n",
947 nbs->box[XX][XX]/(grid->ncx*GPU_NSUBCELL_X),
948 nbs->box[YY][YY]/(grid->ncy*GPU_NSUBCELL_Y),
949 nbs->box[ZZ][ZZ]*grid->ncx*grid->ncy/(grid->nc*GPU_NSUBCELL_Z),
950 ba[XX],ba[YY],ba[ZZ],
951 ba[XX]*grid->ncx*GPU_NSUBCELL_X/nbs->box[XX][XX],
952 ba[YY]*grid->ncy*GPU_NSUBCELL_Y/nbs->box[YY][YY],
953 ba[ZZ]*grid->nc*GPU_NSUBCELL_Z/(grid->ncx*grid->ncy*nbs->box[ZZ][ZZ]));
956 /* Potentially sorts atoms on LJ coefficients !=0 and ==0.
957 * Also sets interaction flags.
959 void sort_on_lj(nbnxn_atomdata_t *nbat,int na_c,
960 int a0,int a1,const int *atinfo,
964 int subc,s,a,n1,n2,a_lj_max,i,j;
965 int sort1[NBNXN_NA_SC_MAX/GPU_NSUBCELL];
966 int sort2[NBNXN_NA_SC_MAX/GPU_NSUBCELL];
972 for(s=a0; s<a1; s+=na_c)
974 /* Make lists for this (sub-)cell on atoms with and without LJ */
979 for(a=s; a<min(s+na_c,a1); a++)
981 haveQ = haveQ || GET_CGINFO_HAS_Q(atinfo[order[a]]);
983 if (GET_CGINFO_HAS_VDW(atinfo[order[a]]))
985 sort1[n1++] = order[a];
990 sort2[n2++] = order[a];
994 /* If we don't have atom with LJ, there's nothing to sort */
997 *flags |= NBNXN_CI_DO_LJ(subc);
1001 /* Only sort when strictly necessary. Ordering particles
1002 * Ordering particles can lead to less accurate summation
1003 * due to rounding, both for LJ and Coulomb interactions.
1005 if (2*(a_lj_max - s) >= na_c)
1009 order[a0+i] = sort1[i];
1013 order[a0+n1+j] = sort2[j];
1017 *flags |= NBNXN_CI_HALF_LJ(subc);
1022 *flags |= NBNXN_CI_DO_COUL(subc);
1028 /* Fill a pair search cell with atoms.
1029 * Potentially sorts atoms and sets the interaction flags.
1031 void fill_cell(const nbnxn_search_t nbs,
1033 nbnxn_atomdata_t *nbat,
1037 int sx,int sy, int sz,
1048 sort_on_lj(nbat,grid->na_c,a0,a1,atinfo,nbs->a,
1049 grid->flags+(a0>>grid->na_c_2log)-grid->cell0);
1052 /* Now we have sorted the atoms, set the cell indices */
1053 for(a=a0; a<a1; a++)
1055 nbs->cell[nbs->a[a]] = a;
1058 copy_rvec_to_nbat_real(nbs->a+a0,a1-a0,grid->na_c,x,
1059 nbat->XFormat,nbat->x,a0,
1062 if (nbat->XFormat == nbatX4)
1064 /* Store the bounding boxes as xyz.xyz. */
1065 offset = ((a0 - grid->cell0*grid->na_sc)>>grid->na_c_2log)*NNBSBB_B;
1066 bb_ptr = grid->bb + offset;
1068 #if defined GMX_DOUBLE && defined NBNXN_SEARCH_SSE
1069 if (2*grid->na_cj == grid->na_c)
1071 calc_bounding_box_x_x4_halves(na,nbat->x+X4_IND_A(a0),bb_ptr,
1072 grid->bbj+offset*2);
1077 calc_bounding_box_x_x4(na,nbat->x+X4_IND_A(a0),bb_ptr);
1080 else if (nbat->XFormat == nbatX8)
1082 /* Store the bounding boxes as xyz.xyz. */
1083 offset = ((a0 - grid->cell0*grid->na_sc)>>grid->na_c_2log)*NNBSBB_B;
1084 bb_ptr = grid->bb + offset;
1086 calc_bounding_box_x_x8(na,nbat->x+X8_IND_A(a0),bb_ptr);
1089 else if (!grid->bSimple)
1091 /* Store the bounding boxes in a format convenient
1092 * for SSE calculations: xxxxyyyyzzzz...
1096 ((a0-grid->cell0*grid->na_sc)>>(grid->na_c_2log+STRIDE_8BB_2LOG))*NNBSBB_XXXX +
1097 (((a0-grid->cell0*grid->na_sc)>>grid->na_c_2log) & (STRIDE_8BB-1));
1099 #ifdef NBNXN_SEARCH_SSE_SINGLE
1100 if (nbat->XFormat == nbatXYZQ)
1102 calc_bounding_box_xxxx_sse(na,nbat->x+a0*nbat->xstride,
1108 calc_bounding_box_xxxx(na,nbat->xstride,nbat->x+a0*nbat->xstride,
1113 fprintf(debug,"%2d %2d %2d bb %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f\n",
1115 bb_ptr[0*STRIDE_8BB],bb_ptr[3*STRIDE_8BB],
1116 bb_ptr[1*STRIDE_8BB],bb_ptr[4*STRIDE_8BB],
1117 bb_ptr[2*STRIDE_8BB],bb_ptr[5*STRIDE_8BB]);
1123 /* Store the bounding boxes as xyz.xyz. */
1124 bb_ptr = grid->bb+((a0-grid->cell0*grid->na_sc)>>grid->na_c_2log)*NNBSBB_B;
1126 calc_bounding_box(na,nbat->xstride,nbat->x+a0*nbat->xstride,
1132 bbo = (a0 - grid->cell0*grid->na_sc)/grid->na_c;
1133 fprintf(debug,"%2d %2d %2d bb %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f\n",
1135 (grid->bb+bbo*NNBSBB_B)[BBL_X],
1136 (grid->bb+bbo*NNBSBB_B)[BBU_X],
1137 (grid->bb+bbo*NNBSBB_B)[BBL_Y],
1138 (grid->bb+bbo*NNBSBB_B)[BBU_Y],
1139 (grid->bb+bbo*NNBSBB_B)[BBL_Z],
1140 (grid->bb+bbo*NNBSBB_B)[BBU_Z]);
1145 /* Spatially sort the atoms within one grid column */
1146 static void sort_columns_simple(const nbnxn_search_t nbs,
1152 nbnxn_atomdata_t *nbat,
1153 int cxy_start,int cxy_end,
1157 int cx,cy,cz,ncz,cfilled,c;
1163 fprintf(debug,"cell0 %d sorting columns %d - %d, atoms %d - %d\n",
1164 grid->cell0,cxy_start,cxy_end,a0,a1);
1167 /* Sort the atoms within each x,y column in 3 dimensions */
1168 for(cxy=cxy_start; cxy<cxy_end; cxy++)
1171 cy = cxy - cx*grid->ncy;
1173 na = grid->cxy_na[cxy];
1174 ncz = grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy];
1175 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1177 /* Sort the atoms within each x,y column on z coordinate */
1178 sort_atoms(ZZ,FALSE,
1181 ncz*grid->na_sc*SORT_GRID_OVERSIZE/nbs->box[ZZ][ZZ],
1182 ncz*grid->na_sc*SGSF,sort_work);
1184 /* Fill the ncz cells in this column */
1185 cfilled = grid->cxy_ind[cxy];
1186 for(cz=0; cz<ncz; cz++)
1188 c = grid->cxy_ind[cxy] + cz ;
1190 ash_c = ash + cz*grid->na_sc;
1191 na_c = min(grid->na_sc,na-(ash_c-ash));
1193 fill_cell(nbs,grid,nbat,
1194 ash_c,ash_c+na_c,atinfo,x,
1195 grid->na_sc*cx + (dd_zone >> 2),
1196 grid->na_sc*cy + (dd_zone & 3),
1200 /* This copy to bbcz is not really necessary.
1201 * But it allows to use the same grid search code
1202 * for the simple and supersub cell setups.
1208 grid->bbcz[c*NNBSBB_D ] = grid->bb[cfilled*NNBSBB_B+2];
1209 grid->bbcz[c*NNBSBB_D+1] = grid->bb[cfilled*NNBSBB_B+6];
1212 /* Set the unused atom indices to -1 */
1213 for(ind=na; ind<ncz*grid->na_sc; ind++)
1215 nbs->a[ash+ind] = -1;
1220 /* Spatially sort the atoms within one grid column */
1221 static void sort_columns_supersub(const nbnxn_search_t nbs,
1227 nbnxn_atomdata_t *nbat,
1228 int cxy_start,int cxy_end,
1232 int cx,cy,cz=-1,c=-1,ncz;
1233 int na,ash,na_c,ind,a;
1234 int subdiv_z,sub_z,na_z,ash_z;
1235 int subdiv_y,sub_y,na_y,ash_y;
1236 int subdiv_x,sub_x,na_x,ash_x;
1238 /* cppcheck-suppress unassignedVariable */
1239 float bb_work_array[NNBSBB_B+3],*bb_work_align;
1241 bb_work_align = (float *)(((size_t)(bb_work_array+3)) & (~((size_t)15)));
1245 fprintf(debug,"cell0 %d sorting columns %d - %d, atoms %d - %d\n",
1246 grid->cell0,cxy_start,cxy_end,a0,a1);
1249 subdiv_x = grid->na_c;
1250 subdiv_y = GPU_NSUBCELL_X*subdiv_x;
1251 subdiv_z = GPU_NSUBCELL_Y*subdiv_y;
1253 /* Sort the atoms within each x,y column in 3 dimensions */
1254 for(cxy=cxy_start; cxy<cxy_end; cxy++)
1257 cy = cxy - cx*grid->ncy;
1259 na = grid->cxy_na[cxy];
1260 ncz = grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy];
1261 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1263 /* Sort the atoms within each x,y column on z coordinate */
1264 sort_atoms(ZZ,FALSE,
1267 ncz*grid->na_sc*SORT_GRID_OVERSIZE/nbs->box[ZZ][ZZ],
1268 ncz*grid->na_sc*SGSF,sort_work);
1270 /* This loop goes over the supercells and subcells along z at once */
1271 for(sub_z=0; sub_z<ncz*GPU_NSUBCELL_Z; sub_z++)
1273 ash_z = ash + sub_z*subdiv_z;
1274 na_z = min(subdiv_z,na-(ash_z-ash));
1276 /* We have already sorted on z */
1278 if (sub_z % GPU_NSUBCELL_Z == 0)
1280 cz = sub_z/GPU_NSUBCELL_Z;
1281 c = grid->cxy_ind[cxy] + cz ;
1283 /* The number of atoms in this supercell */
1284 na_c = min(grid->na_sc,na-(ash_z-ash));
1286 grid->nsubc[c] = min(GPU_NSUBCELL,(na_c+grid->na_c-1)/grid->na_c);
1288 /* Store the z-boundaries of the super cell */
1289 grid->bbcz[c*NNBSBB_D ] = x[nbs->a[ash_z]][ZZ];
1290 grid->bbcz[c*NNBSBB_D+1] = x[nbs->a[ash_z+na_c-1]][ZZ];
1293 #if GPU_NSUBCELL_Y > 1
1294 /* Sort the atoms along y */
1295 sort_atoms(YY,(sub_z & 1),
1296 nbs->a+ash_z,na_z,x,
1297 grid->c0[YY]+cy*grid->sy,grid->inv_sy,
1298 subdiv_y*SGSF,sort_work);
1301 for(sub_y=0; sub_y<GPU_NSUBCELL_Y; sub_y++)
1303 ash_y = ash_z + sub_y*subdiv_y;
1304 na_y = min(subdiv_y,na-(ash_y-ash));
1306 #if GPU_NSUBCELL_X > 1
1307 /* Sort the atoms along x */
1308 sort_atoms(XX,((cz*GPU_NSUBCELL_Y + sub_y) & 1),
1309 nbs->a+ash_y,na_y,x,
1310 grid->c0[XX]+cx*grid->sx,grid->inv_sx,
1311 subdiv_x*SGSF,sort_work);
1314 for(sub_x=0; sub_x<GPU_NSUBCELL_X; sub_x++)
1316 ash_x = ash_y + sub_x*subdiv_x;
1317 na_x = min(subdiv_x,na-(ash_x-ash));
1319 fill_cell(nbs,grid,nbat,
1320 ash_x,ash_x+na_x,atinfo,x,
1321 grid->na_c*(cx*GPU_NSUBCELL_X+sub_x) + (dd_zone >> 2),
1322 grid->na_c*(cy*GPU_NSUBCELL_Y+sub_y) + (dd_zone & 3),
1329 /* Set the unused atom indices to -1 */
1330 for(ind=na; ind<ncz*grid->na_sc; ind++)
1332 nbs->a[ash+ind] = -1;
1337 /* Determine in which grid column atoms should go */
1338 static void calc_column_indices(nbnxn_grid_t *grid,
1340 rvec *x,const int *move,
1341 int thread,int nthread,
1348 /* We add one extra cell for particles which moved during DD */
1349 for(i=0; i<grid->ncx*grid->ncy+1; i++)
1354 n0 = a0 + (int)((thread+0)*(a1 - a0))/nthread;
1355 n1 = a0 + (int)((thread+1)*(a1 - a0))/nthread;
1356 for(i=n0; i<n1; i++)
1358 if (move == NULL || move[i] >= 0)
1360 /* We need to be careful with rounding,
1361 * particles might be a few bits outside the local box.
1362 * The int cast takes care of the lower bound,
1363 * we need to explicitly take care of the upper bound.
1365 cx = (int)((x[i][XX] - grid->c0[XX])*grid->inv_sx);
1366 if (cx == grid->ncx)
1370 cy = (int)((x[i][YY] - grid->c0[YY])*grid->inv_sy);
1371 if (cy == grid->ncy)
1375 /* For the moment cell contains only the, grid local,
1376 * x and y indices, not z.
1378 cell[i] = cx*grid->ncy + cy;
1380 #ifdef DEBUG_NBNXN_GRIDDING
1381 if (cell[i] < 0 || cell[i] >= grid->ncx*grid->ncy)
1384 "grid cell cx %d cy %d out of range (max %d %d)",
1385 cx,cy,grid->ncx,grid->ncy);
1391 /* Put this moved particle after the end of the grid,
1392 * so we can process it later without using conditionals.
1394 cell[i] = grid->ncx*grid->ncy;
1401 /* Determine in which grid cells the atoms should go */
1402 static void calc_cell_indices(const nbnxn_search_t nbs,
1409 nbnxn_atomdata_t *nbat)
1412 int cx,cy,cxy,ncz_max,ncz;
1414 int *cxy_na,cxy_na_i;
1416 nthread = gmx_omp_nthreads_get(emntPairsearch);
1418 #pragma omp parallel for num_threads(nthread) schedule(static)
1419 for(thread=0; thread<nthread; thread++)
1421 calc_column_indices(grid,a0,a1,x,move,thread,nthread,
1422 nbs->cell,nbs->work[thread].cxy_na);
1425 /* Make the cell index as a function of x and y */
1428 grid->cxy_ind[0] = 0;
1429 for(i=0; i<grid->ncx*grid->ncy+1; i++)
1431 /* We set ncz_max at the beginning of the loop iso at the end
1432 * to skip i=grid->ncx*grid->ncy which are moved particles
1433 * that do not need to be ordered on the grid.
1439 cxy_na_i = nbs->work[0].cxy_na[i];
1440 for(thread=1; thread<nthread; thread++)
1442 cxy_na_i += nbs->work[thread].cxy_na[i];
1444 ncz = (cxy_na_i + grid->na_sc - 1)/grid->na_sc;
1445 if (nbat->XFormat == nbatX8)
1447 /* Make the number of cell a multiple of 2 */
1448 ncz = (ncz + 1) & ~1;
1450 grid->cxy_ind[i+1] = grid->cxy_ind[i] + ncz;
1451 /* Clear cxy_na, so we can reuse the array below */
1452 grid->cxy_na[i] = 0;
1454 grid->nc = grid->cxy_ind[grid->ncx*grid->ncy] - grid->cxy_ind[0];
1456 nbat->natoms = (grid->cell0 + grid->nc)*grid->na_sc;
1460 fprintf(debug,"ns na_sc %d na_c %d super-cells: %d x %d y %d z %.1f maxz %d\n",
1461 grid->na_sc,grid->na_c,grid->nc,
1462 grid->ncx,grid->ncy,grid->nc/((double)(grid->ncx*grid->ncy)),
1467 for(cy=0; cy<grid->ncy; cy++)
1469 for(cx=0; cx<grid->ncx; cx++)
1471 fprintf(debug," %2d",grid->cxy_ind[i+1]-grid->cxy_ind[i]);
1474 fprintf(debug,"\n");
1479 /* Make sure the work array for sorting is large enough */
1480 if (ncz_max*grid->na_sc*SGSF > nbs->work[0].sort_work_nalloc)
1482 for(thread=0; thread<nbs->nthread_max; thread++)
1484 nbs->work[thread].sort_work_nalloc =
1485 over_alloc_large(ncz_max*grid->na_sc*SGSF);
1486 srenew(nbs->work[thread].sort_work,
1487 nbs->work[thread].sort_work_nalloc);
1491 /* Now we know the dimensions we can fill the grid.
1492 * This is the first, unsorted fill. We sort the columns after this.
1494 for(i=a0; i<a1; i++)
1496 /* At this point nbs->cell contains the local grid x,y indices */
1498 nbs->a[(grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc + grid->cxy_na[cxy]++] = i;
1501 /* Set the cell indices for the moved particles */
1502 n0 = grid->nc*grid->na_sc;
1503 n1 = grid->nc*grid->na_sc+grid->cxy_na[grid->ncx*grid->ncy];
1504 for(i=n0; i<n1; i++)
1506 nbs->cell[nbs->a[i]] = i;
1509 /* Sort the super-cell columns along z into the sub-cells. */
1510 #pragma omp parallel for num_threads(nbs->nthread_max) schedule(static)
1511 for(thread=0; thread<nbs->nthread_max; thread++)
1515 sort_columns_simple(nbs,dd_zone,grid,a0,a1,atinfo,x,nbat,
1516 ((thread+0)*grid->ncx*grid->ncy)/nthread,
1517 ((thread+1)*grid->ncx*grid->ncy)/nthread,
1518 nbs->work[thread].sort_work);
1522 sort_columns_supersub(nbs,dd_zone,grid,a0,a1,atinfo,x,nbat,
1523 ((thread+0)*grid->ncx*grid->ncy)/nthread,
1524 ((thread+1)*grid->ncx*grid->ncy)/nthread,
1525 nbs->work[thread].sort_work);
1529 #ifdef NBNXN_SEARCH_SSE
1530 if (grid->bSimple && nbat->XFormat == nbatX8)
1532 combine_bounding_box_pairs(grid,grid->bb);
1538 grid->nsubc_tot = 0;
1539 for(i=0; i<grid->nc; i++)
1541 grid->nsubc_tot += grid->nsubc[i];
1549 print_bbsizes_simple(debug,nbs,grid);
1553 fprintf(debug,"ns non-zero sub-cells: %d average atoms %.2f\n",
1554 grid->nsubc_tot,(a1-a0)/(double)grid->nsubc_tot);
1556 print_bbsizes_supersub(debug,nbs,grid);
1561 static void init_grid_flags(nbnxn_cellblock_flags *flags,
1562 const nbnxn_grid_t *grid)
1566 flags->ncb = (grid->nc + NBNXN_CELLBLOCK_SIZE - 1)/NBNXN_CELLBLOCK_SIZE;
1567 if (flags->ncb > flags->flag_nalloc)
1569 flags->flag_nalloc = over_alloc_large(flags->ncb);
1570 srenew(flags->flag,flags->flag_nalloc);
1572 for(cb=0; cb<flags->ncb; cb++)
1574 flags->flag[cb] = 0;
1580 /* Sets up a grid and puts the atoms on the grid.
1581 * This function only operates on one domain of the domain decompostion.
1582 * Note that without domain decomposition there is only one domain.
1584 void nbnxn_put_on_grid(nbnxn_search_t nbs,
1585 int ePBC,matrix box,
1587 rvec corner0,rvec corner1,
1592 int nmoved,int *move,
1594 nbnxn_atomdata_t *nbat)
1598 int nc_max_grid,nc_max;
1600 grid = &nbs->grid[dd_zone];
1602 nbs_cycle_start(&nbs->cc[enbsCCgrid]);
1604 grid->bSimple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
1606 grid->na_c = nbnxn_kernel_to_ci_size(nb_kernel_type);
1607 grid->na_cj = nbnxn_kernel_to_cj_size(nb_kernel_type);
1608 grid->na_sc = (grid->bSimple ? 1 : GPU_NSUBCELL)*grid->na_c;
1609 grid->na_c_2log = get_2log(grid->na_c);
1611 nbat->na_c = grid->na_c;
1620 (nbs->grid[dd_zone-1].cell0 + nbs->grid[dd_zone-1].nc)*
1621 nbs->grid[dd_zone-1].na_sc/grid->na_sc;
1629 copy_mat(box,nbs->box);
1631 if (atom_density >= 0)
1633 grid->atom_density = atom_density;
1637 grid->atom_density = grid_atom_density(n-nmoved,corner0,corner1);
1642 nbs->natoms_local = a1 - nmoved;
1643 /* We assume that nbnxn_put_on_grid is called first
1644 * for the local atoms (dd_zone=0).
1646 nbs->natoms_nonlocal = a1 - nmoved;
1650 nbs->natoms_nonlocal = max(nbs->natoms_nonlocal,a1);
1653 nc_max_grid = set_grid_size_xy(nbs,grid,n-nmoved,corner0,corner1,
1654 nbs->grid[0].atom_density,
1657 nc_max = grid->cell0 + nc_max_grid;
1659 if (a1 > nbs->cell_nalloc)
1661 nbs->cell_nalloc = over_alloc_large(a1);
1662 srenew(nbs->cell,nbs->cell_nalloc);
1665 /* To avoid conditionals we store the moved particles at the end of a,
1666 * make sure we have enough space.
1668 if (nc_max*grid->na_sc + nmoved > nbs->a_nalloc)
1670 nbs->a_nalloc = over_alloc_large(nc_max*grid->na_sc + nmoved);
1671 srenew(nbs->a,nbs->a_nalloc);
1674 if (nc_max*grid->na_sc > nbat->nalloc)
1676 nbnxn_atomdata_realloc(nbat,nc_max*grid->na_sc);
1679 calc_cell_indices(nbs,dd_zone,grid,a0,a1,atinfo,x,move,nbat);
1683 nbat->natoms_local = nbat->natoms;
1686 init_grid_flags(&grid->cellblock_flags,grid);
1688 nbs_cycle_stop(&nbs->cc[enbsCCgrid]);
1691 /* Calls nbnxn_put_on_grid for all non-local domains */
1692 void nbnxn_put_on_grid_nonlocal(nbnxn_search_t nbs,
1693 const gmx_domdec_zones_t *zones,
1697 nbnxn_atomdata_t *nbat)
1702 for(zone=1; zone<zones->n; zone++)
1704 for(d=0; d<DIM; d++)
1706 c0[d] = zones->size[zone].bb_x0[d];
1707 c1[d] = zones->size[zone].bb_x1[d];
1710 nbnxn_put_on_grid(nbs,nbs->ePBC,NULL,
1712 zones->cg_range[zone],
1713 zones->cg_range[zone+1],
1723 /* Add simple grid type information to the local super/sub grid */
1724 void nbnxn_grid_add_simple(nbnxn_search_t nbs,
1725 nbnxn_atomdata_t *nbat)
1731 grid = &nbs->grid[0];
1735 gmx_incons("nbnxn_grid_simple called with a simple grid");
1738 ncd = grid->na_sc/NBNXN_CPU_CLUSTER_I_SIZE;
1740 if (grid->nc*ncd > grid->nc_nalloc_simple)
1742 grid->nc_nalloc_simple = over_alloc_large(grid->nc*ncd);
1743 srenew(grid->bbcz_simple,grid->nc_nalloc_simple*NNBSBB_D);
1744 srenew(grid->bb_simple,grid->nc_nalloc_simple*NNBSBB_B);
1745 srenew(grid->flags_simple,grid->nc_nalloc_simple);
1748 sfree_aligned(grid->bbj);
1749 snew_aligned(grid->bbj,grid->nc_nalloc_simple/2,16);
1753 bbcz = grid->bbcz_simple;
1754 bb = grid->bb_simple;
1756 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntPairsearch)) schedule(static)
1757 for(sc=0; sc<grid->nc; sc++)
1761 for(c=0; c<ncd; c++)
1765 na = NBNXN_CPU_CLUSTER_I_SIZE;
1767 nbat->type[tx*NBNXN_CPU_CLUSTER_I_SIZE+na-1] == nbat->ntype-1)
1774 switch (nbat->XFormat)
1777 /* PACK_X4==NBNXN_CPU_CLUSTER_I_SIZE, so this is simple */
1778 calc_bounding_box_x_x4(na,nbat->x+tx*STRIDE_P4,
1782 /* PACK_X8>NBNXN_CPU_CLUSTER_I_SIZE, more complicated */
1783 calc_bounding_box_x_x8(na,nbat->x+X8_IND_A(tx*NBNXN_CPU_CLUSTER_I_SIZE),
1787 calc_bounding_box(na,nbat->xstride,
1788 nbat->x+tx*NBNXN_CPU_CLUSTER_I_SIZE*nbat->xstride,
1792 bbcz[tx*NNBSBB_D+0] = bb[tx*NNBSBB_B +ZZ];
1793 bbcz[tx*NNBSBB_D+1] = bb[tx*NNBSBB_B+NNBSBB_C+ZZ];
1795 /* No interaction optimization yet here */
1796 grid->flags_simple[tx] = NBNXN_CI_DO_LJ(0) | NBNXN_CI_DO_COUL(0);
1800 grid->flags_simple[tx] = 0;
1805 #ifdef NBNXN_SEARCH_SSE
1806 if (grid->bSimple && nbat->XFormat == nbatX8)
1808 combine_bounding_box_pairs(grid,grid->bb_simple);
1813 void nbnxn_get_ncells(nbnxn_search_t nbs,int *ncx,int *ncy)
1815 *ncx = nbs->grid[0].ncx;
1816 *ncy = nbs->grid[0].ncy;
1819 void nbnxn_get_atomorder(nbnxn_search_t nbs,int **a,int *n)
1821 const nbnxn_grid_t *grid;
1823 grid = &nbs->grid[0];
1825 /* Return the atom order for the home cell (index 0) */
1828 *n = grid->cxy_ind[grid->ncx*grid->ncy]*grid->na_sc;
1831 void nbnxn_set_atomorder(nbnxn_search_t nbs)
1834 int ao,cx,cy,cxy,cz,j;
1836 /* Set the atom order for the home cell (index 0) */
1837 grid = &nbs->grid[0];
1840 for(cx=0; cx<grid->ncx; cx++)
1842 for(cy=0; cy<grid->ncy; cy++)
1844 cxy = cx*grid->ncy + cy;
1845 j = grid->cxy_ind[cxy]*grid->na_sc;
1846 for(cz=0; cz<grid->cxy_na[cxy]; cz++)
1857 /* Determines the cell range along one dimension that
1858 * the bounding box b0 - b1 sees.
1860 static void get_cell_range(real b0,real b1,
1861 int nc,real c0,real s,real invs,
1862 real d2,real r2,int *cf,int *cl)
1864 *cf = max((int)((b0 - c0)*invs),0);
1866 while (*cf > 0 && d2 + sqr((b0 - c0) - (*cf-1+1)*s) < r2)
1871 *cl = min((int)((b1 - c0)*invs),nc-1);
1872 while (*cl < nc-1 && d2 + sqr((*cl+1)*s - (b1 - c0)) < r2)
1878 /* Reference code calculating the distance^2 between two bounding boxes */
1879 static float box_dist2(float bx0,float bx1,float by0,
1880 float by1,float bz0,float bz1,
1888 dl = bx0 - bb[BBU_X];
1889 dh = bb[BBL_X] - bx1;
1894 dl = by0 - bb[BBU_Y];
1895 dh = bb[BBL_Y] - by1;
1900 dl = bz0 - bb[BBU_Z];
1901 dh = bb[BBL_Z] - bz1;
1909 /* Plain C code calculating the distance^2 between two bounding boxes */
1910 static float subc_bb_dist2(int si,const float *bb_i_ci,
1911 int csj,const float *bb_j_all)
1913 const float *bb_i,*bb_j;
1917 bb_i = bb_i_ci + si*NNBSBB_B;
1918 bb_j = bb_j_all + csj*NNBSBB_B;
1922 dl = bb_i[BBL_X] - bb_j[BBU_X];
1923 dh = bb_j[BBL_X] - bb_i[BBU_X];
1928 dl = bb_i[BBL_Y] - bb_j[BBU_Y];
1929 dh = bb_j[BBL_Y] - bb_i[BBU_Y];
1934 dl = bb_i[BBL_Z] - bb_j[BBU_Z];
1935 dh = bb_j[BBL_Z] - bb_i[BBU_Z];
1943 #ifdef NBNXN_SEARCH_SSE
1945 /* SSE code for bb distance for bb format xyz0 */
1946 static float subc_bb_dist2_sse(int na_c,
1947 int si,const float *bb_i_ci,
1948 int csj,const float *bb_j_all)
1950 const float *bb_i,*bb_j;
1952 __m128 bb_i_SSE0,bb_i_SSE1;
1953 __m128 bb_j_SSE0,bb_j_SSE1;
1959 #ifndef GMX_X86_SSE4_1
1960 float d2_array[7],*d2_align;
1962 d2_align = (float *)(((size_t)(d2_array+3)) & (~((size_t)15)));
1967 bb_i = bb_i_ci + si*NNBSBB_B;
1968 bb_j = bb_j_all + csj*NNBSBB_B;
1970 bb_i_SSE0 = _mm_load_ps(bb_i);
1971 bb_i_SSE1 = _mm_load_ps(bb_i+NNBSBB_C);
1972 bb_j_SSE0 = _mm_load_ps(bb_j);
1973 bb_j_SSE1 = _mm_load_ps(bb_j+NNBSBB_C);
1975 dl_SSE = _mm_sub_ps(bb_i_SSE0,bb_j_SSE1);
1976 dh_SSE = _mm_sub_ps(bb_j_SSE0,bb_i_SSE1);
1978 dm_SSE = _mm_max_ps(dl_SSE,dh_SSE);
1979 dm0_SSE = _mm_max_ps(dm_SSE,_mm_setzero_ps());
1980 #ifndef GMX_X86_SSE4_1
1981 d2_SSE = _mm_mul_ps(dm0_SSE,dm0_SSE);
1983 _mm_store_ps(d2_align,d2_SSE);
1985 return d2_align[0] + d2_align[1] + d2_align[2];
1987 /* SSE4.1 dot product of components 0,1,2 */
1988 d2_SSE = _mm_dp_ps(dm0_SSE,dm0_SSE,0x71);
1990 _mm_store_ss(&d2,d2_SSE);
1996 /* Calculate bb bounding distances of bb_i[si,...,si+3] and store them in d2 */
1997 #define SUBC_BB_DIST2_SSE_XXXX_INNER(si,bb_i,d2) \
2001 __m128 dx_0,dy_0,dz_0; \
2002 __m128 dx_1,dy_1,dz_1; \
2005 __m128 m0x,m0y,m0z; \
2007 __m128 d2x,d2y,d2z; \
2010 shi = si*NNBSBB_D*DIM; \
2012 xi_l = _mm_load_ps(bb_i+shi+0*STRIDE_8BB); \
2013 yi_l = _mm_load_ps(bb_i+shi+1*STRIDE_8BB); \
2014 zi_l = _mm_load_ps(bb_i+shi+2*STRIDE_8BB); \
2015 xi_h = _mm_load_ps(bb_i+shi+3*STRIDE_8BB); \
2016 yi_h = _mm_load_ps(bb_i+shi+4*STRIDE_8BB); \
2017 zi_h = _mm_load_ps(bb_i+shi+5*STRIDE_8BB); \
2019 dx_0 = _mm_sub_ps(xi_l,xj_h); \
2020 dy_0 = _mm_sub_ps(yi_l,yj_h); \
2021 dz_0 = _mm_sub_ps(zi_l,zj_h); \
2023 dx_1 = _mm_sub_ps(xj_l,xi_h); \
2024 dy_1 = _mm_sub_ps(yj_l,yi_h); \
2025 dz_1 = _mm_sub_ps(zj_l,zi_h); \
2027 mx = _mm_max_ps(dx_0,dx_1); \
2028 my = _mm_max_ps(dy_0,dy_1); \
2029 mz = _mm_max_ps(dz_0,dz_1); \
2031 m0x = _mm_max_ps(mx,zero); \
2032 m0y = _mm_max_ps(my,zero); \
2033 m0z = _mm_max_ps(mz,zero); \
2035 d2x = _mm_mul_ps(m0x,m0x); \
2036 d2y = _mm_mul_ps(m0y,m0y); \
2037 d2z = _mm_mul_ps(m0z,m0z); \
2039 d2s = _mm_add_ps(d2x,d2y); \
2040 d2t = _mm_add_ps(d2s,d2z); \
2042 _mm_store_ps(d2+si,d2t); \
2045 /* SSE code for nsi bb distances for bb format xxxxyyyyzzzz */
2046 static void subc_bb_dist2_sse_xxxx(const float *bb_j,
2047 int nsi,const float *bb_i,
2050 __m128 xj_l,yj_l,zj_l;
2051 __m128 xj_h,yj_h,zj_h;
2052 __m128 xi_l,yi_l,zi_l;
2053 __m128 xi_h,yi_h,zi_h;
2057 zero = _mm_setzero_ps();
2059 xj_l = _mm_set1_ps(bb_j[0*STRIDE_8BB]);
2060 yj_l = _mm_set1_ps(bb_j[1*STRIDE_8BB]);
2061 zj_l = _mm_set1_ps(bb_j[2*STRIDE_8BB]);
2062 xj_h = _mm_set1_ps(bb_j[3*STRIDE_8BB]);
2063 yj_h = _mm_set1_ps(bb_j[4*STRIDE_8BB]);
2064 zj_h = _mm_set1_ps(bb_j[5*STRIDE_8BB]);
2066 /* Here we "loop" over si (0,STRIDE_8BB) from 0 to nsi with step STRIDE_8BB.
2067 * But as we know the number of iterations is 1 or 2, we unroll manually.
2069 SUBC_BB_DIST2_SSE_XXXX_INNER(0,bb_i,d2);
2070 if (STRIDE_8BB < nsi)
2072 SUBC_BB_DIST2_SSE_XXXX_INNER(STRIDE_8BB,bb_i,d2);
2076 #endif /* NBNXN_SEARCH_SSE */
2078 /* Plain C function which determines if any atom pair between two cells
2079 * is within distance sqrt(rl2).
2081 static gmx_bool subc_in_range_x(int na_c,
2082 int si,const real *x_i,
2083 int csj,int stride,const real *x_j,
2089 for(i=0; i<na_c; i++)
2091 i0 = (si*na_c + i)*DIM;
2092 for(j=0; j<na_c; j++)
2094 j0 = (csj*na_c + j)*stride;
2096 d2 = sqr(x_i[i0 ] - x_j[j0 ]) +
2097 sqr(x_i[i0+1] - x_j[j0+1]) +
2098 sqr(x_i[i0+2] - x_j[j0+2]);
2110 /* SSE function which determines if any atom pair between two cells,
2111 * both with 8 atoms, is within distance sqrt(rl2).
2113 static gmx_bool subc_in_range_sse8(int na_c,
2114 int si,const real *x_i,
2115 int csj,int stride,const real *x_j,
2118 #ifdef NBNXN_SEARCH_SSE_SINGLE
2119 __m128 ix_SSE0,iy_SSE0,iz_SSE0;
2120 __m128 ix_SSE1,iy_SSE1,iz_SSE1;
2127 rc2_SSE = _mm_set1_ps(rl2);
2129 na_c_sse = NBNXN_GPU_CLUSTER_SIZE/STRIDE_8BB;
2130 ix_SSE0 = _mm_load_ps(x_i+(si*na_c_sse*DIM+0)*STRIDE_8BB);
2131 iy_SSE0 = _mm_load_ps(x_i+(si*na_c_sse*DIM+1)*STRIDE_8BB);
2132 iz_SSE0 = _mm_load_ps(x_i+(si*na_c_sse*DIM+2)*STRIDE_8BB);
2133 ix_SSE1 = _mm_load_ps(x_i+(si*na_c_sse*DIM+3)*STRIDE_8BB);
2134 iy_SSE1 = _mm_load_ps(x_i+(si*na_c_sse*DIM+4)*STRIDE_8BB);
2135 iz_SSE1 = _mm_load_ps(x_i+(si*na_c_sse*DIM+5)*STRIDE_8BB);
2137 /* We loop from the outer to the inner particles to maximize
2138 * the chance that we find a pair in range quickly and return.
2144 __m128 jx0_SSE,jy0_SSE,jz0_SSE;
2145 __m128 jx1_SSE,jy1_SSE,jz1_SSE;
2147 __m128 dx_SSE0,dy_SSE0,dz_SSE0;
2148 __m128 dx_SSE1,dy_SSE1,dz_SSE1;
2149 __m128 dx_SSE2,dy_SSE2,dz_SSE2;
2150 __m128 dx_SSE3,dy_SSE3,dz_SSE3;
2161 __m128 wco_any_SSE01,wco_any_SSE23,wco_any_SSE;
2163 jx0_SSE = _mm_load1_ps(x_j+j0*stride+0);
2164 jy0_SSE = _mm_load1_ps(x_j+j0*stride+1);
2165 jz0_SSE = _mm_load1_ps(x_j+j0*stride+2);
2167 jx1_SSE = _mm_load1_ps(x_j+j1*stride+0);
2168 jy1_SSE = _mm_load1_ps(x_j+j1*stride+1);
2169 jz1_SSE = _mm_load1_ps(x_j+j1*stride+2);
2171 /* Calculate distance */
2172 dx_SSE0 = _mm_sub_ps(ix_SSE0,jx0_SSE);
2173 dy_SSE0 = _mm_sub_ps(iy_SSE0,jy0_SSE);
2174 dz_SSE0 = _mm_sub_ps(iz_SSE0,jz0_SSE);
2175 dx_SSE1 = _mm_sub_ps(ix_SSE1,jx0_SSE);
2176 dy_SSE1 = _mm_sub_ps(iy_SSE1,jy0_SSE);
2177 dz_SSE1 = _mm_sub_ps(iz_SSE1,jz0_SSE);
2178 dx_SSE2 = _mm_sub_ps(ix_SSE0,jx1_SSE);
2179 dy_SSE2 = _mm_sub_ps(iy_SSE0,jy1_SSE);
2180 dz_SSE2 = _mm_sub_ps(iz_SSE0,jz1_SSE);
2181 dx_SSE3 = _mm_sub_ps(ix_SSE1,jx1_SSE);
2182 dy_SSE3 = _mm_sub_ps(iy_SSE1,jy1_SSE);
2183 dz_SSE3 = _mm_sub_ps(iz_SSE1,jz1_SSE);
2185 /* rsq = dx*dx+dy*dy+dz*dz */
2186 rsq_SSE0 = gmx_mm_calc_rsq_ps(dx_SSE0,dy_SSE0,dz_SSE0);
2187 rsq_SSE1 = gmx_mm_calc_rsq_ps(dx_SSE1,dy_SSE1,dz_SSE1);
2188 rsq_SSE2 = gmx_mm_calc_rsq_ps(dx_SSE2,dy_SSE2,dz_SSE2);
2189 rsq_SSE3 = gmx_mm_calc_rsq_ps(dx_SSE3,dy_SSE3,dz_SSE3);
2191 wco_SSE0 = _mm_cmplt_ps(rsq_SSE0,rc2_SSE);
2192 wco_SSE1 = _mm_cmplt_ps(rsq_SSE1,rc2_SSE);
2193 wco_SSE2 = _mm_cmplt_ps(rsq_SSE2,rc2_SSE);
2194 wco_SSE3 = _mm_cmplt_ps(rsq_SSE3,rc2_SSE);
2196 wco_any_SSE01 = _mm_or_ps(wco_SSE0,wco_SSE1);
2197 wco_any_SSE23 = _mm_or_ps(wco_SSE2,wco_SSE3);
2198 wco_any_SSE = _mm_or_ps(wco_any_SSE01,wco_any_SSE23);
2200 if (_mm_movemask_ps(wco_any_SSE))
2212 gmx_incons("SSE function called without SSE support");
2218 /* Returns the j sub-cell for index cj_ind */
2219 static int nbl_cj(const nbnxn_pairlist_t *nbl,int cj_ind)
2221 return nbl->cj4[cj_ind>>2].cj[cj_ind & 3];
2224 /* Returns the i-interaction mask of the j sub-cell for index cj_ind */
2225 static unsigned nbl_imask0(const nbnxn_pairlist_t *nbl,int cj_ind)
2227 return nbl->cj4[cj_ind>>2].imei[0].imask;
2230 /* Ensures there is enough space for extra extra exclusion masks */
2231 static void check_excl_space(nbnxn_pairlist_t *nbl,int extra)
2233 if (nbl->nexcl+extra > nbl->excl_nalloc)
2235 nbl->excl_nalloc = over_alloc_small(nbl->nexcl+extra);
2236 nbnxn_realloc_void((void **)&nbl->excl,
2237 nbl->nexcl*sizeof(*nbl->excl),
2238 nbl->excl_nalloc*sizeof(*nbl->excl),
2239 nbl->alloc,nbl->free);
2243 /* Ensures there is enough space for ncell extra j-cells in the list */
2244 static void check_subcell_list_space_simple(nbnxn_pairlist_t *nbl,
2249 cj_max = nbl->ncj + ncell;
2251 if (cj_max > nbl->cj_nalloc)
2253 nbl->cj_nalloc = over_alloc_small(cj_max);
2254 nbnxn_realloc_void((void **)&nbl->cj,
2255 nbl->ncj*sizeof(*nbl->cj),
2256 nbl->cj_nalloc*sizeof(*nbl->cj),
2257 nbl->alloc,nbl->free);
2261 /* Ensures there is enough space for ncell extra j-subcells in the list */
2262 static void check_subcell_list_space_supersub(nbnxn_pairlist_t *nbl,
2265 int ncj4_max,j4,j,w,t;
2268 #define WARP_SIZE 32
2270 /* We can have maximally nsupercell*GPU_NSUBCELL sj lists */
2271 /* We can store 4 j-subcell - i-supercell pairs in one struct.
2272 * since we round down, we need one extra entry.
2274 ncj4_max = ((nbl->work->cj_ind + nsupercell*GPU_NSUBCELL + 4-1) >> 2);
2276 if (ncj4_max > nbl->cj4_nalloc)
2278 nbl->cj4_nalloc = over_alloc_small(ncj4_max);
2279 nbnxn_realloc_void((void **)&nbl->cj4,
2280 nbl->work->cj4_init*sizeof(*nbl->cj4),
2281 nbl->cj4_nalloc*sizeof(*nbl->cj4),
2282 nbl->alloc,nbl->free);
2285 if (ncj4_max > nbl->work->cj4_init)
2287 for(j4=nbl->work->cj4_init; j4<ncj4_max; j4++)
2289 /* No i-subcells and no excl's in the list initially */
2290 for(w=0; w<NWARP; w++)
2292 nbl->cj4[j4].imei[w].imask = 0U;
2293 nbl->cj4[j4].imei[w].excl_ind = 0;
2297 nbl->work->cj4_init = ncj4_max;
2301 /* Set all excl masks for one GPU warp no exclusions */
2302 static void set_no_excls(nbnxn_excl_t *excl)
2306 for(t=0; t<WARP_SIZE; t++)
2308 /* Turn all interaction bits on */
2309 excl->pair[t] = NBNXN_INT_MASK_ALL;
2313 /* Initializes a single nbnxn_pairlist_t data structure */
2314 static void nbnxn_init_pairlist(nbnxn_pairlist_t *nbl,
2316 nbnxn_alloc_t *alloc,
2321 nbl->alloc = nbnxn_alloc_aligned;
2329 nbl->free = nbnxn_free_aligned;
2336 nbl->bSimple = bSimple;
2347 /* We need one element extra in sj, so alloc initially with 1 */
2348 nbl->cj4_nalloc = 0;
2355 nbl->excl_nalloc = 0;
2357 check_excl_space(nbl,1);
2359 set_no_excls(&nbl->excl[0]);
2364 snew_aligned(nbl->work->bb_ci,GPU_NSUBCELL/STRIDE_8BB*NNBSBB_XXXX,16);
2366 snew_aligned(nbl->work->bb_ci,GPU_NSUBCELL*NNBSBB_B,16);
2368 snew_aligned(nbl->work->x_ci,NBNXN_NA_SC_MAX*DIM,16);
2369 #ifdef NBNXN_SEARCH_SSE
2370 snew_aligned(nbl->work->x_ci_x86_simd128,1,16);
2371 #ifdef GMX_X86_AVX_256
2372 snew_aligned(nbl->work->x_ci_x86_simd256,1,32);
2375 snew_aligned(nbl->work->d2,GPU_NSUBCELL,16);
2378 void nbnxn_init_pairlist_set(nbnxn_pairlist_set_t *nbl_list,
2379 gmx_bool bSimple, gmx_bool bCombined,
2380 nbnxn_alloc_t *alloc,
2385 nbl_list->bSimple = bSimple;
2386 nbl_list->bCombined = bCombined;
2388 nbl_list->nnbl = gmx_omp_nthreads_get(emntNonbonded);
2390 if (!nbl_list->bCombined &&
2391 nbl_list->nnbl > NBNXN_CELLBLOCK_MAX_THREADS)
2393 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.",
2394 nbl_list->nnbl,NBNXN_CELLBLOCK_MAX_THREADS,NBNXN_CELLBLOCK_MAX_THREADS);
2397 snew(nbl_list->nbl,nbl_list->nnbl);
2398 /* Execute in order to avoid memory interleaving between threads */
2399 #pragma omp parallel for num_threads(nbl_list->nnbl) schedule(static)
2400 for(i=0; i<nbl_list->nnbl; i++)
2402 /* Allocate the nblist data structure locally on each thread
2403 * to optimize memory access for NUMA architectures.
2405 snew(nbl_list->nbl[i],1);
2407 /* Only list 0 is used on the GPU, use normal allocation for i>0 */
2410 nbnxn_init_pairlist(nbl_list->nbl[i],nbl_list->bSimple,alloc,free);
2414 nbnxn_init_pairlist(nbl_list->nbl[i],nbl_list->bSimple,NULL,NULL);
2419 /* Print statistics of a pair list, used for debug output */
2420 static void print_nblist_statistics_simple(FILE *fp,const nbnxn_pairlist_t *nbl,
2421 const nbnxn_search_t nbs,real rl)
2423 const nbnxn_grid_t *grid;
2428 /* This code only produces correct statistics with domain decomposition */
2429 grid = &nbs->grid[0];
2431 fprintf(fp,"nbl nci %d ncj %d\n",
2433 fprintf(fp,"nbl na_sc %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
2434 nbl->na_sc,rl,nbl->ncj,nbl->ncj/(double)grid->nc,
2435 nbl->ncj/(double)grid->nc*grid->na_sc,
2436 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)));
2438 fprintf(fp,"nbl average j cell list length %.1f\n",
2439 0.25*nbl->ncj/(double)nbl->nci);
2441 for(s=0; s<SHIFTS; s++)
2446 for(i=0; i<nbl->nci; i++)
2448 cs[nbl->ci[i].shift & NBNXN_CI_SHIFT] +=
2449 nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start;
2451 j = nbl->ci[i].cj_ind_start;
2452 while (j < nbl->ci[i].cj_ind_end &&
2453 nbl->cj[j].excl != NBNXN_INT_MASK_ALL)
2459 fprintf(fp,"nbl cell pairs, total: %d excl: %d %.1f%%\n",
2460 nbl->ncj,npexcl,100*npexcl/(double)nbl->ncj);
2461 for(s=0; s<SHIFTS; s++)
2465 fprintf(fp,"nbl shift %2d ncj %3d\n",s,cs[s]);
2470 /* Print statistics of a pair lists, used for debug output */
2471 static void print_nblist_statistics_supersub(FILE *fp,const nbnxn_pairlist_t *nbl,
2472 const nbnxn_search_t nbs,real rl)
2474 const nbnxn_grid_t *grid;
2476 int c[GPU_NSUBCELL+1];
2478 /* This code only produces correct statistics with domain decomposition */
2479 grid = &nbs->grid[0];
2481 fprintf(fp,"nbl nsci %d ncj4 %d nsi %d excl4 %d\n",
2482 nbl->nsci,nbl->ncj4,nbl->nci_tot,nbl->nexcl);
2483 fprintf(fp,"nbl na_c %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
2484 nbl->na_ci,rl,nbl->nci_tot,nbl->nci_tot/(double)grid->nsubc_tot,
2485 nbl->nci_tot/(double)grid->nsubc_tot*grid->na_c,
2486 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)));
2488 fprintf(fp,"nbl average j super cell list length %.1f\n",
2489 0.25*nbl->ncj4/(double)nbl->nsci);
2490 fprintf(fp,"nbl average i sub cell list length %.1f\n",
2491 nbl->nci_tot/(0.25*nbl->ncj4));
2493 for(si=0; si<=GPU_NSUBCELL; si++)
2497 for(i=0; i<nbl->nsci; i++)
2499 for(j4=nbl->sci[i].cj4_ind_start; j4<nbl->sci[i].cj4_ind_end; j4++)
2504 for(si=0; si<GPU_NSUBCELL; si++)
2506 if (nbl->cj4[j4].imei[0].imask & (1U << (j*GPU_NSUBCELL + si)))
2515 for(b=0; b<=GPU_NSUBCELL; b++)
2517 fprintf(fp,"nbl j-list #i-subcell %d %7d %4.1f\n",
2518 b,c[b],100.0*c[b]/(double)(nbl->ncj4*NBNXN_GPU_JGROUP_SIZE));
2522 /* Print the full pair list, used for debug output */
2523 static void print_supersub_nsp(const char *fn,
2524 const nbnxn_pairlist_t *nbl,
2531 sprintf(buf,"%s_%s.xvg",fn,NONLOCAL_I(iloc) ? "nl" : "l");
2532 fp = ffopen(buf,"w");
2534 for(i=0; i<nbl->nci; i++)
2537 for(j4=nbl->sci[i].cj4_ind_start; j4<nbl->sci[i].cj4_ind_end; j4++)
2539 for(p=0; p<NBNXN_GPU_JGROUP_SIZE*GPU_NSUBCELL; p++)
2541 nsp += (nbl->cj4[j4].imei[0].imask >> p) & 1;
2544 fprintf(fp,"%4d %3d %3d\n",
2547 nbl->sci[i].cj4_ind_end-nbl->sci[i].cj4_ind_start);
2553 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp */
2554 static void low_get_nbl_exclusions(nbnxn_pairlist_t *nbl,int cj4,
2555 int warp,nbnxn_excl_t **excl)
2557 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
2559 /* No exclusions set, make a new list entry */
2560 nbl->cj4[cj4].imei[warp].excl_ind = nbl->nexcl;
2562 *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
2563 set_no_excls(*excl);
2567 /* We already have some exclusions, new ones can be added to the list */
2568 *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
2572 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp,
2573 * allocates extra memory, if necessary.
2575 static void get_nbl_exclusions_1(nbnxn_pairlist_t *nbl,int cj4,
2576 int warp,nbnxn_excl_t **excl)
2578 if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
2580 /* We need to make a new list entry, check if we have space */
2581 check_excl_space(nbl,1);
2583 low_get_nbl_exclusions(nbl,cj4,warp,excl);
2586 /* Returns pointers to the exclusion mask for cj4-unit cj4 for both warps,
2587 * allocates extra memory, if necessary.
2589 static void get_nbl_exclusions_2(nbnxn_pairlist_t *nbl,int cj4,
2590 nbnxn_excl_t **excl_w0,
2591 nbnxn_excl_t **excl_w1)
2593 /* Check for space we might need */
2594 check_excl_space(nbl,2);
2596 low_get_nbl_exclusions(nbl,cj4,0,excl_w0);
2597 low_get_nbl_exclusions(nbl,cj4,1,excl_w1);
2600 /* Sets the self exclusions i=j and pair exclusions i>j */
2601 static void set_self_and_newton_excls_supersub(nbnxn_pairlist_t *nbl,
2602 int cj4_ind,int sj_offset,
2605 nbnxn_excl_t *excl[2];
2608 /* Here we only set the set self and double pair exclusions */
2610 get_nbl_exclusions_2(nbl,cj4_ind,&excl[0],&excl[1]);
2612 /* Only minor < major bits set */
2613 for(ej=0; ej<nbl->na_ci; ej++)
2616 for(ei=ej; ei<nbl->na_ci; ei++)
2618 excl[w]->pair[(ej&(4-1))*nbl->na_ci+ei] &=
2619 ~(1U << (sj_offset*GPU_NSUBCELL+si));
2624 /* Returns a diagonal or off-diagonal interaction mask for plain C lists */
2625 static unsigned int get_imask(gmx_bool rdiag,int ci,int cj)
2627 return (rdiag && ci == cj ? NBNXN_INT_MASK_DIAG : NBNXN_INT_MASK_ALL);
2630 #ifdef NBNXN_SEARCH_SSE
2631 /* Returns a diagonal or off-diagonal interaction mask for SIMD128 lists */
2632 static unsigned int get_imask_x86_simd128(gmx_bool rdiag,int ci,int cj)
2634 #ifndef GMX_DOUBLE /* cj-size = 4 */
2635 return (rdiag && ci == cj ? NBNXN_INT_MASK_DIAG : NBNXN_INT_MASK_ALL);
2636 #else /* cj-size = 2 */
2637 return (rdiag && ci*2 == cj ? NBNXN_INT_MASK_DIAG_J2_0 :
2638 (rdiag && ci*2+1 == cj ? NBNXN_INT_MASK_DIAG_J2_1 :
2639 NBNXN_INT_MASK_ALL));
2643 #ifdef GMX_X86_AVX_256
2644 /* Returns a diagonal or off-diagonal interaction mask for SIMD256 lists */
2645 static unsigned int get_imask_x86_simd256(gmx_bool rdiag,int ci,int cj)
2647 #ifndef GMX_DOUBLE /* cj-size = 8 */
2648 return (rdiag && ci == cj*2 ? NBNXN_INT_MASK_DIAG_J8_0 :
2649 (rdiag && ci == cj*2+1 ? NBNXN_INT_MASK_DIAG_J8_1 :
2650 NBNXN_INT_MASK_ALL));
2651 #else /* cj-size = 2 */
2652 return (rdiag && ci == cj ? NBNXN_INT_MASK_DIAG : NBNXN_INT_MASK_ALL);
2656 #endif /* NBNXN_SEARCH_SSE */
2658 /* Plain C code for making a pair list of cell ci vs cell cjf-cjl.
2659 * Checks bounding box distances and possibly atom pair distances.
2661 static void make_cluster_list_simple(const nbnxn_grid_t *gridj,
2662 nbnxn_pairlist_t *nbl,
2663 int ci,int cjf,int cjl,
2664 gmx_bool remove_sub_diag,
2666 real rl2,float rbb2,
2669 const nbnxn_list_work_t *work;
2676 int cjf_gl,cjl_gl,cj;
2680 bb_ci = nbl->work->bb_ci;
2681 x_ci = nbl->work->x_ci;
2684 while (!InRange && cjf <= cjl)
2686 d2 = subc_bb_dist2(0,bb_ci,cjf,gridj->bb);
2689 /* Check if the distance is within the distance where
2690 * we use only the bounding box distance rbb,
2691 * or within the cut-off and there is at least one atom pair
2692 * within the cut-off.
2702 cjf_gl = gridj->cell0 + cjf;
2703 for(i=0; i<NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
2705 for(j=0; j<NBNXN_CPU_CLUSTER_I_SIZE; j++)
2707 InRange = InRange ||
2708 (sqr(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
2709 sqr(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
2710 sqr(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rl2);
2713 *ndistc += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
2726 while (!InRange && cjl > cjf)
2728 d2 = subc_bb_dist2(0,bb_ci,cjl,gridj->bb);
2731 /* Check if the distance is within the distance where
2732 * we use only the bounding box distance rbb,
2733 * or within the cut-off and there is at least one atom pair
2734 * within the cut-off.
2744 cjl_gl = gridj->cell0 + cjl;
2745 for(i=0; i<NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
2747 for(j=0; j<NBNXN_CPU_CLUSTER_I_SIZE; j++)
2749 InRange = InRange ||
2750 (sqr(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
2751 sqr(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
2752 sqr(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rl2);
2755 *ndistc += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
2765 for(cj=cjf; cj<=cjl; cj++)
2767 /* Store cj and the interaction mask */
2768 nbl->cj[nbl->ncj].cj = gridj->cell0 + cj;
2769 nbl->cj[nbl->ncj].excl = get_imask(remove_sub_diag,ci,cj);
2772 /* Increase the closing index in i super-cell list */
2773 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
2777 #ifdef NBNXN_SEARCH_SSE
2778 /* Include make_cluster_list_x86_simd128/256 */
2779 #define GMX_MM128_HERE
2780 #include "gmx_x86_simd_macros.h"
2781 #define STRIDE_S PACK_X4
2782 #include "nbnxn_search_x86_simd.h"
2784 #undef GMX_MM128_HERE
2785 #ifdef GMX_X86_AVX_256
2786 /* Include make_cluster_list_x86_simd128/256 */
2787 #define GMX_MM256_HERE
2788 #include "gmx_x86_simd_macros.h"
2789 #define STRIDE_S GMX_X86_SIMD_WIDTH_HERE
2790 #include "nbnxn_search_x86_simd.h"
2792 #undef GMX_MM256_HERE
2796 /* Plain C or SSE code for making a pair list of super-cell sci vs scj.
2797 * Checks bounding box distances and possibly atom pair distances.
2799 static void make_cluster_list_supersub(const nbnxn_search_t nbs,
2800 const nbnxn_grid_t *gridi,
2801 const nbnxn_grid_t *gridj,
2802 nbnxn_pairlist_t *nbl,
2804 gmx_bool sci_equals_scj,
2805 int stride,const real *x,
2806 real rl2,float rbb2,
2811 int cjo,ci1,ci,cj,cj_gl;
2812 int cj4_ind,cj_offset;
2819 #define PRUNE_LIST_CPU_ONE
2820 #ifdef PRUNE_LIST_CPU_ONE
2824 d2l = nbl->work->d2;
2826 bb_ci = nbl->work->bb_ci;
2827 x_ci = nbl->work->x_ci;
2831 for(cjo=0; cjo<gridj->nsubc[scj]; cjo++)
2833 cj4_ind = (nbl->work->cj_ind >> 2);
2834 cj_offset = nbl->work->cj_ind - cj4_ind*NBNXN_GPU_JGROUP_SIZE;
2835 cj4 = &nbl->cj4[cj4_ind];
2837 cj = scj*GPU_NSUBCELL + cjo;
2839 cj_gl = gridj->cell0*GPU_NSUBCELL + cj;
2841 /* Initialize this j-subcell i-subcell list */
2842 cj4->cj[cj_offset] = cj_gl;
2851 ci1 = gridi->nsubc[sci];
2855 /* Determine all ci1 bb distances in one call with SSE */
2856 subc_bb_dist2_sse_xxxx(gridj->bb+(cj>>STRIDE_8BB_2LOG)*NNBSBB_XXXX+(cj & (STRIDE_8BB-1)),
2862 /* We use a fixed upper-bound instead of ci1 to help optimization */
2863 for(ci=0; ci<GPU_NSUBCELL; ci++)
2870 #ifndef NBNXN_BBXXXX
2871 /* Determine the bb distance between ci and cj */
2872 d2l[ci] = subc_bb_dist2(ci,bb_ci,cj,gridj->bb);
2877 #ifdef PRUNE_LIST_CPU_ALL
2878 /* Check if the distance is within the distance where
2879 * we use only the bounding box distance rbb,
2880 * or within the cut-off and there is at least one atom pair
2881 * within the cut-off. This check is very costly.
2883 *ndistc += na_c*na_c;
2885 (d2 < rl2 && subc_in_range_x(na_c,ci,x_ci,cj_gl,stride,x,rl2)))
2887 /* Check if the distance between the two bounding boxes
2888 * in within the pair-list cut-off.
2893 /* Flag this i-subcell to be taken into account */
2894 imask |= (1U << (cj_offset*GPU_NSUBCELL+ci));
2896 #ifdef PRUNE_LIST_CPU_ONE
2904 #ifdef PRUNE_LIST_CPU_ONE
2905 /* If we only found 1 pair, check if any atoms are actually
2906 * within the cut-off, so we could get rid of it.
2908 if (npair == 1 && d2l[ci_last] >= rbb2)
2910 /* Avoid using function pointers here, as it's slower */
2912 #ifdef NBNXN_8BB_SSE
2917 (na_c,ci_last,x_ci,cj_gl,stride,x,rl2))
2919 imask &= ~(1U << (cj_offset*GPU_NSUBCELL+ci_last));
2927 /* We have a useful sj entry, close it now */
2929 /* Set the exclucions for the ci== sj entry.
2930 * Here we don't bother to check if this entry is actually flagged,
2931 * as it will nearly always be in the list.
2935 set_self_and_newton_excls_supersub(nbl,cj4_ind,cj_offset,cjo);
2938 /* Copy the cluster interaction mask to the list */
2939 for(w=0; w<NWARP; w++)
2941 cj4->imei[w].imask |= imask;
2944 nbl->work->cj_ind++;
2946 /* Keep the count */
2947 nbl->nci_tot += npair;
2949 /* Increase the closing index in i super-cell list */
2950 nbl->sci[nbl->nsci].cj4_ind_end = ((nbl->work->cj_ind+4-1)>>2);
2955 /* Set all atom-pair exclusions from the topology stored in excl
2956 * as masks in the pair-list for simple list i-entry nbl_ci
2958 static void set_ci_top_excls(const nbnxn_search_t nbs,
2959 nbnxn_pairlist_t *nbl,
2960 gmx_bool diagRemoved,
2963 const nbnxn_ci_t *nbl_ci,
2964 const t_blocka *excl)
2968 int cj_ind_first,cj_ind_last;
2969 int cj_first,cj_last;
2971 int i,ai,aj,si,eind,ge,se;
2972 int found,cj_ind_0,cj_ind_1,cj_ind_m;
2976 nbnxn_excl_t *nbl_excl;
2977 int inner_i,inner_e;
2981 if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
2989 cj_ind_first = nbl_ci->cj_ind_start;
2990 cj_ind_last = nbl->ncj - 1;
2992 cj_first = nbl->cj[cj_ind_first].cj;
2993 cj_last = nbl->cj[cj_ind_last].cj;
2995 /* Determine how many contiguous j-cells we have starting
2996 * from the first i-cell. This number can be used to directly
2997 * calculate j-cell indices for excluded atoms.
3000 if (na_ci_2log == na_cj_2log)
3002 while (cj_ind_first + ndirect <= cj_ind_last &&
3003 nbl->cj[cj_ind_first+ndirect].cj == ci + ndirect)
3008 #ifdef NBNXN_SEARCH_SSE
3011 while (cj_ind_first + ndirect <= cj_ind_last &&
3012 nbl->cj[cj_ind_first+ndirect].cj == ci_to_cj(na_cj_2log,ci) + ndirect)
3019 /* Loop over the atoms in the i super-cell */
3020 for(i=0; i<nbl->na_sc; i++)
3022 ai = nbs->a[ci*nbl->na_sc+i];
3025 si = (i>>na_ci_2log);
3027 /* Loop over the topology-based exclusions for this i-atom */
3028 for(eind=excl->index[ai]; eind<excl->index[ai+1]; eind++)
3034 /* The self exclusion are already set, save some time */
3040 /* Without shifts we only calculate interactions j>i
3041 * for one-way pair-lists.
3043 if (diagRemoved && ge <= ci*nbl->na_sc + i)
3048 se = (ge >> na_cj_2log);
3050 /* Could the cluster se be in our list? */
3051 if (se >= cj_first && se <= cj_last)
3053 if (se < cj_first + ndirect)
3055 /* We can calculate cj_ind directly from se */
3056 found = cj_ind_first + se - cj_first;
3060 /* Search for se using bisection */
3062 cj_ind_0 = cj_ind_first + ndirect;
3063 cj_ind_1 = cj_ind_last + 1;
3064 while (found == -1 && cj_ind_0 < cj_ind_1)
3066 cj_ind_m = (cj_ind_0 + cj_ind_1)>>1;
3068 cj_m = nbl->cj[cj_ind_m].cj;
3076 cj_ind_1 = cj_ind_m;
3080 cj_ind_0 = cj_ind_m + 1;
3087 inner_i = i - (si << na_ci_2log);
3088 inner_e = ge - (se << na_cj_2log);
3090 nbl->cj[found].excl &= ~(1U<<((inner_i<<na_cj_2log) + inner_e));
3098 /* Set all atom-pair exclusions from the topology stored in excl
3099 * as masks in the pair-list for i-super-cell entry nbl_sci
3101 static void set_sci_top_excls(const nbnxn_search_t nbs,
3102 nbnxn_pairlist_t *nbl,
3103 gmx_bool diagRemoved,
3105 const nbnxn_sci_t *nbl_sci,
3106 const t_blocka *excl)
3111 int cj_ind_first,cj_ind_last;
3112 int cj_first,cj_last;
3114 int i,ai,aj,si,eind,ge,se;
3115 int found,cj_ind_0,cj_ind_1,cj_ind_m;
3119 nbnxn_excl_t *nbl_excl;
3120 int inner_i,inner_e,w;
3126 if (nbl_sci->cj4_ind_end == nbl_sci->cj4_ind_start)
3134 cj_ind_first = nbl_sci->cj4_ind_start*NBNXN_GPU_JGROUP_SIZE;
3135 cj_ind_last = nbl->work->cj_ind - 1;
3137 cj_first = nbl->cj4[nbl_sci->cj4_ind_start].cj[0];
3138 cj_last = nbl_cj(nbl,cj_ind_last);
3140 /* Determine how many contiguous j-clusters we have starting
3141 * from the first i-cluster. This number can be used to directly
3142 * calculate j-cluster indices for excluded atoms.
3145 while (cj_ind_first + ndirect <= cj_ind_last &&
3146 nbl_cj(nbl,cj_ind_first+ndirect) == sci*GPU_NSUBCELL + ndirect)
3151 /* Loop over the atoms in the i super-cell */
3152 for(i=0; i<nbl->na_sc; i++)
3154 ai = nbs->a[sci*nbl->na_sc+i];
3157 si = (i>>na_c_2log);
3159 /* Loop over the topology-based exclusions for this i-atom */
3160 for(eind=excl->index[ai]; eind<excl->index[ai+1]; eind++)
3166 /* The self exclusion are already set, save some time */
3172 /* Without shifts we only calculate interactions j>i
3173 * for one-way pair-lists.
3175 if (diagRemoved && ge <= sci*nbl->na_sc + i)
3181 /* Could the cluster se be in our list? */
3182 if (se >= cj_first && se <= cj_last)
3184 if (se < cj_first + ndirect)
3186 /* We can calculate cj_ind directly from se */
3187 found = cj_ind_first + se - cj_first;
3191 /* Search for se using bisection */
3193 cj_ind_0 = cj_ind_first + ndirect;
3194 cj_ind_1 = cj_ind_last + 1;
3195 while (found == -1 && cj_ind_0 < cj_ind_1)
3197 cj_ind_m = (cj_ind_0 + cj_ind_1)>>1;
3199 cj_m = nbl_cj(nbl,cj_ind_m);
3207 cj_ind_1 = cj_ind_m;
3211 cj_ind_0 = cj_ind_m + 1;
3218 inner_i = i - si*na_c;
3219 inner_e = ge - se*na_c;
3221 /* Macro for getting the index of atom a within a cluster */
3222 #define AMODI(a) ((a) & (NBNXN_CPU_CLUSTER_I_SIZE - 1))
3223 /* Macro for converting an atom number to a cluster number */
3224 #define A2CI(a) ((a) >> NBNXN_CPU_CLUSTER_I_SIZE_2LOG)
3226 if (nbl_imask0(nbl,found) & (1U << (AMODI(found)*GPU_NSUBCELL + si)))
3230 get_nbl_exclusions_1(nbl,A2CI(found),w,&nbl_excl);
3232 nbl_excl->pair[AMODI(inner_e)*nbl->na_ci+inner_i] &=
3233 ~(1U << (AMODI(found)*GPU_NSUBCELL + si));
3245 /* Reallocate the simple ci list for at least n entries */
3246 static void nb_realloc_ci(nbnxn_pairlist_t *nbl,int n)
3248 nbl->ci_nalloc = over_alloc_small(n);
3249 nbnxn_realloc_void((void **)&nbl->ci,
3250 nbl->nci*sizeof(*nbl->ci),
3251 nbl->ci_nalloc*sizeof(*nbl->ci),
3252 nbl->alloc,nbl->free);
3255 /* Reallocate the super-cell sci list for at least n entries */
3256 static void nb_realloc_sci(nbnxn_pairlist_t *nbl,int n)
3258 nbl->sci_nalloc = over_alloc_small(n);
3259 nbnxn_realloc_void((void **)&nbl->sci,
3260 nbl->nsci*sizeof(*nbl->sci),
3261 nbl->sci_nalloc*sizeof(*nbl->sci),
3262 nbl->alloc,nbl->free);
3265 /* Make a new ci entry at index nbl->nci */
3266 static void new_ci_entry(nbnxn_pairlist_t *nbl,int ci,int shift,int flags,
3267 nbnxn_list_work_t *work)
3269 if (nbl->nci + 1 > nbl->ci_nalloc)
3271 nb_realloc_ci(nbl,nbl->nci+1);
3273 nbl->ci[nbl->nci].ci = ci;
3274 nbl->ci[nbl->nci].shift = shift;
3275 /* Store the interaction flags along with the shift */
3276 nbl->ci[nbl->nci].shift |= flags;
3277 nbl->ci[nbl->nci].cj_ind_start = nbl->ncj;
3278 nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
3281 /* Make a new sci entry at index nbl->nsci */
3282 static void new_sci_entry(nbnxn_pairlist_t *nbl,int sci,int shift,int flags,
3283 nbnxn_list_work_t *work)
3285 if (nbl->nsci + 1 > nbl->sci_nalloc)
3287 nb_realloc_sci(nbl,nbl->nsci+1);
3289 nbl->sci[nbl->nsci].sci = sci;
3290 nbl->sci[nbl->nsci].shift = shift;
3291 nbl->sci[nbl->nsci].cj4_ind_start = nbl->ncj4;
3292 nbl->sci[nbl->nsci].cj4_ind_end = nbl->ncj4;
3295 /* Sort the simple j-list cj on exclusions.
3296 * Entries with exclusions will all be sorted to the beginning of the list.
3298 static void sort_cj_excl(nbnxn_cj_t *cj,int ncj,
3299 nbnxn_list_work_t *work)
3303 if (ncj > work->cj_nalloc)
3305 work->cj_nalloc = over_alloc_large(ncj);
3306 srenew(work->cj,work->cj_nalloc);
3309 /* Make a list of the j-cells involving exclusions */
3311 for(j=0; j<ncj; j++)
3313 if (cj[j].excl != NBNXN_INT_MASK_ALL)
3315 work->cj[jnew++] = cj[j];
3318 /* Check if there are exclusions at all or not just the first entry */
3319 if (!((jnew == 0) ||
3320 (jnew == 1 && cj[0].excl != NBNXN_INT_MASK_ALL)))
3322 for(j=0; j<ncj; j++)
3324 if (cj[j].excl == NBNXN_INT_MASK_ALL)
3326 work->cj[jnew++] = cj[j];
3329 for(j=0; j<ncj; j++)
3331 cj[j] = work->cj[j];
3336 /* Close this simple list i entry */
3337 static void close_ci_entry_simple(nbnxn_pairlist_t *nbl)
3341 /* All content of the new ci entry have already been filled correctly,
3342 * we only need to increase the count here (for non empty lists).
3344 jlen = nbl->ci[nbl->nci].cj_ind_end - nbl->ci[nbl->nci].cj_ind_start;
3347 sort_cj_excl(nbl->cj+nbl->ci[nbl->nci].cj_ind_start,jlen,nbl->work);
3349 if (nbl->ci[nbl->nci].shift & NBNXN_CI_HALF_LJ(0))
3351 nbl->work->ncj_hlj += jlen;
3353 else if (!(nbl->ci[nbl->nci].shift & NBNXN_CI_DO_COUL(0)))
3355 nbl->work->ncj_noq += jlen;
3362 /* Split sci entry for load balancing on the GPU.
3363 * As we only now the current count on our own thread,
3364 * we will need to estimate the current total amount of i-entries.
3365 * As the lists get concatenated later, this estimate depends
3366 * both on nthread and our own thread index thread.
3368 static void split_sci_entry(nbnxn_pairlist_t *nbl,
3369 int nsp_max_av,gmx_bool progBal,int nc_bal,
3370 int thread,int nthread)
3374 int cj4_start,cj4_end,j4len,cj4;
3376 int nsp,nsp_sci,nsp_cj4,nsp_cj4_e,nsp_cj4_p;
3379 /* Estimate the total numbers of ci's of the nblist combined
3380 * over all threads using the target number of ci's.
3382 nsci_est = nc_bal*thread/nthread + nbl->nsci;
3385 /* The first ci blocks should be larger, to avoid overhead.
3386 * The last ci blocks should be smaller, to improve load balancing.
3389 nsp_max_av*nc_bal*3/(2*(nsci_est - 1 + nc_bal)));
3393 nsp_max = nsp_max_av;
3396 cj4_start = nbl->sci[nbl->nsci-1].cj4_ind_start;
3397 cj4_end = nbl->sci[nbl->nsci-1].cj4_ind_end;
3398 j4len = cj4_end - cj4_start;
3400 if (j4len > 1 && j4len*GPU_NSUBCELL*NBNXN_GPU_JGROUP_SIZE > nsp_max)
3402 /* Remove the last ci entry and process the cj4's again */
3411 while (cj4 < cj4_end)
3413 nsp_cj4_p = nsp_cj4;
3415 for(p=0; p<GPU_NSUBCELL*NBNXN_GPU_JGROUP_SIZE; p++)
3417 nsp_cj4 += (nbl->cj4[cj4].imei[0].imask >> p) & 1;
3421 if (nsp > nsp_max && nsp > nsp_cj4)
3423 nbl->sci[sci].cj4_ind_end = cj4;
3426 if (nbl->nsci+1 > nbl->sci_nalloc)
3428 nb_realloc_sci(nbl,nbl->nsci+1);
3430 nbl->sci[sci].sci = nbl->sci[nbl->nsci-1].sci;
3431 nbl->sci[sci].shift = nbl->sci[nbl->nsci-1].shift;
3432 nbl->sci[sci].cj4_ind_start = cj4;
3433 nsp_sci = nsp - nsp_cj4;
3434 nsp_cj4_e = nsp_cj4_p;
3441 /* Put the remaining cj4's in a new ci entry */
3442 nbl->sci[sci].cj4_ind_end = cj4_end;
3444 /* Possibly balance out the last two ci's
3445 * by moving the last cj4 of the second last ci.
3447 if (nsp_sci - nsp_cj4_e >= nsp + nsp_cj4_e)
3449 nbl->sci[sci-1].cj4_ind_end--;
3450 nbl->sci[sci].cj4_ind_start--;
3458 /* Clost this super/sub list i entry */
3459 static void close_ci_entry_supersub(nbnxn_pairlist_t *nbl,
3461 gmx_bool progBal,int nc_bal,
3462 int thread,int nthread)
3467 /* All content of the new ci entry have already been filled correctly,
3468 * we only need to increase the count here (for non empty lists).
3470 j4len = nbl->sci[nbl->nsci].cj4_ind_end - nbl->sci[nbl->nsci].cj4_ind_start;
3473 /* We can only have complete blocks of 4 j-entries in a list,
3474 * so round the count up before closing.
3476 nbl->ncj4 = ((nbl->work->cj_ind + 4-1) >> 2);
3477 nbl->work->cj_ind = nbl->ncj4*NBNXN_GPU_JGROUP_SIZE;
3483 split_sci_entry(nbl,nsp_max_av,progBal,nc_bal,thread,nthread);
3488 /* Syncs the working array before adding another grid pair to the list */
3489 static void sync_work(nbnxn_pairlist_t *nbl)
3493 nbl->work->cj_ind = nbl->ncj4*NBNXN_GPU_JGROUP_SIZE;
3494 nbl->work->cj4_init = nbl->ncj4;
3498 /* Clears an nbnxn_pairlist_t data structure */
3499 static void clear_pairlist(nbnxn_pairlist_t *nbl)
3508 nbl->work->ncj_noq = 0;
3509 nbl->work->ncj_hlj = 0;
3512 /* Sets a simple list i-cell bounding box, including PBC shift */
3513 static void set_icell_bb_simple(const float *bb,int ci,
3514 real shx,real shy,real shz,
3520 bb_ci[BBL_X] = bb[ia+BBL_X] + shx;
3521 bb_ci[BBL_Y] = bb[ia+BBL_Y] + shy;
3522 bb_ci[BBL_Z] = bb[ia+BBL_Z] + shz;
3523 bb_ci[BBU_X] = bb[ia+BBU_X] + shx;
3524 bb_ci[BBU_Y] = bb[ia+BBU_Y] + shy;
3525 bb_ci[BBU_Z] = bb[ia+BBU_Z] + shz;
3528 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
3529 static void set_icell_bb_supersub(const float *bb,int ci,
3530 real shx,real shy,real shz,
3536 ia = ci*(GPU_NSUBCELL>>STRIDE_8BB_2LOG)*NNBSBB_XXXX;
3537 for(m=0; m<(GPU_NSUBCELL>>STRIDE_8BB_2LOG)*NNBSBB_XXXX; m+=NNBSBB_XXXX)
3539 for(i=0; i<STRIDE_8BB; i++)
3541 bb_ci[m+0*STRIDE_8BB+i] = bb[ia+m+0*STRIDE_8BB+i] + shx;
3542 bb_ci[m+1*STRIDE_8BB+i] = bb[ia+m+1*STRIDE_8BB+i] + shy;
3543 bb_ci[m+2*STRIDE_8BB+i] = bb[ia+m+2*STRIDE_8BB+i] + shz;
3544 bb_ci[m+3*STRIDE_8BB+i] = bb[ia+m+3*STRIDE_8BB+i] + shx;
3545 bb_ci[m+4*STRIDE_8BB+i] = bb[ia+m+4*STRIDE_8BB+i] + shy;
3546 bb_ci[m+5*STRIDE_8BB+i] = bb[ia+m+5*STRIDE_8BB+i] + shz;
3550 ia = ci*GPU_NSUBCELL*NNBSBB_B;
3551 for(i=0; i<GPU_NSUBCELL*NNBSBB_B; i+=NNBSBB_B)
3553 bb_ci[BBL_X] = bb[ia+BBL_X] + shx;
3554 bb_ci[BBL_Y] = bb[ia+BBL_Y] + shy;
3555 bb_ci[BBL_Z] = bb[ia+BBL_Z] + shz;
3556 bb_ci[BBU_X] = bb[ia+BBU_X] + shx;
3557 bb_ci[BBU_Y] = bb[ia+BBU_Y] + shy;
3558 bb_ci[BBU_Z] = bb[ia+BBU_Z] + shz;
3563 /* Copies PBC shifted i-cell atom coordinates x,y,z to working array */
3564 static void icell_set_x_simple(int ci,
3565 real shx,real shy,real shz,
3567 int stride,const real *x,
3568 nbnxn_list_work_t *work)
3572 ia = ci*NBNXN_CPU_CLUSTER_I_SIZE;
3574 for(i=0; i<NBNXN_CPU_CLUSTER_I_SIZE; i++)
3576 work->x_ci[i*STRIDE_XYZ+XX] = x[(ia+i)*stride+XX] + shx;
3577 work->x_ci[i*STRIDE_XYZ+YY] = x[(ia+i)*stride+YY] + shy;
3578 work->x_ci[i*STRIDE_XYZ+ZZ] = x[(ia+i)*stride+ZZ] + shz;
3582 /* Copies PBC shifted super-cell atom coordinates x,y,z to working array */
3583 static void icell_set_x_supersub(int ci,
3584 real shx,real shy,real shz,
3586 int stride,const real *x,
3587 nbnxn_list_work_t *work)
3594 ia = ci*GPU_NSUBCELL*na_c;
3595 for(i=0; i<GPU_NSUBCELL*na_c; i++)
3597 x_ci[i*DIM + XX] = x[(ia+i)*stride + XX] + shx;
3598 x_ci[i*DIM + YY] = x[(ia+i)*stride + YY] + shy;
3599 x_ci[i*DIM + ZZ] = x[(ia+i)*stride + ZZ] + shz;
3603 #ifdef NBNXN_SEARCH_SSE
3604 /* Copies PBC shifted super-cell packed atom coordinates to working array */
3605 static void icell_set_x_supersub_sse8(int ci,
3606 real shx,real shy,real shz,
3608 int stride,const real *x,
3609 nbnxn_list_work_t *work)
3616 for(si=0; si<GPU_NSUBCELL; si++)
3618 for(i=0; i<na_c; i+=STRIDE_8BB)
3621 ia = ci*GPU_NSUBCELL*na_c + io;
3622 for(j=0; j<STRIDE_8BB; j++)
3624 x_ci[io*DIM + j + XX*STRIDE_8BB] = x[(ia+j)*stride+XX] + shx;
3625 x_ci[io*DIM + j + YY*STRIDE_8BB] = x[(ia+j)*stride+YY] + shy;
3626 x_ci[io*DIM + j + ZZ*STRIDE_8BB] = x[(ia+j)*stride+ZZ] + shz;
3633 static real nbnxn_rlist_inc_nonloc_fac = 0.6;
3635 /* Due to the cluster size the effective pair-list is longer than
3636 * that of a simple atom pair-list. This function gives the extra distance.
3638 real nbnxn_get_rlist_effective_inc(int cluster_size,real atom_density)
3640 return ((0.5 + nbnxn_rlist_inc_nonloc_fac)*sqr(((cluster_size) - 1.0)/(cluster_size))*pow((cluster_size)/(atom_density),1.0/3.0));
3643 /* Estimates the interaction volume^2 for non-local interactions */
3644 static real nonlocal_vol2(const gmx_domdec_zones_t *zones,rvec ls,real r)
3653 /* Here we simply add up the volumes of 1, 2 or 3 1D decomposition
3654 * not home interaction volume^2. As these volumes are not additive,
3655 * this is an overestimate, but it would only be significant in the limit
3656 * of small cells, where we anyhow need to split the lists into
3657 * as small parts as possible.
3660 for(z=0; z<zones->n; z++)
3662 if (zones->shift[z][XX] + zones->shift[z][YY] + zones->shift[z][ZZ] == 1)
3667 for(d=0; d<DIM; d++)
3669 if (zones->shift[z][d] == 0)
3673 za *= zones->size[z].x1[d] - zones->size[z].x0[d];
3677 /* 4 octants of a sphere */
3678 vold_est = 0.25*M_PI*r*r*r*r;
3679 /* 4 quarter pie slices on the edges */
3680 vold_est += 4*cl*M_PI/6.0*r*r*r;
3681 /* One rectangular volume on a face */
3682 vold_est += ca*0.5*r*r;
3684 vol2_est_tot += vold_est*za;
3688 return vol2_est_tot;
3691 /* Estimates the average size of a full j-list for super/sub setup */
3692 static int get_nsubpair_max(const nbnxn_search_t nbs,
3695 int min_ci_balanced)
3697 const nbnxn_grid_t *grid;
3699 real xy_diag2,r_eff_sup,vol_est,nsp_est,nsp_est_nl;
3702 grid = &nbs->grid[0];
3704 ls[XX] = (grid->c1[XX] - grid->c0[XX])/(grid->ncx*GPU_NSUBCELL_X);
3705 ls[YY] = (grid->c1[YY] - grid->c0[YY])/(grid->ncy*GPU_NSUBCELL_Y);
3706 ls[ZZ] = (grid->c1[ZZ] - grid->c0[ZZ])*grid->ncx*grid->ncy/(grid->nc*GPU_NSUBCELL_Z);
3708 /* The average squared length of the diagonal of a sub cell */
3709 xy_diag2 = ls[XX]*ls[XX] + ls[YY]*ls[YY] + ls[ZZ]*ls[ZZ];
3711 /* The formulas below are a heuristic estimate of the average nsj per si*/
3712 r_eff_sup = rlist + nbnxn_rlist_inc_nonloc_fac*sqr((grid->na_c - 1.0)/grid->na_c)*sqrt(xy_diag2/3);
3714 if (!nbs->DomDec || nbs->zones->n == 1)
3721 sqr(grid->atom_density/grid->na_c)*
3722 nonlocal_vol2(nbs->zones,ls,r_eff_sup);
3727 /* Sub-cell interacts with itself */
3728 vol_est = ls[XX]*ls[YY]*ls[ZZ];
3729 /* 6/2 rectangular volume on the faces */
3730 vol_est += (ls[XX]*ls[YY] + ls[XX]*ls[ZZ] + ls[YY]*ls[ZZ])*r_eff_sup;
3731 /* 12/2 quarter pie slices on the edges */
3732 vol_est += 2*(ls[XX] + ls[YY] + ls[ZZ])*0.25*M_PI*sqr(r_eff_sup);
3733 /* 4 octants of a sphere */
3734 vol_est += 0.5*4.0/3.0*M_PI*pow(r_eff_sup,3);
3736 nsp_est = grid->nsubc_tot*vol_est*grid->atom_density/grid->na_c;
3738 /* Subtract the non-local pair count */
3739 nsp_est -= nsp_est_nl;
3743 fprintf(debug,"nsp_est local %5.1f non-local %5.1f\n",
3744 nsp_est,nsp_est_nl);
3749 nsp_est = nsp_est_nl;
3752 if (min_ci_balanced <= 0 || grid->nc >= min_ci_balanced || grid->nc == 0)
3754 /* We don't need to worry */
3759 /* Thus the (average) maximum j-list size should be as follows */
3760 nsubpair_max = max(1,(int)(nsp_est/min_ci_balanced+0.5));
3762 /* Since the target value is a maximum (this avoid high outliers,
3763 * which lead to load imbalance), not average, we get more lists
3764 * than we ask for (to compensate we need to add GPU_NSUBCELL*4/4).
3765 * But more importantly, the optimal GPU performance moves
3766 * to lower number of block for very small blocks.
3767 * To compensate we add the maximum pair count per cj4.
3769 nsubpair_max += GPU_NSUBCELL*NBNXN_CPU_CLUSTER_I_SIZE;
3774 fprintf(debug,"nbl nsp estimate %.1f, nsubpair_max %d\n",
3775 nsp_est,nsubpair_max);
3778 return nsubpair_max;
3781 /* Debug list print function */
3782 static void print_nblist_ci_cj(FILE *fp,const nbnxn_pairlist_t *nbl)
3786 for(i=0; i<nbl->nci; i++)
3788 fprintf(fp,"ci %4d shift %2d ncj %3d\n",
3789 nbl->ci[i].ci,nbl->ci[i].shift,
3790 nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start);
3792 for(j=nbl->ci[i].cj_ind_start; j<nbl->ci[i].cj_ind_end; j++)
3794 fprintf(fp," cj %5d imask %x\n",
3801 /* Debug list print function */
3802 static void print_nblist_sci_cj(FILE *fp,const nbnxn_pairlist_t *nbl)
3806 for(i=0; i<nbl->nsci; i++)
3808 fprintf(fp,"ci %4d shift %2d ncj4 %2d\n",
3809 nbl->sci[i].sci,nbl->sci[i].shift,
3810 nbl->sci[i].cj4_ind_end - nbl->sci[i].cj4_ind_start);
3812 for(j4=nbl->sci[i].cj4_ind_start; j4<nbl->sci[i].cj4_ind_end; j4++)
3816 fprintf(fp," sj %5d imask %x\n",
3818 nbl->cj4[j4].imei[0].imask);
3824 /* Combine pair lists *nbl generated on multiple threads nblc */
3825 static void combine_nblists(int nnbl,nbnxn_pairlist_t **nbl,
3826 nbnxn_pairlist_t *nblc)
3828 int nsci,ncj4,nexcl;
3833 gmx_incons("combine_nblists does not support simple lists");
3838 nexcl = nblc->nexcl;
3839 for(i=0; i<nnbl; i++)
3841 nsci += nbl[i]->nsci;
3842 ncj4 += nbl[i]->ncj4;
3843 nexcl += nbl[i]->nexcl;
3846 if (nsci > nblc->sci_nalloc)
3848 nb_realloc_sci(nblc,nsci);
3850 if (ncj4 > nblc->cj4_nalloc)
3852 nblc->cj4_nalloc = over_alloc_small(ncj4);
3853 nbnxn_realloc_void((void **)&nblc->cj4,
3854 nblc->ncj4*sizeof(*nblc->cj4),
3855 nblc->cj4_nalloc*sizeof(*nblc->cj4),
3856 nblc->alloc,nblc->free);
3858 if (nexcl > nblc->excl_nalloc)
3860 nblc->excl_nalloc = over_alloc_small(nexcl);
3861 nbnxn_realloc_void((void **)&nblc->excl,
3862 nblc->nexcl*sizeof(*nblc->excl),
3863 nblc->excl_nalloc*sizeof(*nblc->excl),
3864 nblc->alloc,nblc->free);
3867 /* Each thread should copy its own data to the combined arrays,
3868 * as otherwise data will go back and forth between different caches.
3870 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntPairsearch)) schedule(static)
3871 for(n=0; n<nnbl; n++)
3878 const nbnxn_pairlist_t *nbli;
3880 /* Determine the offset in the combined data for our thread */
3881 sci_offset = nblc->nsci;
3882 cj4_offset = nblc->ncj4;
3883 ci_offset = nblc->nci_tot;
3884 excl_offset = nblc->nexcl;
3888 sci_offset += nbl[i]->nsci;
3889 cj4_offset += nbl[i]->ncj4;
3890 ci_offset += nbl[i]->nci_tot;
3891 excl_offset += nbl[i]->nexcl;
3896 for(i=0; i<nbli->nsci; i++)
3898 nblc->sci[sci_offset+i] = nbli->sci[i];
3899 nblc->sci[sci_offset+i].cj4_ind_start += cj4_offset;
3900 nblc->sci[sci_offset+i].cj4_ind_end += cj4_offset;
3903 for(j4=0; j4<nbli->ncj4; j4++)
3905 nblc->cj4[cj4_offset+j4] = nbli->cj4[j4];
3906 nblc->cj4[cj4_offset+j4].imei[0].excl_ind += excl_offset;
3907 nblc->cj4[cj4_offset+j4].imei[1].excl_ind += excl_offset;
3910 for(j4=0; j4<nbli->nexcl; j4++)
3912 nblc->excl[excl_offset+j4] = nbli->excl[j4];
3916 for(n=0; n<nnbl; n++)
3918 nblc->nsci += nbl[n]->nsci;
3919 nblc->ncj4 += nbl[n]->ncj4;
3920 nblc->nci_tot += nbl[n]->nci_tot;
3921 nblc->nexcl += nbl[n]->nexcl;
3925 /* Returns the next ci to be processes by our thread */
3926 static gmx_bool next_ci(const nbnxn_grid_t *grid,
3928 int nth,int ci_block,
3929 int *ci_x,int *ci_y,
3935 if (*ci_b == ci_block)
3937 /* Jump to the next block assigned to this task */
3938 *ci += (nth - 1)*ci_block;
3942 if (*ci >= grid->nc*conv)
3947 while (*ci >= grid->cxy_ind[*ci_x*grid->ncy + *ci_y + 1]*conv)
3950 if (*ci_y == grid->ncy)
3960 /* Returns the distance^2 for which we put cell pairs in the list
3961 * without checking atom pair distances. This is usually < rlist^2.
3963 static float boundingbox_only_distance2(const nbnxn_grid_t *gridi,
3964 const nbnxn_grid_t *gridj,
3968 /* If the distance between two sub-cell bounding boxes is less
3969 * than this distance, do not check the distance between
3970 * all particle pairs in the sub-cell, since then it is likely
3971 * that the box pair has atom pairs within the cut-off.
3972 * We use the nblist cut-off minus 0.5 times the average x/y diagonal
3973 * spacing of the sub-cells. Around 40% of the checked pairs are pruned.
3974 * Using more than 0.5 gains at most 0.5%.
3975 * If forces are calculated more than twice, the performance gain
3976 * in the force calculation outweighs the cost of checking.
3977 * Note that with subcell lists, the atom-pair distance check
3978 * is only performed when only 1 out of 8 sub-cells in within range,
3979 * this is because the GPU is much faster than the cpu.
3984 bbx = 0.5*(gridi->sx + gridj->sx);
3985 bby = 0.5*(gridi->sy + gridj->sy);
3988 bbx /= GPU_NSUBCELL_X;
3989 bby /= GPU_NSUBCELL_Y;
3992 rbb2 = sqr(max(0,rlist - 0.5*sqrt(bbx*bbx + bby*bby)));
3997 return (float)((1+GMX_FLOAT_EPS)*rbb2);
4001 static int get_ci_block_size(const nbnxn_grid_t *gridi,
4002 gmx_bool bDomDec, int nth,
4003 gmx_bool *bFBufferFlag)
4005 const int ci_block_enum = 5;
4006 const int ci_block_denom = 11;
4007 const int ci_block_min_atoms = 16;
4010 /* Here we decide how to distribute the blocks over the threads.
4011 * We use prime numbers to try to avoid that the grid size becomes
4012 * a multiple of the number of threads, which would lead to some
4013 * threads getting "inner" pairs and others getting boundary pairs,
4014 * which in turns will lead to load imbalance between threads.
4015 * Set the block size as 5/11/ntask times the average number of cells
4016 * in a y,z slab. This should ensure a quite uniform distribution
4017 * of the grid parts of the different thread along all three grid
4018 * zone boundaries with 3D domain decomposition. At the same time
4019 * the blocks will not become too small.
4021 ci_block = (gridi->nc*ci_block_enum)/(ci_block_denom*gridi->ncx*nth);
4023 /* Ensure the blocks are not too small: avoids cache invalidation */
4024 if (ci_block*gridi->na_sc < ci_block_min_atoms)
4026 ci_block = (ci_block_min_atoms + gridi->na_sc - 1)/gridi->na_sc;
4029 /* Without domain decomposition
4030 * or with less than 3 blocks per task, divide in nth blocks.
4032 if (!bDomDec || ci_block*3*nth > gridi->nc)
4034 ci_block = (gridi->nc + nth - 1)/nth;
4035 /* With non-interleaved blocks it makes sense to flag which
4036 * part of the force output thread buffer we access.
4037 * We use bit flags, so we have to check if it fits.
4039 *bFBufferFlag = (nth > 1 && nth <= sizeof(unsigned int)*8);
4043 *bFBufferFlag = FALSE;
4049 /* Generates the part of pair-list nbl assigned to our thread */
4050 static void nbnxn_make_pairlist_part(const nbnxn_search_t nbs,
4051 const nbnxn_grid_t *gridi,
4052 const nbnxn_grid_t *gridj,
4053 nbnxn_search_work_t *work,
4054 const nbnxn_atomdata_t *nbat,
4055 const t_blocka *excl,
4059 gmx_bool bFBufferFlag,
4062 int min_ci_balanced,
4064 nbnxn_pairlist_t *nbl)
4071 int ci_b,ci,ci_x,ci_y,ci_xy,cj;
4078 const float *bb_i,*bbcz_i,*bbcz_j;
4080 real bx0,bx1,by0,by1,bz0,bz1;
4082 real d2cx,d2z,d2z_cx,d2z_cy,d2zx,d2zxy,d2xy;
4083 int cxf,cxl,cyf,cyf_x,cyl;
4088 int gridj_flag_shift=0,cj_offset=0;
4089 unsigned *gridj_flag=NULL;
4090 int ncj_old_i,ncj_old_j;
4092 nbs_cycle_start(&work->cc[enbsCCsearch]);
4094 if (gridj->bSimple != nbl->bSimple)
4096 gmx_incons("Grid incompatible with pair-list");
4100 nbl->na_sc = gridj->na_sc;
4101 nbl->na_ci = gridj->na_c;
4102 nbl->na_cj = nbnxn_kernel_to_cj_size(nb_kernel_type);
4103 na_cj_2log = get_2log(nbl->na_cj);
4109 init_grid_flags(&work->gridi_flags,gridi);
4110 init_grid_flags(&work->gridj_flags,gridj);
4112 /* To flag j-blocks for gridj, we need to convert j-clusters to flag blocks */
4113 gridj_flag_shift = 0;
4114 while ((nbl->na_cj<<gridj_flag_shift) < NBNXN_CELLBLOCK_SIZE*nbl->na_ci)
4118 /* We will subtract the cell offset, which is not a multiple of the block size */
4119 cj_offset = ci_to_cj(get_2log(nbl->na_cj),gridj->cell0);
4121 gridj_flag = work->gridj_flags.flag;
4124 copy_mat(nbs->box,box);
4126 rl2 = nbl->rlist*nbl->rlist;
4128 rbb2 = boundingbox_only_distance2(gridi,gridj,nbl->rlist,nbl->bSimple);
4132 fprintf(debug,"nbl bounding box only distance %f\n",sqrt(rbb2));
4135 /* Set the shift range */
4136 for(d=0; d<DIM; d++)
4138 /* Check if we need periodicity shifts.
4139 * Without PBC or with domain decomposition we don't need them.
4141 if (d >= ePBC2npbcdim(nbs->ePBC) || nbs->dd_dim[d])
4148 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
4159 if (nbl->bSimple && !gridi->bSimple)
4161 conv_i = gridi->na_sc/gridj->na_sc;
4162 bb_i = gridi->bb_simple;
4163 bbcz_i = gridi->bbcz_simple;
4164 flags_i = gridi->flags_simple;
4170 bbcz_i = gridi->bbcz;
4171 flags_i = gridi->flags;
4173 cell0_i = gridi->cell0*conv_i;
4175 bbcz_j = gridj->bbcz;
4179 /* Blocks of the conversion factor - 1 give a large repeat count
4180 * combined with a small block size. This should result in good
4181 * load balancing for both small and large domains.
4183 ci_block = conv_i - 1;
4187 fprintf(debug,"nbl nc_i %d col.av. %.1f ci_block %d\n",
4188 gridi->nc,gridi->nc/(double)(gridi->ncx*gridi->ncy),ci_block);
4194 /* Initially ci_b and ci to 1 before where we want them to start,
4195 * as they will both be incremented in next_ci.
4198 ci = th*ci_block - 1;
4201 while (next_ci(gridi,conv_i,nth,ci_block,&ci_x,&ci_y,&ci_b,&ci))
4203 if (nbl->bSimple && flags_i[ci] == 0)
4208 ncj_old_i = nbl->ncj;
4211 if (gridj != gridi && shp[XX] == 0)
4215 bx1 = bb_i[ci*NNBSBB_B+NNBSBB_C+XX];
4219 bx1 = gridi->c0[XX] + (ci_x+1)*gridi->sx;
4221 if (bx1 < gridj->c0[XX])
4223 d2cx = sqr(gridj->c0[XX] - bx1);
4232 ci_xy = ci_x*gridi->ncy + ci_y;
4234 /* Loop over shift vectors in three dimensions */
4235 for (tz=-shp[ZZ]; tz<=shp[ZZ]; tz++)
4237 shz = tz*box[ZZ][ZZ];
4239 bz0 = bbcz_i[ci*NNBSBB_D ] + shz;
4240 bz1 = bbcz_i[ci*NNBSBB_D+1] + shz;
4252 d2z = sqr(bz0 - box[ZZ][ZZ]);
4255 d2z_cx = d2z + d2cx;
4263 bz1/((real)(gridi->cxy_ind[ci_xy+1] - gridi->cxy_ind[ci_xy]));
4268 /* The check with bz1_frac close to or larger than 1 comes later */
4270 for (ty=-shp[YY]; ty<=shp[YY]; ty++)
4272 shy = ty*box[YY][YY] + tz*box[ZZ][YY];
4276 by0 = bb_i[ci*NNBSBB_B +YY] + shy;
4277 by1 = bb_i[ci*NNBSBB_B+NNBSBB_C+YY] + shy;
4281 by0 = gridi->c0[YY] + (ci_y )*gridi->sy + shy;
4282 by1 = gridi->c0[YY] + (ci_y+1)*gridi->sy + shy;
4285 get_cell_range(by0,by1,
4286 gridj->ncy,gridj->c0[YY],gridj->sy,gridj->inv_sy,
4296 if (by1 < gridj->c0[YY])
4298 d2z_cy += sqr(gridj->c0[YY] - by1);
4300 else if (by0 > gridj->c1[YY])
4302 d2z_cy += sqr(by0 - gridj->c1[YY]);
4305 for (tx=-shp[XX]; tx<=shp[XX]; tx++)
4307 shift = XYZ2IS(tx,ty,tz);
4309 #ifdef NBNXN_SHIFT_BACKWARD
4310 if (gridi == gridj && shift > CENTRAL)
4316 shx = tx*box[XX][XX] + ty*box[YY][XX] + tz*box[ZZ][XX];
4320 bx0 = bb_i[ci*NNBSBB_B +XX] + shx;
4321 bx1 = bb_i[ci*NNBSBB_B+NNBSBB_C+XX] + shx;
4325 bx0 = gridi->c0[XX] + (ci_x )*gridi->sx + shx;
4326 bx1 = gridi->c0[XX] + (ci_x+1)*gridi->sx + shx;
4329 get_cell_range(bx0,bx1,
4330 gridj->ncx,gridj->c0[XX],gridj->sx,gridj->inv_sx,
4341 new_ci_entry(nbl,cell0_i+ci,shift,flags_i[ci],
4346 new_sci_entry(nbl,cell0_i+ci,shift,flags_i[ci],
4350 #ifndef NBNXN_SHIFT_BACKWARD
4353 if (shift == CENTRAL && gridi == gridj &&
4357 /* Leave the pairs with i > j.
4358 * x is the major index, so skip half of it.
4365 set_icell_bb_simple(bb_i,ci,shx,shy,shz,
4370 set_icell_bb_supersub(bb_i,ci,shx,shy,shz,
4374 nbs->icell_set_x(cell0_i+ci,shx,shy,shz,
4375 gridi->na_c,nbat->xstride,nbat->x,
4378 for(cx=cxf; cx<=cxl; cx++)
4381 if (gridj->c0[XX] + cx*gridj->sx > bx1)
4383 d2zx += sqr(gridj->c0[XX] + cx*gridj->sx - bx1);
4385 else if (gridj->c0[XX] + (cx+1)*gridj->sx < bx0)
4387 d2zx += sqr(gridj->c0[XX] + (cx+1)*gridj->sx - bx0);
4390 #ifndef NBNXN_SHIFT_BACKWARD
4391 if (gridi == gridj &&
4392 cx == 0 && cyf < ci_y)
4394 if (gridi == gridj &&
4395 cx == 0 && shift == CENTRAL && cyf < ci_y)
4398 /* Leave the pairs with i > j.
4399 * Skip half of y when i and j have the same x.
4408 for(cy=cyf_x; cy<=cyl; cy++)
4410 c0 = gridj->cxy_ind[cx*gridj->ncy+cy];
4411 c1 = gridj->cxy_ind[cx*gridj->ncy+cy+1];
4412 #ifdef NBNXN_SHIFT_BACKWARD
4413 if (gridi == gridj &&
4414 shift == CENTRAL && c0 < ci)
4421 if (gridj->c0[YY] + cy*gridj->sy > by1)
4423 d2zxy += sqr(gridj->c0[YY] + cy*gridj->sy - by1);
4425 else if (gridj->c0[YY] + (cy+1)*gridj->sy < by0)
4427 d2zxy += sqr(gridj->c0[YY] + (cy+1)*gridj->sy - by0);
4429 if (c1 > c0 && d2zxy < rl2)
4431 cs = c0 + (int)(bz1_frac*(c1 - c0));
4439 /* Find the lowest cell that can possibly
4444 (bbcz_j[cf*NNBSBB_D+1] >= bz0 ||
4445 d2xy + sqr(bbcz_j[cf*NNBSBB_D+1] - bz0) < rl2))
4450 /* Find the highest cell that can possibly
4455 (bbcz_j[cl*NNBSBB_D] <= bz1 ||
4456 d2xy + sqr(bbcz_j[cl*NNBSBB_D] - bz1) < rl2))
4461 #ifdef NBNXN_REFCODE
4463 /* Simple reference code, for debugging,
4464 * overrides the more complex code above.
4469 for(k=c0; k<c1; k++)
4471 if (box_dist2(bx0,bx1,by0,by1,bz0,bz1,
4472 bb+k*NNBSBB_B) < rl2 &&
4477 if (box_dist2(bx0,bx1,by0,by1,bz0,bz1,
4478 bb+k*NNBSBB_B) < rl2 &&
4489 /* We want each atom/cell pair only once,
4490 * only use cj >= ci.
4492 #ifndef NBNXN_SHIFT_BACKWARD
4495 if (shift == CENTRAL)
4504 /* For f buffer flags with simple lists */
4505 ncj_old_j = nbl->ncj;
4507 switch (nb_kernel_type)
4510 check_subcell_list_space_simple(nbl,cl-cf+1);
4512 make_cluster_list_simple(gridj,
4514 (gridi == gridj && shift == CENTRAL),
4519 #ifdef NBNXN_SEARCH_SSE
4520 case nbk4xN_X86_SIMD128:
4521 check_subcell_list_space_simple(nbl,ci_to_cj(na_cj_2log,cl-cf)+2);
4522 make_cluster_list_x86_simd128(gridj,
4524 (gridi == gridj && shift == CENTRAL),
4529 #ifdef GMX_X86_AVX_256
4530 case nbk4xN_X86_SIMD256:
4531 check_subcell_list_space_simple(nbl,ci_to_cj(na_cj_2log,cl-cf)+2);
4532 make_cluster_list_x86_simd256(gridj,
4534 (gridi == gridj && shift == CENTRAL),
4541 case nbk8x8x8_PlainC:
4543 check_subcell_list_space_supersub(nbl,cl-cf+1);
4544 for(cj=cf; cj<=cl; cj++)
4546 make_cluster_list_supersub(nbs,gridi,gridj,
4548 (gridi == gridj && shift == CENTRAL && ci == cj),
4549 nbat->xstride,nbat->x,
4555 ncpcheck += cl - cf + 1;
4557 if (bFBufferFlag && nbl->ncj > ncj_old_j)
4561 cbf = (nbl->cj[ncj_old_j].cj - cj_offset) >> gridj_flag_shift;
4562 cbl = (nbl->cj[nbl->ncj-1].cj - cj_offset) >> gridj_flag_shift;
4563 for(cb=cbf; cb<=cbl; cb++)
4565 gridj_flag[cb] = 1U<<th;
4573 /* Set the exclusions for this ci list */
4576 set_ci_top_excls(nbs,
4578 shift == CENTRAL && gridi == gridj,
4581 &(nbl->ci[nbl->nci]),
4586 set_sci_top_excls(nbs,
4588 shift == CENTRAL && gridi == gridj,
4590 &(nbl->sci[nbl->nsci]),
4594 /* Close this ci list */
4597 close_ci_entry_simple(nbl);
4601 close_ci_entry_supersub(nbl,
4603 progBal,min_ci_balanced,
4610 if (bFBufferFlag && nbl->ncj > ncj_old_i)
4612 work->gridi_flags.flag[ci>>NBNXN_CELLBLOCK_SIZE_2LOG] = 1U<<th;
4616 work->ndistc = ndistc;
4618 nbs_cycle_stop(&work->cc[enbsCCsearch]);
4622 fprintf(debug,"number of distance checks %d\n",ndistc);
4623 fprintf(debug,"ncpcheck %s %d\n",gridi==gridj ? "local" : "non-local",
4628 print_nblist_statistics_simple(debug,nbl,nbs,rlist);
4632 print_nblist_statistics_supersub(debug,nbl,nbs,rlist);
4638 static void reduce_cellblock_flags(const nbnxn_search_t nbs,
4640 const nbnxn_grid_t *gridi,
4641 const nbnxn_grid_t *gridj)
4644 const unsigned *flag;
4646 if (gridi->cellblock_flags.bUse)
4648 for(nbl=0; nbl<nnbl; nbl++)
4650 flag = nbs->work[nbl].gridi_flags.flag;
4652 for(cb=0; cb<gridi->cellblock_flags.ncb; cb++)
4654 gridi->cellblock_flags.flag[cb] |= flag[cb];
4658 if (gridj->cellblock_flags.bUse)
4660 for(nbl=0; nbl<nnbl; nbl++)
4662 flag = nbs->work[nbl].gridj_flags.flag;
4664 for(cb=0; cb<gridj->cellblock_flags.ncb; cb++)
4666 gridj->cellblock_flags.flag[cb] |= flag[cb];
4672 static void print_reduction_cost(const nbnxn_grid_t *grids,int ngrid,int nnbl)
4675 const nbnxn_grid_t *grid;
4677 for(g=0; g<ngrid; g++)
4682 if (grid->cellblock_flags.bUse)
4685 for(cb=0; cb<grid->cellblock_flags.ncb; cb++)
4687 for(nbl=0; nbl<nnbl; nbl++)
4689 if (grid->cellblock_flags.flag[cb] == 1)
4693 else if (grid->cellblock_flags.flag[cb] & (1U<<nbl))
4702 c = nnbl*grid->cellblock_flags.ncb;
4704 fprintf(debug,"nbnxn reduction buffers, grid %d: %d flag %d only buf. 0: %4.2f av. reduction: %4.2f\n",
4705 g,nnbl,grid->cellblock_flags.bUse,
4706 c0/(double)(grid->cellblock_flags.ncb),
4707 c/(double)(grid->cellblock_flags.ncb));
4711 /* Make a local or non-local pair-list, depending on iloc */
4712 void nbnxn_make_pairlist(const nbnxn_search_t nbs,
4713 const nbnxn_atomdata_t *nbat,
4714 const t_blocka *excl,
4716 int min_ci_balanced,
4717 nbnxn_pairlist_set_t *nbl_list,
4722 nbnxn_grid_t *gridi,*gridj;
4723 int nzi,zi,zj0,zj1,zj;
4727 nbnxn_pairlist_t **nbl;
4729 gmx_bool CombineNBLists,bFBufferFlag;
4730 int np_tot,np_noq,np_hlj,nap;
4732 nnbl = nbl_list->nnbl;
4733 nbl = nbl_list->nbl;
4734 CombineNBLists = nbl_list->bCombined;
4738 fprintf(debug,"ns making %d nblists\n", nnbl);
4741 if (nbl_list->bSimple)
4743 switch (nb_kernel_type)
4745 #ifdef NBNXN_SEARCH_SSE
4746 case nbk4xN_X86_SIMD128:
4747 nbs->icell_set_x = icell_set_x_x86_simd128;
4749 #ifdef GMX_X86_AVX_256
4750 case nbk4xN_X86_SIMD256:
4751 nbs->icell_set_x = icell_set_x_x86_simd256;
4756 nbs->icell_set_x = icell_set_x_simple;
4762 #ifdef NBNXN_SEARCH_SSE
4763 nbs->icell_set_x = icell_set_x_supersub_sse8;
4765 nbs->icell_set_x = icell_set_x_supersub;
4771 /* Only zone (grid) 0 vs 0 */
4778 nzi = nbs->zones->nizone;
4781 if (!nbl_list->bSimple && min_ci_balanced > 0)
4783 nsubpair_max = get_nsubpair_max(nbs,iloc,rlist,min_ci_balanced);
4790 /* Clear all pair-lists */
4791 for(th=0; th<nnbl; th++)
4793 clear_pairlist(nbl[th]);
4796 for(zi=0; zi<nzi; zi++)
4798 gridi = &nbs->grid[zi];
4800 if (NONLOCAL_I(iloc))
4802 zj0 = nbs->zones->izone[zi].j0;
4803 zj1 = nbs->zones->izone[zi].j1;
4809 for(zj=zj0; zj<zj1; zj++)
4811 gridj = &nbs->grid[zj];
4815 fprintf(debug,"ns search grid %d vs %d\n",zi,zj);
4818 nbs_cycle_start(&nbs->cc[enbsCCsearch]);
4820 if (nbl[0]->bSimple && !gridi->bSimple)
4822 /* Hybrid list, determine blocking later */
4824 bFBufferFlag = FALSE;
4828 ci_block = get_ci_block_size(gridi,nbs->DomDec,nnbl,
4832 bFBufferFlag = FALSE;
4837 fprintf(debug,"grid %d %d F buffer flags %d\n",
4838 zi,zj,bFBufferFlag);
4841 #pragma omp parallel for num_threads(nnbl) schedule(static)
4842 for(th=0; th<nnbl; th++)
4844 if (CombineNBLists && th > 0)
4846 clear_pairlist(nbl[th]);
4849 /* Divide the i super cell equally over the nblists */
4850 nbnxn_make_pairlist_part(nbs,gridi,gridj,
4851 &nbs->work[th],nbat,excl,
4857 (LOCAL_I(iloc) || nbs->zones->n <= 2),
4862 nbs_cycle_stop(&nbs->cc[enbsCCsearch]);
4867 for(th=0; th<nnbl; th++)
4869 inc_nrnb(nrnb,eNR_NBNXN_DIST2,nbs->work[th].ndistc);
4871 if (nbl_list->bSimple)
4873 np_tot += nbl[th]->ncj;
4874 np_noq += nbl[th]->work->ncj_noq;
4875 np_hlj += nbl[th]->work->ncj_hlj;
4879 /* This count ignores potential subsequent pair pruning */
4880 np_tot += nbl[th]->nci_tot;
4883 nap = nbl[0]->na_ci*nbl[0]->na_cj;
4884 nbl_list->natpair_ljq = (np_tot - np_noq)*nap - np_hlj*nap/2;
4885 nbl_list->natpair_lj = np_noq*nap;
4886 nbl_list->natpair_q = np_hlj*nap/2;
4888 if (CombineNBLists && nnbl > 1)
4890 nbs_cycle_start(&nbs->cc[enbsCCcombine]);
4892 combine_nblists(nnbl-1,nbl+1,nbl[0]);
4894 nbs_cycle_stop(&nbs->cc[enbsCCcombine]);
4899 reduce_cellblock_flags(nbs,nnbl,gridi,gridj);
4903 gridi->cellblock_flags.bUse = FALSE;
4904 gridj->cellblock_flags.bUse = FALSE;
4910 print_supersub_nsp("nsubpair",nbl[0],iloc);
4913 /* Special performance logging stuff (env.var. GMX_NBNXN_CYCLE) */
4916 nbs->search_count++;
4918 if (nbs->print_cycles &&
4919 (!nbs->DomDec || (nbs->DomDec && !LOCAL_I(iloc))) &&
4920 nbs->search_count % 100 == 0)
4922 nbs_cycle_print(stderr,nbs);
4925 if (debug && (CombineNBLists && nnbl > 1))
4927 if (nbl[0]->bSimple)
4929 print_nblist_statistics_simple(debug,nbl[0],nbs,rlist);
4933 print_nblist_statistics_supersub(debug,nbl[0],nbs,rlist);
4939 if (nbl[0]->bSimple)
4941 print_nblist_ci_cj(debug,nbl[0]);
4945 print_nblist_sci_cj(debug,nbl[0]);
4948 print_reduction_cost(nbs->grid,nbs->ngrid,nnbl);