Merge branch release-4-6
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_search.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2012, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  */
32
33 #ifdef HAVE_CONFIG_H
34 #include <config.h>
35 #endif
36
37 #include <math.h>
38 #include <string.h>
39 #include "sysstuff.h"
40 #include "smalloc.h"
41 #include "macros.h"
42 #include "maths.h"
43 #include "vec.h"
44 #include "pbc.h"
45 #include "nbnxn_consts.h"
46 /* nbnxn_internal.h included gmx_simd_macros.h */
47 #include "nbnxn_internal.h"
48 #ifdef GMX_NBNXN_SIMD
49 #include "gmx_simd_vec.h"
50 #endif
51 #include "nbnxn_atomdata.h"
52 #include "nbnxn_search.h"
53 #include "gmx_cyclecounter.h"
54 #include "gmxfio.h"
55 #include "gmx_omp_nthreads.h"
56 #include "nrnb.h"
57
58
59 #ifdef NBNXN_SEARCH_BB_SSE
60 /* We use SSE or AVX-128bit for bounding box calculations */
61
62 #ifndef GMX_DOUBLE
63 /* Single precision BBs + coordinates, we can also load coordinates using SSE */
64 #define NBNXN_SEARCH_SSE_SINGLE
65 #endif
66
67 /* Include basic SSE2 stuff */
68 #include <emmintrin.h>
69
70 #if defined NBNXN_SEARCH_SSE_SINGLE && (GPU_NSUBCELL == 4 || GPU_NSUBCELL == 8)
71 /* Store bounding boxes with x, y and z coordinates in packs of 4 */
72 #define NBNXN_PBB_SSE
73 #endif
74
75 /* The width of SSE/AVX128 with single precision for bounding boxes with GPU.
76  * Here AVX-256 turns out to be slightly slower than AVX-128.
77  */
78 #define STRIDE_PBB        4
79 #define STRIDE_PBB_2LOG   2
80
81 #endif /* NBNXN_SEARCH_BB_SSE */
82
83 #ifdef GMX_NBNXN_SIMD
84
85 /* The functions below are macros as they are performance sensitive */
86
87 /* 4x4 list, pack=4: no complex conversion required */
88 /* i-cluster to j-cluster conversion */
89 #define CI_TO_CJ_J4(ci)   (ci)
90 /* cluster index to coordinate array index conversion */
91 #define X_IND_CI_J4(ci)  ((ci)*STRIDE_P4)
92 #define X_IND_CJ_J4(cj)  ((cj)*STRIDE_P4)
93
94 /* 4x2 list, pack=4: j-cluster size is half the packing width */
95 /* i-cluster to j-cluster conversion */
96 #define CI_TO_CJ_J2(ci)  ((ci)<<1)
97 /* cluster index to coordinate array index conversion */
98 #define X_IND_CI_J2(ci)  ((ci)*STRIDE_P4)
99 #define X_IND_CJ_J2(cj)  (((cj)>>1)*STRIDE_P4 + ((cj) & 1)*(PACK_X4>>1))
100
101 /* 4x8 list, pack=8: i-cluster size is half the packing width */
102 /* i-cluster to j-cluster conversion */
103 #define CI_TO_CJ_J8(ci)  ((ci)>>1)
104 /* cluster index to coordinate array index conversion */
105 #define X_IND_CI_J8(ci)  (((ci)>>1)*STRIDE_P8 + ((ci) & 1)*(PACK_X8>>1))
106 #define X_IND_CJ_J8(cj)  ((cj)*STRIDE_P8)
107
108 /* The j-cluster size is matched to the SIMD width */
109 #if GMX_SIMD_WIDTH_HERE == 2
110 #define CI_TO_CJ_SIMD_4XN(ci)  CI_TO_CJ_J2(ci)
111 #define X_IND_CI_SIMD_4XN(ci)  X_IND_CI_J2(ci)
112 #define X_IND_CJ_SIMD_4XN(cj)  X_IND_CJ_J2(cj)
113 #else
114 #if GMX_SIMD_WIDTH_HERE == 4
115 #define CI_TO_CJ_SIMD_4XN(ci)  CI_TO_CJ_J4(ci)
116 #define X_IND_CI_SIMD_4XN(ci)  X_IND_CI_J4(ci)
117 #define X_IND_CJ_SIMD_4XN(cj)  X_IND_CJ_J4(cj)
118 #else
119 #if GMX_SIMD_WIDTH_HERE == 8
120 #define CI_TO_CJ_SIMD_4XN(ci)  CI_TO_CJ_J8(ci)
121 #define X_IND_CI_SIMD_4XN(ci)  X_IND_CI_J8(ci)
122 #define X_IND_CJ_SIMD_4XN(cj)  X_IND_CJ_J8(cj)
123 /* Half SIMD with j-cluster size */
124 #define CI_TO_CJ_SIMD_2XNN(ci) CI_TO_CJ_J4(ci)
125 #define X_IND_CI_SIMD_2XNN(ci) X_IND_CI_J4(ci)
126 #define X_IND_CJ_SIMD_2XNN(cj) X_IND_CJ_J4(cj)
127 #else
128 #if GMX_SIMD_WIDTH_HERE == 16
129 #define CI_TO_CJ_SIMD_2XNN(ci) CI_TO_CJ_J8(ci)
130 #define X_IND_CI_SIMD_2XNN(ci) X_IND_CI_J8(ci)
131 #define X_IND_CJ_SIMD_2XNN(cj) X_IND_CJ_J8(cj)
132 #else
133 #error "unsupported GMX_NBNXN_SIMD_WIDTH"
134 #endif
135 #endif
136 #endif
137 #endif
138
139 #endif /* GMX_NBNXN_SIMD */
140
141
142 #ifdef NBNXN_SEARCH_BB_SSE
143 /* Store bounding boxes corners as quadruplets: xxxxyyyyzzzz */
144 #define NBNXN_BBXXXX
145 /* Size of bounding box corners quadruplet */
146 #define NNBSBB_XXXX      (NNBSBB_D*DIM*STRIDE_PBB)
147 #endif
148
149 /* We shift the i-particles backward for PBC.
150  * This leads to more conditionals than shifting forward.
151  * We do this to get more balanced pair lists.
152  */
153 #define NBNXN_SHIFT_BACKWARD
154
155
156 /* This define is a lazy way to avoid interdependence of the grid
157  * and searching data structures.
158  */
159 #define NBNXN_NA_SC_MAX (GPU_NSUBCELL*NBNXN_GPU_CLUSTER_SIZE)
160
161
162 static void nbs_cycle_clear(nbnxn_cycle_t *cc)
163 {
164     int i;
165
166     for (i = 0; i < enbsCCnr; i++)
167     {
168         cc[i].count = 0;
169         cc[i].c     = 0;
170     }
171 }
172
173 static double Mcyc_av(const nbnxn_cycle_t *cc)
174 {
175     return (double)cc->c*1e-6/cc->count;
176 }
177
178 static void nbs_cycle_print(FILE *fp, const nbnxn_search_t nbs)
179 {
180     int n;
181     int t;
182
183     fprintf(fp, "\n");
184     fprintf(fp, "ns %4d grid %4.1f search %4.1f red.f %5.3f",
185             nbs->cc[enbsCCgrid].count,
186             Mcyc_av(&nbs->cc[enbsCCgrid]),
187             Mcyc_av(&nbs->cc[enbsCCsearch]),
188             Mcyc_av(&nbs->cc[enbsCCreducef]));
189
190     if (nbs->nthread_max > 1)
191     {
192         if (nbs->cc[enbsCCcombine].count > 0)
193         {
194             fprintf(fp, " comb %5.2f",
195                     Mcyc_av(&nbs->cc[enbsCCcombine]));
196         }
197         fprintf(fp, " s. th");
198         for (t = 0; t < nbs->nthread_max; t++)
199         {
200             fprintf(fp, " %4.1f",
201                     Mcyc_av(&nbs->work[t].cc[enbsCCsearch]));
202         }
203     }
204     fprintf(fp, "\n");
205 }
206
207 static void nbnxn_grid_init(nbnxn_grid_t * grid)
208 {
209     grid->cxy_na      = NULL;
210     grid->cxy_ind     = NULL;
211     grid->cxy_nalloc  = 0;
212     grid->bb          = NULL;
213     grid->bbj         = NULL;
214     grid->nc_nalloc   = 0;
215 }
216
217 static int get_2log(int n)
218 {
219     int log2;
220
221     log2 = 0;
222     while ((1<<log2) < n)
223     {
224         log2++;
225     }
226     if ((1<<log2) != n)
227     {
228         gmx_fatal(FARGS, "nbnxn na_c (%d) is not a power of 2", n);
229     }
230
231     return log2;
232 }
233
234 static int nbnxn_kernel_to_ci_size(int nb_kernel_type)
235 {
236     switch (nb_kernel_type)
237     {
238         case nbnxnk4x4_PlainC:
239         case nbnxnk4xN_SIMD_4xN:
240         case nbnxnk4xN_SIMD_2xNN:
241             return NBNXN_CPU_CLUSTER_I_SIZE;
242         case nbnxnk8x8x8_CUDA:
243         case nbnxnk8x8x8_PlainC:
244             /* The cluster size for super/sub lists is only set here.
245              * Any value should work for the pair-search and atomdata code.
246              * The kernels, of course, might require a particular value.
247              */
248             return NBNXN_GPU_CLUSTER_SIZE;
249         default:
250             gmx_incons("unknown kernel type");
251     }
252
253     return 0;
254 }
255
256 int nbnxn_kernel_to_cj_size(int nb_kernel_type)
257 {
258     int nbnxn_simd_width = 0;
259     int cj_size          = 0;
260
261 #ifdef GMX_NBNXN_SIMD
262     nbnxn_simd_width = GMX_SIMD_WIDTH_HERE;
263 #endif
264
265     switch (nb_kernel_type)
266     {
267         case nbnxnk4x4_PlainC:
268             cj_size = NBNXN_CPU_CLUSTER_I_SIZE;
269             break;
270         case nbnxnk4xN_SIMD_4xN:
271             cj_size = nbnxn_simd_width;
272             break;
273         case nbnxnk4xN_SIMD_2xNN:
274             cj_size = nbnxn_simd_width/2;
275             break;
276         case nbnxnk8x8x8_CUDA:
277         case nbnxnk8x8x8_PlainC:
278             cj_size = nbnxn_kernel_to_ci_size(nb_kernel_type);
279             break;
280         default:
281             gmx_incons("unknown kernel type");
282     }
283
284     return cj_size;
285 }
286
287 static int ci_to_cj(int na_cj_2log, int ci)
288 {
289     switch (na_cj_2log)
290     {
291         case 2: return ci;     break;
292         case 1: return (ci<<1); break;
293         case 3: return (ci>>1); break;
294     }
295
296     return 0;
297 }
298
299 gmx_bool nbnxn_kernel_pairlist_simple(int nb_kernel_type)
300 {
301     if (nb_kernel_type == nbnxnkNotSet)
302     {
303         gmx_fatal(FARGS, "Non-bonded kernel type not set for Verlet-style pair-list.");
304     }
305
306     switch (nb_kernel_type)
307     {
308         case nbnxnk8x8x8_CUDA:
309         case nbnxnk8x8x8_PlainC:
310             return FALSE;
311
312         case nbnxnk4x4_PlainC:
313         case nbnxnk4xN_SIMD_4xN:
314         case nbnxnk4xN_SIMD_2xNN:
315             return TRUE;
316
317         default:
318             gmx_incons("Invalid nonbonded kernel type passed!");
319             return FALSE;
320     }
321 }
322
323 void nbnxn_init_search(nbnxn_search_t    * nbs_ptr,
324                        ivec               *n_dd_cells,
325                        gmx_domdec_zones_t *zones,
326                        int                 nthread_max)
327 {
328     nbnxn_search_t nbs;
329     int            d, g, t;
330
331     snew(nbs, 1);
332     *nbs_ptr = nbs;
333
334     nbs->DomDec = (n_dd_cells != NULL);
335
336     clear_ivec(nbs->dd_dim);
337     nbs->ngrid = 1;
338     if (nbs->DomDec)
339     {
340         nbs->zones = zones;
341
342         for (d = 0; d < DIM; d++)
343         {
344             if ((*n_dd_cells)[d] > 1)
345             {
346                 nbs->dd_dim[d] = 1;
347                 /* Each grid matches a DD zone */
348                 nbs->ngrid *= 2;
349             }
350         }
351     }
352
353     snew(nbs->grid, nbs->ngrid);
354     for (g = 0; g < nbs->ngrid; g++)
355     {
356         nbnxn_grid_init(&nbs->grid[g]);
357     }
358     nbs->cell        = NULL;
359     nbs->cell_nalloc = 0;
360     nbs->a           = NULL;
361     nbs->a_nalloc    = 0;
362
363     nbs->nthread_max = nthread_max;
364
365     /* Initialize the work data structures for each thread */
366     snew(nbs->work, nbs->nthread_max);
367     for (t = 0; t < nbs->nthread_max; t++)
368     {
369         nbs->work[t].cxy_na           = NULL;
370         nbs->work[t].cxy_na_nalloc    = 0;
371         nbs->work[t].sort_work        = NULL;
372         nbs->work[t].sort_work_nalloc = 0;
373     }
374
375     /* Initialize detailed nbsearch cycle counting */
376     nbs->print_cycles = (getenv("GMX_NBNXN_CYCLE") != 0);
377     nbs->search_count = 0;
378     nbs_cycle_clear(nbs->cc);
379     for (t = 0; t < nbs->nthread_max; t++)
380     {
381         nbs_cycle_clear(nbs->work[t].cc);
382     }
383 }
384
385 static real grid_atom_density(int n, rvec corner0, rvec corner1)
386 {
387     rvec size;
388
389     rvec_sub(corner1, corner0, size);
390
391     return n/(size[XX]*size[YY]*size[ZZ]);
392 }
393
394 static int set_grid_size_xy(const nbnxn_search_t nbs,
395                             nbnxn_grid_t *grid,
396                             int dd_zone,
397                             int n, rvec corner0, rvec corner1,
398                             real atom_density)
399 {
400     rvec size;
401     int  na_c;
402     real adens, tlen, tlen_x, tlen_y, nc_max;
403     int  t;
404
405     rvec_sub(corner1, corner0, size);
406
407     if (n > grid->na_sc)
408     {
409         /* target cell length */
410         if (grid->bSimple)
411         {
412             /* To minimize the zero interactions, we should make
413              * the largest of the i/j cell cubic.
414              */
415             na_c = max(grid->na_c, grid->na_cj);
416
417             /* Approximately cubic cells */
418             tlen   = pow(na_c/atom_density, 1.0/3.0);
419             tlen_x = tlen;
420             tlen_y = tlen;
421         }
422         else
423         {
424             /* Approximately cubic sub cells */
425             tlen   = pow(grid->na_c/atom_density, 1.0/3.0);
426             tlen_x = tlen*GPU_NSUBCELL_X;
427             tlen_y = tlen*GPU_NSUBCELL_Y;
428         }
429         /* We round ncx and ncy down, because we get less cell pairs
430          * in the nbsist when the fixed cell dimensions (x,y) are
431          * larger than the variable one (z) than the other way around.
432          */
433         grid->ncx = max(1, (int)(size[XX]/tlen_x));
434         grid->ncy = max(1, (int)(size[YY]/tlen_y));
435     }
436     else
437     {
438         grid->ncx = 1;
439         grid->ncy = 1;
440     }
441
442     grid->sx     = size[XX]/grid->ncx;
443     grid->sy     = size[YY]/grid->ncy;
444     grid->inv_sx = 1/grid->sx;
445     grid->inv_sy = 1/grid->sy;
446
447     if (dd_zone > 0)
448     {
449         /* This is a non-home zone, add an extra row of cells
450          * for particles communicated for bonded interactions.
451          * These can be beyond the cut-off. It doesn't matter where
452          * they end up on the grid, but for performance it's better
453          * if they don't end up in cells that can be within cut-off range.
454          */
455         grid->ncx++;
456         grid->ncy++;
457     }
458
459     /* We need one additional cell entry for particles moved by DD */
460     if (grid->ncx*grid->ncy+1 > grid->cxy_nalloc)
461     {
462         grid->cxy_nalloc = over_alloc_large(grid->ncx*grid->ncy+1);
463         srenew(grid->cxy_na, grid->cxy_nalloc);
464         srenew(grid->cxy_ind, grid->cxy_nalloc+1);
465     }
466     for (t = 0; t < nbs->nthread_max; t++)
467     {
468         if (grid->ncx*grid->ncy+1 > nbs->work[t].cxy_na_nalloc)
469         {
470             nbs->work[t].cxy_na_nalloc = over_alloc_large(grid->ncx*grid->ncy+1);
471             srenew(nbs->work[t].cxy_na, nbs->work[t].cxy_na_nalloc);
472         }
473     }
474
475     /* Worst case scenario of 1 atom in each last cell */
476     if (grid->na_cj <= grid->na_c)
477     {
478         nc_max = n/grid->na_sc + grid->ncx*grid->ncy;
479     }
480     else
481     {
482         nc_max = n/grid->na_sc + grid->ncx*grid->ncy*grid->na_cj/grid->na_c;
483     }
484
485     if (nc_max > grid->nc_nalloc)
486     {
487         grid->nc_nalloc = over_alloc_large(nc_max);
488         srenew(grid->nsubc, grid->nc_nalloc);
489         srenew(grid->bbcz, grid->nc_nalloc*NNBSBB_D);
490
491         sfree_aligned(grid->bb);
492         /* This snew also zeros the contents, this avoid possible
493          * floating exceptions in SSE with the unused bb elements.
494          */
495         if (grid->bSimple)
496         {
497             snew_aligned(grid->bb, grid->nc_nalloc, 16);
498         }
499         else
500         {
501 #ifdef NBNXN_BBXXXX
502             int pbb_nalloc;
503
504             pbb_nalloc = grid->nc_nalloc*GPU_NSUBCELL/STRIDE_PBB*NNBSBB_XXXX;
505             snew_aligned(grid->pbb, pbb_nalloc, 16);
506 #else
507             snew_aligned(grid->bb, grid->nc_nalloc*GPU_NSUBCELL, 16);
508 #endif
509         }
510
511         if (grid->bSimple)
512         {
513             if (grid->na_cj == grid->na_c)
514             {
515                 grid->bbj = grid->bb;
516             }
517             else
518             {
519                 sfree_aligned(grid->bbj);
520                 snew_aligned(grid->bbj, grid->nc_nalloc*grid->na_c/grid->na_cj, 16);
521             }
522         }
523
524         srenew(grid->flags, grid->nc_nalloc);
525     }
526
527     copy_rvec(corner0, grid->c0);
528     copy_rvec(corner1, grid->c1);
529
530     return nc_max;
531 }
532
533 /* We need to sort paricles in grid columns on z-coordinate.
534  * As particle are very often distributed homogeneously, we a sorting
535  * algorithm similar to pigeonhole sort. We multiply the z-coordinate
536  * by a factor, cast to an int and try to store in that hole. If the hole
537  * is full, we move this or another particle. A second pass is needed to make
538  * contiguous elements. SORT_GRID_OVERSIZE is the ratio of holes to particles.
539  * 4 is the optimal value for homogeneous particle distribution and allows
540  * for an O(#particles) sort up till distributions were all particles are
541  * concentrated in 1/4 of the space. No NlogN fallback is implemented,
542  * as it can be expensive to detect imhomogeneous particle distributions.
543  * SGSF is the maximum ratio of holes used, in the worst case all particles
544  * end up in the last hole and we need #particles extra holes at the end.
545  */
546 #define SORT_GRID_OVERSIZE 4
547 #define SGSF (SORT_GRID_OVERSIZE + 1)
548
549 /* Sort particle index a on coordinates x along dim.
550  * Backwards tells if we want decreasing iso increasing coordinates.
551  * h0 is the minimum of the coordinate range.
552  * invh is the 1/length of the sorting range.
553  * n_per_h (>=n) is the expected average number of particles per 1/invh
554  * sort is the sorting work array.
555  * sort should have a size of at least n_per_h*SORT_GRID_OVERSIZE + n,
556  * or easier, allocate at least n*SGSF elements.
557  */
558 static void sort_atoms(int dim, gmx_bool Backwards,
559                        int *a, int n, rvec *x,
560                        real h0, real invh, int n_per_h,
561                        int *sort)
562 {
563     int nsort, i, c;
564     int zi, zim, zi_min, zi_max;
565     int cp, tmp;
566
567     if (n <= 1)
568     {
569         /* Nothing to do */
570         return;
571     }
572
573 #ifndef NDEBUG
574     if (n > n_per_h)
575     {
576         gmx_incons("n > n_per_h");
577     }
578 #endif
579
580     /* Transform the inverse range height into the inverse hole height */
581     invh *= n_per_h*SORT_GRID_OVERSIZE;
582
583     /* Set nsort to the maximum possible number of holes used.
584      * In worst case all n elements end up in the last bin.
585      */
586     nsort = n_per_h*SORT_GRID_OVERSIZE + n;
587
588     /* Determine the index range used, so we can limit it for the second pass */
589     zi_min = INT_MAX;
590     zi_max = -1;
591
592     /* Sort the particles using a simple index sort */
593     for (i = 0; i < n; i++)
594     {
595         /* The cast takes care of float-point rounding effects below zero.
596          * This code assumes particles are less than 1/SORT_GRID_OVERSIZE
597          * times the box height out of the box.
598          */
599         zi = (int)((x[a[i]][dim] - h0)*invh);
600
601 #ifndef NDEBUG
602         /* As we can have rounding effect, we use > iso >= here */
603         if (zi < 0 || zi > n_per_h*SORT_GRID_OVERSIZE)
604         {
605             gmx_fatal(FARGS, "(int)((x[%d][%c]=%f - %f)*%f) = %d, not in 0 - %d*%d\n",
606                       a[i], 'x'+dim, x[a[i]][dim], h0, invh, zi,
607                       n_per_h, SORT_GRID_OVERSIZE);
608         }
609 #endif
610
611         /* Ideally this particle should go in sort cell zi,
612          * but that might already be in use,
613          * in that case find the first empty cell higher up
614          */
615         if (sort[zi] < 0)
616         {
617             sort[zi] = a[i];
618             zi_min   = min(zi_min, zi);
619             zi_max   = max(zi_max, zi);
620         }
621         else
622         {
623             /* We have multiple atoms in the same sorting slot.
624              * Sort on real z for minimal bounding box size.
625              * There is an extra check for identical z to ensure
626              * well-defined output order, independent of input order
627              * to ensure binary reproducibility after restarts.
628              */
629             while (sort[zi] >= 0 && ( x[a[i]][dim] >  x[sort[zi]][dim] ||
630                                       (x[a[i]][dim] == x[sort[zi]][dim] &&
631                                        a[i] > sort[zi])))
632             {
633                 zi++;
634             }
635
636             if (sort[zi] >= 0)
637             {
638                 /* Shift all elements by one slot until we find an empty slot */
639                 cp  = sort[zi];
640                 zim = zi + 1;
641                 while (sort[zim] >= 0)
642                 {
643                     tmp       = sort[zim];
644                     sort[zim] = cp;
645                     cp        = tmp;
646                     zim++;
647                 }
648                 sort[zim] = cp;
649                 zi_max    = max(zi_max, zim);
650             }
651             sort[zi] = a[i];
652             zi_max   = max(zi_max, zi);
653         }
654     }
655
656     c = 0;
657     if (!Backwards)
658     {
659         for (zi = 0; zi < nsort; zi++)
660         {
661             if (sort[zi] >= 0)
662             {
663                 a[c++]   = sort[zi];
664                 sort[zi] = -1;
665             }
666         }
667     }
668     else
669     {
670         for (zi = zi_max; zi >= zi_min; zi--)
671         {
672             if (sort[zi] >= 0)
673             {
674                 a[c++]   = sort[zi];
675                 sort[zi] = -1;
676             }
677         }
678     }
679     if (c < n)
680     {
681         gmx_incons("Lost particles while sorting");
682     }
683 }
684
685 #ifdef GMX_DOUBLE
686 #define R2F_D(x) ((float)((x) >= 0 ? ((1-GMX_FLOAT_EPS)*(x)) : ((1+GMX_FLOAT_EPS)*(x))))
687 #define R2F_U(x) ((float)((x) >= 0 ? ((1+GMX_FLOAT_EPS)*(x)) : ((1-GMX_FLOAT_EPS)*(x))))
688 #else
689 #define R2F_D(x) (x)
690 #define R2F_U(x) (x)
691 #endif
692
693 /* Coordinate order x,y,z, bb order xyz0 */
694 static void calc_bounding_box(int na, int stride, const real *x, nbnxn_bb_t *bb)
695 {
696     int  i, j;
697     real xl, xh, yl, yh, zl, zh;
698
699     i  = 0;
700     xl = x[i+XX];
701     xh = x[i+XX];
702     yl = x[i+YY];
703     yh = x[i+YY];
704     zl = x[i+ZZ];
705     zh = x[i+ZZ];
706     i += stride;
707     for (j = 1; j < na; j++)
708     {
709         xl = min(xl, x[i+XX]);
710         xh = max(xh, x[i+XX]);
711         yl = min(yl, x[i+YY]);
712         yh = max(yh, x[i+YY]);
713         zl = min(zl, x[i+ZZ]);
714         zh = max(zh, x[i+ZZ]);
715         i += stride;
716     }
717     /* Note: possible double to float conversion here */
718     bb->lower[BB_X] = R2F_D(xl);
719     bb->lower[BB_Y] = R2F_D(yl);
720     bb->lower[BB_Z] = R2F_D(zl);
721     bb->upper[BB_X] = R2F_U(xh);
722     bb->upper[BB_Y] = R2F_U(yh);
723     bb->upper[BB_Z] = R2F_U(zh);
724 }
725
726 /* Packed coordinates, bb order xyz0 */
727 static void calc_bounding_box_x_x4(int na, const real *x, nbnxn_bb_t *bb)
728 {
729     int  j;
730     real xl, xh, yl, yh, zl, zh;
731
732     xl = x[XX*PACK_X4];
733     xh = x[XX*PACK_X4];
734     yl = x[YY*PACK_X4];
735     yh = x[YY*PACK_X4];
736     zl = x[ZZ*PACK_X4];
737     zh = x[ZZ*PACK_X4];
738     for (j = 1; j < na; j++)
739     {
740         xl = min(xl, x[j+XX*PACK_X4]);
741         xh = max(xh, x[j+XX*PACK_X4]);
742         yl = min(yl, x[j+YY*PACK_X4]);
743         yh = max(yh, x[j+YY*PACK_X4]);
744         zl = min(zl, x[j+ZZ*PACK_X4]);
745         zh = max(zh, x[j+ZZ*PACK_X4]);
746     }
747     /* Note: possible double to float conversion here */
748     bb->lower[BB_X] = R2F_D(xl);
749     bb->lower[BB_Y] = R2F_D(yl);
750     bb->lower[BB_Z] = R2F_D(zl);
751     bb->upper[BB_X] = R2F_U(xh);
752     bb->upper[BB_Y] = R2F_U(yh);
753     bb->upper[BB_Z] = R2F_U(zh);
754 }
755
756 /* Packed coordinates, bb order xyz0 */
757 static void calc_bounding_box_x_x8(int na, const real *x, nbnxn_bb_t *bb)
758 {
759     int  j;
760     real xl, xh, yl, yh, zl, zh;
761
762     xl = x[XX*PACK_X8];
763     xh = x[XX*PACK_X8];
764     yl = x[YY*PACK_X8];
765     yh = x[YY*PACK_X8];
766     zl = x[ZZ*PACK_X8];
767     zh = x[ZZ*PACK_X8];
768     for (j = 1; j < na; j++)
769     {
770         xl = min(xl, x[j+XX*PACK_X8]);
771         xh = max(xh, x[j+XX*PACK_X8]);
772         yl = min(yl, x[j+YY*PACK_X8]);
773         yh = max(yh, x[j+YY*PACK_X8]);
774         zl = min(zl, x[j+ZZ*PACK_X8]);
775         zh = max(zh, x[j+ZZ*PACK_X8]);
776     }
777     /* Note: possible double to float conversion here */
778     bb->lower[BB_X] = R2F_D(xl);
779     bb->lower[BB_Y] = R2F_D(yl);
780     bb->lower[BB_Z] = R2F_D(zl);
781     bb->upper[BB_X] = R2F_U(xh);
782     bb->upper[BB_Y] = R2F_U(yh);
783     bb->upper[BB_Z] = R2F_U(zh);
784 }
785
786 /* Packed coordinates, bb order xyz0 */
787 static void calc_bounding_box_x_x4_halves(int na, const real *x,
788                                           nbnxn_bb_t *bb, nbnxn_bb_t *bbj)
789 {
790     calc_bounding_box_x_x4(min(na, 2), x, bbj);
791
792     if (na > 2)
793     {
794         calc_bounding_box_x_x4(min(na-2, 2), x+(PACK_X4>>1), bbj+1);
795     }
796     else
797     {
798         /* Set the "empty" bounding box to the same as the first one,
799          * so we don't need to treat special cases in the rest of the code.
800          */
801 #ifdef NBNXN_SEARCH_BB_SSE
802         _mm_store_ps(&bbj[1].lower[0], _mm_load_ps(&bbj[0].lower[0]));
803         _mm_store_ps(&bbj[1].upper[0], _mm_load_ps(&bbj[0].upper[0]));
804 #else
805         bbj[1] = bbj[0];
806 #endif
807     }
808
809 #ifdef NBNXN_SEARCH_BB_SSE
810     _mm_store_ps(&bb->lower[0], _mm_min_ps(_mm_load_ps(&bbj[0].lower[0]),
811                                            _mm_load_ps(&bbj[1].lower[0])));
812     _mm_store_ps(&bb->upper[0], _mm_max_ps(_mm_load_ps(&bbj[0].upper[0]),
813                                            _mm_load_ps(&bbj[1].upper[0])));
814 #else
815     {
816         int i;
817
818         for (i = 0; i < NNBSBB_C; i++)
819         {
820             bb->lower[i] = min(bbj[0].lower[i], bbj[1].lower[i]);
821             bb->upper[i] = max(bbj[0].upper[i], bbj[1].upper[i]);
822         }
823     }
824 #endif
825 }
826
827 #ifdef NBNXN_SEARCH_BB_SSE
828
829 /* Coordinate order xyz, bb order xxxxyyyyzzzz */
830 static void calc_bounding_box_xxxx(int na, int stride, const real *x, float *bb)
831 {
832     int  i, j;
833     real xl, xh, yl, yh, zl, zh;
834
835     i  = 0;
836     xl = x[i+XX];
837     xh = x[i+XX];
838     yl = x[i+YY];
839     yh = x[i+YY];
840     zl = x[i+ZZ];
841     zh = x[i+ZZ];
842     i += stride;
843     for (j = 1; j < na; j++)
844     {
845         xl = min(xl, x[i+XX]);
846         xh = max(xh, x[i+XX]);
847         yl = min(yl, x[i+YY]);
848         yh = max(yh, x[i+YY]);
849         zl = min(zl, x[i+ZZ]);
850         zh = max(zh, x[i+ZZ]);
851         i += stride;
852     }
853     /* Note: possible double to float conversion here */
854     bb[0*STRIDE_PBB] = R2F_D(xl);
855     bb[1*STRIDE_PBB] = R2F_D(yl);
856     bb[2*STRIDE_PBB] = R2F_D(zl);
857     bb[3*STRIDE_PBB] = R2F_U(xh);
858     bb[4*STRIDE_PBB] = R2F_U(yh);
859     bb[5*STRIDE_PBB] = R2F_U(zh);
860 }
861
862 #endif /* NBNXN_SEARCH_BB_SSE */
863
864 #ifdef NBNXN_SEARCH_SSE_SINGLE
865
866 /* Coordinate order xyz?, bb order xyz0 */
867 static void calc_bounding_box_sse(int na, const float *x, nbnxn_bb_t *bb)
868 {
869     __m128 bb_0_SSE, bb_1_SSE;
870     __m128 x_SSE;
871
872     int    i;
873
874     bb_0_SSE = _mm_load_ps(x);
875     bb_1_SSE = bb_0_SSE;
876
877     for (i = 1; i < na; i++)
878     {
879         x_SSE    = _mm_load_ps(x+i*NNBSBB_C);
880         bb_0_SSE = _mm_min_ps(bb_0_SSE, x_SSE);
881         bb_1_SSE = _mm_max_ps(bb_1_SSE, x_SSE);
882     }
883
884     _mm_store_ps(&bb->lower[0], bb_0_SSE);
885     _mm_store_ps(&bb->upper[0], bb_1_SSE);
886 }
887
888 /* Coordinate order xyz?, bb order xxxxyyyyzzzz */
889 static void calc_bounding_box_xxxx_sse(int na, const float *x,
890                                        nbnxn_bb_t *bb_work_aligned,
891                                        real *bb)
892 {
893     calc_bounding_box_sse(na, x, bb_work_aligned);
894
895     bb[0*STRIDE_PBB] = bb_work_aligned->lower[BB_X];
896     bb[1*STRIDE_PBB] = bb_work_aligned->lower[BB_Y];
897     bb[2*STRIDE_PBB] = bb_work_aligned->lower[BB_Z];
898     bb[3*STRIDE_PBB] = bb_work_aligned->upper[BB_X];
899     bb[4*STRIDE_PBB] = bb_work_aligned->upper[BB_Y];
900     bb[5*STRIDE_PBB] = bb_work_aligned->upper[BB_Z];
901 }
902
903 #endif /* NBNXN_SEARCH_SSE_SINGLE */
904
905
906 /* Combines pairs of consecutive bounding boxes */
907 static void combine_bounding_box_pairs(nbnxn_grid_t *grid, const nbnxn_bb_t *bb)
908 {
909     int    i, j, sc2, nc2, c2;
910
911     for (i = 0; i < grid->ncx*grid->ncy; i++)
912     {
913         /* Starting bb in a column is expected to be 2-aligned */
914         sc2 = grid->cxy_ind[i]>>1;
915         /* For odd numbers skip the last bb here */
916         nc2 = (grid->cxy_na[i]+3)>>(2+1);
917         for (c2 = sc2; c2 < sc2+nc2; c2++)
918         {
919 #ifdef NBNXN_SEARCH_BB_SSE
920             __m128 min_SSE, max_SSE;
921
922             min_SSE = _mm_min_ps(_mm_load_ps(&bb[c2*2+0].lower[0]),
923                                  _mm_load_ps(&bb[c2*2+1].lower[0]));
924             max_SSE = _mm_max_ps(_mm_load_ps(&bb[c2*2+0].upper[0]),
925                                  _mm_load_ps(&bb[c2*2+1].upper[0]));
926             _mm_store_ps(&grid->bbj[c2].lower[0], min_SSE);
927             _mm_store_ps(&grid->bbj[c2].upper[0], max_SSE);
928 #else
929             for (j = 0; j < NNBSBB_C; j++)
930             {
931                 grid->bbj[c2].lower[j] = min(bb[c2*2].lower[j],
932                                              bb[c2*2].lower[j]);
933                 grid->bbj[c2].upper[j] = max(bb[c2*2].upper[j],
934                                              bb[c2*2].upper[j]);
935             }
936 #endif
937         }
938         if (((grid->cxy_na[i]+3)>>2) & 1)
939         {
940             /* Copy the last bb for odd bb count in this column */
941             for (j = 0; j < NNBSBB_C; j++)
942             {
943                 grid->bbj[c2].lower[j] = bb[c2*2].lower[j];
944                 grid->bbj[c2].upper[j] = bb[c2*2].upper[j];
945             }
946         }
947     }
948 }
949
950
951 /* Prints the average bb size, used for debug output */
952 static void print_bbsizes_simple(FILE                *fp,
953                                  const nbnxn_search_t nbs,
954                                  const nbnxn_grid_t  *grid)
955 {
956     int  c, d;
957     dvec ba;
958
959     clear_dvec(ba);
960     for (c = 0; c < grid->nc; c++)
961     {
962         for (d = 0; d < DIM; d++)
963         {
964             ba[d] += grid->bb[c].upper[d] - grid->bb[c].lower[d];
965         }
966     }
967     dsvmul(1.0/grid->nc, ba, ba);
968
969     fprintf(fp, "ns bb: %4.2f %4.2f %4.2f  %4.2f %4.2f %4.2f rel %4.2f %4.2f %4.2f\n",
970             nbs->box[XX][XX]/grid->ncx,
971             nbs->box[YY][YY]/grid->ncy,
972             nbs->box[ZZ][ZZ]*grid->ncx*grid->ncy/grid->nc,
973             ba[XX], ba[YY], ba[ZZ],
974             ba[XX]*grid->ncx/nbs->box[XX][XX],
975             ba[YY]*grid->ncy/nbs->box[YY][YY],
976             ba[ZZ]*grid->nc/(grid->ncx*grid->ncy*nbs->box[ZZ][ZZ]));
977 }
978
979 /* Prints the average bb size, used for debug output */
980 static void print_bbsizes_supersub(FILE                *fp,
981                                    const nbnxn_search_t nbs,
982                                    const nbnxn_grid_t  *grid)
983 {
984     int  ns, c, s;
985     dvec ba;
986
987     clear_dvec(ba);
988     ns = 0;
989     for (c = 0; c < grid->nc; c++)
990     {
991 #ifdef NBNXN_BBXXXX
992         for (s = 0; s < grid->nsubc[c]; s += STRIDE_PBB)
993         {
994             int cs_w, i, d;
995
996             cs_w = (c*GPU_NSUBCELL + s)/STRIDE_PBB;
997             for (i = 0; i < STRIDE_PBB; i++)
998             {
999                 for (d = 0; d < DIM; d++)
1000                 {
1001                     ba[d] +=
1002                         grid->pbb[cs_w*NNBSBB_XXXX+(DIM+d)*STRIDE_PBB+i] -
1003                         grid->pbb[cs_w*NNBSBB_XXXX+     d *STRIDE_PBB+i];
1004                 }
1005             }
1006         }
1007 #else
1008         for (s = 0; s < grid->nsubc[c]; s++)
1009         {
1010             int cs, d;
1011
1012             cs = c*GPU_NSUBCELL + s;
1013             for (d = 0; d < DIM; d++)
1014             {
1015                 ba[d] += grid->bb[cs].upper[d] - grid->bb[cs].lower[d];
1016             }
1017         }
1018 #endif
1019         ns += grid->nsubc[c];
1020     }
1021     dsvmul(1.0/ns, ba, ba);
1022
1023     fprintf(fp, "ns bb: %4.2f %4.2f %4.2f  %4.2f %4.2f %4.2f rel %4.2f %4.2f %4.2f\n",
1024             nbs->box[XX][XX]/(grid->ncx*GPU_NSUBCELL_X),
1025             nbs->box[YY][YY]/(grid->ncy*GPU_NSUBCELL_Y),
1026             nbs->box[ZZ][ZZ]*grid->ncx*grid->ncy/(grid->nc*GPU_NSUBCELL_Z),
1027             ba[XX], ba[YY], ba[ZZ],
1028             ba[XX]*grid->ncx*GPU_NSUBCELL_X/nbs->box[XX][XX],
1029             ba[YY]*grid->ncy*GPU_NSUBCELL_Y/nbs->box[YY][YY],
1030             ba[ZZ]*grid->nc*GPU_NSUBCELL_Z/(grid->ncx*grid->ncy*nbs->box[ZZ][ZZ]));
1031 }
1032
1033 /* Potentially sorts atoms on LJ coefficients !=0 and ==0.
1034  * Also sets interaction flags.
1035  */
1036 void sort_on_lj(int na_c,
1037                 int a0, int a1, const int *atinfo,
1038                 int *order,
1039                 int *flags)
1040 {
1041     int      subc, s, a, n1, n2, a_lj_max, i, j;
1042     int      sort1[NBNXN_NA_SC_MAX/GPU_NSUBCELL];
1043     int      sort2[NBNXN_NA_SC_MAX/GPU_NSUBCELL];
1044     gmx_bool haveQ;
1045
1046     *flags = 0;
1047
1048     subc = 0;
1049     for (s = a0; s < a1; s += na_c)
1050     {
1051         /* Make lists for this (sub-)cell on atoms with and without LJ */
1052         n1       = 0;
1053         n2       = 0;
1054         haveQ    = FALSE;
1055         a_lj_max = -1;
1056         for (a = s; a < min(s+na_c, a1); a++)
1057         {
1058             haveQ = haveQ || GET_CGINFO_HAS_Q(atinfo[order[a]]);
1059
1060             if (GET_CGINFO_HAS_VDW(atinfo[order[a]]))
1061             {
1062                 sort1[n1++] = order[a];
1063                 a_lj_max    = a;
1064             }
1065             else
1066             {
1067                 sort2[n2++] = order[a];
1068             }
1069         }
1070
1071         /* If we don't have atom with LJ, there's nothing to sort */
1072         if (n1 > 0)
1073         {
1074             *flags |= NBNXN_CI_DO_LJ(subc);
1075
1076             if (2*n1 <= na_c)
1077             {
1078                 /* Only sort when strictly necessary. Ordering particles
1079                  * Ordering particles can lead to less accurate summation
1080                  * due to rounding, both for LJ and Coulomb interactions.
1081                  */
1082                 if (2*(a_lj_max - s) >= na_c)
1083                 {
1084                     for (i = 0; i < n1; i++)
1085                     {
1086                         order[a0+i] = sort1[i];
1087                     }
1088                     for (j = 0; j < n2; j++)
1089                     {
1090                         order[a0+n1+j] = sort2[j];
1091                     }
1092                 }
1093
1094                 *flags |= NBNXN_CI_HALF_LJ(subc);
1095             }
1096         }
1097         if (haveQ)
1098         {
1099             *flags |= NBNXN_CI_DO_COUL(subc);
1100         }
1101         subc++;
1102     }
1103 }
1104
1105 /* Fill a pair search cell with atoms.
1106  * Potentially sorts atoms and sets the interaction flags.
1107  */
1108 void fill_cell(const nbnxn_search_t nbs,
1109                nbnxn_grid_t *grid,
1110                nbnxn_atomdata_t *nbat,
1111                int a0, int a1,
1112                const int *atinfo,
1113                rvec *x,
1114                int sx, int sy, int sz,
1115                nbnxn_bb_t *bb_work_aligned)
1116 {
1117     int        na, a;
1118     size_t     offset;
1119     nbnxn_bb_t *bb_ptr;
1120 #ifdef NBNXN_BBXXXX
1121     float      *pbb_ptr;
1122 #endif
1123
1124     na = a1 - a0;
1125
1126     if (grid->bSimple)
1127     {
1128         sort_on_lj(grid->na_c, a0, a1, atinfo, nbs->a,
1129                    grid->flags+(a0>>grid->na_c_2log)-grid->cell0);
1130     }
1131
1132     /* Now we have sorted the atoms, set the cell indices */
1133     for (a = a0; a < a1; a++)
1134     {
1135         nbs->cell[nbs->a[a]] = a;
1136     }
1137
1138     copy_rvec_to_nbat_real(nbs->a+a0, a1-a0, grid->na_c, x,
1139                            nbat->XFormat, nbat->x, a0,
1140                            sx, sy, sz);
1141
1142     if (nbat->XFormat == nbatX4)
1143     {
1144         /* Store the bounding boxes as xyz.xyz. */
1145         offset = (a0 - grid->cell0*grid->na_sc) >> grid->na_c_2log;
1146         bb_ptr = grid->bb + offset;
1147
1148 #if defined GMX_NBNXN_SIMD && GMX_SIMD_WIDTH_HERE == 2
1149         if (2*grid->na_cj == grid->na_c)
1150         {
1151             calc_bounding_box_x_x4_halves(na, nbat->x+X4_IND_A(a0), bb_ptr,
1152                                           grid->bbj+offset*2);
1153         }
1154         else
1155 #endif
1156         {
1157             calc_bounding_box_x_x4(na, nbat->x+X4_IND_A(a0), bb_ptr);
1158         }
1159     }
1160     else if (nbat->XFormat == nbatX8)
1161     {
1162         /* Store the bounding boxes as xyz.xyz. */
1163         offset = (a0 - grid->cell0*grid->na_sc) >> grid->na_c_2log;
1164         bb_ptr = grid->bb + offset;
1165
1166         calc_bounding_box_x_x8(na, nbat->x+X8_IND_A(a0), bb_ptr);
1167     }
1168 #ifdef NBNXN_BBXXXX
1169     else if (!grid->bSimple)
1170     {
1171         /* Store the bounding boxes in a format convenient
1172          * for SSE calculations: xxxxyyyyzzzz...
1173          */
1174         pbb_ptr =
1175             grid->pbb +
1176             ((a0-grid->cell0*grid->na_sc)>>(grid->na_c_2log+STRIDE_PBB_2LOG))*NNBSBB_XXXX +
1177             (((a0-grid->cell0*grid->na_sc)>>grid->na_c_2log) & (STRIDE_PBB-1));
1178
1179 #ifdef NBNXN_SEARCH_SSE_SINGLE
1180         if (nbat->XFormat == nbatXYZQ)
1181         {
1182             calc_bounding_box_xxxx_sse(na, nbat->x+a0*nbat->xstride,
1183                                        bb_work_aligned, pbb_ptr);
1184         }
1185         else
1186 #endif
1187         {
1188             calc_bounding_box_xxxx(na, nbat->xstride, nbat->x+a0*nbat->xstride,
1189                                    pbb_ptr);
1190         }
1191         if (gmx_debug_at)
1192         {
1193             fprintf(debug, "%2d %2d %2d bb %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f\n",
1194                     sx, sy, sz,
1195                     pbb_ptr[0*STRIDE_PBB], pbb_ptr[3*STRIDE_PBB],
1196                     pbb_ptr[1*STRIDE_PBB], pbb_ptr[4*STRIDE_PBB],
1197                     pbb_ptr[2*STRIDE_PBB], pbb_ptr[5*STRIDE_PBB]);
1198         }
1199     }
1200 #endif
1201     else
1202     {
1203         /* Store the bounding boxes as xyz.xyz. */
1204         bb_ptr = grid->bb+((a0-grid->cell0*grid->na_sc)>>grid->na_c_2log);
1205
1206         calc_bounding_box(na, nbat->xstride, nbat->x+a0*nbat->xstride,
1207                           bb_ptr);
1208
1209         if (gmx_debug_at)
1210         {
1211             int bbo;
1212             bbo = (a0 - grid->cell0*grid->na_sc)/grid->na_c;
1213             fprintf(debug, "%2d %2d %2d bb %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f\n",
1214                     sx, sy, sz,
1215                     grid->bb[bbo].lower[BB_X],
1216                     grid->bb[bbo].lower[BB_Y],
1217                     grid->bb[bbo].lower[BB_Z],
1218                     grid->bb[bbo].upper[BB_X],
1219                     grid->bb[bbo].upper[BB_Y],
1220                     grid->bb[bbo].upper[BB_Z]);
1221         }
1222     }
1223 }
1224
1225 /* Spatially sort the atoms within one grid column */
1226 static void sort_columns_simple(const nbnxn_search_t nbs,
1227                                 int dd_zone,
1228                                 nbnxn_grid_t *grid,
1229                                 int a0, int a1,
1230                                 const int *atinfo,
1231                                 rvec *x,
1232                                 nbnxn_atomdata_t *nbat,
1233                                 int cxy_start, int cxy_end,
1234                                 int *sort_work)
1235 {
1236     int  cxy;
1237     int  cx, cy, cz, ncz, cfilled, c;
1238     int  na, ash, ind, a;
1239     int  na_c, ash_c;
1240
1241     if (debug)
1242     {
1243         fprintf(debug, "cell0 %d sorting columns %d - %d, atoms %d - %d\n",
1244                 grid->cell0, cxy_start, cxy_end, a0, a1);
1245     }
1246
1247     /* Sort the atoms within each x,y column in 3 dimensions */
1248     for (cxy = cxy_start; cxy < cxy_end; cxy++)
1249     {
1250         cx = cxy/grid->ncy;
1251         cy = cxy - cx*grid->ncy;
1252
1253         na  = grid->cxy_na[cxy];
1254         ncz = grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy];
1255         ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1256
1257         /* Sort the atoms within each x,y column on z coordinate */
1258         sort_atoms(ZZ, FALSE,
1259                    nbs->a+ash, na, x,
1260                    grid->c0[ZZ],
1261                    1.0/nbs->box[ZZ][ZZ], ncz*grid->na_sc,
1262                    sort_work);
1263
1264         /* Fill the ncz cells in this column */
1265         cfilled = grid->cxy_ind[cxy];
1266         for (cz = 0; cz < ncz; cz++)
1267         {
1268             c  = grid->cxy_ind[cxy] + cz;
1269
1270             ash_c = ash + cz*grid->na_sc;
1271             na_c  = min(grid->na_sc, na-(ash_c-ash));
1272
1273             fill_cell(nbs, grid, nbat,
1274                       ash_c, ash_c+na_c, atinfo, x,
1275                       grid->na_sc*cx + (dd_zone >> 2),
1276                       grid->na_sc*cy + (dd_zone & 3),
1277                       grid->na_sc*cz,
1278                       NULL);
1279
1280             /* This copy to bbcz is not really necessary.
1281              * But it allows to use the same grid search code
1282              * for the simple and supersub cell setups.
1283              */
1284             if (na_c > 0)
1285             {
1286                 cfilled = c;
1287             }
1288             grid->bbcz[c*NNBSBB_D  ] = grid->bb[cfilled].lower[BB_Z];
1289             grid->bbcz[c*NNBSBB_D+1] = grid->bb[cfilled].upper[BB_Z];
1290         }
1291
1292         /* Set the unused atom indices to -1 */
1293         for (ind = na; ind < ncz*grid->na_sc; ind++)
1294         {
1295             nbs->a[ash+ind] = -1;
1296         }
1297     }
1298 }
1299
1300 /* Spatially sort the atoms within one grid column */
1301 static void sort_columns_supersub(const nbnxn_search_t nbs,
1302                                   int dd_zone,
1303                                   nbnxn_grid_t *grid,
1304                                   int a0, int a1,
1305                                   const int *atinfo,
1306                                   rvec *x,
1307                                   nbnxn_atomdata_t *nbat,
1308                                   int cxy_start, int cxy_end,
1309                                   int *sort_work)
1310 {
1311     int  cxy;
1312     int  cx, cy, cz = -1, c = -1, ncz;
1313     int  na, ash, na_c, ind, a;
1314     int  subdiv_z, sub_z, na_z, ash_z;
1315     int  subdiv_y, sub_y, na_y, ash_y;
1316     int  subdiv_x, sub_x, na_x, ash_x;
1317
1318     /* cppcheck-suppress unassignedVariable */
1319     nbnxn_bb_t bb_work_array[2], *bb_work_aligned;
1320
1321     bb_work_aligned = (nbnxn_bb_t *)(((size_t)(bb_work_array+1)) & (~((size_t)15)));
1322
1323     if (debug)
1324     {
1325         fprintf(debug, "cell0 %d sorting columns %d - %d, atoms %d - %d\n",
1326                 grid->cell0, cxy_start, cxy_end, a0, a1);
1327     }
1328
1329     subdiv_x = grid->na_c;
1330     subdiv_y = GPU_NSUBCELL_X*subdiv_x;
1331     subdiv_z = GPU_NSUBCELL_Y*subdiv_y;
1332
1333     /* Sort the atoms within each x,y column in 3 dimensions */
1334     for (cxy = cxy_start; cxy < cxy_end; cxy++)
1335     {
1336         cx = cxy/grid->ncy;
1337         cy = cxy - cx*grid->ncy;
1338
1339         na  = grid->cxy_na[cxy];
1340         ncz = grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy];
1341         ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1342
1343         /* Sort the atoms within each x,y column on z coordinate */
1344         sort_atoms(ZZ, FALSE,
1345                    nbs->a+ash, na, x,
1346                    grid->c0[ZZ],
1347                    1.0/nbs->box[ZZ][ZZ], ncz*grid->na_sc,
1348                    sort_work);
1349
1350         /* This loop goes over the supercells and subcells along z at once */
1351         for (sub_z = 0; sub_z < ncz*GPU_NSUBCELL_Z; sub_z++)
1352         {
1353             ash_z = ash + sub_z*subdiv_z;
1354             na_z  = min(subdiv_z, na-(ash_z-ash));
1355
1356             /* We have already sorted on z */
1357
1358             if (sub_z % GPU_NSUBCELL_Z == 0)
1359             {
1360                 cz = sub_z/GPU_NSUBCELL_Z;
1361                 c  = grid->cxy_ind[cxy] + cz;
1362
1363                 /* The number of atoms in this supercell */
1364                 na_c = min(grid->na_sc, na-(ash_z-ash));
1365
1366                 grid->nsubc[c] = min(GPU_NSUBCELL, (na_c+grid->na_c-1)/grid->na_c);
1367
1368                 /* Store the z-boundaries of the super cell */
1369                 grid->bbcz[c*NNBSBB_D  ] = x[nbs->a[ash_z]][ZZ];
1370                 grid->bbcz[c*NNBSBB_D+1] = x[nbs->a[ash_z+na_c-1]][ZZ];
1371             }
1372
1373 #if GPU_NSUBCELL_Y > 1
1374             /* Sort the atoms along y */
1375             sort_atoms(YY, (sub_z & 1),
1376                        nbs->a+ash_z, na_z, x,
1377                        grid->c0[YY]+cy*grid->sy,
1378                        grid->inv_sy, subdiv_z,
1379                        sort_work);
1380 #endif
1381
1382             for (sub_y = 0; sub_y < GPU_NSUBCELL_Y; sub_y++)
1383             {
1384                 ash_y = ash_z + sub_y*subdiv_y;
1385                 na_y  = min(subdiv_y, na-(ash_y-ash));
1386
1387 #if GPU_NSUBCELL_X > 1
1388                 /* Sort the atoms along x */
1389                 sort_atoms(XX, ((cz*GPU_NSUBCELL_Y + sub_y) & 1),
1390                            nbs->a+ash_y, na_y, x,
1391                            grid->c0[XX]+cx*grid->sx,
1392                            grid->inv_sx, subdiv_y,
1393                            sort_work);
1394 #endif
1395
1396                 for (sub_x = 0; sub_x < GPU_NSUBCELL_X; sub_x++)
1397                 {
1398                     ash_x = ash_y + sub_x*subdiv_x;
1399                     na_x  = min(subdiv_x, na-(ash_x-ash));
1400
1401                     fill_cell(nbs, grid, nbat,
1402                               ash_x, ash_x+na_x, atinfo, x,
1403                               grid->na_c*(cx*GPU_NSUBCELL_X+sub_x) + (dd_zone >> 2),
1404                               grid->na_c*(cy*GPU_NSUBCELL_Y+sub_y) + (dd_zone & 3),
1405                               grid->na_c*sub_z,
1406                               bb_work_aligned);
1407                 }
1408             }
1409         }
1410
1411         /* Set the unused atom indices to -1 */
1412         for (ind = na; ind < ncz*grid->na_sc; ind++)
1413         {
1414             nbs->a[ash+ind] = -1;
1415         }
1416     }
1417 }
1418
1419 /* Determine in which grid column atoms should go */
1420 static void calc_column_indices(nbnxn_grid_t *grid,
1421                                 int a0, int a1,
1422                                 rvec *x,
1423                                 int dd_zone, const int *move,
1424                                 int thread, int nthread,
1425                                 int *cell,
1426                                 int *cxy_na)
1427 {
1428     int  n0, n1, i;
1429     int  cx, cy;
1430
1431     /* We add one extra cell for particles which moved during DD */
1432     for (i = 0; i < grid->ncx*grid->ncy+1; i++)
1433     {
1434         cxy_na[i] = 0;
1435     }
1436
1437     n0 = a0 + (int)((thread+0)*(a1 - a0))/nthread;
1438     n1 = a0 + (int)((thread+1)*(a1 - a0))/nthread;
1439     if (dd_zone == 0)
1440     {
1441         /* Home zone */
1442         for (i = n0; i < n1; i++)
1443         {
1444             if (move == NULL || move[i] >= 0)
1445             {
1446                 /* We need to be careful with rounding,
1447                  * particles might be a few bits outside the local zone.
1448                  * The int cast takes care of the lower bound,
1449                  * we will explicitly take care of the upper bound.
1450                  */
1451                 cx = (int)((x[i][XX] - grid->c0[XX])*grid->inv_sx);
1452                 cy = (int)((x[i][YY] - grid->c0[YY])*grid->inv_sy);
1453
1454 #ifndef NDEBUG
1455                 if (cx < 0 || cx > grid->ncx ||
1456                     cy < 0 || cy > grid->ncy)
1457                 {
1458                     gmx_fatal(FARGS,
1459                               "grid cell cx %d cy %d out of range (max %d %d)\n"
1460                               "atom %f %f %f, grid->c0 %f %f",
1461                               cx, cy, grid->ncx, grid->ncy,
1462                               x[i][XX], x[i][YY], x[i][ZZ], grid->c0[XX], grid->c0[YY]);
1463                 }
1464 #endif
1465                 /* Take care of potential rouding issues */
1466                 cx = min(cx, grid->ncx - 1);
1467                 cy = min(cy, grid->ncy - 1);
1468
1469                 /* For the moment cell will contain only the, grid local,
1470                  * x and y indices, not z.
1471                  */
1472                 cell[i] = cx*grid->ncy + cy;
1473             }
1474             else
1475             {
1476                 /* Put this moved particle after the end of the grid,
1477                  * so we can process it later without using conditionals.
1478                  */
1479                 cell[i] = grid->ncx*grid->ncy;
1480             }
1481
1482             cxy_na[cell[i]]++;
1483         }
1484     }
1485     else
1486     {
1487         /* Non-home zone */
1488         for (i = n0; i < n1; i++)
1489         {
1490             cx = (int)((x[i][XX] - grid->c0[XX])*grid->inv_sx);
1491             cy = (int)((x[i][YY] - grid->c0[YY])*grid->inv_sy);
1492
1493             /* For non-home zones there could be particles outside
1494              * the non-bonded cut-off range, which have been communicated
1495              * for bonded interactions only. For the result it doesn't
1496              * matter where these end up on the grid. For performance
1497              * we put them in an extra row at the border.
1498              */
1499             cx = max(cx, 0);
1500             cx = min(cx, grid->ncx - 1);
1501             cy = max(cy, 0);
1502             cy = min(cy, grid->ncy - 1);
1503
1504             /* For the moment cell will contain only the, grid local,
1505              * x and y indices, not z.
1506              */
1507             cell[i] = cx*grid->ncy + cy;
1508
1509             cxy_na[cell[i]]++;
1510         }
1511     }
1512 }
1513
1514 /* Determine in which grid cells the atoms should go */
1515 static void calc_cell_indices(const nbnxn_search_t nbs,
1516                               int dd_zone,
1517                               nbnxn_grid_t *grid,
1518                               int a0, int a1,
1519                               const int *atinfo,
1520                               rvec *x,
1521                               const int *move,
1522                               nbnxn_atomdata_t *nbat)
1523 {
1524     int   n0, n1, i;
1525     int   cx, cy, cxy, ncz_max, ncz;
1526     int   nthread, thread;
1527     int  *cxy_na, cxy_na_i;
1528
1529     nthread = gmx_omp_nthreads_get(emntPairsearch);
1530
1531 #pragma omp parallel for num_threads(nthread) schedule(static)
1532     for (thread = 0; thread < nthread; thread++)
1533     {
1534         calc_column_indices(grid, a0, a1, x, dd_zone, move, thread, nthread,
1535                             nbs->cell, nbs->work[thread].cxy_na);
1536     }
1537
1538     /* Make the cell index as a function of x and y */
1539     ncz_max          = 0;
1540     ncz              = 0;
1541     grid->cxy_ind[0] = 0;
1542     for (i = 0; i < grid->ncx*grid->ncy+1; i++)
1543     {
1544         /* We set ncz_max at the beginning of the loop iso at the end
1545          * to skip i=grid->ncx*grid->ncy which are moved particles
1546          * that do not need to be ordered on the grid.
1547          */
1548         if (ncz > ncz_max)
1549         {
1550             ncz_max = ncz;
1551         }
1552         cxy_na_i = nbs->work[0].cxy_na[i];
1553         for (thread = 1; thread < nthread; thread++)
1554         {
1555             cxy_na_i += nbs->work[thread].cxy_na[i];
1556         }
1557         ncz = (cxy_na_i + grid->na_sc - 1)/grid->na_sc;
1558         if (nbat->XFormat == nbatX8)
1559         {
1560             /* Make the number of cell a multiple of 2 */
1561             ncz = (ncz + 1) & ~1;
1562         }
1563         grid->cxy_ind[i+1] = grid->cxy_ind[i] + ncz;
1564         /* Clear cxy_na, so we can reuse the array below */
1565         grid->cxy_na[i] = 0;
1566     }
1567     grid->nc = grid->cxy_ind[grid->ncx*grid->ncy] - grid->cxy_ind[0];
1568
1569     nbat->natoms = (grid->cell0 + grid->nc)*grid->na_sc;
1570
1571     if (debug)
1572     {
1573         fprintf(debug, "ns na_sc %d na_c %d super-cells: %d x %d y %d z %.1f maxz %d\n",
1574                 grid->na_sc, grid->na_c, grid->nc,
1575                 grid->ncx, grid->ncy, grid->nc/((double)(grid->ncx*grid->ncy)),
1576                 ncz_max);
1577         if (gmx_debug_at)
1578         {
1579             i = 0;
1580             for (cy = 0; cy < grid->ncy; cy++)
1581             {
1582                 for (cx = 0; cx < grid->ncx; cx++)
1583                 {
1584                     fprintf(debug, " %2d", grid->cxy_ind[i+1]-grid->cxy_ind[i]);
1585                     i++;
1586                 }
1587                 fprintf(debug, "\n");
1588             }
1589         }
1590     }
1591
1592     /* Make sure the work array for sorting is large enough */
1593     if (ncz_max*grid->na_sc*SGSF > nbs->work[0].sort_work_nalloc)
1594     {
1595         for (thread = 0; thread < nbs->nthread_max; thread++)
1596         {
1597             nbs->work[thread].sort_work_nalloc =
1598                 over_alloc_large(ncz_max*grid->na_sc*SGSF);
1599             srenew(nbs->work[thread].sort_work,
1600                    nbs->work[thread].sort_work_nalloc);
1601             /* When not in use, all elements should be -1 */
1602             for (i = 0; i < nbs->work[thread].sort_work_nalloc; i++)
1603             {
1604                 nbs->work[thread].sort_work[i] = -1;
1605             }
1606         }
1607     }
1608
1609     /* Now we know the dimensions we can fill the grid.
1610      * This is the first, unsorted fill. We sort the columns after this.
1611      */
1612     for (i = a0; i < a1; i++)
1613     {
1614         /* At this point nbs->cell contains the local grid x,y indices */
1615         cxy = nbs->cell[i];
1616         nbs->a[(grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc + grid->cxy_na[cxy]++] = i;
1617     }
1618
1619     if (dd_zone == 0)
1620     {
1621         /* Set the cell indices for the moved particles */
1622         n0 = grid->nc*grid->na_sc;
1623         n1 = grid->nc*grid->na_sc+grid->cxy_na[grid->ncx*grid->ncy];
1624         if (dd_zone == 0)
1625         {
1626             for (i = n0; i < n1; i++)
1627             {
1628                 nbs->cell[nbs->a[i]] = i;
1629             }
1630         }
1631     }
1632
1633     /* Sort the super-cell columns along z into the sub-cells. */
1634 #pragma omp parallel for num_threads(nbs->nthread_max) schedule(static)
1635     for (thread = 0; thread < nbs->nthread_max; thread++)
1636     {
1637         if (grid->bSimple)
1638         {
1639             sort_columns_simple(nbs, dd_zone, grid, a0, a1, atinfo, x, nbat,
1640                                 ((thread+0)*grid->ncx*grid->ncy)/nthread,
1641                                 ((thread+1)*grid->ncx*grid->ncy)/nthread,
1642                                 nbs->work[thread].sort_work);
1643         }
1644         else
1645         {
1646             sort_columns_supersub(nbs, dd_zone, grid, a0, a1, atinfo, x, nbat,
1647                                   ((thread+0)*grid->ncx*grid->ncy)/nthread,
1648                                   ((thread+1)*grid->ncx*grid->ncy)/nthread,
1649                                   nbs->work[thread].sort_work);
1650         }
1651     }
1652
1653     if (grid->bSimple && nbat->XFormat == nbatX8)
1654     {
1655         combine_bounding_box_pairs(grid, grid->bb);
1656     }
1657
1658     if (!grid->bSimple)
1659     {
1660         grid->nsubc_tot = 0;
1661         for (i = 0; i < grid->nc; i++)
1662         {
1663             grid->nsubc_tot += grid->nsubc[i];
1664         }
1665     }
1666
1667     if (debug)
1668     {
1669         if (grid->bSimple)
1670         {
1671             print_bbsizes_simple(debug, nbs, grid);
1672         }
1673         else
1674         {
1675             fprintf(debug, "ns non-zero sub-cells: %d average atoms %.2f\n",
1676                     grid->nsubc_tot, (a1-a0)/(double)grid->nsubc_tot);
1677
1678             print_bbsizes_supersub(debug, nbs, grid);
1679         }
1680     }
1681 }
1682
1683 static void init_buffer_flags(nbnxn_buffer_flags_t *flags,
1684                               int                   natoms)
1685 {
1686     int b;
1687
1688     flags->nflag = (natoms + NBNXN_BUFFERFLAG_SIZE - 1)/NBNXN_BUFFERFLAG_SIZE;
1689     if (flags->nflag > flags->flag_nalloc)
1690     {
1691         flags->flag_nalloc = over_alloc_large(flags->nflag);
1692         srenew(flags->flag, flags->flag_nalloc);
1693     }
1694     for (b = 0; b < flags->nflag; b++)
1695     {
1696         flags->flag[b] = 0;
1697     }
1698 }
1699
1700 /* Sets up a grid and puts the atoms on the grid.
1701  * This function only operates on one domain of the domain decompostion.
1702  * Note that without domain decomposition there is only one domain.
1703  */
1704 void nbnxn_put_on_grid(nbnxn_search_t nbs,
1705                        int ePBC, matrix box,
1706                        int dd_zone,
1707                        rvec corner0, rvec corner1,
1708                        int a0, int a1,
1709                        real atom_density,
1710                        const int *atinfo,
1711                        rvec *x,
1712                        int nmoved, int *move,
1713                        int nb_kernel_type,
1714                        nbnxn_atomdata_t *nbat)
1715 {
1716     nbnxn_grid_t *grid;
1717     int           n;
1718     int           nc_max_grid, nc_max;
1719
1720     grid = &nbs->grid[dd_zone];
1721
1722     nbs_cycle_start(&nbs->cc[enbsCCgrid]);
1723
1724     grid->bSimple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
1725
1726     grid->na_c      = nbnxn_kernel_to_ci_size(nb_kernel_type);
1727     grid->na_cj     = nbnxn_kernel_to_cj_size(nb_kernel_type);
1728     grid->na_sc     = (grid->bSimple ? 1 : GPU_NSUBCELL)*grid->na_c;
1729     grid->na_c_2log = get_2log(grid->na_c);
1730
1731     nbat->na_c = grid->na_c;
1732
1733     if (dd_zone == 0)
1734     {
1735         grid->cell0 = 0;
1736     }
1737     else
1738     {
1739         grid->cell0 =
1740             (nbs->grid[dd_zone-1].cell0 + nbs->grid[dd_zone-1].nc)*
1741             nbs->grid[dd_zone-1].na_sc/grid->na_sc;
1742     }
1743
1744     n = a1 - a0;
1745
1746     if (dd_zone == 0)
1747     {
1748         nbs->ePBC = ePBC;
1749         copy_mat(box, nbs->box);
1750
1751         if (atom_density >= 0)
1752         {
1753             grid->atom_density = atom_density;
1754         }
1755         else
1756         {
1757             grid->atom_density = grid_atom_density(n-nmoved, corner0, corner1);
1758         }
1759
1760         grid->cell0 = 0;
1761
1762         nbs->natoms_local    = a1 - nmoved;
1763         /* We assume that nbnxn_put_on_grid is called first
1764          * for the local atoms (dd_zone=0).
1765          */
1766         nbs->natoms_nonlocal = a1 - nmoved;
1767     }
1768     else
1769     {
1770         nbs->natoms_nonlocal = max(nbs->natoms_nonlocal, a1);
1771     }
1772
1773     nc_max_grid = set_grid_size_xy(nbs, grid,
1774                                    dd_zone, n-nmoved, corner0, corner1,
1775                                    nbs->grid[0].atom_density);
1776
1777     nc_max = grid->cell0 + nc_max_grid;
1778
1779     if (a1 > nbs->cell_nalloc)
1780     {
1781         nbs->cell_nalloc = over_alloc_large(a1);
1782         srenew(nbs->cell, nbs->cell_nalloc);
1783     }
1784
1785     /* To avoid conditionals we store the moved particles at the end of a,
1786      * make sure we have enough space.
1787      */
1788     if (nc_max*grid->na_sc + nmoved > nbs->a_nalloc)
1789     {
1790         nbs->a_nalloc = over_alloc_large(nc_max*grid->na_sc + nmoved);
1791         srenew(nbs->a, nbs->a_nalloc);
1792     }
1793
1794     /* We need padding up to a multiple of the buffer flag size: simply add */
1795     if (nc_max*grid->na_sc + NBNXN_BUFFERFLAG_SIZE > nbat->nalloc)
1796     {
1797         nbnxn_atomdata_realloc(nbat, nc_max*grid->na_sc+NBNXN_BUFFERFLAG_SIZE);
1798     }
1799
1800     calc_cell_indices(nbs, dd_zone, grid, a0, a1, atinfo, x, move, nbat);
1801
1802     if (dd_zone == 0)
1803     {
1804         nbat->natoms_local = nbat->natoms;
1805     }
1806
1807     nbs_cycle_stop(&nbs->cc[enbsCCgrid]);
1808 }
1809
1810 /* Calls nbnxn_put_on_grid for all non-local domains */
1811 void nbnxn_put_on_grid_nonlocal(nbnxn_search_t            nbs,
1812                                 const gmx_domdec_zones_t *zones,
1813                                 const int                *atinfo,
1814                                 rvec                     *x,
1815                                 int                       nb_kernel_type,
1816                                 nbnxn_atomdata_t         *nbat)
1817 {
1818     int  zone, d;
1819     rvec c0, c1;
1820
1821     for (zone = 1; zone < zones->n; zone++)
1822     {
1823         for (d = 0; d < DIM; d++)
1824         {
1825             c0[d] = zones->size[zone].bb_x0[d];
1826             c1[d] = zones->size[zone].bb_x1[d];
1827         }
1828
1829         nbnxn_put_on_grid(nbs, nbs->ePBC, NULL,
1830                           zone, c0, c1,
1831                           zones->cg_range[zone],
1832                           zones->cg_range[zone+1],
1833                           -1,
1834                           atinfo,
1835                           x,
1836                           0, NULL,
1837                           nb_kernel_type,
1838                           nbat);
1839     }
1840 }
1841
1842 /* Add simple grid type information to the local super/sub grid */
1843 void nbnxn_grid_add_simple(nbnxn_search_t    nbs,
1844                            nbnxn_atomdata_t *nbat)
1845 {
1846     nbnxn_grid_t *grid;
1847     float        *bbcz;
1848     nbnxn_bb_t   *bb;
1849     int           ncd, sc;
1850
1851     grid = &nbs->grid[0];
1852
1853     if (grid->bSimple)
1854     {
1855         gmx_incons("nbnxn_grid_simple called with a simple grid");
1856     }
1857
1858     ncd = grid->na_sc/NBNXN_CPU_CLUSTER_I_SIZE;
1859
1860     if (grid->nc*ncd > grid->nc_nalloc_simple)
1861     {
1862         grid->nc_nalloc_simple = over_alloc_large(grid->nc*ncd);
1863         srenew(grid->bbcz_simple, grid->nc_nalloc_simple*NNBSBB_D);
1864         srenew(grid->bb_simple, grid->nc_nalloc_simple);
1865         srenew(grid->flags_simple, grid->nc_nalloc_simple);
1866         if (nbat->XFormat)
1867         {
1868             sfree_aligned(grid->bbj);
1869             snew_aligned(grid->bbj, grid->nc_nalloc_simple/2, 16);
1870         }
1871     }
1872
1873     bbcz = grid->bbcz_simple;
1874     bb   = grid->bb_simple;
1875
1876 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntPairsearch)) schedule(static)
1877     for (sc = 0; sc < grid->nc; sc++)
1878     {
1879         int c, tx, na;
1880
1881         for (c = 0; c < ncd; c++)
1882         {
1883             tx = sc*ncd + c;
1884
1885             na = NBNXN_CPU_CLUSTER_I_SIZE;
1886             while (na > 0 &&
1887                    nbat->type[tx*NBNXN_CPU_CLUSTER_I_SIZE+na-1] == nbat->ntype-1)
1888             {
1889                 na--;
1890             }
1891
1892             if (na > 0)
1893             {
1894                 switch (nbat->XFormat)
1895                 {
1896                     case nbatX4:
1897                         /* PACK_X4==NBNXN_CPU_CLUSTER_I_SIZE, so this is simple */
1898                         calc_bounding_box_x_x4(na, nbat->x+tx*STRIDE_P4,
1899                                                bb+tx);
1900                         break;
1901                     case nbatX8:
1902                         /* PACK_X8>NBNXN_CPU_CLUSTER_I_SIZE, more complicated */
1903                         calc_bounding_box_x_x8(na, nbat->x+X8_IND_A(tx*NBNXN_CPU_CLUSTER_I_SIZE),
1904                                                bb+tx);
1905                         break;
1906                     default:
1907                         calc_bounding_box(na, nbat->xstride,
1908                                           nbat->x+tx*NBNXN_CPU_CLUSTER_I_SIZE*nbat->xstride,
1909                                           bb+tx);
1910                         break;
1911                 }
1912                 bbcz[tx*NNBSBB_D+0] = bb[tx].lower[BB_Z];
1913                 bbcz[tx*NNBSBB_D+1] = bb[tx].upper[BB_Z];
1914
1915                 /* No interaction optimization yet here */
1916                 grid->flags_simple[tx] = NBNXN_CI_DO_LJ(0) | NBNXN_CI_DO_COUL(0);
1917             }
1918             else
1919             {
1920                 grid->flags_simple[tx] = 0;
1921             }
1922         }
1923     }
1924
1925     if (grid->bSimple && nbat->XFormat == nbatX8)
1926     {
1927         combine_bounding_box_pairs(grid, grid->bb_simple);
1928     }
1929 }
1930
1931 void nbnxn_get_ncells(nbnxn_search_t nbs, int *ncx, int *ncy)
1932 {
1933     *ncx = nbs->grid[0].ncx;
1934     *ncy = nbs->grid[0].ncy;
1935 }
1936
1937 void nbnxn_get_atomorder(nbnxn_search_t nbs, int **a, int *n)
1938 {
1939     const nbnxn_grid_t *grid;
1940
1941     grid = &nbs->grid[0];
1942
1943     /* Return the atom order for the home cell (index 0) */
1944     *a  = nbs->a;
1945
1946     *n = grid->cxy_ind[grid->ncx*grid->ncy]*grid->na_sc;
1947 }
1948
1949 void nbnxn_set_atomorder(nbnxn_search_t nbs)
1950 {
1951     nbnxn_grid_t *grid;
1952     int           ao, cx, cy, cxy, cz, j;
1953
1954     /* Set the atom order for the home cell (index 0) */
1955     grid = &nbs->grid[0];
1956
1957     ao = 0;
1958     for (cx = 0; cx < grid->ncx; cx++)
1959     {
1960         for (cy = 0; cy < grid->ncy; cy++)
1961         {
1962             cxy = cx*grid->ncy + cy;
1963             j   = grid->cxy_ind[cxy]*grid->na_sc;
1964             for (cz = 0; cz < grid->cxy_na[cxy]; cz++)
1965             {
1966                 nbs->a[j]     = ao;
1967                 nbs->cell[ao] = j;
1968                 ao++;
1969                 j++;
1970             }
1971         }
1972     }
1973 }
1974
1975 /* Determines the cell range along one dimension that
1976  * the bounding box b0 - b1 sees.
1977  */
1978 static void get_cell_range(real b0, real b1,
1979                            int nc, real c0, real s, real invs,
1980                            real d2, real r2, int *cf, int *cl)
1981 {
1982     *cf = max((int)((b0 - c0)*invs), 0);
1983
1984     while (*cf > 0 && d2 + sqr((b0 - c0) - (*cf-1+1)*s) < r2)
1985     {
1986         (*cf)--;
1987     }
1988
1989     *cl = min((int)((b1 - c0)*invs), nc-1);
1990     while (*cl < nc-1 && d2 + sqr((*cl+1)*s - (b1 - c0)) < r2)
1991     {
1992         (*cl)++;
1993     }
1994 }
1995
1996 /* Reference code calculating the distance^2 between two bounding boxes */
1997 static float box_dist2(float bx0, float bx1, float by0,
1998                        float by1, float bz0, float bz1,
1999                        const nbnxn_bb_t *bb)
2000 {
2001     float d2;
2002     float dl, dh, dm, dm0;
2003
2004     d2 = 0;
2005
2006     dl  = bx0 - bb->upper[BB_X];
2007     dh  = bb->lower[BB_X] - bx1;
2008     dm  = max(dl, dh);
2009     dm0 = max(dm, 0);
2010     d2 += dm0*dm0;
2011
2012     dl  = by0 - bb->upper[BB_Y];
2013     dh  = bb->lower[BB_Y] - by1;
2014     dm  = max(dl, dh);
2015     dm0 = max(dm, 0);
2016     d2 += dm0*dm0;
2017
2018     dl  = bz0 - bb->upper[BB_Z];
2019     dh  = bb->lower[BB_Z] - bz1;
2020     dm  = max(dl, dh);
2021     dm0 = max(dm, 0);
2022     d2 += dm0*dm0;
2023
2024     return d2;
2025 }
2026
2027 /* Plain C code calculating the distance^2 between two bounding boxes */
2028 static float subc_bb_dist2(int si, const nbnxn_bb_t *bb_i_ci,
2029                            int csj, const nbnxn_bb_t *bb_j_all)
2030 {
2031     const nbnxn_bb_t *bb_i, *bb_j;
2032     float             d2;
2033     float             dl, dh, dm, dm0;
2034
2035     bb_i = bb_i_ci  +  si;
2036     bb_j = bb_j_all + csj;
2037
2038     d2 = 0;
2039
2040     dl  = bb_i->lower[BB_X] - bb_j->upper[BB_X];
2041     dh  = bb_j->lower[BB_X] - bb_i->upper[BB_X];
2042     dm  = max(dl, dh);
2043     dm0 = max(dm, 0);
2044     d2 += dm0*dm0;
2045
2046     dl  = bb_i->lower[BB_Y] - bb_j->upper[BB_Y];
2047     dh  = bb_j->lower[BB_Y] - bb_i->upper[BB_Y];
2048     dm  = max(dl, dh);
2049     dm0 = max(dm, 0);
2050     d2 += dm0*dm0;
2051
2052     dl  = bb_i->lower[BB_Z] - bb_j->upper[BB_Z];
2053     dh  = bb_j->lower[BB_Z] - bb_i->upper[BB_Z];
2054     dm  = max(dl, dh);
2055     dm0 = max(dm, 0);
2056     d2 += dm0*dm0;
2057
2058     return d2;
2059 }
2060
2061 #ifdef NBNXN_SEARCH_BB_SSE
2062
2063 /* SSE code for bb distance for bb format xyz0 */
2064 static float subc_bb_dist2_sse(int si, const nbnxn_bb_t *bb_i_ci,
2065                                int csj, const nbnxn_bb_t *bb_j_all)
2066 {
2067     __m128       bb_i_SSE0, bb_i_SSE1;
2068     __m128       bb_j_SSE0, bb_j_SSE1;
2069     __m128       dl_SSE;
2070     __m128       dh_SSE;
2071     __m128       dm_SSE;
2072     __m128       dm0_SSE;
2073     __m128       d2_SSE;
2074 #ifndef GMX_X86_SSE4_1
2075     float        d2_array[7], *d2_align;
2076
2077     d2_align = (float *)(((size_t)(d2_array+3)) & (~((size_t)15)));
2078 #else
2079     float d2;
2080 #endif
2081
2082     bb_i_SSE0 = _mm_load_ps(&bb_i_ci[si].lower[0]);
2083     bb_i_SSE1 = _mm_load_ps(&bb_i_ci[si].upper[0]);
2084     bb_j_SSE0 = _mm_load_ps(&bb_j_all[csj].lower[0]);
2085     bb_j_SSE1 = _mm_load_ps(&bb_j_all[csj].upper[0]);
2086
2087     dl_SSE    = _mm_sub_ps(bb_i_SSE0, bb_j_SSE1);
2088     dh_SSE    = _mm_sub_ps(bb_j_SSE0, bb_i_SSE1);
2089
2090     dm_SSE    = _mm_max_ps(dl_SSE, dh_SSE);
2091     dm0_SSE   = _mm_max_ps(dm_SSE, _mm_setzero_ps());
2092 #ifndef GMX_X86_SSE4_1
2093     d2_SSE    = _mm_mul_ps(dm0_SSE, dm0_SSE);
2094
2095     _mm_store_ps(d2_align, d2_SSE);
2096
2097     return d2_align[0] + d2_align[1] + d2_align[2];
2098 #else
2099     /* SSE4.1 dot product of components 0,1,2 */
2100     d2_SSE    = _mm_dp_ps(dm0_SSE, dm0_SSE, 0x71);
2101
2102     _mm_store_ss(&d2, d2_SSE);
2103
2104     return d2;
2105 #endif
2106 }
2107
2108 /* Calculate bb bounding distances of bb_i[si,...,si+3] and store them in d2 */
2109 #define SUBC_BB_DIST2_SSE_XXXX_INNER(si, bb_i, d2) \
2110     {                                                \
2111         int    shi;                                  \
2112                                                  \
2113         __m128 dx_0, dy_0, dz_0;                       \
2114         __m128 dx_1, dy_1, dz_1;                       \
2115                                                  \
2116         __m128 mx, my, mz;                             \
2117         __m128 m0x, m0y, m0z;                          \
2118                                                  \
2119         __m128 d2x, d2y, d2z;                          \
2120         __m128 d2s, d2t;                              \
2121                                                  \
2122         shi = si*NNBSBB_D*DIM;                       \
2123                                                  \
2124         xi_l = _mm_load_ps(bb_i+shi+0*STRIDE_PBB);   \
2125         yi_l = _mm_load_ps(bb_i+shi+1*STRIDE_PBB);   \
2126         zi_l = _mm_load_ps(bb_i+shi+2*STRIDE_PBB);   \
2127         xi_h = _mm_load_ps(bb_i+shi+3*STRIDE_PBB);   \
2128         yi_h = _mm_load_ps(bb_i+shi+4*STRIDE_PBB);   \
2129         zi_h = _mm_load_ps(bb_i+shi+5*STRIDE_PBB);   \
2130                                                  \
2131         dx_0 = _mm_sub_ps(xi_l, xj_h);                \
2132         dy_0 = _mm_sub_ps(yi_l, yj_h);                \
2133         dz_0 = _mm_sub_ps(zi_l, zj_h);                \
2134                                                  \
2135         dx_1 = _mm_sub_ps(xj_l, xi_h);                \
2136         dy_1 = _mm_sub_ps(yj_l, yi_h);                \
2137         dz_1 = _mm_sub_ps(zj_l, zi_h);                \
2138                                                  \
2139         mx   = _mm_max_ps(dx_0, dx_1);                \
2140         my   = _mm_max_ps(dy_0, dy_1);                \
2141         mz   = _mm_max_ps(dz_0, dz_1);                \
2142                                                  \
2143         m0x  = _mm_max_ps(mx, zero);                  \
2144         m0y  = _mm_max_ps(my, zero);                  \
2145         m0z  = _mm_max_ps(mz, zero);                  \
2146                                                  \
2147         d2x  = _mm_mul_ps(m0x, m0x);                  \
2148         d2y  = _mm_mul_ps(m0y, m0y);                  \
2149         d2z  = _mm_mul_ps(m0z, m0z);                  \
2150                                                  \
2151         d2s  = _mm_add_ps(d2x, d2y);                  \
2152         d2t  = _mm_add_ps(d2s, d2z);                  \
2153                                                  \
2154         _mm_store_ps(d2+si, d2t);                     \
2155     }
2156
2157 /* SSE code for nsi bb distances for bb format xxxxyyyyzzzz */
2158 static void subc_bb_dist2_sse_xxxx(const float *bb_j,
2159                                    int nsi, const float *bb_i,
2160                                    float *d2)
2161 {
2162     __m128 xj_l, yj_l, zj_l;
2163     __m128 xj_h, yj_h, zj_h;
2164     __m128 xi_l, yi_l, zi_l;
2165     __m128 xi_h, yi_h, zi_h;
2166
2167     __m128 zero;
2168
2169     zero = _mm_setzero_ps();
2170
2171     xj_l = _mm_set1_ps(bb_j[0*STRIDE_PBB]);
2172     yj_l = _mm_set1_ps(bb_j[1*STRIDE_PBB]);
2173     zj_l = _mm_set1_ps(bb_j[2*STRIDE_PBB]);
2174     xj_h = _mm_set1_ps(bb_j[3*STRIDE_PBB]);
2175     yj_h = _mm_set1_ps(bb_j[4*STRIDE_PBB]);
2176     zj_h = _mm_set1_ps(bb_j[5*STRIDE_PBB]);
2177
2178     /* Here we "loop" over si (0,STRIDE_PBB) from 0 to nsi with step STRIDE_PBB.
2179      * But as we know the number of iterations is 1 or 2, we unroll manually.
2180      */
2181     SUBC_BB_DIST2_SSE_XXXX_INNER(0, bb_i, d2);
2182     if (STRIDE_PBB < nsi)
2183     {
2184         SUBC_BB_DIST2_SSE_XXXX_INNER(STRIDE_PBB, bb_i, d2);
2185     }
2186 }
2187
2188 #endif /* NBNXN_SEARCH_BB_SSE */
2189
2190 /* Plain C function which determines if any atom pair between two cells
2191  * is within distance sqrt(rl2).
2192  */
2193 static gmx_bool subc_in_range_x(int na_c,
2194                                 int si, const real *x_i,
2195                                 int csj, int stride, const real *x_j,
2196                                 real rl2)
2197 {
2198     int  i, j, i0, j0;
2199     real d2;
2200
2201     for (i = 0; i < na_c; i++)
2202     {
2203         i0 = (si*na_c + i)*DIM;
2204         for (j = 0; j < na_c; j++)
2205         {
2206             j0 = (csj*na_c + j)*stride;
2207
2208             d2 = sqr(x_i[i0  ] - x_j[j0  ]) +
2209                 sqr(x_i[i0+1] - x_j[j0+1]) +
2210                 sqr(x_i[i0+2] - x_j[j0+2]);
2211
2212             if (d2 < rl2)
2213             {
2214                 return TRUE;
2215             }
2216         }
2217     }
2218
2219     return FALSE;
2220 }
2221
2222 #ifdef NBNXN_SEARCH_SSE_SINGLE
2223 /* When we make seperate single/double precision SIMD vector operation
2224  * include files, this function should be moved there (also using FMA).
2225  */
2226 static inline __m128
2227 gmx_mm_calc_rsq_ps(__m128 x, __m128 y, __m128 z)
2228 {
2229     return _mm_add_ps( _mm_add_ps( _mm_mul_ps(x, x), _mm_mul_ps(y, y) ), _mm_mul_ps(z, z) );
2230 }
2231 #endif
2232
2233 /* SSE function which determines if any atom pair between two cells,
2234  * both with 8 atoms, is within distance sqrt(rl2).
2235  * Not performance critical, so only uses plain SSE.
2236  */
2237 static gmx_bool subc_in_range_sse8(int na_c,
2238                                    int si, const real *x_i,
2239                                    int csj, int stride, const real *x_j,
2240                                    real rl2)
2241 {
2242 #ifdef NBNXN_SEARCH_SSE_SINGLE
2243     __m128 ix_SSE0, iy_SSE0, iz_SSE0;
2244     __m128 ix_SSE1, iy_SSE1, iz_SSE1;
2245
2246     __m128 rc2_SSE;
2247
2248     int    na_c_sse;
2249     int    j0, j1;
2250
2251     rc2_SSE   = _mm_set1_ps(rl2);
2252
2253     na_c_sse = NBNXN_GPU_CLUSTER_SIZE/STRIDE_PBB;
2254     ix_SSE0  = _mm_load_ps(x_i+(si*na_c_sse*DIM+0)*STRIDE_PBB);
2255     iy_SSE0  = _mm_load_ps(x_i+(si*na_c_sse*DIM+1)*STRIDE_PBB);
2256     iz_SSE0  = _mm_load_ps(x_i+(si*na_c_sse*DIM+2)*STRIDE_PBB);
2257     ix_SSE1  = _mm_load_ps(x_i+(si*na_c_sse*DIM+3)*STRIDE_PBB);
2258     iy_SSE1  = _mm_load_ps(x_i+(si*na_c_sse*DIM+4)*STRIDE_PBB);
2259     iz_SSE1  = _mm_load_ps(x_i+(si*na_c_sse*DIM+5)*STRIDE_PBB);
2260
2261     /* We loop from the outer to the inner particles to maximize
2262      * the chance that we find a pair in range quickly and return.
2263      */
2264     j0 = csj*na_c;
2265     j1 = j0 + na_c - 1;
2266     while (j0 < j1)
2267     {
2268         __m128 jx0_SSE, jy0_SSE, jz0_SSE;
2269         __m128 jx1_SSE, jy1_SSE, jz1_SSE;
2270
2271         __m128 dx_SSE0, dy_SSE0, dz_SSE0;
2272         __m128 dx_SSE1, dy_SSE1, dz_SSE1;
2273         __m128 dx_SSE2, dy_SSE2, dz_SSE2;
2274         __m128 dx_SSE3, dy_SSE3, dz_SSE3;
2275
2276         __m128 rsq_SSE0;
2277         __m128 rsq_SSE1;
2278         __m128 rsq_SSE2;
2279         __m128 rsq_SSE3;
2280
2281         __m128 wco_SSE0;
2282         __m128 wco_SSE1;
2283         __m128 wco_SSE2;
2284         __m128 wco_SSE3;
2285         __m128 wco_any_SSE01, wco_any_SSE23, wco_any_SSE;
2286
2287         jx0_SSE = _mm_load1_ps(x_j+j0*stride+0);
2288         jy0_SSE = _mm_load1_ps(x_j+j0*stride+1);
2289         jz0_SSE = _mm_load1_ps(x_j+j0*stride+2);
2290
2291         jx1_SSE = _mm_load1_ps(x_j+j1*stride+0);
2292         jy1_SSE = _mm_load1_ps(x_j+j1*stride+1);
2293         jz1_SSE = _mm_load1_ps(x_j+j1*stride+2);
2294
2295         /* Calculate distance */
2296         dx_SSE0            = _mm_sub_ps(ix_SSE0, jx0_SSE);
2297         dy_SSE0            = _mm_sub_ps(iy_SSE0, jy0_SSE);
2298         dz_SSE0            = _mm_sub_ps(iz_SSE0, jz0_SSE);
2299         dx_SSE1            = _mm_sub_ps(ix_SSE1, jx0_SSE);
2300         dy_SSE1            = _mm_sub_ps(iy_SSE1, jy0_SSE);
2301         dz_SSE1            = _mm_sub_ps(iz_SSE1, jz0_SSE);
2302         dx_SSE2            = _mm_sub_ps(ix_SSE0, jx1_SSE);
2303         dy_SSE2            = _mm_sub_ps(iy_SSE0, jy1_SSE);
2304         dz_SSE2            = _mm_sub_ps(iz_SSE0, jz1_SSE);
2305         dx_SSE3            = _mm_sub_ps(ix_SSE1, jx1_SSE);
2306         dy_SSE3            = _mm_sub_ps(iy_SSE1, jy1_SSE);
2307         dz_SSE3            = _mm_sub_ps(iz_SSE1, jz1_SSE);
2308
2309         /* rsq = dx*dx+dy*dy+dz*dz */
2310         rsq_SSE0           = gmx_mm_calc_rsq_ps(dx_SSE0, dy_SSE0, dz_SSE0);
2311         rsq_SSE1           = gmx_mm_calc_rsq_ps(dx_SSE1, dy_SSE1, dz_SSE1);
2312         rsq_SSE2           = gmx_mm_calc_rsq_ps(dx_SSE2, dy_SSE2, dz_SSE2);
2313         rsq_SSE3           = gmx_mm_calc_rsq_ps(dx_SSE3, dy_SSE3, dz_SSE3);
2314
2315         wco_SSE0           = _mm_cmplt_ps(rsq_SSE0, rc2_SSE);
2316         wco_SSE1           = _mm_cmplt_ps(rsq_SSE1, rc2_SSE);
2317         wco_SSE2           = _mm_cmplt_ps(rsq_SSE2, rc2_SSE);
2318         wco_SSE3           = _mm_cmplt_ps(rsq_SSE3, rc2_SSE);
2319
2320         wco_any_SSE01      = _mm_or_ps(wco_SSE0, wco_SSE1);
2321         wco_any_SSE23      = _mm_or_ps(wco_SSE2, wco_SSE3);
2322         wco_any_SSE        = _mm_or_ps(wco_any_SSE01, wco_any_SSE23);
2323
2324         if (_mm_movemask_ps(wco_any_SSE))
2325         {
2326             return TRUE;
2327         }
2328
2329         j0++;
2330         j1--;
2331     }
2332     return FALSE;
2333
2334 #else
2335     /* No SSE */
2336     gmx_incons("SSE function called without SSE support");
2337
2338     return TRUE;
2339 #endif
2340 }
2341
2342 /* Returns the j sub-cell for index cj_ind */
2343 static int nbl_cj(const nbnxn_pairlist_t *nbl, int cj_ind)
2344 {
2345     return nbl->cj4[cj_ind >> NBNXN_GPU_JGROUP_SIZE_2LOG].cj[cj_ind & (NBNXN_GPU_JGROUP_SIZE - 1)];
2346 }
2347
2348 /* Returns the i-interaction mask of the j sub-cell for index cj_ind */
2349 static unsigned nbl_imask0(const nbnxn_pairlist_t *nbl, int cj_ind)
2350 {
2351     return nbl->cj4[cj_ind >> NBNXN_GPU_JGROUP_SIZE_2LOG].imei[0].imask;
2352 }
2353
2354 /* Ensures there is enough space for extra extra exclusion masks */
2355 static void check_excl_space(nbnxn_pairlist_t *nbl, int extra)
2356 {
2357     if (nbl->nexcl+extra > nbl->excl_nalloc)
2358     {
2359         nbl->excl_nalloc = over_alloc_small(nbl->nexcl+extra);
2360         nbnxn_realloc_void((void **)&nbl->excl,
2361                            nbl->nexcl*sizeof(*nbl->excl),
2362                            nbl->excl_nalloc*sizeof(*nbl->excl),
2363                            nbl->alloc, nbl->free);
2364     }
2365 }
2366
2367 /* Ensures there is enough space for ncell extra j-cells in the list */
2368 static void check_subcell_list_space_simple(nbnxn_pairlist_t *nbl,
2369                                             int               ncell)
2370 {
2371     int cj_max;
2372
2373     cj_max = nbl->ncj + ncell;
2374
2375     if (cj_max > nbl->cj_nalloc)
2376     {
2377         nbl->cj_nalloc = over_alloc_small(cj_max);
2378         nbnxn_realloc_void((void **)&nbl->cj,
2379                            nbl->ncj*sizeof(*nbl->cj),
2380                            nbl->cj_nalloc*sizeof(*nbl->cj),
2381                            nbl->alloc, nbl->free);
2382     }
2383 }
2384
2385 /* Ensures there is enough space for ncell extra j-subcells in the list */
2386 static void check_subcell_list_space_supersub(nbnxn_pairlist_t *nbl,
2387                                               int               nsupercell)
2388 {
2389     int ncj4_max, j4, j, w, t;
2390
2391 #define NWARP       2
2392 #define WARP_SIZE  32
2393
2394     /* We can have maximally nsupercell*GPU_NSUBCELL sj lists */
2395     /* We can store 4 j-subcell - i-supercell pairs in one struct.
2396      * since we round down, we need one extra entry.
2397      */
2398     ncj4_max = ((nbl->work->cj_ind + nsupercell*GPU_NSUBCELL + NBNXN_GPU_JGROUP_SIZE - 1) >> NBNXN_GPU_JGROUP_SIZE_2LOG);
2399
2400     if (ncj4_max > nbl->cj4_nalloc)
2401     {
2402         nbl->cj4_nalloc = over_alloc_small(ncj4_max);
2403         nbnxn_realloc_void((void **)&nbl->cj4,
2404                            nbl->work->cj4_init*sizeof(*nbl->cj4),
2405                            nbl->cj4_nalloc*sizeof(*nbl->cj4),
2406                            nbl->alloc, nbl->free);
2407     }
2408
2409     if (ncj4_max > nbl->work->cj4_init)
2410     {
2411         for (j4 = nbl->work->cj4_init; j4 < ncj4_max; j4++)
2412         {
2413             /* No i-subcells and no excl's in the list initially */
2414             for (w = 0; w < NWARP; w++)
2415             {
2416                 nbl->cj4[j4].imei[w].imask    = 0U;
2417                 nbl->cj4[j4].imei[w].excl_ind = 0;
2418
2419             }
2420         }
2421         nbl->work->cj4_init = ncj4_max;
2422     }
2423 }
2424
2425 /* Set all excl masks for one GPU warp no exclusions */
2426 static void set_no_excls(nbnxn_excl_t *excl)
2427 {
2428     int t;
2429
2430     for (t = 0; t < WARP_SIZE; t++)
2431     {
2432         /* Turn all interaction bits on */
2433         excl->pair[t] = NBNXN_INTERACTION_MASK_ALL;
2434     }
2435 }
2436
2437 /* Initializes a single nbnxn_pairlist_t data structure */
2438 static void nbnxn_init_pairlist(nbnxn_pairlist_t *nbl,
2439                                 gmx_bool          bSimple,
2440                                 nbnxn_alloc_t    *alloc,
2441                                 nbnxn_free_t     *free)
2442 {
2443     if (alloc == NULL)
2444     {
2445         nbl->alloc = nbnxn_alloc_aligned;
2446     }
2447     else
2448     {
2449         nbl->alloc = alloc;
2450     }
2451     if (free == NULL)
2452     {
2453         nbl->free = nbnxn_free_aligned;
2454     }
2455     else
2456     {
2457         nbl->free = free;
2458     }
2459
2460     nbl->bSimple     = bSimple;
2461     nbl->na_sc       = 0;
2462     nbl->na_ci       = 0;
2463     nbl->na_cj       = 0;
2464     nbl->nci         = 0;
2465     nbl->ci          = NULL;
2466     nbl->ci_nalloc   = 0;
2467     nbl->ncj         = 0;
2468     nbl->cj          = NULL;
2469     nbl->cj_nalloc   = 0;
2470     nbl->ncj4        = 0;
2471     /* We need one element extra in sj, so alloc initially with 1 */
2472     nbl->cj4_nalloc  = 0;
2473     nbl->cj4         = NULL;
2474     nbl->nci_tot     = 0;
2475
2476     if (!nbl->bSimple)
2477     {
2478         nbl->excl        = NULL;
2479         nbl->excl_nalloc = 0;
2480         nbl->nexcl       = 0;
2481         check_excl_space(nbl, 1);
2482         nbl->nexcl       = 1;
2483         set_no_excls(&nbl->excl[0]);
2484     }
2485
2486     snew(nbl->work, 1);
2487     if (nbl->bSimple)
2488     {
2489         snew_aligned(nbl->work->bb_ci, 1, NBNXN_MEM_ALIGN);
2490     }
2491     else
2492     {
2493 #ifdef NBNXN_BBXXXX
2494         snew_aligned(nbl->work->pbb_ci, GPU_NSUBCELL/STRIDE_PBB*NNBSBB_XXXX, NBNXN_MEM_ALIGN);
2495 #else
2496         snew_aligned(nbl->work->bb_ci, GPU_NSUBCELL, NBNXN_MEM_ALIGN);
2497 #endif
2498     }
2499     snew_aligned(nbl->work->x_ci, NBNXN_NA_SC_MAX*DIM, NBNXN_MEM_ALIGN);
2500 #ifdef GMX_NBNXN_SIMD
2501     snew_aligned(nbl->work->x_ci_simd_4xn, 1, NBNXN_MEM_ALIGN);
2502     snew_aligned(nbl->work->x_ci_simd_2xnn, 1, NBNXN_MEM_ALIGN);
2503 #endif
2504     snew_aligned(nbl->work->d2, GPU_NSUBCELL, NBNXN_MEM_ALIGN);
2505
2506     nbl->work->sort            = NULL;
2507     nbl->work->sort_nalloc     = 0;
2508     nbl->work->sci_sort        = NULL;
2509     nbl->work->sci_sort_nalloc = 0;
2510 }
2511
2512 void nbnxn_init_pairlist_set(nbnxn_pairlist_set_t *nbl_list,
2513                              gmx_bool bSimple, gmx_bool bCombined,
2514                              nbnxn_alloc_t *alloc,
2515                              nbnxn_free_t  *free)
2516 {
2517     int i;
2518
2519     nbl_list->bSimple   = bSimple;
2520     nbl_list->bCombined = bCombined;
2521
2522     nbl_list->nnbl = gmx_omp_nthreads_get(emntNonbonded);
2523
2524     if (!nbl_list->bCombined &&
2525         nbl_list->nnbl > NBNXN_BUFFERFLAG_MAX_THREADS)
2526     {
2527         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.",
2528                   nbl_list->nnbl, NBNXN_BUFFERFLAG_MAX_THREADS, NBNXN_BUFFERFLAG_MAX_THREADS);
2529     }
2530
2531     snew(nbl_list->nbl, nbl_list->nnbl);
2532     /* Execute in order to avoid memory interleaving between threads */
2533 #pragma omp parallel for num_threads(nbl_list->nnbl) schedule(static)
2534     for (i = 0; i < nbl_list->nnbl; i++)
2535     {
2536         /* Allocate the nblist data structure locally on each thread
2537          * to optimize memory access for NUMA architectures.
2538          */
2539         snew(nbl_list->nbl[i], 1);
2540
2541         /* Only list 0 is used on the GPU, use normal allocation for i>0 */
2542         if (i == 0)
2543         {
2544             nbnxn_init_pairlist(nbl_list->nbl[i], nbl_list->bSimple, alloc, free);
2545         }
2546         else
2547         {
2548             nbnxn_init_pairlist(nbl_list->nbl[i], nbl_list->bSimple, NULL, NULL);
2549         }
2550     }
2551 }
2552
2553 /* Print statistics of a pair list, used for debug output */
2554 static void print_nblist_statistics_simple(FILE *fp, const nbnxn_pairlist_t *nbl,
2555                                            const nbnxn_search_t nbs, real rl)
2556 {
2557     const nbnxn_grid_t *grid;
2558     int                 cs[SHIFTS];
2559     int                 s, i, j;
2560     int                 npexcl;
2561
2562     /* This code only produces correct statistics with domain decomposition */
2563     grid = &nbs->grid[0];
2564
2565     fprintf(fp, "nbl nci %d ncj %d\n",
2566             nbl->nci, nbl->ncj);
2567     fprintf(fp, "nbl na_sc %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
2568             nbl->na_sc, rl, nbl->ncj, nbl->ncj/(double)grid->nc,
2569             nbl->ncj/(double)grid->nc*grid->na_sc,
2570             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)));
2571
2572     fprintf(fp, "nbl average j cell list length %.1f\n",
2573             0.25*nbl->ncj/(double)nbl->nci);
2574
2575     for (s = 0; s < SHIFTS; s++)
2576     {
2577         cs[s] = 0;
2578     }
2579     npexcl = 0;
2580     for (i = 0; i < nbl->nci; i++)
2581     {
2582         cs[nbl->ci[i].shift & NBNXN_CI_SHIFT] +=
2583             nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start;
2584
2585         j = nbl->ci[i].cj_ind_start;
2586         while (j < nbl->ci[i].cj_ind_end &&
2587                nbl->cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
2588         {
2589             npexcl++;
2590             j++;
2591         }
2592     }
2593     fprintf(fp, "nbl cell pairs, total: %d excl: %d %.1f%%\n",
2594             nbl->ncj, npexcl, 100*npexcl/(double)nbl->ncj);
2595     for (s = 0; s < SHIFTS; s++)
2596     {
2597         if (cs[s] > 0)
2598         {
2599             fprintf(fp, "nbl shift %2d ncj %3d\n", s, cs[s]);
2600         }
2601     }
2602 }
2603
2604 /* Print statistics of a pair lists, used for debug output */
2605 static void print_nblist_statistics_supersub(FILE *fp, const nbnxn_pairlist_t *nbl,
2606                                              const nbnxn_search_t nbs, real rl)
2607 {
2608     const nbnxn_grid_t *grid;
2609     int                 i, j4, j, si, b;
2610     int                 c[GPU_NSUBCELL+1];
2611
2612     /* This code only produces correct statistics with domain decomposition */
2613     grid = &nbs->grid[0];
2614
2615     fprintf(fp, "nbl nsci %d ncj4 %d nsi %d excl4 %d\n",
2616             nbl->nsci, nbl->ncj4, nbl->nci_tot, nbl->nexcl);
2617     fprintf(fp, "nbl na_c %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
2618             nbl->na_ci, rl, nbl->nci_tot, nbl->nci_tot/(double)grid->nsubc_tot,
2619             nbl->nci_tot/(double)grid->nsubc_tot*grid->na_c,
2620             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)));
2621
2622     fprintf(fp, "nbl average j super cell list length %.1f\n",
2623             0.25*nbl->ncj4/(double)nbl->nsci);
2624     fprintf(fp, "nbl average i sub cell list length %.1f\n",
2625             nbl->nci_tot/((double)nbl->ncj4));
2626
2627     for (si = 0; si <= GPU_NSUBCELL; si++)
2628     {
2629         c[si] = 0;
2630     }
2631     for (i = 0; i < nbl->nsci; i++)
2632     {
2633         for (j4 = nbl->sci[i].cj4_ind_start; j4 < nbl->sci[i].cj4_ind_end; j4++)
2634         {
2635             for (j = 0; j < NBNXN_GPU_JGROUP_SIZE; j++)
2636             {
2637                 b = 0;
2638                 for (si = 0; si < GPU_NSUBCELL; si++)
2639                 {
2640                     if (nbl->cj4[j4].imei[0].imask & (1U << (j*GPU_NSUBCELL + si)))
2641                     {
2642                         b++;
2643                     }
2644                 }
2645                 c[b]++;
2646             }
2647         }
2648     }
2649     for (b = 0; b <= GPU_NSUBCELL; b++)
2650     {
2651         fprintf(fp, "nbl j-list #i-subcell %d %7d %4.1f\n",
2652                 b, c[b], 100.0*c[b]/(double)(nbl->ncj4*NBNXN_GPU_JGROUP_SIZE));
2653     }
2654 }
2655
2656 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp */
2657 static void low_get_nbl_exclusions(nbnxn_pairlist_t *nbl, int cj4,
2658                                    int warp, nbnxn_excl_t **excl)
2659 {
2660     if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
2661     {
2662         /* No exclusions set, make a new list entry */
2663         nbl->cj4[cj4].imei[warp].excl_ind = nbl->nexcl;
2664         nbl->nexcl++;
2665         *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
2666         set_no_excls(*excl);
2667     }
2668     else
2669     {
2670         /* We already have some exclusions, new ones can be added to the list */
2671         *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
2672     }
2673 }
2674
2675 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp,
2676  * allocates extra memory, if necessary.
2677  */
2678 static void get_nbl_exclusions_1(nbnxn_pairlist_t *nbl, int cj4,
2679                                  int warp, nbnxn_excl_t **excl)
2680 {
2681     if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
2682     {
2683         /* We need to make a new list entry, check if we have space */
2684         check_excl_space(nbl, 1);
2685     }
2686     low_get_nbl_exclusions(nbl, cj4, warp, excl);
2687 }
2688
2689 /* Returns pointers to the exclusion mask for cj4-unit cj4 for both warps,
2690  * allocates extra memory, if necessary.
2691  */
2692 static void get_nbl_exclusions_2(nbnxn_pairlist_t *nbl, int cj4,
2693                                  nbnxn_excl_t **excl_w0,
2694                                  nbnxn_excl_t **excl_w1)
2695 {
2696     /* Check for space we might need */
2697     check_excl_space(nbl, 2);
2698
2699     low_get_nbl_exclusions(nbl, cj4, 0, excl_w0);
2700     low_get_nbl_exclusions(nbl, cj4, 1, excl_w1);
2701 }
2702
2703 /* Sets the self exclusions i=j and pair exclusions i>j */
2704 static void set_self_and_newton_excls_supersub(nbnxn_pairlist_t *nbl,
2705                                                int cj4_ind, int sj_offset,
2706                                                int si)
2707 {
2708     nbnxn_excl_t *excl[2];
2709     int           ei, ej, w;
2710
2711     /* Here we only set the set self and double pair exclusions */
2712
2713     get_nbl_exclusions_2(nbl, cj4_ind, &excl[0], &excl[1]);
2714
2715     /* Only minor < major bits set */
2716     for (ej = 0; ej < nbl->na_ci; ej++)
2717     {
2718         w = (ej>>2);
2719         for (ei = ej; ei < nbl->na_ci; ei++)
2720         {
2721             excl[w]->pair[(ej & (NBNXN_GPU_JGROUP_SIZE-1))*nbl->na_ci + ei] &=
2722                 ~(1U << (sj_offset*GPU_NSUBCELL + si));
2723         }
2724     }
2725 }
2726
2727 /* Returns a diagonal or off-diagonal interaction mask for plain C lists */
2728 static unsigned int get_imask(gmx_bool rdiag, int ci, int cj)
2729 {
2730     return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
2731 }
2732
2733 /* Returns a diagonal or off-diagonal interaction mask for cj-size=2 */
2734 static unsigned int get_imask_simd_j2(gmx_bool rdiag, int ci, int cj)
2735 {
2736     return (rdiag && ci*2 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_0 :
2737             (rdiag && ci*2+1 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_1 :
2738              NBNXN_INTERACTION_MASK_ALL));
2739 }
2740
2741 /* Returns a diagonal or off-diagonal interaction mask for cj-size=4 */
2742 static unsigned int get_imask_simd_j4(gmx_bool rdiag, int ci, int cj)
2743 {
2744     return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
2745 }
2746
2747 /* Returns a diagonal or off-diagonal interaction mask for cj-size=8 */
2748 static unsigned int get_imask_simd_j8(gmx_bool rdiag, int ci, int cj)
2749 {
2750     return (rdiag && ci == cj*2 ? NBNXN_INTERACTION_MASK_DIAG_J8_0 :
2751             (rdiag && ci == cj*2+1 ? NBNXN_INTERACTION_MASK_DIAG_J8_1 :
2752              NBNXN_INTERACTION_MASK_ALL));
2753 }
2754
2755 #ifdef GMX_NBNXN_SIMD
2756 #if GMX_SIMD_WIDTH_HERE == 2
2757 #define get_imask_simd_4xn  get_imask_simd_j2
2758 #endif
2759 #if GMX_SIMD_WIDTH_HERE == 4
2760 #define get_imask_simd_4xn  get_imask_simd_j4
2761 #endif
2762 #if GMX_SIMD_WIDTH_HERE == 8
2763 #define get_imask_simd_4xn  get_imask_simd_j8
2764 #define get_imask_simd_2xnn get_imask_simd_j4
2765 #endif
2766 #if GMX_SIMD_WIDTH_HERE == 16
2767 #define get_imask_simd_2xnn get_imask_simd_j8
2768 #endif
2769 #endif
2770
2771 /* Plain C code for making a pair list of cell ci vs cell cjf-cjl.
2772  * Checks bounding box distances and possibly atom pair distances.
2773  */
2774 static void make_cluster_list_simple(const nbnxn_grid_t *gridj,
2775                                      nbnxn_pairlist_t *nbl,
2776                                      int ci, int cjf, int cjl,
2777                                      gmx_bool remove_sub_diag,
2778                                      const real *x_j,
2779                                      real rl2, float rbb2,
2780                                      int *ndistc)
2781 {
2782     const nbnxn_list_work_t *work;
2783
2784     const nbnxn_bb_t        *bb_ci;
2785     const real              *x_ci;
2786
2787     gmx_bool                 InRange;
2788     real                     d2;
2789     int                      cjf_gl, cjl_gl, cj;
2790
2791     work = nbl->work;
2792
2793     bb_ci = nbl->work->bb_ci;
2794     x_ci  = nbl->work->x_ci;
2795
2796     InRange = FALSE;
2797     while (!InRange && cjf <= cjl)
2798     {
2799         d2       = subc_bb_dist2(0, bb_ci, cjf, gridj->bb);
2800         *ndistc += 2;
2801
2802         /* Check if the distance is within the distance where
2803          * we use only the bounding box distance rbb,
2804          * or within the cut-off and there is at least one atom pair
2805          * within the cut-off.
2806          */
2807         if (d2 < rbb2)
2808         {
2809             InRange = TRUE;
2810         }
2811         else if (d2 < rl2)
2812         {
2813             int i, j;
2814
2815             cjf_gl = gridj->cell0 + cjf;
2816             for (i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
2817             {
2818                 for (j = 0; j < NBNXN_CPU_CLUSTER_I_SIZE; j++)
2819                 {
2820                     InRange = InRange ||
2821                         (sqr(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
2822                          sqr(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
2823                          sqr(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rl2);
2824                 }
2825             }
2826             *ndistc += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
2827         }
2828         if (!InRange)
2829         {
2830             cjf++;
2831         }
2832     }
2833     if (!InRange)
2834     {
2835         return;
2836     }
2837
2838     InRange = FALSE;
2839     while (!InRange && cjl > cjf)
2840     {
2841         d2       = subc_bb_dist2(0, bb_ci, cjl, gridj->bb);
2842         *ndistc += 2;
2843
2844         /* Check if the distance is within the distance where
2845          * we use only the bounding box distance rbb,
2846          * or within the cut-off and there is at least one atom pair
2847          * within the cut-off.
2848          */
2849         if (d2 < rbb2)
2850         {
2851             InRange = TRUE;
2852         }
2853         else if (d2 < rl2)
2854         {
2855             int i, j;
2856
2857             cjl_gl = gridj->cell0 + cjl;
2858             for (i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
2859             {
2860                 for (j = 0; j < NBNXN_CPU_CLUSTER_I_SIZE; j++)
2861                 {
2862                     InRange = InRange ||
2863                         (sqr(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
2864                          sqr(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
2865                          sqr(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rl2);
2866                 }
2867             }
2868             *ndistc += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
2869         }
2870         if (!InRange)
2871         {
2872             cjl--;
2873         }
2874     }
2875
2876     if (cjf <= cjl)
2877     {
2878         for (cj = cjf; cj <= cjl; cj++)
2879         {
2880             /* Store cj and the interaction mask */
2881             nbl->cj[nbl->ncj].cj   = gridj->cell0 + cj;
2882             nbl->cj[nbl->ncj].excl = get_imask(remove_sub_diag, ci, cj);
2883             nbl->ncj++;
2884         }
2885         /* Increase the closing index in i super-cell list */
2886         nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
2887     }
2888 }
2889
2890 #ifdef GMX_NBNXN_SIMD_4XN
2891 #include "nbnxn_search_simd_4xn.h"
2892 #endif
2893 #ifdef GMX_NBNXN_SIMD_2XNN
2894 #include "nbnxn_search_simd_2xnn.h"
2895 #endif
2896
2897 /* Plain C or SSE code for making a pair list of super-cell sci vs scj.
2898  * Checks bounding box distances and possibly atom pair distances.
2899  */
2900 static void make_cluster_list_supersub(const nbnxn_grid_t *gridi,
2901                                        const nbnxn_grid_t *gridj,
2902                                        nbnxn_pairlist_t *nbl,
2903                                        int sci, int scj,
2904                                        gmx_bool sci_equals_scj,
2905                                        int stride, const real *x,
2906                                        real rl2, float rbb2,
2907                                        int *ndistc)
2908 {
2909     int          na_c;
2910     int          npair;
2911     int          cjo, ci1, ci, cj, cj_gl;
2912     int          cj4_ind, cj_offset;
2913     unsigned     imask;
2914     nbnxn_cj4_t *cj4;
2915 #ifdef NBNXN_BBXXXX
2916     const float      *pbb_ci;
2917 #else
2918     const nbnxn_bb_t *bb_ci;
2919 #endif
2920     const real  *x_ci;
2921     float       *d2l, d2;
2922     int          w;
2923 #define PRUNE_LIST_CPU_ONE
2924 #ifdef PRUNE_LIST_CPU_ONE
2925     int  ci_last = -1;
2926 #endif
2927
2928     d2l = nbl->work->d2;
2929
2930 #ifdef NBNXN_BBXXXX
2931     pbb_ci = nbl->work->pbb_ci;
2932 #else
2933     bb_ci  = nbl->work->bb_ci;
2934 #endif
2935     x_ci   = nbl->work->x_ci;
2936
2937     na_c = gridj->na_c;
2938
2939     for (cjo = 0; cjo < gridj->nsubc[scj]; cjo++)
2940     {
2941         cj4_ind   = (nbl->work->cj_ind >> NBNXN_GPU_JGROUP_SIZE_2LOG);
2942         cj_offset = nbl->work->cj_ind - cj4_ind*NBNXN_GPU_JGROUP_SIZE;
2943         cj4       = &nbl->cj4[cj4_ind];
2944
2945         cj = scj*GPU_NSUBCELL + cjo;
2946
2947         cj_gl = gridj->cell0*GPU_NSUBCELL + cj;
2948
2949         /* Initialize this j-subcell i-subcell list */
2950         cj4->cj[cj_offset] = cj_gl;
2951         imask              = 0;
2952
2953         if (sci_equals_scj)
2954         {
2955             ci1 = cjo + 1;
2956         }
2957         else
2958         {
2959             ci1 = gridi->nsubc[sci];
2960         }
2961
2962 #ifdef NBNXN_BBXXXX
2963         /* Determine all ci1 bb distances in one call with SSE */
2964         subc_bb_dist2_sse_xxxx(gridj->pbb+(cj>>STRIDE_PBB_2LOG)*NNBSBB_XXXX+(cj & (STRIDE_PBB-1)),
2965                                ci1, pbb_ci, d2l);
2966         *ndistc += na_c*2;
2967 #endif
2968
2969         npair = 0;
2970         /* We use a fixed upper-bound instead of ci1 to help optimization */
2971         for (ci = 0; ci < GPU_NSUBCELL; ci++)
2972         {
2973             if (ci == ci1)
2974             {
2975                 break;
2976             }
2977
2978 #ifndef NBNXN_BBXXXX
2979             /* Determine the bb distance between ci and cj */
2980             d2l[ci]  = subc_bb_dist2(ci, bb_ci, cj, gridj->bb);
2981             *ndistc += 2;
2982 #endif
2983             d2 = d2l[ci];
2984
2985 #ifdef PRUNE_LIST_CPU_ALL
2986             /* Check if the distance is within the distance where
2987              * we use only the bounding box distance rbb,
2988              * or within the cut-off and there is at least one atom pair
2989              * within the cut-off. This check is very costly.
2990              */
2991             *ndistc += na_c*na_c;
2992             if (d2 < rbb2 ||
2993                 (d2 < rl2 &&
2994 #ifdef NBNXN_PBB_SSE
2995                  subc_in_range_sse8
2996 #else
2997                  subc_in_range_x
2998 #endif
2999                      (na_c, ci, x_ci, cj_gl, stride, x, rl2)))
3000 #else
3001             /* Check if the distance between the two bounding boxes
3002              * in within the pair-list cut-off.
3003              */
3004             if (d2 < rl2)
3005 #endif
3006             {
3007                 /* Flag this i-subcell to be taken into account */
3008                 imask |= (1U << (cj_offset*GPU_NSUBCELL+ci));
3009
3010 #ifdef PRUNE_LIST_CPU_ONE
3011                 ci_last = ci;
3012 #endif
3013
3014                 npair++;
3015             }
3016         }
3017
3018 #ifdef PRUNE_LIST_CPU_ONE
3019         /* If we only found 1 pair, check if any atoms are actually
3020          * within the cut-off, so we could get rid of it.
3021          */
3022         if (npair == 1 && d2l[ci_last] >= rbb2)
3023         {
3024             /* Avoid using function pointers here, as it's slower */
3025             if (
3026 #ifdef NBNXN_PBB_SSE
3027                 !subc_in_range_sse8
3028 #else
3029                 !subc_in_range_x
3030 #endif
3031                     (na_c, ci_last, x_ci, cj_gl, stride, x, rl2))
3032             {
3033                 imask &= ~(1U << (cj_offset*GPU_NSUBCELL+ci_last));
3034                 npair--;
3035             }
3036         }
3037 #endif
3038
3039         if (npair > 0)
3040         {
3041             /* We have a useful sj entry, close it now */
3042
3043             /* Set the exclucions for the ci== sj entry.
3044              * Here we don't bother to check if this entry is actually flagged,
3045              * as it will nearly always be in the list.
3046              */
3047             if (sci_equals_scj)
3048             {
3049                 set_self_and_newton_excls_supersub(nbl, cj4_ind, cj_offset, cjo);
3050             }
3051
3052             /* Copy the cluster interaction mask to the list */
3053             for (w = 0; w < NWARP; w++)
3054             {
3055                 cj4->imei[w].imask |= imask;
3056             }
3057
3058             nbl->work->cj_ind++;
3059
3060             /* Keep the count */
3061             nbl->nci_tot += npair;
3062
3063             /* Increase the closing index in i super-cell list */
3064             nbl->sci[nbl->nsci].cj4_ind_end =
3065                 ((nbl->work->cj_ind+NBNXN_GPU_JGROUP_SIZE-1) >> NBNXN_GPU_JGROUP_SIZE_2LOG);
3066         }
3067     }
3068 }
3069
3070 /* Set all atom-pair exclusions from the topology stored in excl
3071  * as masks in the pair-list for simple list i-entry nbl_ci
3072  */
3073 static void set_ci_top_excls(const nbnxn_search_t nbs,
3074                              nbnxn_pairlist_t    *nbl,
3075                              gmx_bool             diagRemoved,
3076                              int                  na_ci_2log,
3077                              int                  na_cj_2log,
3078                              const nbnxn_ci_t    *nbl_ci,
3079                              const t_blocka      *excl)
3080 {
3081     const int    *cell;
3082     int           ci;
3083     int           cj_ind_first, cj_ind_last;
3084     int           cj_first, cj_last;
3085     int           ndirect;
3086     int           i, ai, aj, si, eind, ge, se;
3087     int           found, cj_ind_0, cj_ind_1, cj_ind_m;
3088     int           cj_m;
3089     gmx_bool      Found_si;
3090     int           si_ind;
3091     nbnxn_excl_t *nbl_excl;
3092     int           inner_i, inner_e;
3093
3094     cell = nbs->cell;
3095
3096     if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
3097     {
3098         /* Empty list */
3099         return;
3100     }
3101
3102     ci = nbl_ci->ci;
3103
3104     cj_ind_first = nbl_ci->cj_ind_start;
3105     cj_ind_last  = nbl->ncj - 1;
3106
3107     cj_first = nbl->cj[cj_ind_first].cj;
3108     cj_last  = nbl->cj[cj_ind_last].cj;
3109
3110     /* Determine how many contiguous j-cells we have starting
3111      * from the first i-cell. This number can be used to directly
3112      * calculate j-cell indices for excluded atoms.
3113      */
3114     ndirect = 0;
3115     if (na_ci_2log == na_cj_2log)
3116     {
3117         while (cj_ind_first + ndirect <= cj_ind_last &&
3118                nbl->cj[cj_ind_first+ndirect].cj == ci + ndirect)
3119         {
3120             ndirect++;
3121         }
3122     }
3123 #ifdef NBNXN_SEARCH_BB_SSE
3124     else
3125     {
3126         while (cj_ind_first + ndirect <= cj_ind_last &&
3127                nbl->cj[cj_ind_first+ndirect].cj == ci_to_cj(na_cj_2log, ci) + ndirect)
3128         {
3129             ndirect++;
3130         }
3131     }
3132 #endif
3133
3134     /* Loop over the atoms in the i super-cell */
3135     for (i = 0; i < nbl->na_sc; i++)
3136     {
3137         ai = nbs->a[ci*nbl->na_sc+i];
3138         if (ai >= 0)
3139         {
3140             si  = (i>>na_ci_2log);
3141
3142             /* Loop over the topology-based exclusions for this i-atom */
3143             for (eind = excl->index[ai]; eind < excl->index[ai+1]; eind++)
3144             {
3145                 aj = excl->a[eind];
3146
3147                 if (aj == ai)
3148                 {
3149                     /* The self exclusion are already set, save some time */
3150                     continue;
3151                 }
3152
3153                 ge = cell[aj];
3154
3155                 /* Without shifts we only calculate interactions j>i
3156                  * for one-way pair-lists.
3157                  */
3158                 if (diagRemoved && ge <= ci*nbl->na_sc + i)
3159                 {
3160                     continue;
3161                 }
3162
3163                 se = (ge >> na_cj_2log);
3164
3165                 /* Could the cluster se be in our list? */
3166                 if (se >= cj_first && se <= cj_last)
3167                 {
3168                     if (se < cj_first + ndirect)
3169                     {
3170                         /* We can calculate cj_ind directly from se */
3171                         found = cj_ind_first + se - cj_first;
3172                     }
3173                     else
3174                     {
3175                         /* Search for se using bisection */
3176                         found    = -1;
3177                         cj_ind_0 = cj_ind_first + ndirect;
3178                         cj_ind_1 = cj_ind_last + 1;
3179                         while (found == -1 && cj_ind_0 < cj_ind_1)
3180                         {
3181                             cj_ind_m = (cj_ind_0 + cj_ind_1)>>1;
3182
3183                             cj_m = nbl->cj[cj_ind_m].cj;
3184
3185                             if (se == cj_m)
3186                             {
3187                                 found = cj_ind_m;
3188                             }
3189                             else if (se < cj_m)
3190                             {
3191                                 cj_ind_1 = cj_ind_m;
3192                             }
3193                             else
3194                             {
3195                                 cj_ind_0 = cj_ind_m + 1;
3196                             }
3197                         }
3198                     }
3199
3200                     if (found >= 0)
3201                     {
3202                         inner_i = i  - (si << na_ci_2log);
3203                         inner_e = ge - (se << na_cj_2log);
3204
3205                         nbl->cj[found].excl &= ~(1U<<((inner_i<<na_cj_2log) + inner_e));
3206                     }
3207                 }
3208             }
3209         }
3210     }
3211 }
3212
3213 /* Set all atom-pair exclusions from the topology stored in excl
3214  * as masks in the pair-list for i-super-cell entry nbl_sci
3215  */
3216 static void set_sci_top_excls(const nbnxn_search_t nbs,
3217                               nbnxn_pairlist_t    *nbl,
3218                               gmx_bool             diagRemoved,
3219                               int                  na_c_2log,
3220                               const nbnxn_sci_t   *nbl_sci,
3221                               const t_blocka      *excl)
3222 {
3223     const int    *cell;
3224     int           na_c;
3225     int           sci;
3226     int           cj_ind_first, cj_ind_last;
3227     int           cj_first, cj_last;
3228     int           ndirect;
3229     int           i, ai, aj, si, eind, ge, se;
3230     int           found, cj_ind_0, cj_ind_1, cj_ind_m;
3231     int           cj_m;
3232     gmx_bool      Found_si;
3233     int           si_ind;
3234     nbnxn_excl_t *nbl_excl;
3235     int           inner_i, inner_e, w;
3236
3237     cell = nbs->cell;
3238
3239     na_c = nbl->na_ci;
3240
3241     if (nbl_sci->cj4_ind_end == nbl_sci->cj4_ind_start)
3242     {
3243         /* Empty list */
3244         return;
3245     }
3246
3247     sci = nbl_sci->sci;
3248
3249     cj_ind_first = nbl_sci->cj4_ind_start*NBNXN_GPU_JGROUP_SIZE;
3250     cj_ind_last  = nbl->work->cj_ind - 1;
3251
3252     cj_first = nbl->cj4[nbl_sci->cj4_ind_start].cj[0];
3253     cj_last  = nbl_cj(nbl, cj_ind_last);
3254
3255     /* Determine how many contiguous j-clusters we have starting
3256      * from the first i-cluster. This number can be used to directly
3257      * calculate j-cluster indices for excluded atoms.
3258      */
3259     ndirect = 0;
3260     while (cj_ind_first + ndirect <= cj_ind_last &&
3261            nbl_cj(nbl, cj_ind_first+ndirect) == sci*GPU_NSUBCELL + ndirect)
3262     {
3263         ndirect++;
3264     }
3265
3266     /* Loop over the atoms in the i super-cell */
3267     for (i = 0; i < nbl->na_sc; i++)
3268     {
3269         ai = nbs->a[sci*nbl->na_sc+i];
3270         if (ai >= 0)
3271         {
3272             si  = (i>>na_c_2log);
3273
3274             /* Loop over the topology-based exclusions for this i-atom */
3275             for (eind = excl->index[ai]; eind < excl->index[ai+1]; eind++)
3276             {
3277                 aj = excl->a[eind];
3278
3279                 if (aj == ai)
3280                 {
3281                     /* The self exclusion are already set, save some time */
3282                     continue;
3283                 }
3284
3285                 ge = cell[aj];
3286
3287                 /* Without shifts we only calculate interactions j>i
3288                  * for one-way pair-lists.
3289                  */
3290                 if (diagRemoved && ge <= sci*nbl->na_sc + i)
3291                 {
3292                     continue;
3293                 }
3294
3295                 se = ge>>na_c_2log;
3296                 /* Could the cluster se be in our list? */
3297                 if (se >= cj_first && se <= cj_last)
3298                 {
3299                     if (se < cj_first + ndirect)
3300                     {
3301                         /* We can calculate cj_ind directly from se */
3302                         found = cj_ind_first + se - cj_first;
3303                     }
3304                     else
3305                     {
3306                         /* Search for se using bisection */
3307                         found    = -1;
3308                         cj_ind_0 = cj_ind_first + ndirect;
3309                         cj_ind_1 = cj_ind_last + 1;
3310                         while (found == -1 && cj_ind_0 < cj_ind_1)
3311                         {
3312                             cj_ind_m = (cj_ind_0 + cj_ind_1)>>1;
3313
3314                             cj_m = nbl_cj(nbl, cj_ind_m);
3315
3316                             if (se == cj_m)
3317                             {
3318                                 found = cj_ind_m;
3319                             }
3320                             else if (se < cj_m)
3321                             {
3322                                 cj_ind_1 = cj_ind_m;
3323                             }
3324                             else
3325                             {
3326                                 cj_ind_0 = cj_ind_m + 1;
3327                             }
3328                         }
3329                     }
3330
3331                     if (found >= 0)
3332                     {
3333                         inner_i = i  - si*na_c;
3334                         inner_e = ge - se*na_c;
3335
3336 /* Macro for getting the index of atom a within a cluster */
3337 #define AMODCJ4(a)  ((a) & (NBNXN_GPU_JGROUP_SIZE - 1))
3338 /* Macro for converting an atom number to a cluster number */
3339 #define A2CJ4(a)    ((a) >> NBNXN_GPU_JGROUP_SIZE_2LOG)
3340 /* Macro for getting the index of an i-atom within a warp */
3341 #define AMODWI(a)   ((a) & (NBNXN_GPU_CLUSTER_SIZE/2 - 1))
3342
3343                         if (nbl_imask0(nbl, found) & (1U << (AMODCJ4(found)*GPU_NSUBCELL + si)))
3344                         {
3345                             w       = (inner_e >> 2);
3346
3347                             get_nbl_exclusions_1(nbl, A2CJ4(found), w, &nbl_excl);
3348
3349                             nbl_excl->pair[AMODWI(inner_e)*nbl->na_ci+inner_i] &=
3350                                 ~(1U << (AMODCJ4(found)*GPU_NSUBCELL + si));
3351                         }
3352
3353 #undef AMODCJ4
3354 #undef A2CJ4
3355 #undef AMODWI
3356                     }
3357                 }
3358             }
3359         }
3360     }
3361 }
3362
3363 /* Reallocate the simple ci list for at least n entries */
3364 static void nb_realloc_ci(nbnxn_pairlist_t *nbl, int n)
3365 {
3366     nbl->ci_nalloc = over_alloc_small(n);
3367     nbnxn_realloc_void((void **)&nbl->ci,
3368                        nbl->nci*sizeof(*nbl->ci),
3369                        nbl->ci_nalloc*sizeof(*nbl->ci),
3370                        nbl->alloc, nbl->free);
3371 }
3372
3373 /* Reallocate the super-cell sci list for at least n entries */
3374 static void nb_realloc_sci(nbnxn_pairlist_t *nbl, int n)
3375 {
3376     nbl->sci_nalloc = over_alloc_small(n);
3377     nbnxn_realloc_void((void **)&nbl->sci,
3378                        nbl->nsci*sizeof(*nbl->sci),
3379                        nbl->sci_nalloc*sizeof(*nbl->sci),
3380                        nbl->alloc, nbl->free);
3381 }
3382
3383 /* Make a new ci entry at index nbl->nci */
3384 static void new_ci_entry(nbnxn_pairlist_t *nbl, int ci, int shift, int flags)
3385 {
3386     if (nbl->nci + 1 > nbl->ci_nalloc)
3387     {
3388         nb_realloc_ci(nbl, nbl->nci+1);
3389     }
3390     nbl->ci[nbl->nci].ci            = ci;
3391     nbl->ci[nbl->nci].shift         = shift;
3392     /* Store the interaction flags along with the shift */
3393     nbl->ci[nbl->nci].shift        |= flags;
3394     nbl->ci[nbl->nci].cj_ind_start  = nbl->ncj;
3395     nbl->ci[nbl->nci].cj_ind_end    = nbl->ncj;
3396 }
3397
3398 /* Make a new sci entry at index nbl->nsci */
3399 static void new_sci_entry(nbnxn_pairlist_t *nbl, int sci, int shift)
3400 {
3401     if (nbl->nsci + 1 > nbl->sci_nalloc)
3402     {
3403         nb_realloc_sci(nbl, nbl->nsci+1);
3404     }
3405     nbl->sci[nbl->nsci].sci           = sci;
3406     nbl->sci[nbl->nsci].shift         = shift;
3407     nbl->sci[nbl->nsci].cj4_ind_start = nbl->ncj4;
3408     nbl->sci[nbl->nsci].cj4_ind_end   = nbl->ncj4;
3409 }
3410
3411 /* Sort the simple j-list cj on exclusions.
3412  * Entries with exclusions will all be sorted to the beginning of the list.
3413  */
3414 static void sort_cj_excl(nbnxn_cj_t *cj, int ncj,
3415                          nbnxn_list_work_t *work)
3416 {
3417     int jnew, j;
3418
3419     if (ncj > work->cj_nalloc)
3420     {
3421         work->cj_nalloc = over_alloc_large(ncj);
3422         srenew(work->cj, work->cj_nalloc);
3423     }
3424
3425     /* Make a list of the j-cells involving exclusions */
3426     jnew = 0;
3427     for (j = 0; j < ncj; j++)
3428     {
3429         if (cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
3430         {
3431             work->cj[jnew++] = cj[j];
3432         }
3433     }
3434     /* Check if there are exclusions at all or not just the first entry */
3435     if (!((jnew == 0) ||
3436           (jnew == 1 && cj[0].excl != NBNXN_INTERACTION_MASK_ALL)))
3437     {
3438         for (j = 0; j < ncj; j++)
3439         {
3440             if (cj[j].excl == NBNXN_INTERACTION_MASK_ALL)
3441             {
3442                 work->cj[jnew++] = cj[j];
3443             }
3444         }
3445         for (j = 0; j < ncj; j++)
3446         {
3447             cj[j] = work->cj[j];
3448         }
3449     }
3450 }
3451
3452 /* Close this simple list i entry */
3453 static void close_ci_entry_simple(nbnxn_pairlist_t *nbl)
3454 {
3455     int jlen;
3456
3457     /* All content of the new ci entry have already been filled correctly,
3458      * we only need to increase the count here (for non empty lists).
3459      */
3460     jlen = nbl->ci[nbl->nci].cj_ind_end - nbl->ci[nbl->nci].cj_ind_start;
3461     if (jlen > 0)
3462     {
3463         sort_cj_excl(nbl->cj+nbl->ci[nbl->nci].cj_ind_start, jlen, nbl->work);
3464
3465         /* The counts below are used for non-bonded pair/flop counts
3466          * and should therefore match the available kernel setups.
3467          */
3468         if (!(nbl->ci[nbl->nci].shift & NBNXN_CI_DO_COUL(0)))
3469         {
3470             nbl->work->ncj_noq += jlen;
3471         }
3472         else if ((nbl->ci[nbl->nci].shift & NBNXN_CI_HALF_LJ(0)) ||
3473                  !(nbl->ci[nbl->nci].shift & NBNXN_CI_DO_LJ(0)))
3474         {
3475             nbl->work->ncj_hlj += jlen;
3476         }
3477
3478         nbl->nci++;
3479     }
3480 }
3481
3482 /* Split sci entry for load balancing on the GPU.
3483  * Splitting ensures we have enough lists to fully utilize the whole GPU.
3484  * With progBal we generate progressively smaller lists, which improves
3485  * load balancing. As we only know the current count on our own thread,
3486  * we will need to estimate the current total amount of i-entries.
3487  * As the lists get concatenated later, this estimate depends
3488  * both on nthread and our own thread index.
3489  */
3490 static void split_sci_entry(nbnxn_pairlist_t *nbl,
3491                             int nsp_max_av, gmx_bool progBal, int nc_bal,
3492                             int thread, int nthread)
3493 {
3494     int nsci_est;
3495     int nsp_max;
3496     int cj4_start, cj4_end, j4len, cj4;
3497     int sci;
3498     int nsp, nsp_sci, nsp_cj4, nsp_cj4_e, nsp_cj4_p;
3499     int p;
3500
3501     if (progBal)
3502     {
3503         /* Estimate the total numbers of ci's of the nblist combined
3504          * over all threads using the target number of ci's.
3505          */
3506         nsci_est = nc_bal*thread/nthread + nbl->nsci;
3507
3508         /* The first ci blocks should be larger, to avoid overhead.
3509          * The last ci blocks should be smaller, to improve load balancing.
3510          */
3511         nsp_max = max(1,
3512                       nsp_max_av*nc_bal*3/(2*(nsci_est - 1 + nc_bal)));
3513     }
3514     else
3515     {
3516         nsp_max = nsp_max_av;
3517     }
3518
3519     cj4_start = nbl->sci[nbl->nsci-1].cj4_ind_start;
3520     cj4_end   = nbl->sci[nbl->nsci-1].cj4_ind_end;
3521     j4len     = cj4_end - cj4_start;
3522
3523     if (j4len > 1 && j4len*GPU_NSUBCELL*NBNXN_GPU_JGROUP_SIZE > nsp_max)
3524     {
3525         /* Remove the last ci entry and process the cj4's again */
3526         nbl->nsci -= 1;
3527
3528         sci        = nbl->nsci;
3529         nsp        = 0;
3530         nsp_sci    = 0;
3531         nsp_cj4_e  = 0;
3532         nsp_cj4    = 0;
3533         for (cj4 = cj4_start; cj4 < cj4_end; cj4++)
3534         {
3535             nsp_cj4_p = nsp_cj4;
3536             /* Count the number of cluster pairs in this cj4 group */
3537             nsp_cj4   = 0;
3538             for (p = 0; p < GPU_NSUBCELL*NBNXN_GPU_JGROUP_SIZE; p++)
3539             {
3540                 nsp_cj4 += (nbl->cj4[cj4].imei[0].imask >> p) & 1;
3541             }
3542
3543             if (nsp_cj4 > 0 && nsp + nsp_cj4 > nsp_max)
3544             {
3545                 /* Split the list at cj4 */
3546                 nbl->sci[sci].cj4_ind_end = cj4;
3547                 /* Create a new sci entry */
3548                 sci++;
3549                 nbl->nsci++;
3550                 if (nbl->nsci+1 > nbl->sci_nalloc)
3551                 {
3552                     nb_realloc_sci(nbl, nbl->nsci+1);
3553                 }
3554                 nbl->sci[sci].sci           = nbl->sci[nbl->nsci-1].sci;
3555                 nbl->sci[sci].shift         = nbl->sci[nbl->nsci-1].shift;
3556                 nbl->sci[sci].cj4_ind_start = cj4;
3557                 nsp_sci                     = nsp;
3558                 nsp_cj4_e                   = nsp_cj4_p;
3559                 nsp                         = 0;
3560             }
3561             nsp += nsp_cj4;
3562         }
3563
3564         /* Put the remaining cj4's in the last sci entry */
3565         nbl->sci[sci].cj4_ind_end = cj4_end;
3566
3567         /* Possibly balance out the last two sci's
3568          * by moving the last cj4 of the second last sci.
3569          */
3570         if (nsp_sci - nsp_cj4_e >= nsp + nsp_cj4_e)
3571         {
3572             nbl->sci[sci-1].cj4_ind_end--;
3573             nbl->sci[sci].cj4_ind_start--;
3574         }
3575
3576         nbl->nsci++;
3577     }
3578 }
3579
3580 /* Clost this super/sub list i entry */
3581 static void close_ci_entry_supersub(nbnxn_pairlist_t *nbl,
3582                                     int nsp_max_av,
3583                                     gmx_bool progBal, int nc_bal,
3584                                     int thread, int nthread)
3585 {
3586     int j4len, tlen;
3587     int nb, b;
3588
3589     /* All content of the new ci entry have already been filled correctly,
3590      * we only need to increase the count here (for non empty lists).
3591      */
3592     j4len = nbl->sci[nbl->nsci].cj4_ind_end - nbl->sci[nbl->nsci].cj4_ind_start;
3593     if (j4len > 0)
3594     {
3595         /* We can only have complete blocks of 4 j-entries in a list,
3596          * so round the count up before closing.
3597          */
3598         nbl->ncj4         = ((nbl->work->cj_ind + NBNXN_GPU_JGROUP_SIZE - 1) >> NBNXN_GPU_JGROUP_SIZE_2LOG);
3599         nbl->work->cj_ind = nbl->ncj4*NBNXN_GPU_JGROUP_SIZE;
3600
3601         nbl->nsci++;
3602
3603         if (nsp_max_av > 0)
3604         {
3605             /* Measure the size of the new entry and potentially split it */
3606             split_sci_entry(nbl, nsp_max_av, progBal, nc_bal, thread, nthread);
3607         }
3608     }
3609 }
3610
3611 /* Syncs the working array before adding another grid pair to the list */
3612 static void sync_work(nbnxn_pairlist_t *nbl)
3613 {
3614     if (!nbl->bSimple)
3615     {
3616         nbl->work->cj_ind   = nbl->ncj4*NBNXN_GPU_JGROUP_SIZE;
3617         nbl->work->cj4_init = nbl->ncj4;
3618     }
3619 }
3620
3621 /* Clears an nbnxn_pairlist_t data structure */
3622 static void clear_pairlist(nbnxn_pairlist_t *nbl)
3623 {
3624     nbl->nci           = 0;
3625     nbl->nsci          = 0;
3626     nbl->ncj           = 0;
3627     nbl->ncj4          = 0;
3628     nbl->nci_tot       = 0;
3629     nbl->nexcl         = 1;
3630
3631     nbl->work->ncj_noq = 0;
3632     nbl->work->ncj_hlj = 0;
3633 }
3634
3635 /* Sets a simple list i-cell bounding box, including PBC shift */
3636 static gmx_inline void set_icell_bb_simple(const nbnxn_bb_t *bb, int ci,
3637                                            real shx, real shy, real shz,
3638                                            nbnxn_bb_t *bb_ci)
3639 {
3640     bb_ci->lower[BB_X] = bb[ci].lower[BB_X] + shx;
3641     bb_ci->lower[BB_Y] = bb[ci].lower[BB_Y] + shy;
3642     bb_ci->lower[BB_Z] = bb[ci].lower[BB_Z] + shz;
3643     bb_ci->upper[BB_X] = bb[ci].upper[BB_X] + shx;
3644     bb_ci->upper[BB_Y] = bb[ci].upper[BB_Y] + shy;
3645     bb_ci->upper[BB_Z] = bb[ci].upper[BB_Z] + shz;
3646 }
3647
3648 #ifdef NBNXN_BBXXXX
3649 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
3650 static void set_icell_bbxxxx_supersub(const float *bb, int ci,
3651                                       real shx, real shy, real shz,
3652                                       float *bb_ci)
3653 {
3654     int ia, m, i;
3655
3656     ia = ci*(GPU_NSUBCELL>>STRIDE_PBB_2LOG)*NNBSBB_XXXX;
3657     for (m = 0; m < (GPU_NSUBCELL>>STRIDE_PBB_2LOG)*NNBSBB_XXXX; m += NNBSBB_XXXX)
3658     {
3659         for (i = 0; i < STRIDE_PBB; i++)
3660         {
3661             bb_ci[m+0*STRIDE_PBB+i] = bb[ia+m+0*STRIDE_PBB+i] + shx;
3662             bb_ci[m+1*STRIDE_PBB+i] = bb[ia+m+1*STRIDE_PBB+i] + shy;
3663             bb_ci[m+2*STRIDE_PBB+i] = bb[ia+m+2*STRIDE_PBB+i] + shz;
3664             bb_ci[m+3*STRIDE_PBB+i] = bb[ia+m+3*STRIDE_PBB+i] + shx;
3665             bb_ci[m+4*STRIDE_PBB+i] = bb[ia+m+4*STRIDE_PBB+i] + shy;
3666             bb_ci[m+5*STRIDE_PBB+i] = bb[ia+m+5*STRIDE_PBB+i] + shz;
3667         }
3668     }
3669 }
3670 #endif
3671
3672 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
3673 static void set_icell_bb_supersub(const nbnxn_bb_t *bb, int ci,
3674                                   real shx, real shy, real shz,
3675                                   nbnxn_bb_t *bb_ci)
3676 {
3677     int i;
3678
3679     for (i = 0; i < GPU_NSUBCELL; i++)
3680     {
3681         set_icell_bb_simple(bb, ci*GPU_NSUBCELL+i,
3682                             shx, shy, shz,
3683                             &bb_ci[i]);
3684     }
3685 }
3686
3687 /* Copies PBC shifted i-cell atom coordinates x,y,z to working array */
3688 static void icell_set_x_simple(int ci,
3689                                real shx, real shy, real shz,
3690                                int gmx_unused na_c,
3691                                int stride, const real *x,
3692                                nbnxn_list_work_t *work)
3693 {
3694     int  ia, i;
3695
3696     ia = ci*NBNXN_CPU_CLUSTER_I_SIZE;
3697
3698     for (i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE; i++)
3699     {
3700         work->x_ci[i*STRIDE_XYZ+XX] = x[(ia+i)*stride+XX] + shx;
3701         work->x_ci[i*STRIDE_XYZ+YY] = x[(ia+i)*stride+YY] + shy;
3702         work->x_ci[i*STRIDE_XYZ+ZZ] = x[(ia+i)*stride+ZZ] + shz;
3703     }
3704 }
3705
3706 /* Copies PBC shifted super-cell atom coordinates x,y,z to working array */
3707 static void icell_set_x_supersub(int ci,
3708                                  real shx, real shy, real shz,
3709                                  int na_c,
3710                                  int stride, const real *x,
3711                                  nbnxn_list_work_t *work)
3712 {
3713     int  ia, i;
3714     real *x_ci;
3715
3716     x_ci = work->x_ci;
3717
3718     ia = ci*GPU_NSUBCELL*na_c;
3719     for (i = 0; i < GPU_NSUBCELL*na_c; i++)
3720     {
3721         x_ci[i*DIM + XX] = x[(ia+i)*stride + XX] + shx;
3722         x_ci[i*DIM + YY] = x[(ia+i)*stride + YY] + shy;
3723         x_ci[i*DIM + ZZ] = x[(ia+i)*stride + ZZ] + shz;
3724     }
3725 }
3726
3727 #ifdef NBNXN_SEARCH_BB_SSE
3728 /* Copies PBC shifted super-cell packed atom coordinates to working array */
3729 static void icell_set_x_supersub_sse8(int ci,
3730                                       real shx, real shy, real shz,
3731                                       int na_c,
3732                                       int stride, const real *x,
3733                                       nbnxn_list_work_t *work)
3734 {
3735     int  si, io, ia, i, j;
3736     real *x_ci;
3737
3738     x_ci = work->x_ci;
3739
3740     for (si = 0; si < GPU_NSUBCELL; si++)
3741     {
3742         for (i = 0; i < na_c; i += STRIDE_PBB)
3743         {
3744             io = si*na_c + i;
3745             ia = ci*GPU_NSUBCELL*na_c + io;
3746             for (j = 0; j < STRIDE_PBB; j++)
3747             {
3748                 x_ci[io*DIM + j + XX*STRIDE_PBB] = x[(ia+j)*stride+XX] + shx;
3749                 x_ci[io*DIM + j + YY*STRIDE_PBB] = x[(ia+j)*stride+YY] + shy;
3750                 x_ci[io*DIM + j + ZZ*STRIDE_PBB] = x[(ia+j)*stride+ZZ] + shz;
3751             }
3752         }
3753     }
3754 }
3755 #endif
3756
3757 static real nbnxn_rlist_inc_nonloc_fac = 0.6;
3758
3759 /* Due to the cluster size the effective pair-list is longer than
3760  * that of a simple atom pair-list. This function gives the extra distance.
3761  */
3762 real nbnxn_get_rlist_effective_inc(int cluster_size, real atom_density)
3763 {
3764     return ((0.5 + nbnxn_rlist_inc_nonloc_fac)*sqr(((cluster_size) - 1.0)/(cluster_size))*pow((cluster_size)/(atom_density), 1.0/3.0));
3765 }
3766
3767 /* Estimates the interaction volume^2 for non-local interactions */
3768 static real nonlocal_vol2(const gmx_domdec_zones_t *zones, rvec ls, real r)
3769 {
3770     int  z, d;
3771     real cl, ca, za;
3772     real vold_est;
3773     real vol2_est_tot;
3774
3775     vol2_est_tot = 0;
3776
3777     /* Here we simply add up the volumes of 1, 2 or 3 1D decomposition
3778      * not home interaction volume^2. As these volumes are not additive,
3779      * this is an overestimate, but it would only be significant in the limit
3780      * of small cells, where we anyhow need to split the lists into
3781      * as small parts as possible.
3782      */
3783
3784     for (z = 0; z < zones->n; z++)
3785     {
3786         if (zones->shift[z][XX] + zones->shift[z][YY] + zones->shift[z][ZZ] == 1)
3787         {
3788             cl = 0;
3789             ca = 1;
3790             za = 1;
3791             for (d = 0; d < DIM; d++)
3792             {
3793                 if (zones->shift[z][d] == 0)
3794                 {
3795                     cl += 0.5*ls[d];
3796                     ca *= ls[d];
3797                     za *= zones->size[z].x1[d] - zones->size[z].x0[d];
3798                 }
3799             }
3800
3801             /* 4 octants of a sphere */
3802             vold_est  = 0.25*M_PI*r*r*r*r;
3803             /* 4 quarter pie slices on the edges */
3804             vold_est += 4*cl*M_PI/6.0*r*r*r;
3805             /* One rectangular volume on a face */
3806             vold_est += ca*0.5*r*r;
3807
3808             vol2_est_tot += vold_est*za;
3809         }
3810     }
3811
3812     return vol2_est_tot;
3813 }
3814
3815 /* Estimates the average size of a full j-list for super/sub setup */
3816 static int get_nsubpair_max(const nbnxn_search_t nbs,
3817                             int                  iloc,
3818                             real                 rlist,
3819                             int                  min_ci_balanced)
3820 {
3821     const nbnxn_grid_t *grid;
3822     rvec ls;
3823     real xy_diag2, r_eff_sup, vol_est, nsp_est, nsp_est_nl;
3824     int  nsubpair_max;
3825
3826     grid = &nbs->grid[0];
3827
3828     ls[XX] = (grid->c1[XX] - grid->c0[XX])/(grid->ncx*GPU_NSUBCELL_X);
3829     ls[YY] = (grid->c1[YY] - grid->c0[YY])/(grid->ncy*GPU_NSUBCELL_Y);
3830     ls[ZZ] = (grid->c1[ZZ] - grid->c0[ZZ])*grid->ncx*grid->ncy/(grid->nc*GPU_NSUBCELL_Z);
3831
3832     /* The average squared length of the diagonal of a sub cell */
3833     xy_diag2 = ls[XX]*ls[XX] + ls[YY]*ls[YY] + ls[ZZ]*ls[ZZ];
3834
3835     /* The formulas below are a heuristic estimate of the average nsj per si*/
3836     r_eff_sup = rlist + nbnxn_rlist_inc_nonloc_fac*sqr((grid->na_c - 1.0)/grid->na_c)*sqrt(xy_diag2/3);
3837
3838     if (!nbs->DomDec || nbs->zones->n == 1)
3839     {
3840         nsp_est_nl = 0;
3841     }
3842     else
3843     {
3844         nsp_est_nl =
3845             sqr(grid->atom_density/grid->na_c)*
3846             nonlocal_vol2(nbs->zones, ls, r_eff_sup);
3847     }
3848
3849     if (LOCAL_I(iloc))
3850     {
3851         /* Sub-cell interacts with itself */
3852         vol_est  = ls[XX]*ls[YY]*ls[ZZ];
3853         /* 6/2 rectangular volume on the faces */
3854         vol_est += (ls[XX]*ls[YY] + ls[XX]*ls[ZZ] + ls[YY]*ls[ZZ])*r_eff_sup;
3855         /* 12/2 quarter pie slices on the edges */
3856         vol_est += 2*(ls[XX] + ls[YY] + ls[ZZ])*0.25*M_PI*sqr(r_eff_sup);
3857         /* 4 octants of a sphere */
3858         vol_est += 0.5*4.0/3.0*M_PI*pow(r_eff_sup, 3);
3859
3860         nsp_est = grid->nsubc_tot*vol_est*grid->atom_density/grid->na_c;
3861
3862         /* Subtract the non-local pair count */
3863         nsp_est -= nsp_est_nl;
3864
3865         if (debug)
3866         {
3867             fprintf(debug, "nsp_est local %5.1f non-local %5.1f\n",
3868                     nsp_est, nsp_est_nl);
3869         }
3870     }
3871     else
3872     {
3873         nsp_est = nsp_est_nl;
3874     }
3875
3876     if (min_ci_balanced <= 0 || grid->nc >= min_ci_balanced || grid->nc == 0)
3877     {
3878         /* We don't need to worry */
3879         nsubpair_max = -1;
3880     }
3881     else
3882     {
3883         /* Thus the (average) maximum j-list size should be as follows */
3884         nsubpair_max = max(1, (int)(nsp_est/min_ci_balanced+0.5));
3885
3886         /* Since the target value is a maximum (this avoids high outliers,
3887          * which lead to load imbalance), not average, we add half the
3888          * number of pairs in a cj4 block to get the average about right.
3889          */
3890         nsubpair_max += GPU_NSUBCELL*NBNXN_GPU_JGROUP_SIZE/2;
3891     }
3892
3893     if (debug)
3894     {
3895         fprintf(debug, "nbl nsp estimate %.1f, nsubpair_max %d\n",
3896                 nsp_est, nsubpair_max);
3897     }
3898
3899     return nsubpair_max;
3900 }
3901
3902 /* Debug list print function */
3903 static void print_nblist_ci_cj(FILE *fp, const nbnxn_pairlist_t *nbl)
3904 {
3905     int i, j;
3906
3907     for (i = 0; i < nbl->nci; i++)
3908     {
3909         fprintf(fp, "ci %4d  shift %2d  ncj %3d\n",
3910                 nbl->ci[i].ci, nbl->ci[i].shift,
3911                 nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start);
3912
3913         for (j = nbl->ci[i].cj_ind_start; j < nbl->ci[i].cj_ind_end; j++)
3914         {
3915             fprintf(fp, "  cj %5d  imask %x\n",
3916                     nbl->cj[j].cj,
3917                     nbl->cj[j].excl);
3918         }
3919     }
3920 }
3921
3922 /* Debug list print function */
3923 static void print_nblist_sci_cj(FILE *fp, const nbnxn_pairlist_t *nbl)
3924 {
3925     int i, j4, j, ncp, si;
3926
3927     for (i = 0; i < nbl->nsci; i++)
3928     {
3929         fprintf(fp, "ci %4d  shift %2d  ncj4 %2d\n",
3930                 nbl->sci[i].sci, nbl->sci[i].shift,
3931                 nbl->sci[i].cj4_ind_end - nbl->sci[i].cj4_ind_start);
3932
3933         ncp = 0;
3934         for (j4 = nbl->sci[i].cj4_ind_start; j4 < nbl->sci[i].cj4_ind_end; j4++)
3935         {
3936             for (j = 0; j < NBNXN_GPU_JGROUP_SIZE; j++)
3937             {
3938                 fprintf(fp, "  sj %5d  imask %x\n",
3939                         nbl->cj4[j4].cj[j],
3940                         nbl->cj4[j4].imei[0].imask);
3941                 for (si = 0; si < GPU_NSUBCELL; si++)
3942                 {
3943                     if (nbl->cj4[j4].imei[0].imask & (1U << (j*GPU_NSUBCELL + si)))
3944                     {
3945                         ncp++;
3946                     }
3947                 }
3948             }
3949         }
3950         fprintf(fp, "ci %4d  shift %2d  ncj4 %2d ncp %3d\n",
3951                 nbl->sci[i].sci, nbl->sci[i].shift,
3952                 nbl->sci[i].cj4_ind_end - nbl->sci[i].cj4_ind_start,
3953                 ncp);
3954     }
3955 }
3956
3957 /* Combine pair lists *nbl generated on multiple threads nblc */
3958 static void combine_nblists(int nnbl, nbnxn_pairlist_t **nbl,
3959                             nbnxn_pairlist_t *nblc)
3960 {
3961     int nsci, ncj4, nexcl;
3962     int n, i;
3963
3964     if (nblc->bSimple)
3965     {
3966         gmx_incons("combine_nblists does not support simple lists");
3967     }
3968
3969     nsci  = nblc->nsci;
3970     ncj4  = nblc->ncj4;
3971     nexcl = nblc->nexcl;
3972     for (i = 0; i < nnbl; i++)
3973     {
3974         nsci  += nbl[i]->nsci;
3975         ncj4  += nbl[i]->ncj4;
3976         nexcl += nbl[i]->nexcl;
3977     }
3978
3979     if (nsci > nblc->sci_nalloc)
3980     {
3981         nb_realloc_sci(nblc, nsci);
3982     }
3983     if (ncj4 > nblc->cj4_nalloc)
3984     {
3985         nblc->cj4_nalloc = over_alloc_small(ncj4);
3986         nbnxn_realloc_void((void **)&nblc->cj4,
3987                            nblc->ncj4*sizeof(*nblc->cj4),
3988                            nblc->cj4_nalloc*sizeof(*nblc->cj4),
3989                            nblc->alloc, nblc->free);
3990     }
3991     if (nexcl > nblc->excl_nalloc)
3992     {
3993         nblc->excl_nalloc = over_alloc_small(nexcl);
3994         nbnxn_realloc_void((void **)&nblc->excl,
3995                            nblc->nexcl*sizeof(*nblc->excl),
3996                            nblc->excl_nalloc*sizeof(*nblc->excl),
3997                            nblc->alloc, nblc->free);
3998     }
3999
4000     /* Each thread should copy its own data to the combined arrays,
4001      * as otherwise data will go back and forth between different caches.
4002      */
4003 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntPairsearch)) schedule(static)
4004     for (n = 0; n < nnbl; n++)
4005     {
4006         int sci_offset;
4007         int cj4_offset;
4008         int ci_offset;
4009         int excl_offset;
4010         int i, j4;
4011         const nbnxn_pairlist_t *nbli;
4012
4013         /* Determine the offset in the combined data for our thread */
4014         sci_offset  = nblc->nsci;
4015         cj4_offset  = nblc->ncj4;
4016         ci_offset   = nblc->nci_tot;
4017         excl_offset = nblc->nexcl;
4018
4019         for (i = 0; i < n; i++)
4020         {
4021             sci_offset  += nbl[i]->nsci;
4022             cj4_offset  += nbl[i]->ncj4;
4023             ci_offset   += nbl[i]->nci_tot;
4024             excl_offset += nbl[i]->nexcl;
4025         }
4026
4027         nbli = nbl[n];
4028
4029         for (i = 0; i < nbli->nsci; i++)
4030         {
4031             nblc->sci[sci_offset+i]                = nbli->sci[i];
4032             nblc->sci[sci_offset+i].cj4_ind_start += cj4_offset;
4033             nblc->sci[sci_offset+i].cj4_ind_end   += cj4_offset;
4034         }
4035
4036         for (j4 = 0; j4 < nbli->ncj4; j4++)
4037         {
4038             nblc->cj4[cj4_offset+j4]                   = nbli->cj4[j4];
4039             nblc->cj4[cj4_offset+j4].imei[0].excl_ind += excl_offset;
4040             nblc->cj4[cj4_offset+j4].imei[1].excl_ind += excl_offset;
4041         }
4042
4043         for (j4 = 0; j4 < nbli->nexcl; j4++)
4044         {
4045             nblc->excl[excl_offset+j4] = nbli->excl[j4];
4046         }
4047     }
4048
4049     for (n = 0; n < nnbl; n++)
4050     {
4051         nblc->nsci    += nbl[n]->nsci;
4052         nblc->ncj4    += nbl[n]->ncj4;
4053         nblc->nci_tot += nbl[n]->nci_tot;
4054         nblc->nexcl   += nbl[n]->nexcl;
4055     }
4056 }
4057
4058 /* Returns the next ci to be processes by our thread */
4059 static gmx_bool next_ci(const nbnxn_grid_t *grid,
4060                         int conv,
4061                         int nth, int ci_block,
4062                         int *ci_x, int *ci_y,
4063                         int *ci_b, int *ci)
4064 {
4065     (*ci_b)++;
4066     (*ci)++;
4067
4068     if (*ci_b == ci_block)
4069     {
4070         /* Jump to the next block assigned to this task */
4071         *ci   += (nth - 1)*ci_block;
4072         *ci_b  = 0;
4073     }
4074
4075     if (*ci >= grid->nc*conv)
4076     {
4077         return FALSE;
4078     }
4079
4080     while (*ci >= grid->cxy_ind[*ci_x*grid->ncy + *ci_y + 1]*conv)
4081     {
4082         *ci_y += 1;
4083         if (*ci_y == grid->ncy)
4084         {
4085             *ci_x += 1;
4086             *ci_y  = 0;
4087         }
4088     }
4089
4090     return TRUE;
4091 }
4092
4093 /* Returns the distance^2 for which we put cell pairs in the list
4094  * without checking atom pair distances. This is usually < rlist^2.
4095  */
4096 static float boundingbox_only_distance2(const nbnxn_grid_t *gridi,
4097                                         const nbnxn_grid_t *gridj,
4098                                         real                rlist,
4099                                         gmx_bool            simple)
4100 {
4101     /* If the distance between two sub-cell bounding boxes is less
4102      * than this distance, do not check the distance between
4103      * all particle pairs in the sub-cell, since then it is likely
4104      * that the box pair has atom pairs within the cut-off.
4105      * We use the nblist cut-off minus 0.5 times the average x/y diagonal
4106      * spacing of the sub-cells. Around 40% of the checked pairs are pruned.
4107      * Using more than 0.5 gains at most 0.5%.
4108      * If forces are calculated more than twice, the performance gain
4109      * in the force calculation outweighs the cost of checking.
4110      * Note that with subcell lists, the atom-pair distance check
4111      * is only performed when only 1 out of 8 sub-cells in within range,
4112      * this is because the GPU is much faster than the cpu.
4113      */
4114     real bbx, bby;
4115     real rbb2;
4116
4117     bbx = 0.5*(gridi->sx + gridj->sx);
4118     bby = 0.5*(gridi->sy + gridj->sy);
4119     if (!simple)
4120     {
4121         bbx /= GPU_NSUBCELL_X;
4122         bby /= GPU_NSUBCELL_Y;
4123     }
4124
4125     rbb2 = sqr(max(0, rlist - 0.5*sqrt(bbx*bbx + bby*bby)));
4126
4127 #ifndef GMX_DOUBLE
4128     return rbb2;
4129 #else
4130     return (float)((1+GMX_FLOAT_EPS)*rbb2);
4131 #endif
4132 }
4133
4134 static int get_ci_block_size(const nbnxn_grid_t *gridi,
4135                              gmx_bool bDomDec, int nth)
4136 {
4137     const int ci_block_enum      = 5;
4138     const int ci_block_denom     = 11;
4139     const int ci_block_min_atoms = 16;
4140     int ci_block;
4141
4142     /* Here we decide how to distribute the blocks over the threads.
4143      * We use prime numbers to try to avoid that the grid size becomes
4144      * a multiple of the number of threads, which would lead to some
4145      * threads getting "inner" pairs and others getting boundary pairs,
4146      * which in turns will lead to load imbalance between threads.
4147      * Set the block size as 5/11/ntask times the average number of cells
4148      * in a y,z slab. This should ensure a quite uniform distribution
4149      * of the grid parts of the different thread along all three grid
4150      * zone boundaries with 3D domain decomposition. At the same time
4151      * the blocks will not become too small.
4152      */
4153     ci_block = (gridi->nc*ci_block_enum)/(ci_block_denom*gridi->ncx*nth);
4154
4155     /* Ensure the blocks are not too small: avoids cache invalidation */
4156     if (ci_block*gridi->na_sc < ci_block_min_atoms)
4157     {
4158         ci_block = (ci_block_min_atoms + gridi->na_sc - 1)/gridi->na_sc;
4159     }
4160
4161     /* Without domain decomposition
4162      * or with less than 3 blocks per task, divide in nth blocks.
4163      */
4164     if (!bDomDec || ci_block*3*nth > gridi->nc)
4165     {
4166         ci_block = (gridi->nc + nth - 1)/nth;
4167     }
4168
4169     return ci_block;
4170 }
4171
4172 /* Generates the part of pair-list nbl assigned to our thread */
4173 static void nbnxn_make_pairlist_part(const nbnxn_search_t nbs,
4174                                      const nbnxn_grid_t *gridi,
4175                                      const nbnxn_grid_t *gridj,
4176                                      nbnxn_search_work_t *work,
4177                                      const nbnxn_atomdata_t *nbat,
4178                                      const t_blocka *excl,
4179                                      real rlist,
4180                                      int nb_kernel_type,
4181                                      int ci_block,
4182                                      gmx_bool bFBufferFlag,
4183                                      int nsubpair_max,
4184                                      gmx_bool progBal,
4185                                      int min_ci_balanced,
4186                                      int th, int nth,
4187                                      nbnxn_pairlist_t *nbl)
4188 {
4189     int  na_cj_2log;
4190     matrix box;
4191     real rl2;
4192     float rbb2;
4193     int  d;
4194     int  ci_b, ci, ci_x, ci_y, ci_xy, cj;
4195     ivec shp;
4196     int  tx, ty, tz;
4197     int  shift;
4198     gmx_bool bMakeList;
4199     real shx, shy, shz;
4200     int  conv_i, cell0_i;
4201     const nbnxn_bb_t *bb_i=NULL;
4202 #ifdef NBNXN_BBXXXX
4203     const float *pbb_i=NULL;
4204 #endif
4205     const float *bbcz_i, *bbcz_j;
4206     const int *flags_i;
4207     real bx0, bx1, by0, by1, bz0, bz1;
4208     real bz1_frac;
4209     real d2cx, d2z, d2z_cx, d2z_cy, d2zx, d2zxy, d2xy;
4210     int  cxf, cxl, cyf, cyf_x, cyl;
4211     int  cx, cy;
4212     int  c0, c1, cs, cf, cl;
4213     int  ndistc;
4214     int  ncpcheck;
4215     int  gridi_flag_shift = 0, gridj_flag_shift = 0;
4216     unsigned *gridj_flag  = NULL;
4217     int  ncj_old_i, ncj_old_j;
4218
4219     nbs_cycle_start(&work->cc[enbsCCsearch]);
4220
4221     if (gridj->bSimple != nbl->bSimple)
4222     {
4223         gmx_incons("Grid incompatible with pair-list");
4224     }
4225
4226     sync_work(nbl);
4227     nbl->na_sc = gridj->na_sc;
4228     nbl->na_ci = gridj->na_c;
4229     nbl->na_cj = nbnxn_kernel_to_cj_size(nb_kernel_type);
4230     na_cj_2log = get_2log(nbl->na_cj);
4231
4232     nbl->rlist  = rlist;
4233
4234     if (bFBufferFlag)
4235     {
4236         /* Determine conversion of clusters to flag blocks */
4237         gridi_flag_shift = 0;
4238         while ((nbl->na_ci<<gridi_flag_shift) < NBNXN_BUFFERFLAG_SIZE)
4239         {
4240             gridi_flag_shift++;
4241         }
4242         gridj_flag_shift = 0;
4243         while ((nbl->na_cj<<gridj_flag_shift) < NBNXN_BUFFERFLAG_SIZE)
4244         {
4245             gridj_flag_shift++;
4246         }
4247
4248         gridj_flag = work->buffer_flags.flag;
4249     }
4250
4251     copy_mat(nbs->box, box);
4252
4253     rl2 = nbl->rlist*nbl->rlist;
4254
4255     rbb2 = boundingbox_only_distance2(gridi, gridj, nbl->rlist, nbl->bSimple);
4256
4257     if (debug)
4258     {
4259         fprintf(debug, "nbl bounding box only distance %f\n", sqrt(rbb2));
4260     }
4261
4262     /* Set the shift range */
4263     for (d = 0; d < DIM; d++)
4264     {
4265         /* Check if we need periodicity shifts.
4266          * Without PBC or with domain decomposition we don't need them.
4267          */
4268         if (d >= ePBC2npbcdim(nbs->ePBC) || nbs->dd_dim[d])
4269         {
4270             shp[d] = 0;
4271         }
4272         else
4273         {
4274             if (d == XX &&
4275                 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
4276             {
4277                 shp[d] = 2;
4278             }
4279             else
4280             {
4281                 shp[d] = 1;
4282             }
4283         }
4284     }
4285
4286     if (nbl->bSimple && !gridi->bSimple)
4287     {
4288         conv_i  = gridi->na_sc/gridj->na_sc;
4289         bb_i    = gridi->bb_simple;
4290         bbcz_i  = gridi->bbcz_simple;
4291         flags_i = gridi->flags_simple;
4292     }
4293     else
4294     {
4295         conv_i  = 1;
4296 #ifdef NBNXN_BBXXXX
4297         if (gridi->bSimple)
4298         {
4299             bb_i  = gridi->bb;
4300         }
4301         else
4302         {
4303             pbb_i = gridi->pbb;
4304         }
4305 #else
4306         /* We use the normal bounding box format for both grid types */
4307         bb_i  = gridi->bb;
4308 #endif
4309         bbcz_i  = gridi->bbcz;
4310         flags_i = gridi->flags;
4311     }
4312     cell0_i = gridi->cell0*conv_i;
4313
4314     bbcz_j = gridj->bbcz;
4315
4316     if (conv_i != 1)
4317     {
4318         /* Blocks of the conversion factor - 1 give a large repeat count
4319          * combined with a small block size. This should result in good
4320          * load balancing for both small and large domains.
4321          */
4322         ci_block = conv_i - 1;
4323     }
4324     if (debug)
4325     {
4326         fprintf(debug, "nbl nc_i %d col.av. %.1f ci_block %d\n",
4327                 gridi->nc, gridi->nc/(double)(gridi->ncx*gridi->ncy), ci_block);
4328     }
4329
4330     ndistc   = 0;
4331     ncpcheck = 0;
4332
4333     /* Initially ci_b and ci to 1 before where we want them to start,
4334      * as they will both be incremented in next_ci.
4335      */
4336     ci_b = -1;
4337     ci   = th*ci_block - 1;
4338     ci_x = 0;
4339     ci_y = 0;
4340     while (next_ci(gridi, conv_i, nth, ci_block, &ci_x, &ci_y, &ci_b, &ci))
4341     {
4342         if (nbl->bSimple && flags_i[ci] == 0)
4343         {
4344             continue;
4345         }
4346
4347         ncj_old_i = nbl->ncj;
4348
4349         d2cx = 0;
4350         if (gridj != gridi && shp[XX] == 0)
4351         {
4352             if (nbl->bSimple)
4353             {
4354                 bx1 = bb_i[ci].upper[BB_X];
4355             }
4356             else
4357             {
4358                 bx1 = gridi->c0[XX] + (ci_x+1)*gridi->sx;
4359             }
4360             if (bx1 < gridj->c0[XX])
4361             {
4362                 d2cx = sqr(gridj->c0[XX] - bx1);
4363
4364                 if (d2cx >= rl2)
4365                 {
4366                     continue;
4367                 }
4368             }
4369         }
4370
4371         ci_xy = ci_x*gridi->ncy + ci_y;
4372
4373         /* Loop over shift vectors in three dimensions */
4374         for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
4375         {
4376             shz = tz*box[ZZ][ZZ];
4377
4378             bz0 = bbcz_i[ci*NNBSBB_D  ] + shz;
4379             bz1 = bbcz_i[ci*NNBSBB_D+1] + shz;
4380
4381             if (tz == 0)
4382             {
4383                 d2z = 0;
4384             }
4385             else if (tz < 0)
4386             {
4387                 d2z = sqr(bz1);
4388             }
4389             else
4390             {
4391                 d2z = sqr(bz0 - box[ZZ][ZZ]);
4392             }
4393
4394             d2z_cx = d2z + d2cx;
4395
4396             if (d2z_cx >= rl2)
4397             {
4398                 continue;
4399             }
4400
4401             bz1_frac =
4402                 bz1/((real)(gridi->cxy_ind[ci_xy+1] - gridi->cxy_ind[ci_xy]));
4403             if (bz1_frac < 0)
4404             {
4405                 bz1_frac = 0;
4406             }
4407             /* The check with bz1_frac close to or larger than 1 comes later */
4408
4409             for (ty = -shp[YY]; ty <= shp[YY]; ty++)
4410             {
4411                 shy = ty*box[YY][YY] + tz*box[ZZ][YY];
4412
4413                 if (nbl->bSimple)
4414                 {
4415                     by0 = bb_i[ci].lower[BB_Y] + shy;
4416                     by1 = bb_i[ci].upper[BB_Y] + shy;
4417                 }
4418                 else
4419                 {
4420                     by0 = gridi->c0[YY] + (ci_y  )*gridi->sy + shy;
4421                     by1 = gridi->c0[YY] + (ci_y+1)*gridi->sy + shy;
4422                 }
4423
4424                 get_cell_range(by0, by1,
4425                                gridj->ncy, gridj->c0[YY], gridj->sy, gridj->inv_sy,
4426                                d2z_cx, rl2,
4427                                &cyf, &cyl);
4428
4429                 if (cyf > cyl)
4430                 {
4431                     continue;
4432                 }
4433
4434                 d2z_cy = d2z;
4435                 if (by1 < gridj->c0[YY])
4436                 {
4437                     d2z_cy += sqr(gridj->c0[YY] - by1);
4438                 }
4439                 else if (by0 > gridj->c1[YY])
4440                 {
4441                     d2z_cy += sqr(by0 - gridj->c1[YY]);
4442                 }
4443
4444                 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
4445                 {
4446                     shift = XYZ2IS(tx, ty, tz);
4447
4448 #ifdef NBNXN_SHIFT_BACKWARD
4449                     if (gridi == gridj && shift > CENTRAL)
4450                     {
4451                         continue;
4452                     }
4453 #endif
4454
4455                     shx = tx*box[XX][XX] + ty*box[YY][XX] + tz*box[ZZ][XX];
4456
4457                     if (nbl->bSimple)
4458                     {
4459                         bx0 = bb_i[ci].lower[BB_X] + shx;
4460                         bx1 = bb_i[ci].upper[BB_X] + shx;
4461                     }
4462                     else
4463                     {
4464                         bx0 = gridi->c0[XX] + (ci_x  )*gridi->sx + shx;
4465                         bx1 = gridi->c0[XX] + (ci_x+1)*gridi->sx + shx;
4466                     }
4467
4468                     get_cell_range(bx0, bx1,
4469                                    gridj->ncx, gridj->c0[XX], gridj->sx, gridj->inv_sx,
4470                                    d2z_cy, rl2,
4471                                    &cxf, &cxl);
4472
4473                     if (cxf > cxl)
4474                     {
4475                         continue;
4476                     }
4477
4478                     if (nbl->bSimple)
4479                     {
4480                         new_ci_entry(nbl, cell0_i+ci, shift, flags_i[ci]);
4481                     }
4482                     else
4483                     {
4484                         new_sci_entry(nbl, cell0_i+ci, shift);
4485                     }
4486
4487 #ifndef NBNXN_SHIFT_BACKWARD
4488                     if (cxf < ci_x)
4489 #else
4490                     if (shift == CENTRAL && gridi == gridj &&
4491                         cxf < ci_x)
4492 #endif
4493                     {
4494                         /* Leave the pairs with i > j.
4495                          * x is the major index, so skip half of it.
4496                          */
4497                         cxf = ci_x;
4498                     }
4499
4500                     if (nbl->bSimple)
4501                     {
4502                         set_icell_bb_simple(bb_i, ci, shx, shy, shz,
4503                                             nbl->work->bb_ci);
4504                     }
4505                     else
4506                     {
4507 #ifdef NBNXN_BBXXXX
4508                         set_icell_bbxxxx_supersub(pbb_i, ci, shx, shy, shz,
4509                                                   nbl->work->pbb_ci);
4510 #else
4511                         set_icell_bb_supersub(bb_i, ci, shx, shy, shz,
4512                                               nbl->work->bb_ci);
4513 #endif
4514                     }
4515
4516                     nbs->icell_set_x(cell0_i+ci, shx, shy, shz,
4517                                      gridi->na_c, nbat->xstride, nbat->x,
4518                                      nbl->work);
4519
4520                     for (cx = cxf; cx <= cxl; cx++)
4521                     {
4522                         d2zx = d2z;
4523                         if (gridj->c0[XX] + cx*gridj->sx > bx1)
4524                         {
4525                             d2zx += sqr(gridj->c0[XX] + cx*gridj->sx - bx1);
4526                         }
4527                         else if (gridj->c0[XX] + (cx+1)*gridj->sx < bx0)
4528                         {
4529                             d2zx += sqr(gridj->c0[XX] + (cx+1)*gridj->sx - bx0);
4530                         }
4531
4532 #ifndef NBNXN_SHIFT_BACKWARD
4533                         if (gridi == gridj &&
4534                             cx == 0 && cyf < ci_y)
4535 #else
4536                         if (gridi == gridj &&
4537                             cx == 0 && shift == CENTRAL && cyf < ci_y)
4538 #endif
4539                         {
4540                             /* Leave the pairs with i > j.
4541                              * Skip half of y when i and j have the same x.
4542                              */
4543                             cyf_x = ci_y;
4544                         }
4545                         else
4546                         {
4547                             cyf_x = cyf;
4548                         }
4549
4550                         for (cy = cyf_x; cy <= cyl; cy++)
4551                         {
4552                             c0 = gridj->cxy_ind[cx*gridj->ncy+cy];
4553                             c1 = gridj->cxy_ind[cx*gridj->ncy+cy+1];
4554 #ifdef NBNXN_SHIFT_BACKWARD
4555                             if (gridi == gridj &&
4556                                 shift == CENTRAL && c0 < ci)
4557                             {
4558                                 c0 = ci;
4559                             }
4560 #endif
4561
4562                             d2zxy = d2zx;
4563                             if (gridj->c0[YY] + cy*gridj->sy > by1)
4564                             {
4565                                 d2zxy += sqr(gridj->c0[YY] + cy*gridj->sy - by1);
4566                             }
4567                             else if (gridj->c0[YY] + (cy+1)*gridj->sy < by0)
4568                             {
4569                                 d2zxy += sqr(gridj->c0[YY] + (cy+1)*gridj->sy - by0);
4570                             }
4571                             if (c1 > c0 && d2zxy < rl2)
4572                             {
4573                                 cs = c0 + (int)(bz1_frac*(c1 - c0));
4574                                 if (cs >= c1)
4575                                 {
4576                                     cs = c1 - 1;
4577                                 }
4578
4579                                 d2xy = d2zxy - d2z;
4580
4581                                 /* Find the lowest cell that can possibly
4582                                  * be within range.
4583                                  */
4584                                 cf = cs;
4585                                 while (cf > c0 &&
4586                                        (bbcz_j[cf*NNBSBB_D+1] >= bz0 ||
4587                                         d2xy + sqr(bbcz_j[cf*NNBSBB_D+1] - bz0) < rl2))
4588                                 {
4589                                     cf--;
4590                                 }
4591
4592                                 /* Find the highest cell that can possibly
4593                                  * be within range.
4594                                  */
4595                                 cl = cs;
4596                                 while (cl < c1-1 &&
4597                                        (bbcz_j[cl*NNBSBB_D] <= bz1 ||
4598                                         d2xy + sqr(bbcz_j[cl*NNBSBB_D] - bz1) < rl2))
4599                                 {
4600                                     cl++;
4601                                 }
4602
4603 #ifdef NBNXN_REFCODE
4604                                 {
4605                                     /* Simple reference code, for debugging,
4606                                      * overrides the more complex code above.
4607                                      */
4608                                     int k;
4609                                     cf = c1;
4610                                     cl = -1;
4611                                     for (k = c0; k < c1; k++)
4612                                     {
4613                                         if (box_dist2(bx0, bx1, by0, by1, bz0, bz1, bb+k) < rl2 &&
4614                                             k < cf)
4615                                         {
4616                                             cf = k;
4617                                         }
4618                                         if (box_dist2(bx0, bx1, by0, by1, bz0, bz1, bb+k) < rl2 &&
4619                                             k > cl)
4620                                         {
4621                                             cl = k;
4622                                         }
4623                                     }
4624                                 }
4625 #endif
4626
4627                                 if (gridi == gridj)
4628                                 {
4629                                     /* We want each atom/cell pair only once,
4630                                      * only use cj >= ci.
4631                                      */
4632 #ifndef NBNXN_SHIFT_BACKWARD
4633                                     cf = max(cf, ci);
4634 #else
4635                                     if (shift == CENTRAL)
4636                                     {
4637                                         cf = max(cf, ci);
4638                                     }
4639 #endif
4640                                 }
4641
4642                                 if (cf <= cl)
4643                                 {
4644                                     /* For f buffer flags with simple lists */
4645                                     ncj_old_j = nbl->ncj;
4646
4647                                     switch (nb_kernel_type)
4648                                     {
4649                                         case nbnxnk4x4_PlainC:
4650                                             check_subcell_list_space_simple(nbl, cl-cf+1);
4651
4652                                             make_cluster_list_simple(gridj,
4653                                                                      nbl, ci, cf, cl,
4654                                                                      (gridi == gridj && shift == CENTRAL),
4655                                                                      nbat->x,
4656                                                                      rl2, rbb2,
4657                                                                      &ndistc);
4658                                             break;
4659 #ifdef GMX_NBNXN_SIMD_4XN
4660                                         case nbnxnk4xN_SIMD_4xN:
4661                                             check_subcell_list_space_simple(nbl, ci_to_cj(na_cj_2log, cl-cf)+2);
4662                                             make_cluster_list_simd_4xn(gridj,
4663                                                                        nbl, ci, cf, cl,
4664                                                                        (gridi == gridj && shift == CENTRAL),
4665                                                                        nbat->x,
4666                                                                        rl2, rbb2,
4667                                                                        &ndistc);
4668                                             break;
4669 #endif
4670 #ifdef GMX_NBNXN_SIMD_2XNN
4671                                         case nbnxnk4xN_SIMD_2xNN:
4672                                             check_subcell_list_space_simple(nbl, ci_to_cj(na_cj_2log, cl-cf)+2);
4673                                             make_cluster_list_simd_2xnn(gridj,
4674                                                                         nbl, ci, cf, cl,
4675                                                                         (gridi == gridj && shift == CENTRAL),
4676                                                                         nbat->x,
4677                                                                         rl2, rbb2,
4678                                                                         &ndistc);
4679                                             break;
4680 #endif
4681                                         case nbnxnk8x8x8_PlainC:
4682                                         case nbnxnk8x8x8_CUDA:
4683                                             check_subcell_list_space_supersub(nbl, cl-cf+1);
4684                                             for (cj = cf; cj <= cl; cj++)
4685                                             {
4686                                                 make_cluster_list_supersub(gridi, gridj,
4687                                                                            nbl, ci, cj,
4688                                                                            (gridi == gridj && shift == CENTRAL && ci == cj),
4689                                                                            nbat->xstride, nbat->x,
4690                                                                            rl2, rbb2,
4691                                                                            &ndistc);
4692                                             }
4693                                             break;
4694                                     }
4695                                     ncpcheck += cl - cf + 1;
4696
4697                                     if (bFBufferFlag && nbl->ncj > ncj_old_j)
4698                                     {
4699                                         int cbf, cbl, cb;
4700
4701                                         cbf = nbl->cj[ncj_old_j].cj >> gridj_flag_shift;
4702                                         cbl = nbl->cj[nbl->ncj-1].cj >> gridj_flag_shift;
4703                                         for (cb = cbf; cb <= cbl; cb++)
4704                                         {
4705                                             gridj_flag[cb] = 1U<<th;
4706                                         }
4707                                     }
4708                                 }
4709                             }
4710                         }
4711                     }
4712
4713                     /* Set the exclusions for this ci list */
4714                     if (nbl->bSimple)
4715                     {
4716                         set_ci_top_excls(nbs,
4717                                          nbl,
4718                                          shift == CENTRAL && gridi == gridj,
4719                                          gridj->na_c_2log,
4720                                          na_cj_2log,
4721                                          &(nbl->ci[nbl->nci]),
4722                                          excl);
4723                     }
4724                     else
4725                     {
4726                         set_sci_top_excls(nbs,
4727                                           nbl,
4728                                           shift == CENTRAL && gridi == gridj,
4729                                           gridj->na_c_2log,
4730                                           &(nbl->sci[nbl->nsci]),
4731                                           excl);
4732                     }
4733
4734                     /* Close this ci list */
4735                     if (nbl->bSimple)
4736                     {
4737                         close_ci_entry_simple(nbl);
4738                     }
4739                     else
4740                     {
4741                         close_ci_entry_supersub(nbl,
4742                                                 nsubpair_max,
4743                                                 progBal, min_ci_balanced,
4744                                                 th, nth);
4745                     }
4746                 }
4747             }
4748         }
4749
4750         if (bFBufferFlag && nbl->ncj > ncj_old_i)
4751         {
4752             work->buffer_flags.flag[(gridi->cell0+ci)>>gridi_flag_shift] = 1U<<th;
4753         }
4754     }
4755
4756     work->ndistc = ndistc;
4757
4758     nbs_cycle_stop(&work->cc[enbsCCsearch]);
4759
4760     if (debug)
4761     {
4762         fprintf(debug, "number of distance checks %d\n", ndistc);
4763         fprintf(debug, "ncpcheck %s %d\n", gridi == gridj ? "local" : "non-local",
4764                 ncpcheck);
4765
4766         if (nbl->bSimple)
4767         {
4768             print_nblist_statistics_simple(debug, nbl, nbs, rlist);
4769         }
4770         else
4771         {
4772             print_nblist_statistics_supersub(debug, nbl, nbs, rlist);
4773         }
4774
4775     }
4776 }
4777
4778 static void reduce_buffer_flags(const nbnxn_search_t        nbs,
4779                                 int                         nsrc,
4780                                 const nbnxn_buffer_flags_t *dest)
4781 {
4782     int s, b;
4783     const unsigned *flag;
4784
4785     for (s = 0; s < nsrc; s++)
4786     {
4787         flag = nbs->work[s].buffer_flags.flag;
4788
4789         for (b = 0; b < dest->nflag; b++)
4790         {
4791             dest->flag[b] |= flag[b];
4792         }
4793     }
4794 }
4795
4796 static void print_reduction_cost(const nbnxn_buffer_flags_t *flags, int nout)
4797 {
4798     int nelem, nkeep, ncopy, nred, b, c, out;
4799
4800     nelem = 0;
4801     nkeep = 0;
4802     ncopy = 0;
4803     nred  = 0;
4804     for (b = 0; b < flags->nflag; b++)
4805     {
4806         if (flags->flag[b] == 1)
4807         {
4808             /* Only flag 0 is set, no copy of reduction required */
4809             nelem++;
4810             nkeep++;
4811         }
4812         else if (flags->flag[b] > 0)
4813         {
4814             c = 0;
4815             for (out = 0; out < nout; out++)
4816             {
4817                 if (flags->flag[b] & (1U<<out))
4818                 {
4819                     c++;
4820                 }
4821             }
4822             nelem += c;
4823             if (c == 1)
4824             {
4825                 ncopy++;
4826             }
4827             else
4828             {
4829                 nred += c;
4830             }
4831         }
4832     }
4833
4834     fprintf(debug, "nbnxn reduction: #flag %d #list %d elem %4.2f, keep %4.2f copy %4.2f red %4.2f\n",
4835             flags->nflag, nout,
4836             nelem/(double)(flags->nflag),
4837             nkeep/(double)(flags->nflag),
4838             ncopy/(double)(flags->nflag),
4839             nred/(double)(flags->nflag));
4840 }
4841
4842 /* Perform a count (linear) sort to sort the smaller lists to the end.
4843  * This avoids load imbalance on the GPU, as large lists will be
4844  * scheduled and executed first and the smaller lists later.
4845  * Load balancing between multi-processors only happens at the end
4846  * and there smaller lists lead to more effective load balancing.
4847  * The sorting is done on the cj4 count, not on the actual pair counts.
4848  * Not only does this make the sort faster, but it also results in
4849  * better load balancing than using a list sorted on exact load.
4850  * This function swaps the pointer in the pair list to avoid a copy operation.
4851  */
4852 static void sort_sci(nbnxn_pairlist_t *nbl)
4853 {
4854     nbnxn_list_work_t *work;
4855     int                m, i, s, s0, s1;
4856     nbnxn_sci_t       *sci_sort;
4857
4858     if (nbl->ncj4 <= nbl->nsci)
4859     {
4860         /* nsci = 0 or all sci have size 1, sorting won't change the order */
4861         return;
4862     }
4863
4864     work = nbl->work;
4865
4866     /* We will distinguish differences up to double the average */
4867     m = (2*nbl->ncj4)/nbl->nsci;
4868
4869     if (m + 1 > work->sort_nalloc)
4870     {
4871         work->sort_nalloc = over_alloc_large(m + 1);
4872         srenew(work->sort, work->sort_nalloc);
4873     }
4874
4875     if (work->sci_sort_nalloc != nbl->sci_nalloc)
4876     {
4877         work->sci_sort_nalloc = nbl->sci_nalloc;
4878         nbnxn_realloc_void((void **)&work->sci_sort,
4879                            0,
4880                            work->sci_sort_nalloc*sizeof(*work->sci_sort),
4881                            nbl->alloc, nbl->free);
4882     }
4883
4884     /* Count the entries of each size */
4885     for (i = 0; i <= m; i++)
4886     {
4887         work->sort[i] = 0;
4888     }
4889     for (s = 0; s < nbl->nsci; s++)
4890     {
4891         i = min(m, nbl->sci[s].cj4_ind_end - nbl->sci[s].cj4_ind_start);
4892         work->sort[i]++;
4893     }
4894     /* Calculate the offset for each count */
4895     s0            = work->sort[m];
4896     work->sort[m] = 0;
4897     for (i = m - 1; i >= 0; i--)
4898     {
4899         s1            = work->sort[i];
4900         work->sort[i] = work->sort[i + 1] + s0;
4901         s0            = s1;
4902     }
4903
4904     /* Sort entries directly into place */
4905     sci_sort = work->sci_sort;
4906     for (s = 0; s < nbl->nsci; s++)
4907     {
4908         i = min(m, nbl->sci[s].cj4_ind_end - nbl->sci[s].cj4_ind_start);
4909         sci_sort[work->sort[i]++] = nbl->sci[s];
4910     }
4911
4912     /* Swap the sci pointers so we use the new, sorted list */
4913     work->sci_sort = nbl->sci;
4914     nbl->sci       = sci_sort;
4915 }
4916
4917 /* Make a local or non-local pair-list, depending on iloc */
4918 void nbnxn_make_pairlist(const nbnxn_search_t  nbs,
4919                          nbnxn_atomdata_t     *nbat,
4920                          const t_blocka       *excl,
4921                          real                  rlist,
4922                          int                   min_ci_balanced,
4923                          nbnxn_pairlist_set_t *nbl_list,
4924                          int                   iloc,
4925                          int                   nb_kernel_type,
4926                          t_nrnb               *nrnb)
4927 {
4928     nbnxn_grid_t *gridi, *gridj;
4929     gmx_bool bGPUCPU;
4930     int nzi, zi, zj0, zj1, zj;
4931     int nsubpair_max;
4932     int th;
4933     int nnbl;
4934     nbnxn_pairlist_t **nbl;
4935     int ci_block;
4936     gmx_bool CombineNBLists;
4937     gmx_bool progBal;
4938     int np_tot, np_noq, np_hlj, nap;
4939
4940     /* Check if we are running hybrid GPU + CPU nbnxn mode */
4941     bGPUCPU = (!nbs->grid[0].bSimple && nbl_list->bSimple);
4942
4943     nnbl            = nbl_list->nnbl;
4944     nbl             = nbl_list->nbl;
4945     CombineNBLists  = nbl_list->bCombined;
4946
4947     if (debug)
4948     {
4949         fprintf(debug, "ns making %d nblists\n", nnbl);
4950     }
4951
4952     nbat->bUseBufferFlags = (nbat->nout > 1);
4953     /* We should re-init the flags before making the first list */
4954     if (nbat->bUseBufferFlags && (LOCAL_I(iloc) || bGPUCPU))
4955     {
4956         init_buffer_flags(&nbat->buffer_flags, nbat->natoms);
4957     }
4958
4959     if (nbl_list->bSimple)
4960     {
4961         switch (nb_kernel_type)
4962         {
4963 #ifdef GMX_NBNXN_SIMD_4XN
4964             case nbnxnk4xN_SIMD_4xN:
4965                 nbs->icell_set_x = icell_set_x_simd_4xn;
4966                 break;
4967 #endif
4968 #ifdef GMX_NBNXN_SIMD_2XNN
4969             case nbnxnk4xN_SIMD_2xNN:
4970                 nbs->icell_set_x = icell_set_x_simd_2xnn;
4971                 break;
4972 #endif
4973             default:
4974                 nbs->icell_set_x = icell_set_x_simple;
4975                 break;
4976         }
4977     }
4978     else
4979     {
4980 #ifdef NBNXN_SEARCH_BB_SSE
4981         nbs->icell_set_x = icell_set_x_supersub_sse8;
4982 #else
4983         nbs->icell_set_x = icell_set_x_supersub;
4984 #endif
4985     }
4986
4987     if (LOCAL_I(iloc))
4988     {
4989         /* Only zone (grid) 0 vs 0 */
4990         nzi = 1;
4991         zj0 = 0;
4992         zj1 = 1;
4993     }
4994     else
4995     {
4996         nzi = nbs->zones->nizone;
4997     }
4998
4999     if (!nbl_list->bSimple && min_ci_balanced > 0)
5000     {
5001         nsubpair_max = get_nsubpair_max(nbs, iloc, rlist, min_ci_balanced);
5002     }
5003     else
5004     {
5005         nsubpair_max = 0;
5006     }
5007
5008     /* Clear all pair-lists */
5009     for (th = 0; th < nnbl; th++)
5010     {
5011         clear_pairlist(nbl[th]);
5012     }
5013
5014     for (zi = 0; zi < nzi; zi++)
5015     {
5016         gridi = &nbs->grid[zi];
5017
5018         if (NONLOCAL_I(iloc))
5019         {
5020             zj0 = nbs->zones->izone[zi].j0;
5021             zj1 = nbs->zones->izone[zi].j1;
5022             if (zi == 0)
5023             {
5024                 zj0++;
5025             }
5026         }
5027         for (zj = zj0; zj < zj1; zj++)
5028         {
5029             gridj = &nbs->grid[zj];
5030
5031             if (debug)
5032             {
5033                 fprintf(debug, "ns search grid %d vs %d\n", zi, zj);
5034             }
5035
5036             nbs_cycle_start(&nbs->cc[enbsCCsearch]);
5037
5038             if (nbl[0]->bSimple && !gridi->bSimple)
5039             {
5040                 /* Hybrid list, determine blocking later */
5041                 ci_block = 0;
5042             }
5043             else
5044             {
5045                 ci_block = get_ci_block_size(gridi, nbs->DomDec, nnbl);
5046             }
5047
5048 #pragma omp parallel for num_threads(nnbl) schedule(static)
5049             for (th = 0; th < nnbl; th++)
5050             {
5051                 /* Re-init the thread-local work flag data before making
5052                  * the first list (not an elegant conditional).
5053                  */
5054                 if (nbat->bUseBufferFlags && ((zi == 0 && zj == 0) ||
5055                                               (bGPUCPU && zi == 0 && zj == 1)))
5056                 {
5057                     init_buffer_flags(&nbs->work[th].buffer_flags, nbat->natoms);
5058                 }
5059
5060                 if (CombineNBLists && th > 0)
5061                 {
5062                     clear_pairlist(nbl[th]);
5063                 }
5064
5065                 /* With GPU: generate progressively smaller lists for
5066                  * load balancing for local only or non-local with 2 zones.
5067                  */
5068                 progBal = (LOCAL_I(iloc) || nbs->zones->n <= 2);
5069
5070                 /* Divide the i super cell equally over the nblists */
5071                 nbnxn_make_pairlist_part(nbs, gridi, gridj,
5072                                          &nbs->work[th], nbat, excl,
5073                                          rlist,
5074                                          nb_kernel_type,
5075                                          ci_block,
5076                                          nbat->bUseBufferFlags,
5077                                          nsubpair_max,
5078                                          progBal, min_ci_balanced,
5079                                          th, nnbl,
5080                                          nbl[th]);
5081             }
5082             nbs_cycle_stop(&nbs->cc[enbsCCsearch]);
5083
5084             np_tot = 0;
5085             np_noq = 0;
5086             np_hlj = 0;
5087             for (th = 0; th < nnbl; th++)
5088             {
5089                 inc_nrnb(nrnb, eNR_NBNXN_DIST2, nbs->work[th].ndistc);
5090
5091                 if (nbl_list->bSimple)
5092                 {
5093                     np_tot += nbl[th]->ncj;
5094                     np_noq += nbl[th]->work->ncj_noq;
5095                     np_hlj += nbl[th]->work->ncj_hlj;
5096                 }
5097                 else
5098                 {
5099                     /* This count ignores potential subsequent pair pruning */
5100                     np_tot += nbl[th]->nci_tot;
5101                 }
5102             }
5103             nap                   = nbl[0]->na_ci*nbl[0]->na_cj;
5104             nbl_list->natpair_ljq = (np_tot - np_noq)*nap - np_hlj*nap/2;
5105             nbl_list->natpair_lj  = np_noq*nap;
5106             nbl_list->natpair_q   = np_hlj*nap/2;
5107
5108             if (CombineNBLists && nnbl > 1)
5109             {
5110                 nbs_cycle_start(&nbs->cc[enbsCCcombine]);
5111
5112                 combine_nblists(nnbl-1, nbl+1, nbl[0]);
5113
5114                 nbs_cycle_stop(&nbs->cc[enbsCCcombine]);
5115             }
5116         }
5117     }
5118
5119     if (!nbl_list->bSimple)
5120     {
5121         /* Sort the entries on size, large ones first */
5122         if (CombineNBLists || nnbl == 1)
5123         {
5124             sort_sci(nbl[0]);
5125         }
5126         else
5127         {
5128 #pragma omp parallel for num_threads(nnbl) schedule(static)
5129             for (th = 0; th < nnbl; th++)
5130             {
5131                 sort_sci(nbl[th]);
5132             }
5133         }
5134     }
5135
5136     if (nbat->bUseBufferFlags)
5137     {
5138         reduce_buffer_flags(nbs, nnbl, &nbat->buffer_flags);
5139     }
5140
5141     /* Special performance logging stuff (env.var. GMX_NBNXN_CYCLE) */
5142     if (LOCAL_I(iloc))
5143     {
5144         nbs->search_count++;
5145     }
5146     if (nbs->print_cycles &&
5147         (!nbs->DomDec || (nbs->DomDec && !LOCAL_I(iloc))) &&
5148         nbs->search_count % 100 == 0)
5149     {
5150         nbs_cycle_print(stderr, nbs);
5151     }
5152
5153     if (debug && (CombineNBLists && nnbl > 1))
5154     {
5155         if (nbl[0]->bSimple)
5156         {
5157             print_nblist_statistics_simple(debug, nbl[0], nbs, rlist);
5158         }
5159         else
5160         {
5161             print_nblist_statistics_supersub(debug, nbl[0], nbs, rlist);
5162         }
5163     }
5164
5165     if (debug)
5166     {
5167         if (gmx_debug_at)
5168         {
5169             if (nbl[0]->bSimple)
5170             {
5171                 print_nblist_ci_cj(debug, nbl[0]);
5172             }
5173             else
5174             {
5175                 print_nblist_sci_cj(debug, nbl[0]);
5176             }
5177         }
5178
5179         if (nbat->bUseBufferFlags)
5180         {
5181             print_reduction_cost(&nbat->buffer_flags, nnbl);
5182         }
5183     }
5184 }