Merge branch release-5-1
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_search.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
36 #include "gmxpre.h"
37
38 #include "nbnxn_search.h"
39
40 #include "config.h"
41
42 #include <assert.h>
43 #include <string.h>
44
45 #include <cmath>
46
47 #include <algorithm>
48
49 #include "gromacs/domdec/domdec_struct.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/utilities.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/mdlib/gmx_omp_nthreads.h"
55 #include "gromacs/mdlib/nb_verlet.h"
56 #include "gromacs/mdlib/nbnxn_atomdata.h"
57 #include "gromacs/mdlib/nbnxn_consts.h"
58 #include "gromacs/mdlib/nbnxn_grid.h"
59 #include "gromacs/mdlib/nbnxn_internal.h"
60 #include "gromacs/mdlib/nbnxn_simd.h"
61 #include "gromacs/mdlib/nbnxn_util.h"
62 #include "gromacs/mdlib/ns.h"
63 #include "gromacs/mdtypes/group.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/pbcutil/ishift.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/simd/simd.h"
68 #include "gromacs/simd/vector_operations.h"
69 #include "gromacs/utility/exceptions.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/smalloc.h"
72
73 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
74
75 #if GMX_SIMD
76
77 /* The functions below are macros as they are performance sensitive */
78
79 /* 4x4 list, pack=4: no complex conversion required */
80 /* i-cluster to j-cluster conversion */
81 #define CI_TO_CJ_J4(ci)   (ci)
82 /* cluster index to coordinate array index conversion */
83 #define X_IND_CI_J4(ci)  ((ci)*STRIDE_P4)
84 #define X_IND_CJ_J4(cj)  ((cj)*STRIDE_P4)
85
86 /* 4x2 list, pack=4: j-cluster size is half the packing width */
87 /* i-cluster to j-cluster conversion */
88 #define CI_TO_CJ_J2(ci)  ((ci)<<1)
89 /* cluster index to coordinate array index conversion */
90 #define X_IND_CI_J2(ci)  ((ci)*STRIDE_P4)
91 #define X_IND_CJ_J2(cj)  (((cj)>>1)*STRIDE_P4 + ((cj) & 1)*(PACK_X4>>1))
92
93 /* 4x8 list, pack=8: i-cluster size is half the packing width */
94 /* i-cluster to j-cluster conversion */
95 #define CI_TO_CJ_J8(ci)  ((ci)>>1)
96 /* cluster index to coordinate array index conversion */
97 #define X_IND_CI_J8(ci)  (((ci)>>1)*STRIDE_P8 + ((ci) & 1)*(PACK_X8>>1))
98 #define X_IND_CJ_J8(cj)  ((cj)*STRIDE_P8)
99
100 /* The j-cluster size is matched to the SIMD width */
101 #if GMX_SIMD_REAL_WIDTH == 2
102 #define CI_TO_CJ_SIMD_4XN(ci)  CI_TO_CJ_J2(ci)
103 #define X_IND_CI_SIMD_4XN(ci)  X_IND_CI_J2(ci)
104 #define X_IND_CJ_SIMD_4XN(cj)  X_IND_CJ_J2(cj)
105 #else
106 #if GMX_SIMD_REAL_WIDTH == 4
107 #define CI_TO_CJ_SIMD_4XN(ci)  CI_TO_CJ_J4(ci)
108 #define X_IND_CI_SIMD_4XN(ci)  X_IND_CI_J4(ci)
109 #define X_IND_CJ_SIMD_4XN(cj)  X_IND_CJ_J4(cj)
110 #else
111 #if GMX_SIMD_REAL_WIDTH == 8
112 #define CI_TO_CJ_SIMD_4XN(ci)  CI_TO_CJ_J8(ci)
113 #define X_IND_CI_SIMD_4XN(ci)  X_IND_CI_J8(ci)
114 #define X_IND_CJ_SIMD_4XN(cj)  X_IND_CJ_J8(cj)
115 /* Half SIMD with j-cluster size */
116 #define CI_TO_CJ_SIMD_2XNN(ci) CI_TO_CJ_J4(ci)
117 #define X_IND_CI_SIMD_2XNN(ci) X_IND_CI_J4(ci)
118 #define X_IND_CJ_SIMD_2XNN(cj) X_IND_CJ_J4(cj)
119 #else
120 #if GMX_SIMD_REAL_WIDTH == 16
121 #define CI_TO_CJ_SIMD_2XNN(ci) CI_TO_CJ_J8(ci)
122 #define X_IND_CI_SIMD_2XNN(ci) X_IND_CI_J8(ci)
123 #define X_IND_CJ_SIMD_2XNN(cj) X_IND_CJ_J8(cj)
124 #else
125 #error "unsupported GMX_SIMD_REAL_WIDTH"
126 #endif
127 #endif
128 #endif
129 #endif
130
131 #endif // GMX_SIMD
132
133
134 /* We shift the i-particles backward for PBC.
135  * This leads to more conditionals than shifting forward.
136  * We do this to get more balanced pair lists.
137  */
138 #define NBNXN_SHIFT_BACKWARD
139
140
141 static void nbs_cycle_clear(nbnxn_cycle_t *cc)
142 {
143     for (int i = 0; i < enbsCCnr; i++)
144     {
145         cc[i].count = 0;
146         cc[i].c     = 0;
147     }
148 }
149
150 static double Mcyc_av(const nbnxn_cycle_t *cc)
151 {
152     return (double)cc->c*1e-6/cc->count;
153 }
154
155 static void nbs_cycle_print(FILE *fp, const nbnxn_search_t nbs)
156 {
157     fprintf(fp, "\n");
158     fprintf(fp, "ns %4d grid %4.1f search %4.1f red.f %5.3f",
159             nbs->cc[enbsCCgrid].count,
160             Mcyc_av(&nbs->cc[enbsCCgrid]),
161             Mcyc_av(&nbs->cc[enbsCCsearch]),
162             Mcyc_av(&nbs->cc[enbsCCreducef]));
163
164     if (nbs->nthread_max > 1)
165     {
166         if (nbs->cc[enbsCCcombine].count > 0)
167         {
168             fprintf(fp, " comb %5.2f",
169                     Mcyc_av(&nbs->cc[enbsCCcombine]));
170         }
171         fprintf(fp, " s. th");
172         for (int t = 0; t < nbs->nthread_max; t++)
173         {
174             fprintf(fp, " %4.1f",
175                     Mcyc_av(&nbs->work[t].cc[enbsCCsearch]));
176         }
177     }
178     fprintf(fp, "\n");
179 }
180
181 static gmx_inline int ci_to_cj(int na_cj_2log, int ci)
182 {
183     switch (na_cj_2log)
184     {
185         case 2: return ci;     break;
186         case 1: return (ci<<1); break;
187         case 3: return (ci>>1); break;
188     }
189
190     return 0;
191 }
192
193 gmx_bool nbnxn_kernel_pairlist_simple(int nb_kernel_type)
194 {
195     if (nb_kernel_type == nbnxnkNotSet)
196     {
197         gmx_fatal(FARGS, "Non-bonded kernel type not set for Verlet-style pair-list.");
198     }
199
200     switch (nb_kernel_type)
201     {
202         case nbnxnk8x8x8_GPU:
203         case nbnxnk8x8x8_PlainC:
204             return FALSE;
205
206         case nbnxnk4x4_PlainC:
207         case nbnxnk4xN_SIMD_4xN:
208         case nbnxnk4xN_SIMD_2xNN:
209             return TRUE;
210
211         default:
212             gmx_incons("Invalid nonbonded kernel type passed!");
213             return FALSE;
214     }
215 }
216
217 /* Initializes a single nbnxn_pairlist_t data structure */
218 static void nbnxn_init_pairlist_fep(t_nblist *nl)
219 {
220     nl->type        = GMX_NBLIST_INTERACTION_FREE_ENERGY;
221     nl->igeometry   = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
222     /* The interaction functions are set in the free energy kernel fuction */
223     nl->ivdw        = -1;
224     nl->ivdwmod     = -1;
225     nl->ielec       = -1;
226     nl->ielecmod    = -1;
227
228     nl->maxnri      = 0;
229     nl->maxnrj      = 0;
230     nl->nri         = 0;
231     nl->nrj         = 0;
232     nl->iinr        = NULL;
233     nl->gid         = NULL;
234     nl->shift       = NULL;
235     nl->jindex      = NULL;
236     nl->jjnr        = NULL;
237     nl->excl_fep    = NULL;
238
239 }
240
241 void nbnxn_init_search(nbnxn_search_t           * nbs_ptr,
242                        ivec                      *n_dd_cells,
243                        struct gmx_domdec_zones_t *zones,
244                        gmx_bool                   bFEP,
245                        int                        nthread_max)
246 {
247     nbnxn_search_t nbs;
248     int            ngrid;
249
250     snew(nbs, 1);
251     *nbs_ptr = nbs;
252
253     nbs->bFEP   = bFEP;
254
255     nbs->DomDec = (n_dd_cells != NULL);
256
257     clear_ivec(nbs->dd_dim);
258     ngrid = 1;
259     if (nbs->DomDec)
260     {
261         nbs->zones = zones;
262
263         for (int d = 0; d < DIM; d++)
264         {
265             if ((*n_dd_cells)[d] > 1)
266             {
267                 nbs->dd_dim[d] = 1;
268                 /* Each grid matches a DD zone */
269                 ngrid *= 2;
270             }
271         }
272     }
273
274     nbnxn_grids_init(nbs, ngrid);
275
276     nbs->cell        = NULL;
277     nbs->cell_nalloc = 0;
278     nbs->a           = NULL;
279     nbs->a_nalloc    = 0;
280
281     nbs->nthread_max = nthread_max;
282
283     /* Initialize the work data structures for each thread */
284     snew(nbs->work, nbs->nthread_max);
285     for (int t = 0; t < nbs->nthread_max; t++)
286     {
287         nbs->work[t].cxy_na           = NULL;
288         nbs->work[t].cxy_na_nalloc    = 0;
289         nbs->work[t].sort_work        = NULL;
290         nbs->work[t].sort_work_nalloc = 0;
291
292         snew(nbs->work[t].nbl_fep, 1);
293         nbnxn_init_pairlist_fep(nbs->work[t].nbl_fep);
294     }
295
296     /* Initialize detailed nbsearch cycle counting */
297     nbs->print_cycles = (getenv("GMX_NBNXN_CYCLE") != 0);
298     nbs->search_count = 0;
299     nbs_cycle_clear(nbs->cc);
300     for (int t = 0; t < nbs->nthread_max; t++)
301     {
302         nbs_cycle_clear(nbs->work[t].cc);
303     }
304 }
305
306 static void init_buffer_flags(nbnxn_buffer_flags_t *flags,
307                               int                   natoms)
308 {
309     flags->nflag = (natoms + NBNXN_BUFFERFLAG_SIZE - 1)/NBNXN_BUFFERFLAG_SIZE;
310     if (flags->nflag > flags->flag_nalloc)
311     {
312         flags->flag_nalloc = over_alloc_large(flags->nflag);
313         srenew(flags->flag, flags->flag_nalloc);
314     }
315     for (int b = 0; b < flags->nflag; b++)
316     {
317         bitmask_clear(&(flags->flag[b]));
318     }
319 }
320
321 /* Determines the cell range along one dimension that
322  * the bounding box b0 - b1 sees.
323  */
324 static void get_cell_range(real b0, real b1,
325                            int nc, real c0, real s, real invs,
326                            real d2, real r2, int *cf, int *cl)
327 {
328     *cf = std::max(static_cast<int>((b0 - c0)*invs), 0);
329
330     while (*cf > 0 && d2 + gmx::square((b0 - c0) - (*cf-1+1)*s) < r2)
331     {
332         (*cf)--;
333     }
334
335     *cl = std::min(static_cast<int>((b1 - c0)*invs), nc-1);
336     while (*cl < nc-1 && d2 + gmx::square((*cl+1)*s - (b1 - c0)) < r2)
337     {
338         (*cl)++;
339     }
340 }
341
342 /* Reference code calculating the distance^2 between two bounding boxes */
343 static float box_dist2(float bx0, float bx1, float by0,
344                        float by1, float bz0, float bz1,
345                        const nbnxn_bb_t *bb)
346 {
347     float d2;
348     float dl, dh, dm, dm0;
349
350     d2 = 0;
351
352     dl  = bx0 - bb->upper[BB_X];
353     dh  = bb->lower[BB_X] - bx1;
354     dm  = std::max(dl, dh);
355     dm0 = std::max(dm, 0.0f);
356     d2 += dm0*dm0;
357
358     dl  = by0 - bb->upper[BB_Y];
359     dh  = bb->lower[BB_Y] - by1;
360     dm  = std::max(dl, dh);
361     dm0 = std::max(dm, 0.0f);
362     d2 += dm0*dm0;
363
364     dl  = bz0 - bb->upper[BB_Z];
365     dh  = bb->lower[BB_Z] - bz1;
366     dm  = std::max(dl, dh);
367     dm0 = std::max(dm, 0.0f);
368     d2 += dm0*dm0;
369
370     return d2;
371 }
372
373 /* Plain C code calculating the distance^2 between two bounding boxes */
374 static float subc_bb_dist2(int si, const nbnxn_bb_t *bb_i_ci,
375                            int csj, const nbnxn_bb_t *bb_j_all)
376 {
377     const nbnxn_bb_t *bb_i, *bb_j;
378     float             d2;
379     float             dl, dh, dm, dm0;
380
381     bb_i = bb_i_ci  +  si;
382     bb_j = bb_j_all + csj;
383
384     d2 = 0;
385
386     dl  = bb_i->lower[BB_X] - bb_j->upper[BB_X];
387     dh  = bb_j->lower[BB_X] - bb_i->upper[BB_X];
388     dm  = std::max(dl, dh);
389     dm0 = std::max(dm, 0.0f);
390     d2 += dm0*dm0;
391
392     dl  = bb_i->lower[BB_Y] - bb_j->upper[BB_Y];
393     dh  = bb_j->lower[BB_Y] - bb_i->upper[BB_Y];
394     dm  = std::max(dl, dh);
395     dm0 = std::max(dm, 0.0f);
396     d2 += dm0*dm0;
397
398     dl  = bb_i->lower[BB_Z] - bb_j->upper[BB_Z];
399     dh  = bb_j->lower[BB_Z] - bb_i->upper[BB_Z];
400     dm  = std::max(dl, dh);
401     dm0 = std::max(dm, 0.0f);
402     d2 += dm0*dm0;
403
404     return d2;
405 }
406
407 #ifdef NBNXN_SEARCH_BB_SIMD4
408
409 /* 4-wide SIMD code for bb distance for bb format xyz0 */
410 static float subc_bb_dist2_simd4(int si, const nbnxn_bb_t *bb_i_ci,
411                                  int csj, const nbnxn_bb_t *bb_j_all)
412 {
413     // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
414     using namespace gmx;
415
416     Simd4Float bb_i_S0, bb_i_S1;
417     Simd4Float bb_j_S0, bb_j_S1;
418     Simd4Float dl_S;
419     Simd4Float dh_S;
420     Simd4Float dm_S;
421     Simd4Float dm0_S;
422
423     bb_i_S0 = load4(&bb_i_ci[si].lower[0]);
424     bb_i_S1 = load4(&bb_i_ci[si].upper[0]);
425     bb_j_S0 = load4(&bb_j_all[csj].lower[0]);
426     bb_j_S1 = load4(&bb_j_all[csj].upper[0]);
427
428     dl_S    = bb_i_S0 - bb_j_S1;
429     dh_S    = bb_j_S0 - bb_i_S1;
430
431     dm_S    = max(dl_S, dh_S);
432     dm0_S   = max(dm_S, simd4SetZeroF());
433
434     return dotProduct(dm0_S, dm0_S);
435 }
436
437 /* Calculate bb bounding distances of bb_i[si,...,si+3] and store them in d2 */
438 #define SUBC_BB_DIST2_SIMD4_XXXX_INNER(si, bb_i, d2) \
439     {                                                \
440         int               shi;                                  \
441                                                  \
442         Simd4Float        dx_0, dy_0, dz_0;                    \
443         Simd4Float        dx_1, dy_1, dz_1;                    \
444                                                  \
445         Simd4Float        mx, my, mz;                          \
446         Simd4Float        m0x, m0y, m0z;                       \
447                                                  \
448         Simd4Float        d2x, d2y, d2z;                       \
449         Simd4Float        d2s, d2t;                            \
450                                                  \
451         shi = si*NNBSBB_D*DIM;                       \
452                                                  \
453         xi_l = load4(bb_i+shi+0*STRIDE_PBB);   \
454         yi_l = load4(bb_i+shi+1*STRIDE_PBB);   \
455         zi_l = load4(bb_i+shi+2*STRIDE_PBB);   \
456         xi_h = load4(bb_i+shi+3*STRIDE_PBB);   \
457         yi_h = load4(bb_i+shi+4*STRIDE_PBB);   \
458         zi_h = load4(bb_i+shi+5*STRIDE_PBB);   \
459                                                  \
460         dx_0 = xi_l - xj_h;                 \
461         dy_0 = yi_l - yj_h;                 \
462         dz_0 = zi_l - zj_h;                 \
463                                                  \
464         dx_1 = xj_l - xi_h;                 \
465         dy_1 = yj_l - yi_h;                 \
466         dz_1 = zj_l - zi_h;                 \
467                                                  \
468         mx   = max(dx_0, dx_1);                 \
469         my   = max(dy_0, dy_1);                 \
470         mz   = max(dz_0, dz_1);                 \
471                                                  \
472         m0x  = max(mx, zero);                   \
473         m0y  = max(my, zero);                   \
474         m0z  = max(mz, zero);                   \
475                                                  \
476         d2x  = m0x * m0x;                   \
477         d2y  = m0y * m0y;                   \
478         d2z  = m0z * m0z;                   \
479                                                  \
480         d2s  = d2x + d2y;                   \
481         d2t  = d2s + d2z;                   \
482                                                  \
483         store4(d2+si, d2t);                      \
484     }
485
486 /* 4-wide SIMD code for nsi bb distances for bb format xxxxyyyyzzzz */
487 static void subc_bb_dist2_simd4_xxxx(const float *bb_j,
488                                      int nsi, const float *bb_i,
489                                      float *d2)
490 {
491     // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
492     using namespace gmx;
493
494     Simd4Float xj_l, yj_l, zj_l;
495     Simd4Float xj_h, yj_h, zj_h;
496     Simd4Float xi_l, yi_l, zi_l;
497     Simd4Float xi_h, yi_h, zi_h;
498
499     Simd4Float zero;
500
501     zero = setZero();
502
503     xj_l = Simd4Float(bb_j[0*STRIDE_PBB]);
504     yj_l = Simd4Float(bb_j[1*STRIDE_PBB]);
505     zj_l = Simd4Float(bb_j[2*STRIDE_PBB]);
506     xj_h = Simd4Float(bb_j[3*STRIDE_PBB]);
507     yj_h = Simd4Float(bb_j[4*STRIDE_PBB]);
508     zj_h = Simd4Float(bb_j[5*STRIDE_PBB]);
509
510     /* Here we "loop" over si (0,STRIDE_PBB) from 0 to nsi with step STRIDE_PBB.
511      * But as we know the number of iterations is 1 or 2, we unroll manually.
512      */
513     SUBC_BB_DIST2_SIMD4_XXXX_INNER(0, bb_i, d2);
514     if (STRIDE_PBB < nsi)
515     {
516         SUBC_BB_DIST2_SIMD4_XXXX_INNER(STRIDE_PBB, bb_i, d2);
517     }
518 }
519
520 #endif /* NBNXN_SEARCH_BB_SIMD4 */
521
522
523 /* Returns if any atom pair from two clusters is within distance sqrt(rl2) */
524 static gmx_inline gmx_bool
525 clusterpair_in_range(const nbnxn_list_work_t *work,
526                      int si,
527                      int csj, int stride, const real *x_j,
528                      real rl2)
529 {
530 #ifndef NBNXN_SEARCH_BB_SIMD4
531
532     /* Plain C version.
533      * All coordinates are stored as xyzxyz...
534      */
535
536     const real *x_i = work->x_ci;
537
538     for (int i = 0; i < nbnxn_gpu_cluster_size; i++)
539     {
540         int i0 = (si*nbnxn_gpu_cluster_size + i)*DIM;
541         for (int j = 0; j < nbnxn_gpu_cluster_size; j++)
542         {
543             int  j0 = (csj*nbnxn_gpu_cluster_size + j)*stride;
544
545             real d2 = gmx::square(x_i[i0  ] - x_j[j0  ]) + gmx::square(x_i[i0+1] - x_j[j0+1]) + gmx::square(x_i[i0+2] - x_j[j0+2]);
546
547             if (d2 < rl2)
548             {
549                 return TRUE;
550             }
551         }
552     }
553
554     return FALSE;
555
556 #else /* !NBNXN_SEARCH_BB_SIMD4 */
557
558     /* 4-wide SIMD version.
559      * A cluster is hard-coded to 8 atoms.
560      * The coordinates x_i are stored as xxxxyyyy..., x_j is stored xyzxyz...
561      * Using 8-wide AVX is not faster on Intel Sandy Bridge.
562      */
563     assert(nbnxn_gpu_cluster_size == 8);
564
565     Simd4Real   rc2_S      = Simd4Real(rl2);
566
567     const real *x_i        = work->x_ci_simd;
568
569     int         dim_stride = nbnxn_gpu_cluster_size*DIM;
570     Simd4Real   ix_S0      = load4(x_i + si*dim_stride + 0*GMX_SIMD4_WIDTH);
571     Simd4Real   iy_S0      = load4(x_i + si*dim_stride + 1*GMX_SIMD4_WIDTH);
572     Simd4Real   iz_S0      = load4(x_i + si*dim_stride + 2*GMX_SIMD4_WIDTH);
573     Simd4Real   ix_S1      = load4(x_i + si*dim_stride + 3*GMX_SIMD4_WIDTH);
574     Simd4Real   iy_S1      = load4(x_i + si*dim_stride + 4*GMX_SIMD4_WIDTH);
575     Simd4Real   iz_S1      = load4(x_i + si*dim_stride + 5*GMX_SIMD4_WIDTH);
576
577     /* We loop from the outer to the inner particles to maximize
578      * the chance that we find a pair in range quickly and return.
579      */
580     int j0 = csj*nbnxn_gpu_cluster_size;
581     int j1 = j0 + nbnxn_gpu_cluster_size - 1;
582     while (j0 < j1)
583     {
584         Simd4Real jx0_S, jy0_S, jz0_S;
585         Simd4Real jx1_S, jy1_S, jz1_S;
586
587         Simd4Real dx_S0, dy_S0, dz_S0;
588         Simd4Real dx_S1, dy_S1, dz_S1;
589         Simd4Real dx_S2, dy_S2, dz_S2;
590         Simd4Real dx_S3, dy_S3, dz_S3;
591
592         Simd4Real rsq_S0;
593         Simd4Real rsq_S1;
594         Simd4Real rsq_S2;
595         Simd4Real rsq_S3;
596
597         Simd4Bool wco_S0;
598         Simd4Bool wco_S1;
599         Simd4Bool wco_S2;
600         Simd4Bool wco_S3;
601         Simd4Bool wco_any_S01, wco_any_S23, wco_any_S;
602
603         jx0_S = Simd4Real(x_j[j0*stride+0]);
604         jy0_S = Simd4Real(x_j[j0*stride+1]);
605         jz0_S = Simd4Real(x_j[j0*stride+2]);
606
607         jx1_S = Simd4Real(x_j[j1*stride+0]);
608         jy1_S = Simd4Real(x_j[j1*stride+1]);
609         jz1_S = Simd4Real(x_j[j1*stride+2]);
610
611         /* Calculate distance */
612         dx_S0            = ix_S0 - jx0_S;
613         dy_S0            = iy_S0 - jy0_S;
614         dz_S0            = iz_S0 - jz0_S;
615         dx_S1            = ix_S1 - jx0_S;
616         dy_S1            = iy_S1 - jy0_S;
617         dz_S1            = iz_S1 - jz0_S;
618         dx_S2            = ix_S0 - jx1_S;
619         dy_S2            = iy_S0 - jy1_S;
620         dz_S2            = iz_S0 - jz1_S;
621         dx_S3            = ix_S1 - jx1_S;
622         dy_S3            = iy_S1 - jy1_S;
623         dz_S3            = iz_S1 - jz1_S;
624
625         /* rsq = dx*dx+dy*dy+dz*dz */
626         rsq_S0           = norm2(dx_S0, dy_S0, dz_S0);
627         rsq_S1           = norm2(dx_S1, dy_S1, dz_S1);
628         rsq_S2           = norm2(dx_S2, dy_S2, dz_S2);
629         rsq_S3           = norm2(dx_S3, dy_S3, dz_S3);
630
631         wco_S0           = (rsq_S0 < rc2_S);
632         wco_S1           = (rsq_S1 < rc2_S);
633         wco_S2           = (rsq_S2 < rc2_S);
634         wco_S3           = (rsq_S3 < rc2_S);
635
636         wco_any_S01      = wco_S0 || wco_S1;
637         wco_any_S23      = wco_S2 || wco_S3;
638         wco_any_S        = wco_any_S01 || wco_any_S23;
639
640         if (anyTrue(wco_any_S))
641         {
642             return TRUE;
643         }
644
645         j0++;
646         j1--;
647     }
648
649     return FALSE;
650
651 #endif /* !NBNXN_SEARCH_BB_SIMD4 */
652 }
653
654 /* Returns the j sub-cell for index cj_ind */
655 static int nbl_cj(const nbnxn_pairlist_t *nbl, int cj_ind)
656 {
657     return nbl->cj4[cj_ind/nbnxn_gpu_jgroup_size].cj[cj_ind & (nbnxn_gpu_jgroup_size - 1)];
658 }
659
660 /* Returns the i-interaction mask of the j sub-cell for index cj_ind */
661 static unsigned int nbl_imask0(const nbnxn_pairlist_t *nbl, int cj_ind)
662 {
663     return nbl->cj4[cj_ind/nbnxn_gpu_jgroup_size].imei[0].imask;
664 }
665
666 /* Ensures there is enough space for extra extra exclusion masks */
667 static void check_excl_space(nbnxn_pairlist_t *nbl, int extra)
668 {
669     if (nbl->nexcl+extra > nbl->excl_nalloc)
670     {
671         nbl->excl_nalloc = over_alloc_small(nbl->nexcl+extra);
672         nbnxn_realloc_void((void **)&nbl->excl,
673                            nbl->nexcl*sizeof(*nbl->excl),
674                            nbl->excl_nalloc*sizeof(*nbl->excl),
675                            nbl->alloc, nbl->free);
676     }
677 }
678
679 /* Ensures there is enough space for ncell extra j-cells in the list */
680 static void check_cell_list_space_simple(nbnxn_pairlist_t *nbl,
681                                          int               ncell)
682 {
683     int cj_max;
684
685     cj_max = nbl->ncj + ncell;
686
687     if (cj_max > nbl->cj_nalloc)
688     {
689         nbl->cj_nalloc = over_alloc_small(cj_max);
690         nbnxn_realloc_void((void **)&nbl->cj,
691                            nbl->ncj*sizeof(*nbl->cj),
692                            nbl->cj_nalloc*sizeof(*nbl->cj),
693                            nbl->alloc, nbl->free);
694     }
695 }
696
697 /* Ensures there is enough space for ncell extra j-clusters in the list */
698 static void check_cell_list_space_supersub(nbnxn_pairlist_t *nbl,
699                                            int               ncell)
700 {
701     int ncj4_max, w;
702
703     /* We can have maximally nsupercell*gpu_ncluster_per_cell sj lists */
704     /* We can store 4 j-subcell - i-supercell pairs in one struct.
705      * since we round down, we need one extra entry.
706      */
707     ncj4_max = ((nbl->work->cj_ind + ncell*gpu_ncluster_per_cell + nbnxn_gpu_jgroup_size - 1)/nbnxn_gpu_jgroup_size);
708
709     if (ncj4_max > nbl->cj4_nalloc)
710     {
711         nbl->cj4_nalloc = over_alloc_small(ncj4_max);
712         nbnxn_realloc_void((void **)&nbl->cj4,
713                            nbl->work->cj4_init*sizeof(*nbl->cj4),
714                            nbl->cj4_nalloc*sizeof(*nbl->cj4),
715                            nbl->alloc, nbl->free);
716     }
717
718     if (ncj4_max > nbl->work->cj4_init)
719     {
720         for (int j4 = nbl->work->cj4_init; j4 < ncj4_max; j4++)
721         {
722             /* No i-subcells and no excl's in the list initially */
723             for (w = 0; w < nbnxn_gpu_clusterpair_split; w++)
724             {
725                 nbl->cj4[j4].imei[w].imask    = 0U;
726                 nbl->cj4[j4].imei[w].excl_ind = 0;
727
728             }
729         }
730         nbl->work->cj4_init = ncj4_max;
731     }
732 }
733
734 /* Set all excl masks for one GPU warp no exclusions */
735 static void set_no_excls(nbnxn_excl_t *excl)
736 {
737     for (int t = 0; t < nbnxn_gpu_excl_size; t++)
738     {
739         /* Turn all interaction bits on */
740         excl->pair[t] = NBNXN_INTERACTION_MASK_ALL;
741     }
742 }
743
744 /* Initializes a single nbnxn_pairlist_t data structure */
745 static void nbnxn_init_pairlist(nbnxn_pairlist_t *nbl,
746                                 gmx_bool          bSimple,
747                                 nbnxn_alloc_t    *alloc,
748                                 nbnxn_free_t     *free)
749 {
750     if (alloc == NULL)
751     {
752         nbl->alloc = nbnxn_alloc_aligned;
753     }
754     else
755     {
756         nbl->alloc = alloc;
757     }
758     if (free == NULL)
759     {
760         nbl->free = nbnxn_free_aligned;
761     }
762     else
763     {
764         nbl->free = free;
765     }
766
767     nbl->bSimple     = bSimple;
768     nbl->na_sc       = 0;
769     nbl->na_ci       = 0;
770     nbl->na_cj       = 0;
771     nbl->nci         = 0;
772     nbl->ci          = NULL;
773     nbl->ci_nalloc   = 0;
774     nbl->nsci        = 0;
775     nbl->sci         = NULL;
776     nbl->sci_nalloc  = 0;
777     nbl->ncj         = 0;
778     nbl->cj          = NULL;
779     nbl->cj_nalloc   = 0;
780     nbl->ncj4        = 0;
781     /* We need one element extra in sj, so alloc initially with 1 */
782     nbl->cj4_nalloc  = 0;
783     nbl->cj4         = NULL;
784     nbl->nci_tot     = 0;
785
786     if (!nbl->bSimple)
787     {
788         GMX_ASSERT(nbnxn_gpu_ncluster_per_supercluster == gpu_ncluster_per_cell, "The search code assumes that the a super-cluster matches a search grid cell");
789
790         GMX_ASSERT(sizeof(nbl->cj4[0].imei[0].imask)*8 >= nbnxn_gpu_jgroup_size*gpu_ncluster_per_cell, "The i super-cluster cluster interaction mask does not contain a sufficient number of bits");
791         GMX_ASSERT(sizeof(nbl->excl[0])*8 >= nbnxn_gpu_jgroup_size*gpu_ncluster_per_cell, "The GPU exclusion mask does not contain a sufficient number of bits");
792
793         nbl->excl        = NULL;
794         nbl->excl_nalloc = 0;
795         nbl->nexcl       = 0;
796         check_excl_space(nbl, 1);
797         nbl->nexcl       = 1;
798         set_no_excls(&nbl->excl[0]);
799     }
800
801     snew(nbl->work, 1);
802     if (nbl->bSimple)
803     {
804         snew_aligned(nbl->work->bb_ci, 1, NBNXN_SEARCH_BB_MEM_ALIGN);
805     }
806     else
807     {
808 #ifdef NBNXN_BBXXXX
809         snew_aligned(nbl->work->pbb_ci, gpu_ncluster_per_cell/STRIDE_PBB*NNBSBB_XXXX, NBNXN_SEARCH_BB_MEM_ALIGN);
810 #else
811         snew_aligned(nbl->work->bb_ci, gpu_ncluster_per_cell, NBNXN_SEARCH_BB_MEM_ALIGN);
812 #endif
813     }
814     int gpu_clusterpair_nc = gpu_ncluster_per_cell*nbnxn_gpu_cluster_size*DIM;
815     snew(nbl->work->x_ci, gpu_clusterpair_nc);
816 #if GMX_SIMD
817     snew_aligned(nbl->work->x_ci_simd,
818                  std::max(NBNXN_CPU_CLUSTER_I_SIZE*DIM*GMX_SIMD_REAL_WIDTH,
819                           gpu_clusterpair_nc),
820                  GMX_SIMD_REAL_WIDTH);
821 #endif
822     snew_aligned(nbl->work->d2, gpu_ncluster_per_cell, NBNXN_SEARCH_BB_MEM_ALIGN);
823
824     nbl->work->sort            = NULL;
825     nbl->work->sort_nalloc     = 0;
826     nbl->work->sci_sort        = NULL;
827     nbl->work->sci_sort_nalloc = 0;
828 }
829
830 void nbnxn_init_pairlist_set(nbnxn_pairlist_set_t *nbl_list,
831                              gmx_bool bSimple, gmx_bool bCombined,
832                              nbnxn_alloc_t *alloc,
833                              nbnxn_free_t  *free)
834 {
835     nbl_list->bSimple   = bSimple;
836     nbl_list->bCombined = bCombined;
837
838     nbl_list->nnbl = gmx_omp_nthreads_get(emntNonbonded);
839
840     if (!nbl_list->bCombined &&
841         nbl_list->nnbl > NBNXN_BUFFERFLAG_MAX_THREADS)
842     {
843         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.",
844                   nbl_list->nnbl, NBNXN_BUFFERFLAG_MAX_THREADS, NBNXN_BUFFERFLAG_MAX_THREADS);
845     }
846
847     snew(nbl_list->nbl, nbl_list->nnbl);
848     snew(nbl_list->nbl_fep, nbl_list->nnbl);
849     /* Execute in order to avoid memory interleaving between threads */
850 #pragma omp parallel for num_threads(nbl_list->nnbl) schedule(static)
851     for (int i = 0; i < nbl_list->nnbl; i++)
852     {
853         try
854         {
855             /* Allocate the nblist data structure locally on each thread
856              * to optimize memory access for NUMA architectures.
857              */
858             snew(nbl_list->nbl[i], 1);
859
860             /* Only list 0 is used on the GPU, use normal allocation for i>0 */
861             if (i == 0)
862             {
863                 nbnxn_init_pairlist(nbl_list->nbl[i], nbl_list->bSimple, alloc, free);
864             }
865             else
866             {
867                 nbnxn_init_pairlist(nbl_list->nbl[i], nbl_list->bSimple, NULL, NULL);
868             }
869
870             snew(nbl_list->nbl_fep[i], 1);
871             nbnxn_init_pairlist_fep(nbl_list->nbl_fep[i]);
872         }
873         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
874     }
875 }
876
877 /* Print statistics of a pair list, used for debug output */
878 static void print_nblist_statistics_simple(FILE *fp, const nbnxn_pairlist_t *nbl,
879                                            const nbnxn_search_t nbs, real rl)
880 {
881     const nbnxn_grid_t *grid;
882     int                 cs[SHIFTS];
883     int                 npexcl;
884
885     /* This code only produces correct statistics with domain decomposition */
886     grid = &nbs->grid[0];
887
888     fprintf(fp, "nbl nci %d ncj %d\n",
889             nbl->nci, nbl->ncj);
890     fprintf(fp, "nbl na_sc %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
891             nbl->na_sc, rl, nbl->ncj, nbl->ncj/(double)grid->nc,
892             nbl->ncj/(double)grid->nc*grid->na_sc,
893             nbl->ncj/(double)grid->nc*grid->na_sc/(0.5*4.0/3.0*M_PI*rl*rl*rl*grid->nc*grid->na_sc/(grid->size[XX]*grid->size[YY]*grid->size[ZZ])));
894
895     fprintf(fp, "nbl average j cell list length %.1f\n",
896             0.25*nbl->ncj/(double)std::max(nbl->nci, 1));
897
898     for (int s = 0; s < SHIFTS; s++)
899     {
900         cs[s] = 0;
901     }
902     npexcl = 0;
903     for (int i = 0; i < nbl->nci; i++)
904     {
905         cs[nbl->ci[i].shift & NBNXN_CI_SHIFT] +=
906             nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start;
907
908         int j = nbl->ci[i].cj_ind_start;
909         while (j < nbl->ci[i].cj_ind_end &&
910                nbl->cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
911         {
912             npexcl++;
913             j++;
914         }
915     }
916     fprintf(fp, "nbl cell pairs, total: %d excl: %d %.1f%%\n",
917             nbl->ncj, npexcl, 100*npexcl/(double)std::max(nbl->ncj, 1));
918     for (int s = 0; s < SHIFTS; s++)
919     {
920         if (cs[s] > 0)
921         {
922             fprintf(fp, "nbl shift %2d ncj %3d\n", s, cs[s]);
923         }
924     }
925 }
926
927 /* Print statistics of a pair lists, used for debug output */
928 static void print_nblist_statistics_supersub(FILE *fp, const nbnxn_pairlist_t *nbl,
929                                              const nbnxn_search_t nbs, real rl)
930 {
931     const nbnxn_grid_t *grid;
932     int                 b;
933     int                 c[gpu_ncluster_per_cell + 1];
934     double              sum_nsp, sum_nsp2;
935     int                 nsp_max;
936
937     /* This code only produces correct statistics with domain decomposition */
938     grid = &nbs->grid[0];
939
940     fprintf(fp, "nbl nsci %d ncj4 %d nsi %d excl4 %d\n",
941             nbl->nsci, nbl->ncj4, nbl->nci_tot, nbl->nexcl);
942     fprintf(fp, "nbl na_c %d rl %g ncp %d per cell %.1f atoms %.1f ratio %.2f\n",
943             nbl->na_ci, rl, nbl->nci_tot, nbl->nci_tot/(double)grid->nsubc_tot,
944             nbl->nci_tot/(double)grid->nsubc_tot*grid->na_c,
945             nbl->nci_tot/(double)grid->nsubc_tot*grid->na_c/(0.5*4.0/3.0*M_PI*rl*rl*rl*grid->nsubc_tot*grid->na_c/(grid->size[XX]*grid->size[YY]*grid->size[ZZ])));
946
947     sum_nsp  = 0;
948     sum_nsp2 = 0;
949     nsp_max  = 0;
950     for (int si = 0; si <= gpu_ncluster_per_cell; si++)
951     {
952         c[si] = 0;
953     }
954     for (int i = 0; i < nbl->nsci; i++)
955     {
956         int nsp;
957
958         nsp = 0;
959         for (int j4 = nbl->sci[i].cj4_ind_start; j4 < nbl->sci[i].cj4_ind_end; j4++)
960         {
961             for (int j = 0; j < nbnxn_gpu_jgroup_size; j++)
962             {
963                 b = 0;
964                 for (int si = 0; si < gpu_ncluster_per_cell; si++)
965                 {
966                     if (nbl->cj4[j4].imei[0].imask & (1U << (j*gpu_ncluster_per_cell + si)))
967                     {
968                         b++;
969                     }
970                 }
971                 nsp += b;
972                 c[b]++;
973             }
974         }
975         sum_nsp  += nsp;
976         sum_nsp2 += nsp*nsp;
977         nsp_max   = std::max(nsp_max, nsp);
978     }
979     if (nbl->nsci > 0)
980     {
981         sum_nsp  /= nbl->nsci;
982         sum_nsp2 /= nbl->nsci;
983     }
984     fprintf(fp, "nbl #cluster-pairs: av %.1f stddev %.1f max %d\n",
985             sum_nsp, std::sqrt(sum_nsp2 - sum_nsp*sum_nsp), nsp_max);
986
987     if (nbl->ncj4 > 0)
988     {
989         for (b = 0; b <= gpu_ncluster_per_cell; b++)
990         {
991             fprintf(fp, "nbl j-list #i-subcell %d %7d %4.1f\n",
992                     b, c[b],
993                     100.0*c[b]/(double)(nbl->ncj4*nbnxn_gpu_jgroup_size));
994         }
995     }
996 }
997
998 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp */
999 static void low_get_nbl_exclusions(nbnxn_pairlist_t *nbl, int cj4,
1000                                    int warp, nbnxn_excl_t **excl)
1001 {
1002     if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
1003     {
1004         /* No exclusions set, make a new list entry */
1005         nbl->cj4[cj4].imei[warp].excl_ind = nbl->nexcl;
1006         nbl->nexcl++;
1007         *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
1008         set_no_excls(*excl);
1009     }
1010     else
1011     {
1012         /* We already have some exclusions, new ones can be added to the list */
1013         *excl = &nbl->excl[nbl->cj4[cj4].imei[warp].excl_ind];
1014     }
1015 }
1016
1017 /* Returns a pointer to the exclusion mask for cj4-unit cj4, warp warp,
1018  * generates a new element and allocates extra memory, if necessary.
1019  */
1020 static void get_nbl_exclusions_1(nbnxn_pairlist_t *nbl, int cj4,
1021                                  int warp, nbnxn_excl_t **excl)
1022 {
1023     if (nbl->cj4[cj4].imei[warp].excl_ind == 0)
1024     {
1025         /* We need to make a new list entry, check if we have space */
1026         check_excl_space(nbl, 1);
1027     }
1028     low_get_nbl_exclusions(nbl, cj4, warp, excl);
1029 }
1030
1031 /* Returns pointers to the exclusion masks for cj4-unit cj4 for both warps,
1032  * generates a new element and allocates extra memory, if necessary.
1033  */
1034 static void get_nbl_exclusions_2(nbnxn_pairlist_t *nbl, int cj4,
1035                                  nbnxn_excl_t **excl_w0,
1036                                  nbnxn_excl_t **excl_w1)
1037 {
1038     /* Check for space we might need */
1039     check_excl_space(nbl, 2);
1040
1041     low_get_nbl_exclusions(nbl, cj4, 0, excl_w0);
1042     low_get_nbl_exclusions(nbl, cj4, 1, excl_w1);
1043 }
1044
1045 /* Sets the self exclusions i=j and pair exclusions i>j */
1046 static void set_self_and_newton_excls_supersub(nbnxn_pairlist_t *nbl,
1047                                                int cj4_ind, int sj_offset,
1048                                                int i_cluster_in_cell)
1049 {
1050     nbnxn_excl_t *excl[nbnxn_gpu_clusterpair_split];
1051
1052     /* Here we only set the set self and double pair exclusions */
1053
1054     assert(nbnxn_gpu_clusterpair_split == 2);
1055
1056     get_nbl_exclusions_2(nbl, cj4_ind, &excl[0], &excl[1]);
1057
1058     /* Only minor < major bits set */
1059     for (int ej = 0; ej < nbl->na_ci; ej++)
1060     {
1061         int w = (ej>>2);
1062         for (int ei = ej; ei < nbl->na_ci; ei++)
1063         {
1064             excl[w]->pair[(ej & (nbnxn_gpu_jgroup_size-1))*nbl->na_ci + ei] &=
1065                 ~(1U << (sj_offset*gpu_ncluster_per_cell + i_cluster_in_cell));
1066         }
1067     }
1068 }
1069
1070 /* Returns a diagonal or off-diagonal interaction mask for plain C lists */
1071 static unsigned int get_imask(gmx_bool rdiag, int ci, int cj)
1072 {
1073     return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
1074 }
1075
1076 /* Returns a diagonal or off-diagonal interaction mask for cj-size=2 */
1077 static unsigned int get_imask_simd_j2(gmx_bool rdiag, int ci, int cj)
1078 {
1079     return (rdiag && ci*2 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_0 :
1080             (rdiag && ci*2+1 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_1 :
1081              NBNXN_INTERACTION_MASK_ALL));
1082 }
1083
1084 /* Returns a diagonal or off-diagonal interaction mask for cj-size=4 */
1085 static unsigned int get_imask_simd_j4(gmx_bool rdiag, int ci, int cj)
1086 {
1087     return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
1088 }
1089
1090 /* Returns a diagonal or off-diagonal interaction mask for cj-size=8 */
1091 static unsigned int get_imask_simd_j8(gmx_bool rdiag, int ci, int cj)
1092 {
1093     return (rdiag && ci == cj*2 ? NBNXN_INTERACTION_MASK_DIAG_J8_0 :
1094             (rdiag && ci == cj*2+1 ? NBNXN_INTERACTION_MASK_DIAG_J8_1 :
1095              NBNXN_INTERACTION_MASK_ALL));
1096 }
1097
1098 #if GMX_SIMD
1099 #if GMX_SIMD_REAL_WIDTH == 2
1100 #define get_imask_simd_4xn  get_imask_simd_j2
1101 #endif
1102 #if GMX_SIMD_REAL_WIDTH == 4
1103 #define get_imask_simd_4xn  get_imask_simd_j4
1104 #endif
1105 #if GMX_SIMD_REAL_WIDTH == 8
1106 #define get_imask_simd_4xn  get_imask_simd_j8
1107 #define get_imask_simd_2xnn get_imask_simd_j4
1108 #endif
1109 #if GMX_SIMD_REAL_WIDTH == 16
1110 #define get_imask_simd_2xnn get_imask_simd_j8
1111 #endif
1112 #endif
1113
1114 /* Plain C code for making a pair list of cell ci vs cell cjf-cjl.
1115  * Checks bounding box distances and possibly atom pair distances.
1116  */
1117 static void make_cluster_list_simple(const nbnxn_grid_t *gridj,
1118                                      nbnxn_pairlist_t *nbl,
1119                                      int ci, int cjf, int cjl,
1120                                      gmx_bool remove_sub_diag,
1121                                      const real *x_j,
1122                                      real rl2, float rbb2,
1123                                      int *ndistc)
1124 {
1125     const nbnxn_bb_t        *bb_ci;
1126     const real              *x_ci;
1127
1128     gmx_bool                 InRange;
1129     real                     d2;
1130     int                      cjf_gl, cjl_gl;
1131
1132     bb_ci = nbl->work->bb_ci;
1133     x_ci  = nbl->work->x_ci;
1134
1135     InRange = FALSE;
1136     while (!InRange && cjf <= cjl)
1137     {
1138         d2       = subc_bb_dist2(0, bb_ci, cjf, gridj->bb);
1139         *ndistc += 2;
1140
1141         /* Check if the distance is within the distance where
1142          * we use only the bounding box distance rbb,
1143          * or within the cut-off and there is at least one atom pair
1144          * within the cut-off.
1145          */
1146         if (d2 < rbb2)
1147         {
1148             InRange = TRUE;
1149         }
1150         else if (d2 < rl2)
1151         {
1152             cjf_gl = gridj->cell0 + cjf;
1153             for (int i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
1154             {
1155                 for (int j = 0; j < NBNXN_CPU_CLUSTER_I_SIZE; j++)
1156                 {
1157                     InRange = InRange ||
1158                         (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
1159                          gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
1160                          gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjf_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rl2);
1161                 }
1162             }
1163             *ndistc += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
1164         }
1165         if (!InRange)
1166         {
1167             cjf++;
1168         }
1169     }
1170     if (!InRange)
1171     {
1172         return;
1173     }
1174
1175     InRange = FALSE;
1176     while (!InRange && cjl > cjf)
1177     {
1178         d2       = subc_bb_dist2(0, bb_ci, cjl, gridj->bb);
1179         *ndistc += 2;
1180
1181         /* Check if the distance is within the distance where
1182          * we use only the bounding box distance rbb,
1183          * or within the cut-off and there is at least one atom pair
1184          * within the cut-off.
1185          */
1186         if (d2 < rbb2)
1187         {
1188             InRange = TRUE;
1189         }
1190         else if (d2 < rl2)
1191         {
1192             cjl_gl = gridj->cell0 + cjl;
1193             for (int i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE && !InRange; i++)
1194             {
1195                 for (int j = 0; j < NBNXN_CPU_CLUSTER_I_SIZE; j++)
1196                 {
1197                     InRange = InRange ||
1198                         (gmx::square(x_ci[i*STRIDE_XYZ+XX] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+XX]) +
1199                          gmx::square(x_ci[i*STRIDE_XYZ+YY] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+YY]) +
1200                          gmx::square(x_ci[i*STRIDE_XYZ+ZZ] - x_j[(cjl_gl*NBNXN_CPU_CLUSTER_I_SIZE+j)*STRIDE_XYZ+ZZ]) < rl2);
1201                 }
1202             }
1203             *ndistc += NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
1204         }
1205         if (!InRange)
1206         {
1207             cjl--;
1208         }
1209     }
1210
1211     if (cjf <= cjl)
1212     {
1213         for (int cj = cjf; cj <= cjl; cj++)
1214         {
1215             /* Store cj and the interaction mask */
1216             nbl->cj[nbl->ncj].cj   = gridj->cell0 + cj;
1217             nbl->cj[nbl->ncj].excl = get_imask(remove_sub_diag, ci, cj);
1218             nbl->ncj++;
1219         }
1220         /* Increase the closing index in i super-cell list */
1221         nbl->ci[nbl->nci].cj_ind_end = nbl->ncj;
1222     }
1223 }
1224
1225 #ifdef GMX_NBNXN_SIMD_4XN
1226 #include "gromacs/mdlib/nbnxn_search_simd_4xn.h"
1227 #endif
1228 #ifdef GMX_NBNXN_SIMD_2XNN
1229 #include "gromacs/mdlib/nbnxn_search_simd_2xnn.h"
1230 #endif
1231
1232 /* Plain C or SIMD4 code for making a pair list of super-cell sci vs scj.
1233  * Checks bounding box distances and possibly atom pair distances.
1234  */
1235 static void make_cluster_list_supersub(const nbnxn_grid_t *gridi,
1236                                        const nbnxn_grid_t *gridj,
1237                                        nbnxn_pairlist_t *nbl,
1238                                        int sci, int scj,
1239                                        gmx_bool sci_equals_scj,
1240                                        int stride, const real *x,
1241                                        real rl2, float rbb2,
1242                                        int *ndistc)
1243 {
1244     nbnxn_list_work_t *work   = nbl->work;
1245
1246 #ifdef NBNXN_BBXXXX
1247     const float       *pbb_ci = work->pbb_ci;
1248 #else
1249     const nbnxn_bb_t  *bb_ci  = work->bb_ci;
1250 #endif
1251
1252     const int na_c = nbnxn_gpu_cluster_size;
1253     assert(na_c == gridi->na_c);
1254     assert(na_c == gridj->na_c);
1255
1256     /* We generate the pairlist mainly based on bounding-box distances
1257      * and do atom pair distance based pruning on the GPU.
1258      * Only if a j-group contains a single cluster-pair, we try to prune
1259      * that pair based on atom distances on the CPU to avoid empty j-groups.
1260      */
1261 #define PRUNE_LIST_CPU_ONE
1262 #ifdef PRUNE_LIST_CPU_ONE
1263     int  ci_last = -1;
1264 #endif
1265
1266     float *d2l = work->d2;
1267
1268     for (int subc = 0; subc < gridj->nsubc[scj]; subc++)
1269     {
1270         int          cj4_ind   = nbl->work->cj_ind/nbnxn_gpu_jgroup_size;
1271         int          cj_offset = nbl->work->cj_ind - cj4_ind*nbnxn_gpu_jgroup_size;
1272         nbnxn_cj4_t *cj4       = &nbl->cj4[cj4_ind];
1273
1274         int          cj        = scj*gpu_ncluster_per_cell + subc;
1275
1276         int          cj_gl     = gridj->cell0*gpu_ncluster_per_cell + cj;
1277
1278         /* Initialize this j-subcell i-subcell list */
1279         cj4->cj[cj_offset] = cj_gl;
1280
1281         int ci1;
1282         if (sci_equals_scj)
1283         {
1284             ci1 = subc + 1;
1285         }
1286         else
1287         {
1288             ci1 = gridi->nsubc[sci];
1289         }
1290
1291 #ifdef NBNXN_BBXXXX
1292         /* Determine all ci1 bb distances in one call with SIMD4 */
1293         subc_bb_dist2_simd4_xxxx(gridj->pbb+(cj>>STRIDE_PBB_2LOG)*NNBSBB_XXXX+(cj & (STRIDE_PBB-1)),
1294                                  ci1, pbb_ci, d2l);
1295         *ndistc += na_c*2;
1296 #endif
1297
1298         int          npair = 0;
1299         unsigned int imask = 0;
1300         /* We use a fixed upper-bound instead of ci1 to help optimization */
1301         for (int ci = 0; ci < gpu_ncluster_per_cell; ci++)
1302         {
1303             if (ci == ci1)
1304             {
1305                 break;
1306             }
1307
1308 #ifndef NBNXN_BBXXXX
1309             /* Determine the bb distance between ci and cj */
1310             d2l[ci]  = subc_bb_dist2(ci, bb_ci, cj, gridj->bb);
1311             *ndistc += 2;
1312 #endif
1313             float d2 = d2l[ci];
1314
1315 #ifdef PRUNE_LIST_CPU_ALL
1316             /* Check if the distance is within the distance where
1317              * we use only the bounding box distance rbb,
1318              * or within the cut-off and there is at least one atom pair
1319              * within the cut-off. This check is very costly.
1320              */
1321             *ndistc += na_c*na_c;
1322             if (d2 < rbb2 ||
1323                 (d2 < rl2 &&
1324                  clusterpair_in_range(work, ci, cj_gl, stride, x, rl2)))
1325 #else
1326             /* Check if the distance between the two bounding boxes
1327              * in within the pair-list cut-off.
1328              */
1329             if (d2 < rl2)
1330 #endif
1331             {
1332                 /* Flag this i-subcell to be taken into account */
1333                 imask |= (1U << (cj_offset*gpu_ncluster_per_cell + ci));
1334
1335 #ifdef PRUNE_LIST_CPU_ONE
1336                 ci_last = ci;
1337 #endif
1338
1339                 npair++;
1340             }
1341         }
1342
1343 #ifdef PRUNE_LIST_CPU_ONE
1344         /* If we only found 1 pair, check if any atoms are actually
1345          * within the cut-off, so we could get rid of it.
1346          */
1347         if (npair == 1 && d2l[ci_last] >= rbb2 &&
1348             !clusterpair_in_range(work, ci_last, cj_gl, stride, x, rl2))
1349         {
1350             imask &= ~(1U << (cj_offset*gpu_ncluster_per_cell + ci_last));
1351             npair--;
1352         }
1353 #endif
1354
1355         if (npair > 0)
1356         {
1357             /* We have a useful sj entry, close it now */
1358
1359             /* Set the exclucions for the ci== sj entry.
1360              * Here we don't bother to check if this entry is actually flagged,
1361              * as it will nearly always be in the list.
1362              */
1363             if (sci_equals_scj)
1364             {
1365                 set_self_and_newton_excls_supersub(nbl, cj4_ind, cj_offset, subc);
1366             }
1367
1368             /* Copy the cluster interaction mask to the list */
1369             for (int w = 0; w < nbnxn_gpu_clusterpair_split; w++)
1370             {
1371                 cj4->imei[w].imask |= imask;
1372             }
1373
1374             nbl->work->cj_ind++;
1375
1376             /* Keep the count */
1377             nbl->nci_tot += npair;
1378
1379             /* Increase the closing index in i super-cell list */
1380             nbl->sci[nbl->nsci].cj4_ind_end =
1381                 (nbl->work->cj_ind + nbnxn_gpu_jgroup_size - 1)/nbnxn_gpu_jgroup_size;
1382         }
1383     }
1384 }
1385
1386 /* Set all atom-pair exclusions from the topology stored in excl
1387  * as masks in the pair-list for simple list i-entry nbl_ci
1388  */
1389 static void set_ci_top_excls(const nbnxn_search_t nbs,
1390                              nbnxn_pairlist_t    *nbl,
1391                              gmx_bool             diagRemoved,
1392                              int                  na_ci_2log,
1393                              int                  na_cj_2log,
1394                              const nbnxn_ci_t    *nbl_ci,
1395                              const t_blocka      *excl)
1396 {
1397     const int    *cell;
1398     int           ci;
1399     int           cj_ind_first, cj_ind_last;
1400     int           cj_first, cj_last;
1401     int           ndirect;
1402     int           ai, aj, si, ge, se;
1403     int           found, cj_ind_0, cj_ind_1, cj_ind_m;
1404     int           cj_m;
1405     int           inner_i, inner_e;
1406
1407     cell = nbs->cell;
1408
1409     if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
1410     {
1411         /* Empty list */
1412         return;
1413     }
1414
1415     ci = nbl_ci->ci;
1416
1417     cj_ind_first = nbl_ci->cj_ind_start;
1418     cj_ind_last  = nbl->ncj - 1;
1419
1420     cj_first = nbl->cj[cj_ind_first].cj;
1421     cj_last  = nbl->cj[cj_ind_last].cj;
1422
1423     /* Determine how many contiguous j-cells we have starting
1424      * from the first i-cell. This number can be used to directly
1425      * calculate j-cell indices for excluded atoms.
1426      */
1427     ndirect = 0;
1428     if (na_ci_2log == na_cj_2log)
1429     {
1430         while (cj_ind_first + ndirect <= cj_ind_last &&
1431                nbl->cj[cj_ind_first+ndirect].cj == ci + ndirect)
1432         {
1433             ndirect++;
1434         }
1435     }
1436 #ifdef NBNXN_SEARCH_BB_SIMD4
1437     else
1438     {
1439         while (cj_ind_first + ndirect <= cj_ind_last &&
1440                nbl->cj[cj_ind_first+ndirect].cj == ci_to_cj(na_cj_2log, ci) + ndirect)
1441         {
1442             ndirect++;
1443         }
1444     }
1445 #endif
1446
1447     /* Loop over the atoms in the i super-cell */
1448     for (int i = 0; i < nbl->na_sc; i++)
1449     {
1450         ai = nbs->a[ci*nbl->na_sc+i];
1451         if (ai >= 0)
1452         {
1453             si  = (i>>na_ci_2log);
1454
1455             /* Loop over the topology-based exclusions for this i-atom */
1456             for (int eind = excl->index[ai]; eind < excl->index[ai+1]; eind++)
1457             {
1458                 aj = excl->a[eind];
1459
1460                 if (aj == ai)
1461                 {
1462                     /* The self exclusion are already set, save some time */
1463                     continue;
1464                 }
1465
1466                 ge = cell[aj];
1467
1468                 /* Without shifts we only calculate interactions j>i
1469                  * for one-way pair-lists.
1470                  */
1471                 if (diagRemoved && ge <= ci*nbl->na_sc + i)
1472                 {
1473                     continue;
1474                 }
1475
1476                 se = (ge >> na_cj_2log);
1477
1478                 /* Could the cluster se be in our list? */
1479                 if (se >= cj_first && se <= cj_last)
1480                 {
1481                     if (se < cj_first + ndirect)
1482                     {
1483                         /* We can calculate cj_ind directly from se */
1484                         found = cj_ind_first + se - cj_first;
1485                     }
1486                     else
1487                     {
1488                         /* Search for se using bisection */
1489                         found    = -1;
1490                         cj_ind_0 = cj_ind_first + ndirect;
1491                         cj_ind_1 = cj_ind_last + 1;
1492                         while (found == -1 && cj_ind_0 < cj_ind_1)
1493                         {
1494                             cj_ind_m = (cj_ind_0 + cj_ind_1)>>1;
1495
1496                             cj_m = nbl->cj[cj_ind_m].cj;
1497
1498                             if (se == cj_m)
1499                             {
1500                                 found = cj_ind_m;
1501                             }
1502                             else if (se < cj_m)
1503                             {
1504                                 cj_ind_1 = cj_ind_m;
1505                             }
1506                             else
1507                             {
1508                                 cj_ind_0 = cj_ind_m + 1;
1509                             }
1510                         }
1511                     }
1512
1513                     if (found >= 0)
1514                     {
1515                         inner_i = i  - (si << na_ci_2log);
1516                         inner_e = ge - (se << na_cj_2log);
1517
1518                         nbl->cj[found].excl &= ~(1U<<((inner_i<<na_cj_2log) + inner_e));
1519                     }
1520                 }
1521             }
1522         }
1523     }
1524 }
1525
1526 /* Add a new i-entry to the FEP list and copy the i-properties */
1527 static gmx_inline void fep_list_new_nri_copy(t_nblist *nlist)
1528 {
1529     /* Add a new i-entry */
1530     nlist->nri++;
1531
1532     assert(nlist->nri < nlist->maxnri);
1533
1534     /* Duplicate the last i-entry, except for jindex, which continues */
1535     nlist->iinr[nlist->nri]   = nlist->iinr[nlist->nri-1];
1536     nlist->shift[nlist->nri]  = nlist->shift[nlist->nri-1];
1537     nlist->gid[nlist->nri]    = nlist->gid[nlist->nri-1];
1538     nlist->jindex[nlist->nri] = nlist->nrj;
1539 }
1540
1541 /* For load balancing of the free-energy lists over threads, we set
1542  * the maximum nrj size of an i-entry to 40. This leads to good
1543  * load balancing in the worst case scenario of a single perturbed
1544  * particle on 16 threads, while not introducing significant overhead.
1545  * Note that half of the perturbed pairs will anyhow end up in very small lists,
1546  * since non perturbed i-particles will see few perturbed j-particles).
1547  */
1548 const int max_nrj_fep = 40;
1549
1550 /* Exclude the perturbed pairs from the Verlet list. This is only done to avoid
1551  * singularities for overlapping particles (0/0), since the charges and
1552  * LJ parameters have been zeroed in the nbnxn data structure.
1553  * Simultaneously make a group pair list for the perturbed pairs.
1554  */
1555 static void make_fep_list(const nbnxn_search_t    nbs,
1556                           const nbnxn_atomdata_t *nbat,
1557                           nbnxn_pairlist_t       *nbl,
1558                           gmx_bool                bDiagRemoved,
1559                           nbnxn_ci_t             *nbl_ci,
1560                           const nbnxn_grid_t     *gridi,
1561                           const nbnxn_grid_t     *gridj,
1562                           t_nblist               *nlist)
1563 {
1564     int      ci, cj_ind_start, cj_ind_end, cja, cjr;
1565     int      nri_max;
1566     int      ngid, gid_i = 0, gid_j, gid;
1567     int      egp_shift, egp_mask;
1568     int      gid_cj = 0;
1569     int      ind_i, ind_j, ai, aj;
1570     int      nri;
1571     gmx_bool bFEP_i, bFEP_i_all;
1572
1573     if (nbl_ci->cj_ind_end == nbl_ci->cj_ind_start)
1574     {
1575         /* Empty list */
1576         return;
1577     }
1578
1579     ci = nbl_ci->ci;
1580
1581     cj_ind_start = nbl_ci->cj_ind_start;
1582     cj_ind_end   = nbl_ci->cj_ind_end;
1583
1584     /* In worst case we have alternating energy groups
1585      * and create #atom-pair lists, which means we need the size
1586      * of a cluster pair (na_ci*na_cj) times the number of cj's.
1587      */
1588     nri_max = nbl->na_ci*nbl->na_cj*(cj_ind_end - cj_ind_start);
1589     if (nlist->nri + nri_max > nlist->maxnri)
1590     {
1591         nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1592         reallocate_nblist(nlist);
1593     }
1594
1595     ngid = nbat->nenergrp;
1596
1597     if (static_cast<std::size_t>(ngid*gridj->na_cj) > sizeof(gid_cj)*8)
1598     {
1599         gmx_fatal(FARGS, "The Verlet scheme with %dx%d kernels and free-energy only supports up to %d energy groups",
1600                   gridi->na_c, gridj->na_cj, (sizeof(gid_cj)*8)/gridj->na_cj);
1601     }
1602
1603     egp_shift = nbat->neg_2log;
1604     egp_mask  = (1<<nbat->neg_2log) - 1;
1605
1606     /* Loop over the atoms in the i sub-cell */
1607     bFEP_i_all = TRUE;
1608     for (int i = 0; i < nbl->na_ci; i++)
1609     {
1610         ind_i = ci*nbl->na_ci + i;
1611         ai    = nbs->a[ind_i];
1612         if (ai >= 0)
1613         {
1614             nri                  = nlist->nri;
1615             nlist->jindex[nri+1] = nlist->jindex[nri];
1616             nlist->iinr[nri]     = ai;
1617             /* The actual energy group pair index is set later */
1618             nlist->gid[nri]      = 0;
1619             nlist->shift[nri]    = nbl_ci->shift & NBNXN_CI_SHIFT;
1620
1621             bFEP_i = gridi->fep[ci - gridi->cell0] & (1 << i);
1622
1623             bFEP_i_all = bFEP_i_all && bFEP_i;
1624
1625             if ((nlist->nrj + cj_ind_end - cj_ind_start)*nbl->na_cj > nlist->maxnrj)
1626             {
1627                 nlist->maxnrj = over_alloc_small((nlist->nrj + cj_ind_end - cj_ind_start)*nbl->na_cj);
1628                 srenew(nlist->jjnr,     nlist->maxnrj);
1629                 srenew(nlist->excl_fep, nlist->maxnrj);
1630             }
1631
1632             if (ngid > 1)
1633             {
1634                 gid_i = (nbat->energrp[ci] >> (egp_shift*i)) & egp_mask;
1635             }
1636
1637             for (int cj_ind = cj_ind_start; cj_ind < cj_ind_end; cj_ind++)
1638             {
1639                 unsigned int fep_cj;
1640
1641                 cja = nbl->cj[cj_ind].cj;
1642
1643                 if (gridj->na_cj == gridj->na_c)
1644                 {
1645                     cjr    = cja - gridj->cell0;
1646                     fep_cj = gridj->fep[cjr];
1647                     if (ngid > 1)
1648                     {
1649                         gid_cj = nbat->energrp[cja];
1650                     }
1651                 }
1652                 else if (2*gridj->na_cj == gridj->na_c)
1653                 {
1654                     cjr    = cja - gridj->cell0*2;
1655                     /* Extract half of the ci fep/energrp mask */
1656                     fep_cj = (gridj->fep[cjr>>1] >> ((cjr&1)*gridj->na_cj)) & ((1<<gridj->na_cj) - 1);
1657                     if (ngid > 1)
1658                     {
1659                         gid_cj = nbat->energrp[cja>>1] >> ((cja&1)*gridj->na_cj*egp_shift) & ((1<<(gridj->na_cj*egp_shift)) - 1);
1660                     }
1661                 }
1662                 else
1663                 {
1664                     cjr    = cja - (gridj->cell0>>1);
1665                     /* Combine two ci fep masks/energrp */
1666                     fep_cj = gridj->fep[cjr*2] + (gridj->fep[cjr*2+1] << gridj->na_c);
1667                     if (ngid > 1)
1668                     {
1669                         gid_cj = nbat->energrp[cja*2] + (nbat->energrp[cja*2+1] << (gridj->na_c*egp_shift));
1670                     }
1671                 }
1672
1673                 if (bFEP_i || fep_cj != 0)
1674                 {
1675                     for (int j = 0; j < nbl->na_cj; j++)
1676                     {
1677                         /* Is this interaction perturbed and not excluded? */
1678                         ind_j = cja*nbl->na_cj + j;
1679                         aj    = nbs->a[ind_j];
1680                         if (aj >= 0 &&
1681                             (bFEP_i || (fep_cj & (1 << j))) &&
1682                             (!bDiagRemoved || ind_j >= ind_i))
1683                         {
1684                             if (ngid > 1)
1685                             {
1686                                 gid_j = (gid_cj >> (j*egp_shift)) & egp_mask;
1687                                 gid   = GID(gid_i, gid_j, ngid);
1688
1689                                 if (nlist->nrj > nlist->jindex[nri] &&
1690                                     nlist->gid[nri] != gid)
1691                                 {
1692                                     /* Energy group pair changed: new list */
1693                                     fep_list_new_nri_copy(nlist);
1694                                     nri = nlist->nri;
1695                                 }
1696                                 nlist->gid[nri] = gid;
1697                             }
1698
1699                             if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1700                             {
1701                                 fep_list_new_nri_copy(nlist);
1702                                 nri = nlist->nri;
1703                             }
1704
1705                             /* Add it to the FEP list */
1706                             nlist->jjnr[nlist->nrj]     = aj;
1707                             nlist->excl_fep[nlist->nrj] = (nbl->cj[cj_ind].excl >> (i*nbl->na_cj + j)) & 1;
1708                             nlist->nrj++;
1709
1710                             /* Exclude it from the normal list.
1711                              * Note that the charge has been set to zero,
1712                              * but we need to avoid 0/0, as perturbed atoms
1713                              * can be on top of each other.
1714                              */
1715                             nbl->cj[cj_ind].excl &= ~(1U << (i*nbl->na_cj + j));
1716                         }
1717                     }
1718                 }
1719             }
1720
1721             if (nlist->nrj > nlist->jindex[nri])
1722             {
1723                 /* Actually add this new, non-empty, list */
1724                 nlist->nri++;
1725                 nlist->jindex[nlist->nri] = nlist->nrj;
1726             }
1727         }
1728     }
1729
1730     if (bFEP_i_all)
1731     {
1732         /* All interactions are perturbed, we can skip this entry */
1733         nbl_ci->cj_ind_end = cj_ind_start;
1734     }
1735 }
1736
1737 /* Return the index of atom a within a cluster */
1738 static gmx_inline int cj_mod_cj4(int cj)
1739 {
1740     return cj & (nbnxn_gpu_jgroup_size - 1);
1741 }
1742
1743 /* Convert a j-cluster to a cj4 group */
1744 static gmx_inline int cj_to_cj4(int cj)
1745 {
1746     return cj/nbnxn_gpu_jgroup_size;
1747 }
1748
1749 /* Return the index of an j-atom within a warp */
1750 static gmx_inline int a_mod_wj(int a)
1751 {
1752     return a & (nbnxn_gpu_cluster_size/nbnxn_gpu_clusterpair_split - 1);
1753 }
1754
1755 /* As make_fep_list above, but for super/sub lists. */
1756 static void make_fep_list_supersub(const nbnxn_search_t    nbs,
1757                                    const nbnxn_atomdata_t *nbat,
1758                                    nbnxn_pairlist_t       *nbl,
1759                                    gmx_bool                bDiagRemoved,
1760                                    const nbnxn_sci_t      *nbl_sci,
1761                                    real                    shx,
1762                                    real                    shy,
1763                                    real                    shz,
1764                                    real                    rlist_fep2,
1765                                    const nbnxn_grid_t     *gridi,
1766                                    const nbnxn_grid_t     *gridj,
1767                                    t_nblist               *nlist)
1768 {
1769     int                sci, cj4_ind_start, cj4_ind_end, cjr;
1770     int                nri_max;
1771     int                c_abs;
1772     int                ind_i, ind_j, ai, aj;
1773     int                nri;
1774     gmx_bool           bFEP_i;
1775     real               xi, yi, zi;
1776     const nbnxn_cj4_t *cj4;
1777
1778     if (nbl_sci->cj4_ind_end == nbl_sci->cj4_ind_start)
1779     {
1780         /* Empty list */
1781         return;
1782     }
1783
1784     sci = nbl_sci->sci;
1785
1786     cj4_ind_start = nbl_sci->cj4_ind_start;
1787     cj4_ind_end   = nbl_sci->cj4_ind_end;
1788
1789     /* Here we process one super-cell, max #atoms na_sc, versus a list
1790      * cj4 entries, each with max nbnxn_gpu_jgroup_size cj's, each
1791      * of size na_cj atoms.
1792      * On the GPU we don't support energy groups (yet).
1793      * So for each of the na_sc i-atoms, we need max one FEP list
1794      * for each max_nrj_fep j-atoms.
1795      */
1796     nri_max = nbl->na_sc*nbl->na_cj*(1 + ((cj4_ind_end - cj4_ind_start)*nbnxn_gpu_jgroup_size)/max_nrj_fep);
1797     if (nlist->nri + nri_max > nlist->maxnri)
1798     {
1799         nlist->maxnri = over_alloc_large(nlist->nri + nri_max);
1800         reallocate_nblist(nlist);
1801     }
1802
1803     /* Loop over the atoms in the i super-cluster */
1804     for (int c = 0; c < gpu_ncluster_per_cell; c++)
1805     {
1806         c_abs = sci*gpu_ncluster_per_cell + c;
1807
1808         for (int i = 0; i < nbl->na_ci; i++)
1809         {
1810             ind_i = c_abs*nbl->na_ci + i;
1811             ai    = nbs->a[ind_i];
1812             if (ai >= 0)
1813             {
1814                 nri                  = nlist->nri;
1815                 nlist->jindex[nri+1] = nlist->jindex[nri];
1816                 nlist->iinr[nri]     = ai;
1817                 /* With GPUs, energy groups are not supported */
1818                 nlist->gid[nri]      = 0;
1819                 nlist->shift[nri]    = nbl_sci->shift & NBNXN_CI_SHIFT;
1820
1821                 bFEP_i = (gridi->fep[c_abs - gridi->cell0*gpu_ncluster_per_cell] & (1 << i));
1822
1823                 xi = nbat->x[ind_i*nbat->xstride+XX] + shx;
1824                 yi = nbat->x[ind_i*nbat->xstride+YY] + shy;
1825                 zi = nbat->x[ind_i*nbat->xstride+ZZ] + shz;
1826
1827                 if ((nlist->nrj + cj4_ind_end - cj4_ind_start)*nbnxn_gpu_jgroup_size*nbl->na_cj > nlist->maxnrj)
1828                 {
1829                     nlist->maxnrj = over_alloc_small((nlist->nrj + cj4_ind_end - cj4_ind_start)*nbnxn_gpu_jgroup_size*nbl->na_cj);
1830                     srenew(nlist->jjnr,     nlist->maxnrj);
1831                     srenew(nlist->excl_fep, nlist->maxnrj);
1832                 }
1833
1834                 for (int cj4_ind = cj4_ind_start; cj4_ind < cj4_ind_end; cj4_ind++)
1835                 {
1836                     cj4 = &nbl->cj4[cj4_ind];
1837
1838                     for (int gcj = 0; gcj < nbnxn_gpu_jgroup_size; gcj++)
1839                     {
1840                         unsigned int fep_cj;
1841
1842                         if ((cj4->imei[0].imask & (1U << (gcj*gpu_ncluster_per_cell + c))) == 0)
1843                         {
1844                             /* Skip this ci for this cj */
1845                             continue;
1846                         }
1847
1848                         cjr = cj4->cj[gcj] - gridj->cell0*gpu_ncluster_per_cell;
1849
1850                         fep_cj = gridj->fep[cjr];
1851
1852                         if (bFEP_i || fep_cj != 0)
1853                         {
1854                             for (int j = 0; j < nbl->na_cj; j++)
1855                             {
1856                                 /* Is this interaction perturbed and not excluded? */
1857                                 ind_j = (gridj->cell0*gpu_ncluster_per_cell + cjr)*nbl->na_cj + j;
1858                                 aj    = nbs->a[ind_j];
1859                                 if (aj >= 0 &&
1860                                     (bFEP_i || (fep_cj & (1 << j))) &&
1861                                     (!bDiagRemoved || ind_j >= ind_i))
1862                                 {
1863                                     nbnxn_excl_t *excl;
1864                                     int           excl_pair;
1865                                     unsigned int  excl_bit;
1866                                     real          dx, dy, dz;
1867
1868                                     get_nbl_exclusions_1(nbl, cj4_ind, j>>2, &excl);
1869
1870                                     excl_pair = a_mod_wj(j)*nbl->na_ci + i;
1871                                     excl_bit  = (1U << (gcj*gpu_ncluster_per_cell + c));
1872
1873                                     dx = nbat->x[ind_j*nbat->xstride+XX] - xi;
1874                                     dy = nbat->x[ind_j*nbat->xstride+YY] - yi;
1875                                     dz = nbat->x[ind_j*nbat->xstride+ZZ] - zi;
1876
1877                                     /* The unpruned GPU list has more than 2/3
1878                                      * of the atom pairs beyond rlist. Using
1879                                      * this list will cause a lot of overhead
1880                                      * in the CPU FEP kernels, especially
1881                                      * relative to the fast GPU kernels.
1882                                      * So we prune the FEP list here.
1883                                      */
1884                                     if (dx*dx + dy*dy + dz*dz < rlist_fep2)
1885                                     {
1886                                         if (nlist->nrj - nlist->jindex[nri] >= max_nrj_fep)
1887                                         {
1888                                             fep_list_new_nri_copy(nlist);
1889                                             nri = nlist->nri;
1890                                         }
1891
1892                                         /* Add it to the FEP list */
1893                                         nlist->jjnr[nlist->nrj]     = aj;
1894                                         nlist->excl_fep[nlist->nrj] = (excl->pair[excl_pair] & excl_bit) ? 1 : 0;
1895                                         nlist->nrj++;
1896                                     }
1897
1898                                     /* Exclude it from the normal list.
1899                                      * Note that the charge and LJ parameters have
1900                                      * been set to zero, but we need to avoid 0/0,
1901                                      * as perturbed atoms can be on top of each other.
1902                                      */
1903                                     excl->pair[excl_pair] &= ~excl_bit;
1904                                 }
1905                             }
1906
1907                             /* Note that we could mask out this pair in imask
1908                              * if all i- and/or all j-particles are perturbed.
1909                              * But since the perturbed pairs on the CPU will
1910                              * take an order of magnitude more time, the GPU
1911                              * will finish before the CPU and there is no gain.
1912                              */
1913                         }
1914                     }
1915                 }
1916
1917                 if (nlist->nrj > nlist->jindex[nri])
1918                 {
1919                     /* Actually add this new, non-empty, list */
1920                     nlist->nri++;
1921                     nlist->jindex[nlist->nri] = nlist->nrj;
1922                 }
1923             }
1924         }
1925     }
1926 }
1927
1928 /* Set all atom-pair exclusions from the topology stored in excl
1929  * as masks in the pair-list for i-super-cell entry nbl_sci
1930  */
1931 static void set_sci_top_excls(const nbnxn_search_t nbs,
1932                               nbnxn_pairlist_t    *nbl,
1933                               gmx_bool             diagRemoved,
1934                               int                  na_c_2log,
1935                               const nbnxn_sci_t   *nbl_sci,
1936                               const t_blocka      *excl)
1937 {
1938     const int    *cell;
1939     int           na_c;
1940     int           sci;
1941     int           cj_ind_first, cj_ind_last;
1942     int           cj_first, cj_last;
1943     int           ndirect;
1944     int           ai, aj, si, ge, se;
1945     int           found, cj_ind_0, cj_ind_1, cj_ind_m;
1946     int           cj_m;
1947     nbnxn_excl_t *nbl_excl;
1948     int           inner_i, inner_e, w;
1949
1950     cell = nbs->cell;
1951
1952     na_c = nbl->na_ci;
1953
1954     if (nbl_sci->cj4_ind_end == nbl_sci->cj4_ind_start)
1955     {
1956         /* Empty list */
1957         return;
1958     }
1959
1960     sci = nbl_sci->sci;
1961
1962     cj_ind_first = nbl_sci->cj4_ind_start*nbnxn_gpu_jgroup_size;
1963     cj_ind_last  = nbl->work->cj_ind - 1;
1964
1965     cj_first = nbl->cj4[nbl_sci->cj4_ind_start].cj[0];
1966     cj_last  = nbl_cj(nbl, cj_ind_last);
1967
1968     /* Determine how many contiguous j-clusters we have starting
1969      * from the first i-cluster. This number can be used to directly
1970      * calculate j-cluster indices for excluded atoms.
1971      */
1972     ndirect = 0;
1973     while (cj_ind_first + ndirect <= cj_ind_last &&
1974            nbl_cj(nbl, cj_ind_first+ndirect) == sci*gpu_ncluster_per_cell + ndirect)
1975     {
1976         ndirect++;
1977     }
1978
1979     /* Loop over the atoms in the i super-cell */
1980     for (int i = 0; i < nbl->na_sc; i++)
1981     {
1982         ai = nbs->a[sci*nbl->na_sc+i];
1983         if (ai >= 0)
1984         {
1985             si  = (i>>na_c_2log);
1986
1987             /* Loop over the topology-based exclusions for this i-atom */
1988             for (int eind = excl->index[ai]; eind < excl->index[ai+1]; eind++)
1989             {
1990                 aj = excl->a[eind];
1991
1992                 if (aj == ai)
1993                 {
1994                     /* The self exclusion are already set, save some time */
1995                     continue;
1996                 }
1997
1998                 ge = cell[aj];
1999
2000                 /* Without shifts we only calculate interactions j>i
2001                  * for one-way pair-lists.
2002                  */
2003                 if (diagRemoved && ge <= sci*nbl->na_sc + i)
2004                 {
2005                     continue;
2006                 }
2007
2008                 se = ge>>na_c_2log;
2009                 /* Could the cluster se be in our list? */
2010                 if (se >= cj_first && se <= cj_last)
2011                 {
2012                     if (se < cj_first + ndirect)
2013                     {
2014                         /* We can calculate cj_ind directly from se */
2015                         found = cj_ind_first + se - cj_first;
2016                     }
2017                     else
2018                     {
2019                         /* Search for se using bisection */
2020                         found    = -1;
2021                         cj_ind_0 = cj_ind_first + ndirect;
2022                         cj_ind_1 = cj_ind_last + 1;
2023                         while (found == -1 && cj_ind_0 < cj_ind_1)
2024                         {
2025                             cj_ind_m = (cj_ind_0 + cj_ind_1)>>1;
2026
2027                             cj_m = nbl_cj(nbl, cj_ind_m);
2028
2029                             if (se == cj_m)
2030                             {
2031                                 found = cj_ind_m;
2032                             }
2033                             else if (se < cj_m)
2034                             {
2035                                 cj_ind_1 = cj_ind_m;
2036                             }
2037                             else
2038                             {
2039                                 cj_ind_0 = cj_ind_m + 1;
2040                             }
2041                         }
2042                     }
2043
2044                     if (found >= 0)
2045                     {
2046                         inner_i = i  - si*na_c;
2047                         inner_e = ge - se*na_c;
2048
2049                         if (nbl_imask0(nbl, found) & (1U << (cj_mod_cj4(found)*gpu_ncluster_per_cell + si)))
2050                         {
2051                             w       = (inner_e >> 2);
2052
2053                             get_nbl_exclusions_1(nbl, cj_to_cj4(found), w, &nbl_excl);
2054
2055                             nbl_excl->pair[a_mod_wj(inner_e)*nbl->na_ci+inner_i] &=
2056                                 ~(1U << (cj_mod_cj4(found)*gpu_ncluster_per_cell + si));
2057                         }
2058                     }
2059                 }
2060             }
2061         }
2062     }
2063 }
2064
2065 /* Reallocate the simple ci list for at least n entries */
2066 static void nb_realloc_ci(nbnxn_pairlist_t *nbl, int n)
2067 {
2068     nbl->ci_nalloc = over_alloc_small(n);
2069     nbnxn_realloc_void((void **)&nbl->ci,
2070                        nbl->nci*sizeof(*nbl->ci),
2071                        nbl->ci_nalloc*sizeof(*nbl->ci),
2072                        nbl->alloc, nbl->free);
2073 }
2074
2075 /* Reallocate the super-cell sci list for at least n entries */
2076 static void nb_realloc_sci(nbnxn_pairlist_t *nbl, int n)
2077 {
2078     nbl->sci_nalloc = over_alloc_small(n);
2079     nbnxn_realloc_void((void **)&nbl->sci,
2080                        nbl->nsci*sizeof(*nbl->sci),
2081                        nbl->sci_nalloc*sizeof(*nbl->sci),
2082                        nbl->alloc, nbl->free);
2083 }
2084
2085 /* Make a new ci entry at index nbl->nci */
2086 static void new_ci_entry(nbnxn_pairlist_t *nbl, int ci, int shift, int flags)
2087 {
2088     if (nbl->nci + 1 > nbl->ci_nalloc)
2089     {
2090         nb_realloc_ci(nbl, nbl->nci+1);
2091     }
2092     nbl->ci[nbl->nci].ci            = ci;
2093     nbl->ci[nbl->nci].shift         = shift;
2094     /* Store the interaction flags along with the shift */
2095     nbl->ci[nbl->nci].shift        |= flags;
2096     nbl->ci[nbl->nci].cj_ind_start  = nbl->ncj;
2097     nbl->ci[nbl->nci].cj_ind_end    = nbl->ncj;
2098 }
2099
2100 /* Make a new sci entry at index nbl->nsci */
2101 static void new_sci_entry(nbnxn_pairlist_t *nbl, int sci, int shift)
2102 {
2103     if (nbl->nsci + 1 > nbl->sci_nalloc)
2104     {
2105         nb_realloc_sci(nbl, nbl->nsci+1);
2106     }
2107     nbl->sci[nbl->nsci].sci           = sci;
2108     nbl->sci[nbl->nsci].shift         = shift;
2109     nbl->sci[nbl->nsci].cj4_ind_start = nbl->ncj4;
2110     nbl->sci[nbl->nsci].cj4_ind_end   = nbl->ncj4;
2111 }
2112
2113 /* Sort the simple j-list cj on exclusions.
2114  * Entries with exclusions will all be sorted to the beginning of the list.
2115  */
2116 static void sort_cj_excl(nbnxn_cj_t *cj, int ncj,
2117                          nbnxn_list_work_t *work)
2118 {
2119     int jnew;
2120
2121     if (ncj > work->cj_nalloc)
2122     {
2123         work->cj_nalloc = over_alloc_large(ncj);
2124         srenew(work->cj, work->cj_nalloc);
2125     }
2126
2127     /* Make a list of the j-cells involving exclusions */
2128     jnew = 0;
2129     for (int j = 0; j < ncj; j++)
2130     {
2131         if (cj[j].excl != NBNXN_INTERACTION_MASK_ALL)
2132         {
2133             work->cj[jnew++] = cj[j];
2134         }
2135     }
2136     /* Check if there are exclusions at all or not just the first entry */
2137     if (!((jnew == 0) ||
2138           (jnew == 1 && cj[0].excl != NBNXN_INTERACTION_MASK_ALL)))
2139     {
2140         for (int j = 0; j < ncj; j++)
2141         {
2142             if (cj[j].excl == NBNXN_INTERACTION_MASK_ALL)
2143             {
2144                 work->cj[jnew++] = cj[j];
2145             }
2146         }
2147         for (int j = 0; j < ncj; j++)
2148         {
2149             cj[j] = work->cj[j];
2150         }
2151     }
2152 }
2153
2154 /* Close this simple list i entry */
2155 static void close_ci_entry_simple(nbnxn_pairlist_t *nbl)
2156 {
2157     int jlen;
2158
2159     /* All content of the new ci entry have already been filled correctly,
2160      * we only need to increase the count here (for non empty lists).
2161      */
2162     jlen = nbl->ci[nbl->nci].cj_ind_end - nbl->ci[nbl->nci].cj_ind_start;
2163     if (jlen > 0)
2164     {
2165         sort_cj_excl(nbl->cj+nbl->ci[nbl->nci].cj_ind_start, jlen, nbl->work);
2166
2167         /* The counts below are used for non-bonded pair/flop counts
2168          * and should therefore match the available kernel setups.
2169          */
2170         if (!(nbl->ci[nbl->nci].shift & NBNXN_CI_DO_COUL(0)))
2171         {
2172             nbl->work->ncj_noq += jlen;
2173         }
2174         else if ((nbl->ci[nbl->nci].shift & NBNXN_CI_HALF_LJ(0)) ||
2175                  !(nbl->ci[nbl->nci].shift & NBNXN_CI_DO_LJ(0)))
2176         {
2177             nbl->work->ncj_hlj += jlen;
2178         }
2179
2180         nbl->nci++;
2181     }
2182 }
2183
2184 /* Split sci entry for load balancing on the GPU.
2185  * Splitting ensures we have enough lists to fully utilize the whole GPU.
2186  * With progBal we generate progressively smaller lists, which improves
2187  * load balancing. As we only know the current count on our own thread,
2188  * we will need to estimate the current total amount of i-entries.
2189  * As the lists get concatenated later, this estimate depends
2190  * both on nthread and our own thread index.
2191  */
2192 static void split_sci_entry(nbnxn_pairlist_t *nbl,
2193                             int nsp_target_av,
2194                             gmx_bool progBal, int nsp_tot_est,
2195                             int thread, int nthread)
2196 {
2197     int nsp_est;
2198     int nsp_max;
2199     int cj4_start, cj4_end, j4len;
2200     int sci;
2201     int nsp, nsp_sci, nsp_cj4, nsp_cj4_e, nsp_cj4_p;
2202
2203     if (progBal)
2204     {
2205         /* Estimate the total numbers of ci's of the nblist combined
2206          * over all threads using the target number of ci's.
2207          */
2208         nsp_est = (nsp_tot_est*thread)/nthread + nbl->nci_tot;
2209
2210         /* The first ci blocks should be larger, to avoid overhead.
2211          * The last ci blocks should be smaller, to improve load balancing.
2212          * The factor 3/2 makes the first block 3/2 times the target average
2213          * and ensures that the total number of blocks end up equal to
2214          * that with of equally sized blocks of size nsp_target_av.
2215          */
2216         nsp_max = nsp_target_av*nsp_tot_est*3/(2*(nsp_est + nsp_tot_est));
2217     }
2218     else
2219     {
2220         nsp_max = nsp_target_av;
2221     }
2222
2223     /* Since nsp_max is a maximum/cut-off (this avoids high outliers,
2224      * which lead to load imbalance), not an average, we add half the
2225      * number of pairs in a cj4 block to get the average about right.
2226      */
2227     nsp_max += gpu_ncluster_per_cell*nbnxn_gpu_jgroup_size/2;
2228
2229     cj4_start = nbl->sci[nbl->nsci-1].cj4_ind_start;
2230     cj4_end   = nbl->sci[nbl->nsci-1].cj4_ind_end;
2231     j4len     = cj4_end - cj4_start;
2232
2233     if (j4len > 1 && j4len*gpu_ncluster_per_cell*nbnxn_gpu_jgroup_size > nsp_max)
2234     {
2235         /* Remove the last ci entry and process the cj4's again */
2236         nbl->nsci -= 1;
2237
2238         sci        = nbl->nsci;
2239         nsp        = 0;
2240         nsp_sci    = 0;
2241         nsp_cj4_e  = 0;
2242         nsp_cj4    = 0;
2243         for (int cj4 = cj4_start; cj4 < cj4_end; cj4++)
2244         {
2245             nsp_cj4_p = nsp_cj4;
2246             /* Count the number of cluster pairs in this cj4 group */
2247             nsp_cj4   = 0;
2248             for (int p = 0; p < gpu_ncluster_per_cell*nbnxn_gpu_jgroup_size; p++)
2249             {
2250                 nsp_cj4 += (nbl->cj4[cj4].imei[0].imask >> p) & 1;
2251             }
2252
2253             /* Check if we should split at this cj4 to get a list of size nsp */
2254             if (nsp > 0 && nsp + nsp_cj4 > nsp_max)
2255             {
2256                 /* Split the list at cj4 */
2257                 nbl->sci[sci].cj4_ind_end = cj4;
2258                 /* Create a new sci entry */
2259                 sci++;
2260                 nbl->nsci++;
2261                 if (nbl->nsci+1 > nbl->sci_nalloc)
2262                 {
2263                     nb_realloc_sci(nbl, nbl->nsci+1);
2264                 }
2265                 nbl->sci[sci].sci           = nbl->sci[nbl->nsci-1].sci;
2266                 nbl->sci[sci].shift         = nbl->sci[nbl->nsci-1].shift;
2267                 nbl->sci[sci].cj4_ind_start = cj4;
2268                 nsp_sci                     = nsp;
2269                 nsp_cj4_e                   = nsp_cj4_p;
2270                 nsp                         = 0;
2271             }
2272             nsp += nsp_cj4;
2273         }
2274
2275         /* Put the remaining cj4's in the last sci entry */
2276         nbl->sci[sci].cj4_ind_end = cj4_end;
2277
2278         /* Possibly balance out the last two sci's
2279          * by moving the last cj4 of the second last sci.
2280          */
2281         if (nsp_sci - nsp_cj4_e >= nsp + nsp_cj4_e)
2282         {
2283             nbl->sci[sci-1].cj4_ind_end--;
2284             nbl->sci[sci].cj4_ind_start--;
2285         }
2286
2287         nbl->nsci++;
2288     }
2289 }
2290
2291 /* Clost this super/sub list i entry */
2292 static void close_ci_entry_supersub(nbnxn_pairlist_t *nbl,
2293                                     int nsp_max_av,
2294                                     gmx_bool progBal, int nsp_tot_est,
2295                                     int thread, int nthread)
2296 {
2297     /* All content of the new ci entry have already been filled correctly,
2298      * we only need to increase the count here (for non empty lists).
2299      */
2300     int j4len = nbl->sci[nbl->nsci].cj4_ind_end - nbl->sci[nbl->nsci].cj4_ind_start;
2301     if (j4len > 0)
2302     {
2303         /* We can only have complete blocks of 4 j-entries in a list,
2304          * so round the count up before closing.
2305          */
2306         nbl->ncj4         = (nbl->work->cj_ind + nbnxn_gpu_jgroup_size - 1)/nbnxn_gpu_jgroup_size;
2307         nbl->work->cj_ind = nbl->ncj4*nbnxn_gpu_jgroup_size;
2308
2309         nbl->nsci++;
2310
2311         if (nsp_max_av > 0)
2312         {
2313             /* Measure the size of the new entry and potentially split it */
2314             split_sci_entry(nbl, nsp_max_av, progBal, nsp_tot_est,
2315                             thread, nthread);
2316         }
2317     }
2318 }
2319
2320 /* Syncs the working array before adding another grid pair to the list */
2321 static void sync_work(nbnxn_pairlist_t *nbl)
2322 {
2323     if (!nbl->bSimple)
2324     {
2325         nbl->work->cj_ind   = nbl->ncj4*nbnxn_gpu_jgroup_size;
2326         nbl->work->cj4_init = nbl->ncj4;
2327     }
2328 }
2329
2330 /* Clears an nbnxn_pairlist_t data structure */
2331 static void clear_pairlist(nbnxn_pairlist_t *nbl)
2332 {
2333     nbl->nci           = 0;
2334     nbl->nsci          = 0;
2335     nbl->ncj           = 0;
2336     nbl->ncj4          = 0;
2337     nbl->nci_tot       = 0;
2338     nbl->nexcl         = 1;
2339
2340     nbl->work->ncj_noq = 0;
2341     nbl->work->ncj_hlj = 0;
2342 }
2343
2344 /* Clears a group scheme pair list */
2345 static void clear_pairlist_fep(t_nblist *nl)
2346 {
2347     nl->nri = 0;
2348     nl->nrj = 0;
2349     if (nl->jindex == NULL)
2350     {
2351         snew(nl->jindex, 1);
2352     }
2353     nl->jindex[0] = 0;
2354 }
2355
2356 /* Sets a simple list i-cell bounding box, including PBC shift */
2357 static gmx_inline void set_icell_bb_simple(const nbnxn_bb_t *bb, int ci,
2358                                            real shx, real shy, real shz,
2359                                            nbnxn_bb_t *bb_ci)
2360 {
2361     bb_ci->lower[BB_X] = bb[ci].lower[BB_X] + shx;
2362     bb_ci->lower[BB_Y] = bb[ci].lower[BB_Y] + shy;
2363     bb_ci->lower[BB_Z] = bb[ci].lower[BB_Z] + shz;
2364     bb_ci->upper[BB_X] = bb[ci].upper[BB_X] + shx;
2365     bb_ci->upper[BB_Y] = bb[ci].upper[BB_Y] + shy;
2366     bb_ci->upper[BB_Z] = bb[ci].upper[BB_Z] + shz;
2367 }
2368
2369 #ifdef NBNXN_BBXXXX
2370 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2371 static void set_icell_bbxxxx_supersub(const float *bb, int ci,
2372                                       real shx, real shy, real shz,
2373                                       float *bb_ci)
2374 {
2375     int ia = ci*(gpu_ncluster_per_cell >> STRIDE_PBB_2LOG)*NNBSBB_XXXX;
2376     for (int m = 0; m < (gpu_ncluster_per_cell >> STRIDE_PBB_2LOG)*NNBSBB_XXXX; m += NNBSBB_XXXX)
2377     {
2378         for (int i = 0; i < STRIDE_PBB; i++)
2379         {
2380             bb_ci[m+0*STRIDE_PBB+i] = bb[ia+m+0*STRIDE_PBB+i] + shx;
2381             bb_ci[m+1*STRIDE_PBB+i] = bb[ia+m+1*STRIDE_PBB+i] + shy;
2382             bb_ci[m+2*STRIDE_PBB+i] = bb[ia+m+2*STRIDE_PBB+i] + shz;
2383             bb_ci[m+3*STRIDE_PBB+i] = bb[ia+m+3*STRIDE_PBB+i] + shx;
2384             bb_ci[m+4*STRIDE_PBB+i] = bb[ia+m+4*STRIDE_PBB+i] + shy;
2385             bb_ci[m+5*STRIDE_PBB+i] = bb[ia+m+5*STRIDE_PBB+i] + shz;
2386         }
2387     }
2388 }
2389 #endif
2390
2391 /* Sets a super-cell and sub cell bounding boxes, including PBC shift */
2392 static void set_icell_bb_supersub(const nbnxn_bb_t *bb, int ci,
2393                                   real shx, real shy, real shz,
2394                                   nbnxn_bb_t *bb_ci)
2395 {
2396     for (int i = 0; i < gpu_ncluster_per_cell; i++)
2397     {
2398         set_icell_bb_simple(bb, ci*gpu_ncluster_per_cell+i,
2399                             shx, shy, shz,
2400                             &bb_ci[i]);
2401     }
2402 }
2403
2404 /* Copies PBC shifted i-cell atom coordinates x,y,z to working array */
2405 static void icell_set_x_simple(int ci,
2406                                real shx, real shy, real shz,
2407                                int stride, const real *x,
2408                                nbnxn_list_work_t *work)
2409 {
2410     int ia = ci*NBNXN_CPU_CLUSTER_I_SIZE;
2411
2412     for (int i = 0; i < NBNXN_CPU_CLUSTER_I_SIZE; i++)
2413     {
2414         work->x_ci[i*STRIDE_XYZ+XX] = x[(ia+i)*stride+XX] + shx;
2415         work->x_ci[i*STRIDE_XYZ+YY] = x[(ia+i)*stride+YY] + shy;
2416         work->x_ci[i*STRIDE_XYZ+ZZ] = x[(ia+i)*stride+ZZ] + shz;
2417     }
2418 }
2419
2420 /* Copies PBC shifted super-cell atom coordinates x,y,z to working array */
2421 static void icell_set_x_supersub(int ci,
2422                                  real shx, real shy, real shz,
2423                                  int stride, const real *x,
2424                                  nbnxn_list_work_t *work)
2425 {
2426 #ifndef NBNXN_SEARCH_BB_SIMD4
2427
2428     real * x_ci = work->x_ci;
2429
2430     int    ia = ci*gpu_ncluster_per_cell*nbnxn_gpu_cluster_size;
2431     for (int i = 0; i < gpu_ncluster_per_cell*nbnxn_gpu_cluster_size; i++)
2432     {
2433         x_ci[i*DIM + XX] = x[(ia+i)*stride + XX] + shx;
2434         x_ci[i*DIM + YY] = x[(ia+i)*stride + YY] + shy;
2435         x_ci[i*DIM + ZZ] = x[(ia+i)*stride + ZZ] + shz;
2436     }
2437
2438 #else /* !NBNXN_SEARCH_BB_SIMD4 */
2439
2440     real * x_ci = work->x_ci_simd;
2441
2442     for (int si = 0; si < gpu_ncluster_per_cell; si++)
2443     {
2444         for (int i = 0; i < nbnxn_gpu_cluster_size; i += GMX_SIMD4_WIDTH)
2445         {
2446             int io = si*nbnxn_gpu_cluster_size + i;
2447             int ia = ci*gpu_ncluster_per_cell*nbnxn_gpu_cluster_size + io;
2448             for (int j = 0; j < GMX_SIMD4_WIDTH; j++)
2449             {
2450                 x_ci[io*DIM + j + XX*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + XX] + shx;
2451                 x_ci[io*DIM + j + YY*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + YY] + shy;
2452                 x_ci[io*DIM + j + ZZ*GMX_SIMD4_WIDTH] = x[(ia + j)*stride + ZZ] + shz;
2453             }
2454         }
2455     }
2456
2457 #endif /* !NBNXN_SEARCH_BB_SIMD4 */
2458 }
2459
2460 static real minimum_subgrid_size_xy(const nbnxn_grid_t *grid)
2461 {
2462     if (grid->bSimple)
2463     {
2464         return std::min(grid->sx, grid->sy);
2465     }
2466     else
2467     {
2468         return std::min(grid->sx/gpu_ncluster_per_cell_x,
2469                         grid->sy/gpu_ncluster_per_cell_y);
2470     }
2471 }
2472
2473 static real effective_buffer_1x1_vs_MxN(const nbnxn_grid_t *gridi,
2474                                         const nbnxn_grid_t *gridj)
2475 {
2476     const real eff_1x1_buffer_fac_overest = 0.1;
2477
2478     /* Determine an atom-pair list cut-off buffer size for atom pairs,
2479      * to be added to rlist (including buffer) used for MxN.
2480      * This is for converting an MxN list to a 1x1 list. This means we can't
2481      * use the normal buffer estimate, as we have an MxN list in which
2482      * some atom pairs beyond rlist are missing. We want to capture
2483      * the beneficial effect of buffering by extra pairs just outside rlist,
2484      * while removing the useless pairs that are further away from rlist.
2485      * (Also the buffer could have been set manually not using the estimate.)
2486      * This buffer size is an overestimate.
2487      * We add 10% of the smallest grid sub-cell dimensions.
2488      * Note that the z-size differs per cell and we don't use this,
2489      * so we overestimate.
2490      * With PME, the 10% value gives a buffer that is somewhat larger
2491      * than the effective buffer with a tolerance of 0.005 kJ/mol/ps.
2492      * Smaller tolerances or using RF lead to a smaller effective buffer,
2493      * so 10% gives a safe overestimate.
2494      */
2495     return eff_1x1_buffer_fac_overest*(minimum_subgrid_size_xy(gridi) +
2496                                        minimum_subgrid_size_xy(gridj));
2497 }
2498
2499 /* Clusters at the cut-off only increase rlist by 60% of their size */
2500 static real nbnxn_rlist_inc_outside_fac = 0.6;
2501
2502 /* Due to the cluster size the effective pair-list is longer than
2503  * that of a simple atom pair-list. This function gives the extra distance.
2504  */
2505 real nbnxn_get_rlist_effective_inc(int cluster_size_j, real atom_density)
2506 {
2507     int  cluster_size_i;
2508     real vol_inc_i, vol_inc_j;
2509
2510     /* We should get this from the setup, but currently it's the same for
2511      * all setups, including GPUs.
2512      */
2513     cluster_size_i = NBNXN_CPU_CLUSTER_I_SIZE;
2514
2515     vol_inc_i = (cluster_size_i - 1)/atom_density;
2516     vol_inc_j = (cluster_size_j - 1)/atom_density;
2517
2518     return nbnxn_rlist_inc_outside_fac*std::cbrt(vol_inc_i + vol_inc_j);
2519 }
2520
2521 /* Estimates the interaction volume^2 for non-local interactions */
2522 static real nonlocal_vol2(const struct gmx_domdec_zones_t *zones, rvec ls, real r)
2523 {
2524     real cl, ca, za;
2525     real vold_est;
2526     real vol2_est_tot;
2527
2528     vol2_est_tot = 0;
2529
2530     /* Here we simply add up the volumes of 1, 2 or 3 1D decomposition
2531      * not home interaction volume^2. As these volumes are not additive,
2532      * this is an overestimate, but it would only be significant in the limit
2533      * of small cells, where we anyhow need to split the lists into
2534      * as small parts as possible.
2535      */
2536
2537     for (int z = 0; z < zones->n; z++)
2538     {
2539         if (zones->shift[z][XX] + zones->shift[z][YY] + zones->shift[z][ZZ] == 1)
2540         {
2541             cl = 0;
2542             ca = 1;
2543             za = 1;
2544             for (int d = 0; d < DIM; d++)
2545             {
2546                 if (zones->shift[z][d] == 0)
2547                 {
2548                     cl += 0.5*ls[d];
2549                     ca *= ls[d];
2550                     za *= zones->size[z].x1[d] - zones->size[z].x0[d];
2551                 }
2552             }
2553
2554             /* 4 octants of a sphere */
2555             vold_est  = 0.25*M_PI*r*r*r*r;
2556             /* 4 quarter pie slices on the edges */
2557             vold_est += 4*cl*M_PI/6.0*r*r*r;
2558             /* One rectangular volume on a face */
2559             vold_est += ca*0.5*r*r;
2560
2561             vol2_est_tot += vold_est*za;
2562         }
2563     }
2564
2565     return vol2_est_tot;
2566 }
2567
2568 /* Estimates the average size of a full j-list for super/sub setup */
2569 static void get_nsubpair_target(const nbnxn_search_t  nbs,
2570                                 int                   iloc,
2571                                 real                  rlist,
2572                                 int                   min_ci_balanced,
2573                                 int                  *nsubpair_target,
2574                                 int                  *nsubpair_tot_est)
2575 {
2576     /* The target value of 36 seems to be the optimum for Kepler.
2577      * Maxwell is less sensitive to the exact value.
2578      */
2579     const int           nsubpair_target_min = 36;
2580     const nbnxn_grid_t *grid;
2581     rvec                ls;
2582     real                xy_diag2, r_eff_sup, vol_est, nsp_est, nsp_est_nl;
2583
2584     grid = &nbs->grid[0];
2585
2586     if (min_ci_balanced <= 0 || grid->nc >= min_ci_balanced || grid->nc == 0)
2587     {
2588         /* We don't need to balance the list sizes */
2589         *nsubpair_target  = 0;
2590         *nsubpair_tot_est = 0;
2591
2592         return;
2593     }
2594
2595     ls[XX] = (grid->c1[XX] - grid->c0[XX])/(grid->ncx*gpu_ncluster_per_cell_x);
2596     ls[YY] = (grid->c1[YY] - grid->c0[YY])/(grid->ncy*gpu_ncluster_per_cell_y);
2597     ls[ZZ] = grid->na_c/(grid->atom_density*ls[XX]*ls[YY]);
2598
2599     /* The average squared length of the diagonal of a sub cell */
2600     xy_diag2 = ls[XX]*ls[XX] + ls[YY]*ls[YY] + ls[ZZ]*ls[ZZ];
2601
2602     /* The formulas below are a heuristic estimate of the average nsj per si*/
2603     r_eff_sup = rlist + nbnxn_rlist_inc_outside_fac*gmx::square((grid->na_c - 1.0)/grid->na_c)*std::sqrt(xy_diag2/3);
2604
2605     if (!nbs->DomDec || nbs->zones->n == 1)
2606     {
2607         nsp_est_nl = 0;
2608     }
2609     else
2610     {
2611         nsp_est_nl =
2612             gmx::square(grid->atom_density/grid->na_c)*
2613             nonlocal_vol2(nbs->zones, ls, r_eff_sup);
2614     }
2615
2616     if (LOCAL_I(iloc))
2617     {
2618         /* Sub-cell interacts with itself */
2619         vol_est  = ls[XX]*ls[YY]*ls[ZZ];
2620         /* 6/2 rectangular volume on the faces */
2621         vol_est += (ls[XX]*ls[YY] + ls[XX]*ls[ZZ] + ls[YY]*ls[ZZ])*r_eff_sup;
2622         /* 12/2 quarter pie slices on the edges */
2623         vol_est += 2*(ls[XX] + ls[YY] + ls[ZZ])*0.25*M_PI*gmx::square(r_eff_sup);
2624         /* 4 octants of a sphere */
2625         vol_est += 0.5*4.0/3.0*M_PI*gmx::power3(r_eff_sup);
2626
2627         /* Estimate the number of cluster pairs as the local number of
2628          * clusters times the volume they interact with times the density.
2629          */
2630         nsp_est = grid->nsubc_tot*vol_est*grid->atom_density/grid->na_c;
2631
2632         /* Subtract the non-local pair count */
2633         nsp_est -= nsp_est_nl;
2634
2635         /* For small cut-offs nsp_est will be an underesimate.
2636          * With DD nsp_est_nl is an overestimate so nsp_est can get negative.
2637          * So to avoid too small or negative nsp_est we set a minimum of
2638          * all cells interacting with all 3^3 direct neighbors (3^3-1)/2+1=14.
2639          * This might be a slight overestimate for small non-periodic groups of
2640          * atoms as will occur for a local domain with DD, but for small
2641          * groups of atoms we'll anyhow be limited by nsubpair_target_min,
2642          * so this overestimation will not matter.
2643          */
2644         nsp_est = std::max(nsp_est, grid->nsubc_tot*static_cast<real>(14));
2645
2646         if (debug)
2647         {
2648             fprintf(debug, "nsp_est local %5.1f non-local %5.1f\n",
2649                     nsp_est, nsp_est_nl);
2650         }
2651     }
2652     else
2653     {
2654         nsp_est = nsp_est_nl;
2655     }
2656
2657     /* Thus the (average) maximum j-list size should be as follows.
2658      * Since there is overhead, we shouldn't make the lists too small
2659      * (and we can't chop up j-groups) so we use a minimum target size of 36.
2660      */
2661     *nsubpair_target  = std::max(nsubpair_target_min,
2662                                  static_cast<int>(nsp_est/min_ci_balanced + 0.5));
2663     *nsubpair_tot_est = static_cast<int>(nsp_est);
2664
2665     if (debug)
2666     {
2667         fprintf(debug, "nbl nsp estimate %.1f, nsubpair_target %d\n",
2668                 nsp_est, *nsubpair_target);
2669     }
2670 }
2671
2672 /* Debug list print function */
2673 static void print_nblist_ci_cj(FILE *fp, const nbnxn_pairlist_t *nbl)
2674 {
2675     for (int i = 0; i < nbl->nci; i++)
2676     {
2677         fprintf(fp, "ci %4d  shift %2d  ncj %3d\n",
2678                 nbl->ci[i].ci, nbl->ci[i].shift,
2679                 nbl->ci[i].cj_ind_end - nbl->ci[i].cj_ind_start);
2680
2681         for (int j = nbl->ci[i].cj_ind_start; j < nbl->ci[i].cj_ind_end; j++)
2682         {
2683             fprintf(fp, "  cj %5d  imask %x\n",
2684                     nbl->cj[j].cj,
2685                     nbl->cj[j].excl);
2686         }
2687     }
2688 }
2689
2690 /* Debug list print function */
2691 static void print_nblist_sci_cj(FILE *fp, const nbnxn_pairlist_t *nbl)
2692 {
2693     for (int i = 0; i < nbl->nsci; i++)
2694     {
2695         fprintf(fp, "ci %4d  shift %2d  ncj4 %2d\n",
2696                 nbl->sci[i].sci, nbl->sci[i].shift,
2697                 nbl->sci[i].cj4_ind_end - nbl->sci[i].cj4_ind_start);
2698
2699         int ncp = 0;
2700         for (int j4 = nbl->sci[i].cj4_ind_start; j4 < nbl->sci[i].cj4_ind_end; j4++)
2701         {
2702             for (int j = 0; j < nbnxn_gpu_jgroup_size; j++)
2703             {
2704                 fprintf(fp, "  sj %5d  imask %x\n",
2705                         nbl->cj4[j4].cj[j],
2706                         nbl->cj4[j4].imei[0].imask);
2707                 for (int si = 0; si < gpu_ncluster_per_cell; si++)
2708                 {
2709                     if (nbl->cj4[j4].imei[0].imask & (1U << (j*gpu_ncluster_per_cell + si)))
2710                     {
2711                         ncp++;
2712                     }
2713                 }
2714             }
2715         }
2716         fprintf(fp, "ci %4d  shift %2d  ncj4 %2d ncp %3d\n",
2717                 nbl->sci[i].sci, nbl->sci[i].shift,
2718                 nbl->sci[i].cj4_ind_end - nbl->sci[i].cj4_ind_start,
2719                 ncp);
2720     }
2721 }
2722
2723 /* Combine pair lists *nbl generated on multiple threads nblc */
2724 static void combine_nblists(int nnbl, nbnxn_pairlist_t **nbl,
2725                             nbnxn_pairlist_t *nblc)
2726 {
2727     int nsci, ncj4, nexcl;
2728
2729     if (nblc->bSimple)
2730     {
2731         gmx_incons("combine_nblists does not support simple lists");
2732     }
2733
2734     nsci  = nblc->nsci;
2735     ncj4  = nblc->ncj4;
2736     nexcl = nblc->nexcl;
2737     for (int i = 0; i < nnbl; i++)
2738     {
2739         nsci  += nbl[i]->nsci;
2740         ncj4  += nbl[i]->ncj4;
2741         nexcl += nbl[i]->nexcl;
2742     }
2743
2744     if (nsci > nblc->sci_nalloc)
2745     {
2746         nb_realloc_sci(nblc, nsci);
2747     }
2748     if (ncj4 > nblc->cj4_nalloc)
2749     {
2750         nblc->cj4_nalloc = over_alloc_small(ncj4);
2751         nbnxn_realloc_void((void **)&nblc->cj4,
2752                            nblc->ncj4*sizeof(*nblc->cj4),
2753                            nblc->cj4_nalloc*sizeof(*nblc->cj4),
2754                            nblc->alloc, nblc->free);
2755     }
2756     if (nexcl > nblc->excl_nalloc)
2757     {
2758         nblc->excl_nalloc = over_alloc_small(nexcl);
2759         nbnxn_realloc_void((void **)&nblc->excl,
2760                            nblc->nexcl*sizeof(*nblc->excl),
2761                            nblc->excl_nalloc*sizeof(*nblc->excl),
2762                            nblc->alloc, nblc->free);
2763     }
2764
2765     /* Each thread should copy its own data to the combined arrays,
2766      * as otherwise data will go back and forth between different caches.
2767      */
2768 #if GMX_OPENMP && !(defined __clang_analyzer__)
2769     // cppcheck-suppress unreadVariable
2770     int nthreads = gmx_omp_nthreads_get(emntPairsearch);
2771 #endif
2772
2773 #pragma omp parallel for num_threads(nthreads) schedule(static)
2774     for (int n = 0; n < nnbl; n++)
2775     {
2776         try
2777         {
2778             int                     sci_offset;
2779             int                     cj4_offset;
2780             int                     excl_offset;
2781             const nbnxn_pairlist_t *nbli;
2782
2783             /* Determine the offset in the combined data for our thread */
2784             sci_offset  = nblc->nsci;
2785             cj4_offset  = nblc->ncj4;
2786             excl_offset = nblc->nexcl;
2787
2788             for (int i = 0; i < n; i++)
2789             {
2790                 sci_offset  += nbl[i]->nsci;
2791                 cj4_offset  += nbl[i]->ncj4;
2792                 excl_offset += nbl[i]->nexcl;
2793             }
2794
2795             nbli = nbl[n];
2796
2797             for (int i = 0; i < nbli->nsci; i++)
2798             {
2799                 nblc->sci[sci_offset+i]                = nbli->sci[i];
2800                 nblc->sci[sci_offset+i].cj4_ind_start += cj4_offset;
2801                 nblc->sci[sci_offset+i].cj4_ind_end   += cj4_offset;
2802             }
2803
2804             for (int j4 = 0; j4 < nbli->ncj4; j4++)
2805             {
2806                 nblc->cj4[cj4_offset+j4]                   = nbli->cj4[j4];
2807                 nblc->cj4[cj4_offset+j4].imei[0].excl_ind += excl_offset;
2808                 nblc->cj4[cj4_offset+j4].imei[1].excl_ind += excl_offset;
2809             }
2810
2811             for (int j4 = 0; j4 < nbli->nexcl; j4++)
2812             {
2813                 nblc->excl[excl_offset+j4] = nbli->excl[j4];
2814             }
2815         }
2816         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2817     }
2818
2819     for (int n = 0; n < nnbl; n++)
2820     {
2821         nblc->nsci    += nbl[n]->nsci;
2822         nblc->ncj4    += nbl[n]->ncj4;
2823         nblc->nci_tot += nbl[n]->nci_tot;
2824         nblc->nexcl   += nbl[n]->nexcl;
2825     }
2826 }
2827
2828 static void balance_fep_lists(const nbnxn_search_t  nbs,
2829                               nbnxn_pairlist_set_t *nbl_lists)
2830 {
2831     int       nnbl;
2832     int       nri_tot, nrj_tot, nrj_target;
2833     int       th_dest;
2834     t_nblist *nbld;
2835
2836     nnbl = nbl_lists->nnbl;
2837
2838     if (nnbl == 1)
2839     {
2840         /* Nothing to balance */
2841         return;
2842     }
2843
2844     /* Count the total i-lists and pairs */
2845     nri_tot = 0;
2846     nrj_tot = 0;
2847     for (int th = 0; th < nnbl; th++)
2848     {
2849         nri_tot += nbl_lists->nbl_fep[th]->nri;
2850         nrj_tot += nbl_lists->nbl_fep[th]->nrj;
2851     }
2852
2853     nrj_target = (nrj_tot + nnbl - 1)/nnbl;
2854
2855     assert(gmx_omp_nthreads_get(emntNonbonded) == nnbl);
2856
2857 #pragma omp parallel for schedule(static) num_threads(nnbl)
2858     for (int th = 0; th < nnbl; th++)
2859     {
2860         try
2861         {
2862             t_nblist *nbl;
2863
2864             nbl = nbs->work[th].nbl_fep;
2865
2866             /* Note that here we allocate for the total size, instead of
2867              * a per-thread esimate (which is hard to obtain).
2868              */
2869             if (nri_tot > nbl->maxnri)
2870             {
2871                 nbl->maxnri = over_alloc_large(nri_tot);
2872                 reallocate_nblist(nbl);
2873             }
2874             if (nri_tot > nbl->maxnri || nrj_tot > nbl->maxnrj)
2875             {
2876                 nbl->maxnrj = over_alloc_small(nrj_tot);
2877                 srenew(nbl->jjnr, nbl->maxnrj);
2878                 srenew(nbl->excl_fep, nbl->maxnrj);
2879             }
2880
2881             clear_pairlist_fep(nbl);
2882         }
2883         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2884     }
2885
2886     /* Loop over the source lists and assign and copy i-entries */
2887     th_dest = 0;
2888     nbld    = nbs->work[th_dest].nbl_fep;
2889     for (int th = 0; th < nnbl; th++)
2890     {
2891         t_nblist *nbls;
2892
2893         nbls = nbl_lists->nbl_fep[th];
2894
2895         for (int i = 0; i < nbls->nri; i++)
2896         {
2897             int nrj;
2898
2899             /* The number of pairs in this i-entry */
2900             nrj = nbls->jindex[i+1] - nbls->jindex[i];
2901
2902             /* Decide if list th_dest is too large and we should procede
2903              * to the next destination list.
2904              */
2905             if (th_dest+1 < nnbl && nbld->nrj > 0 &&
2906                 nbld->nrj + nrj - nrj_target > nrj_target - nbld->nrj)
2907             {
2908                 th_dest++;
2909                 nbld = nbs->work[th_dest].nbl_fep;
2910             }
2911
2912             nbld->iinr[nbld->nri]  = nbls->iinr[i];
2913             nbld->gid[nbld->nri]   = nbls->gid[i];
2914             nbld->shift[nbld->nri] = nbls->shift[i];
2915
2916             for (int j = nbls->jindex[i]; j < nbls->jindex[i+1]; j++)
2917             {
2918                 nbld->jjnr[nbld->nrj]     = nbls->jjnr[j];
2919                 nbld->excl_fep[nbld->nrj] = nbls->excl_fep[j];
2920                 nbld->nrj++;
2921             }
2922             nbld->nri++;
2923             nbld->jindex[nbld->nri] = nbld->nrj;
2924         }
2925     }
2926
2927     /* Swap the list pointers */
2928     for (int th = 0; th < nnbl; th++)
2929     {
2930         t_nblist *nbl_tmp;
2931
2932         nbl_tmp                = nbl_lists->nbl_fep[th];
2933         nbl_lists->nbl_fep[th] = nbs->work[th].nbl_fep;
2934         nbs->work[th].nbl_fep  = nbl_tmp;
2935
2936         if (debug)
2937         {
2938             fprintf(debug, "nbl_fep[%d] nri %4d nrj %4d\n",
2939                     th,
2940                     nbl_lists->nbl_fep[th]->nri,
2941                     nbl_lists->nbl_fep[th]->nrj);
2942         }
2943     }
2944 }
2945
2946 /* Returns the next ci to be processes by our thread */
2947 static gmx_bool next_ci(const nbnxn_grid_t *grid,
2948                         int conv,
2949                         int nth, int ci_block,
2950                         int *ci_x, int *ci_y,
2951                         int *ci_b, int *ci)
2952 {
2953     (*ci_b)++;
2954     (*ci)++;
2955
2956     if (*ci_b == ci_block)
2957     {
2958         /* Jump to the next block assigned to this task */
2959         *ci   += (nth - 1)*ci_block;
2960         *ci_b  = 0;
2961     }
2962
2963     if (*ci >= grid->nc*conv)
2964     {
2965         return FALSE;
2966     }
2967
2968     while (*ci >= grid->cxy_ind[*ci_x*grid->ncy + *ci_y + 1]*conv)
2969     {
2970         *ci_y += 1;
2971         if (*ci_y == grid->ncy)
2972         {
2973             *ci_x += 1;
2974             *ci_y  = 0;
2975         }
2976     }
2977
2978     return TRUE;
2979 }
2980
2981 /* Returns the distance^2 for which we put cell pairs in the list
2982  * without checking atom pair distances. This is usually < rlist^2.
2983  */
2984 static float boundingbox_only_distance2(const nbnxn_grid_t *gridi,
2985                                         const nbnxn_grid_t *gridj,
2986                                         real                rlist,
2987                                         gmx_bool            simple)
2988 {
2989     /* If the distance between two sub-cell bounding boxes is less
2990      * than this distance, do not check the distance between
2991      * all particle pairs in the sub-cell, since then it is likely
2992      * that the box pair has atom pairs within the cut-off.
2993      * We use the nblist cut-off minus 0.5 times the average x/y diagonal
2994      * spacing of the sub-cells. Around 40% of the checked pairs are pruned.
2995      * Using more than 0.5 gains at most 0.5%.
2996      * If forces are calculated more than twice, the performance gain
2997      * in the force calculation outweighs the cost of checking.
2998      * Note that with subcell lists, the atom-pair distance check
2999      * is only performed when only 1 out of 8 sub-cells in within range,
3000      * this is because the GPU is much faster than the cpu.
3001      */
3002     real bbx, bby;
3003     real rbb2;
3004
3005     bbx = 0.5*(gridi->sx + gridj->sx);
3006     bby = 0.5*(gridi->sy + gridj->sy);
3007     if (!simple)
3008     {
3009         bbx /= gpu_ncluster_per_cell_x;
3010         bby /= gpu_ncluster_per_cell_y;
3011     }
3012
3013     rbb2 = std::max(0.0, rlist - 0.5*std::sqrt(bbx*bbx + bby*bby));
3014     rbb2 = rbb2 * rbb2;
3015
3016 #if !GMX_DOUBLE
3017     return rbb2;
3018 #else
3019     return (float)((1+GMX_FLOAT_EPS)*rbb2);
3020 #endif
3021 }
3022
3023 static int get_ci_block_size(const nbnxn_grid_t *gridi,
3024                              gmx_bool bDomDec, int nth)
3025 {
3026     const int ci_block_enum      = 5;
3027     const int ci_block_denom     = 11;
3028     const int ci_block_min_atoms = 16;
3029     int       ci_block;
3030
3031     /* Here we decide how to distribute the blocks over the threads.
3032      * We use prime numbers to try to avoid that the grid size becomes
3033      * a multiple of the number of threads, which would lead to some
3034      * threads getting "inner" pairs and others getting boundary pairs,
3035      * which in turns will lead to load imbalance between threads.
3036      * Set the block size as 5/11/ntask times the average number of cells
3037      * in a y,z slab. This should ensure a quite uniform distribution
3038      * of the grid parts of the different thread along all three grid
3039      * zone boundaries with 3D domain decomposition. At the same time
3040      * the blocks will not become too small.
3041      */
3042     ci_block = (gridi->nc*ci_block_enum)/(ci_block_denom*gridi->ncx*nth);
3043
3044     /* Ensure the blocks are not too small: avoids cache invalidation */
3045     if (ci_block*gridi->na_sc < ci_block_min_atoms)
3046     {
3047         ci_block = (ci_block_min_atoms + gridi->na_sc - 1)/gridi->na_sc;
3048     }
3049
3050     /* Without domain decomposition
3051      * or with less than 3 blocks per task, divide in nth blocks.
3052      */
3053     if (!bDomDec || nth*3*ci_block > gridi->nc)
3054     {
3055         ci_block = (gridi->nc + nth - 1)/nth;
3056     }
3057
3058     if (ci_block > 1 && (nth - 1)*ci_block >= gridi->nc)
3059     {
3060         /* Some threads have no work. Although reducing the block size
3061          * does not decrease the block count on the first few threads,
3062          * with GPUs better mixing of "upper" cells that have more empty
3063          * clusters results in a somewhat lower max load over all threads.
3064          * Without GPUs the regime of so few atoms per thread is less
3065          * performance relevant, but with 8-wide SIMD the same reasoning
3066          * applies, since the pair list uses 4 i-atom "sub-clusters".
3067          */
3068         ci_block--;
3069     }
3070
3071     return ci_block;
3072 }
3073
3074 /* Generates the part of pair-list nbl assigned to our thread */
3075 static void nbnxn_make_pairlist_part(const nbnxn_search_t nbs,
3076                                      const nbnxn_grid_t *gridi,
3077                                      const nbnxn_grid_t *gridj,
3078                                      nbnxn_search_work_t *work,
3079                                      const nbnxn_atomdata_t *nbat,
3080                                      const t_blocka *excl,
3081                                      real rlist,
3082                                      int nb_kernel_type,
3083                                      int ci_block,
3084                                      gmx_bool bFBufferFlag,
3085                                      int nsubpair_max,
3086                                      gmx_bool progBal,
3087                                      int nsubpair_tot_est,
3088                                      int th, int nth,
3089                                      nbnxn_pairlist_t *nbl,
3090                                      t_nblist *nbl_fep)
3091 {
3092     int               na_cj_2log;
3093     matrix            box;
3094     real              rl2, rl_fep2 = 0;
3095     float             rbb2;
3096     int               ci_b, ci, ci_x, ci_y, ci_xy, cj;
3097     ivec              shp;
3098     int               shift;
3099     real              shx, shy, shz;
3100     int               conv_i, cell0_i;
3101     const nbnxn_bb_t *bb_i = NULL;
3102 #ifdef NBNXN_BBXXXX
3103     const float      *pbb_i = NULL;
3104 #endif
3105     const float      *bbcz_i, *bbcz_j;
3106     const int        *flags_i;
3107     real              bx0, bx1, by0, by1, bz0, bz1;
3108     real              bz1_frac;
3109     real              d2cx, d2z, d2z_cx, d2z_cy, d2zx, d2zxy, d2xy;
3110     int               cxf, cxl, cyf, cyf_x, cyl;
3111     int               c0, c1, cs, cf, cl;
3112     int               ndistc;
3113     int               ncpcheck;
3114     int               gridi_flag_shift = 0, gridj_flag_shift = 0;
3115     gmx_bitmask_t    *gridj_flag       = NULL;
3116     int               ncj_old_i, ncj_old_j;
3117
3118     nbs_cycle_start(&work->cc[enbsCCsearch]);
3119
3120     if (gridj->bSimple != nbl->bSimple)
3121     {
3122         gmx_incons("Grid incompatible with pair-list");
3123     }
3124
3125     sync_work(nbl);
3126     nbl->na_sc = gridj->na_sc;
3127     nbl->na_ci = gridj->na_c;
3128     nbl->na_cj = nbnxn_kernel_to_cluster_j_size(nb_kernel_type);
3129     na_cj_2log = get_2log(nbl->na_cj);
3130
3131     nbl->rlist  = rlist;
3132
3133     if (bFBufferFlag)
3134     {
3135         /* Determine conversion of clusters to flag blocks */
3136         gridi_flag_shift = 0;
3137         while ((nbl->na_ci<<gridi_flag_shift) < NBNXN_BUFFERFLAG_SIZE)
3138         {
3139             gridi_flag_shift++;
3140         }
3141         gridj_flag_shift = 0;
3142         while ((nbl->na_cj<<gridj_flag_shift) < NBNXN_BUFFERFLAG_SIZE)
3143         {
3144             gridj_flag_shift++;
3145         }
3146
3147         gridj_flag = work->buffer_flags.flag;
3148     }
3149
3150     copy_mat(nbs->box, box);
3151
3152     rl2 = nbl->rlist*nbl->rlist;
3153
3154     if (nbs->bFEP && !nbl->bSimple)
3155     {
3156         /* Determine an atom-pair list cut-off distance for FEP atom pairs.
3157          * We should not simply use rlist, since then we would not have
3158          * the small, effective buffering of the NxN lists.
3159          * The buffer is on overestimate, but the resulting cost for pairs
3160          * beyond rlist is neglible compared to the FEP pairs within rlist.
3161          */
3162         rl_fep2 = nbl->rlist + effective_buffer_1x1_vs_MxN(gridi, gridj);
3163
3164         if (debug)
3165         {
3166             fprintf(debug, "nbl_fep atom-pair rlist %f\n", rl_fep2);
3167         }
3168         rl_fep2 = rl_fep2*rl_fep2;
3169     }
3170
3171     rbb2 = boundingbox_only_distance2(gridi, gridj, nbl->rlist, nbl->bSimple);
3172
3173     if (debug)
3174     {
3175         fprintf(debug, "nbl bounding box only distance %f\n", std::sqrt(rbb2));
3176     }
3177
3178     /* Set the shift range */
3179     for (int d = 0; d < DIM; d++)
3180     {
3181         /* Check if we need periodicity shifts.
3182          * Without PBC or with domain decomposition we don't need them.
3183          */
3184         if (d >= ePBC2npbcdim(nbs->ePBC) || nbs->dd_dim[d])
3185         {
3186             shp[d] = 0;
3187         }
3188         else
3189         {
3190             if (d == XX &&
3191                 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < std::sqrt(rl2))
3192             {
3193                 shp[d] = 2;
3194             }
3195             else
3196             {
3197                 shp[d] = 1;
3198             }
3199         }
3200     }
3201
3202     if (nbl->bSimple && !gridi->bSimple)
3203     {
3204         conv_i  = gridi->na_sc/gridj->na_sc;
3205         bb_i    = gridi->bb_simple;
3206         bbcz_i  = gridi->bbcz_simple;
3207         flags_i = gridi->flags_simple;
3208     }
3209     else
3210     {
3211         conv_i  = 1;
3212 #ifdef NBNXN_BBXXXX
3213         if (gridi->bSimple)
3214         {
3215             bb_i  = gridi->bb;
3216         }
3217         else
3218         {
3219             pbb_i = gridi->pbb;
3220         }
3221 #else
3222         /* We use the normal bounding box format for both grid types */
3223         bb_i  = gridi->bb;
3224 #endif
3225         bbcz_i  = gridi->bbcz;
3226         flags_i = gridi->flags;
3227     }
3228     cell0_i = gridi->cell0*conv_i;
3229
3230     bbcz_j = gridj->bbcz;
3231
3232     if (conv_i != 1)
3233     {
3234         /* Blocks of the conversion factor - 1 give a large repeat count
3235          * combined with a small block size. This should result in good
3236          * load balancing for both small and large domains.
3237          */
3238         ci_block = conv_i - 1;
3239     }
3240     if (debug)
3241     {
3242         fprintf(debug, "nbl nc_i %d col.av. %.1f ci_block %d\n",
3243                 gridi->nc, gridi->nc/(double)(gridi->ncx*gridi->ncy), ci_block);
3244     }
3245
3246     ndistc   = 0;
3247     ncpcheck = 0;
3248
3249     /* Initially ci_b and ci to 1 before where we want them to start,
3250      * as they will both be incremented in next_ci.
3251      */
3252     ci_b = -1;
3253     ci   = th*ci_block - 1;
3254     ci_x = 0;
3255     ci_y = 0;
3256     while (next_ci(gridi, conv_i, nth, ci_block, &ci_x, &ci_y, &ci_b, &ci))
3257     {
3258         if (nbl->bSimple && flags_i[ci] == 0)
3259         {
3260             continue;
3261         }
3262
3263         ncj_old_i = nbl->ncj;
3264
3265         d2cx = 0;
3266         if (gridj != gridi && shp[XX] == 0)
3267         {
3268             if (nbl->bSimple)
3269             {
3270                 bx1 = bb_i[ci].upper[BB_X];
3271             }
3272             else
3273             {
3274                 bx1 = gridi->c0[XX] + (ci_x+1)*gridi->sx;
3275             }
3276             if (bx1 < gridj->c0[XX])
3277             {
3278                 d2cx = gmx::square(gridj->c0[XX] - bx1);
3279
3280                 if (d2cx >= rl2)
3281                 {
3282                     continue;
3283                 }
3284             }
3285         }
3286
3287         ci_xy = ci_x*gridi->ncy + ci_y;
3288
3289         /* Loop over shift vectors in three dimensions */
3290         for (int tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
3291         {
3292             shz = tz*box[ZZ][ZZ];
3293
3294             bz0 = bbcz_i[ci*NNBSBB_D  ] + shz;
3295             bz1 = bbcz_i[ci*NNBSBB_D+1] + shz;
3296
3297             if (tz == 0)
3298             {
3299                 d2z = 0;
3300             }
3301             else if (tz < 0)
3302             {
3303                 d2z = gmx::square(bz1);
3304             }
3305             else
3306             {
3307                 d2z = gmx::square(bz0 - box[ZZ][ZZ]);
3308             }
3309
3310             d2z_cx = d2z + d2cx;
3311
3312             if (d2z_cx >= rl2)
3313             {
3314                 continue;
3315             }
3316
3317             bz1_frac = bz1/(gridi->cxy_ind[ci_xy+1] - gridi->cxy_ind[ci_xy]);
3318             if (bz1_frac < 0)
3319             {
3320                 bz1_frac = 0;
3321             }
3322             /* The check with bz1_frac close to or larger than 1 comes later */
3323
3324             for (int ty = -shp[YY]; ty <= shp[YY]; ty++)
3325             {
3326                 shy = ty*box[YY][YY] + tz*box[ZZ][YY];
3327
3328                 if (nbl->bSimple)
3329                 {
3330                     by0 = bb_i[ci].lower[BB_Y] + shy;
3331                     by1 = bb_i[ci].upper[BB_Y] + shy;
3332                 }
3333                 else
3334                 {
3335                     by0 = gridi->c0[YY] + (ci_y  )*gridi->sy + shy;
3336                     by1 = gridi->c0[YY] + (ci_y+1)*gridi->sy + shy;
3337                 }
3338
3339                 get_cell_range(by0, by1,
3340                                gridj->ncy, gridj->c0[YY], gridj->sy, gridj->inv_sy,
3341                                d2z_cx, rl2,
3342                                &cyf, &cyl);
3343
3344                 if (cyf > cyl)
3345                 {
3346                     continue;
3347                 }
3348
3349                 d2z_cy = d2z;
3350                 if (by1 < gridj->c0[YY])
3351                 {
3352                     d2z_cy += gmx::square(gridj->c0[YY] - by1);
3353                 }
3354                 else if (by0 > gridj->c1[YY])
3355                 {
3356                     d2z_cy += gmx::square(by0 - gridj->c1[YY]);
3357                 }
3358
3359                 for (int tx = -shp[XX]; tx <= shp[XX]; tx++)
3360                 {
3361                     shift = XYZ2IS(tx, ty, tz);
3362
3363 #ifdef NBNXN_SHIFT_BACKWARD
3364                     if (gridi == gridj && shift > CENTRAL)
3365                     {
3366                         continue;
3367                     }
3368 #endif
3369
3370                     shx = tx*box[XX][XX] + ty*box[YY][XX] + tz*box[ZZ][XX];
3371
3372                     if (nbl->bSimple)
3373                     {
3374                         bx0 = bb_i[ci].lower[BB_X] + shx;
3375                         bx1 = bb_i[ci].upper[BB_X] + shx;
3376                     }
3377                     else
3378                     {
3379                         bx0 = gridi->c0[XX] + (ci_x  )*gridi->sx + shx;
3380                         bx1 = gridi->c0[XX] + (ci_x+1)*gridi->sx + shx;
3381                     }
3382
3383                     get_cell_range(bx0, bx1,
3384                                    gridj->ncx, gridj->c0[XX], gridj->sx, gridj->inv_sx,
3385                                    d2z_cy, rl2,
3386                                    &cxf, &cxl);
3387
3388                     if (cxf > cxl)
3389                     {
3390                         continue;
3391                     }
3392
3393                     if (nbl->bSimple)
3394                     {
3395                         new_ci_entry(nbl, cell0_i+ci, shift, flags_i[ci]);
3396                     }
3397                     else
3398                     {
3399                         new_sci_entry(nbl, cell0_i+ci, shift);
3400                     }
3401
3402 #ifndef NBNXN_SHIFT_BACKWARD
3403                     if (cxf < ci_x)
3404 #else
3405                     if (shift == CENTRAL && gridi == gridj &&
3406                         cxf < ci_x)
3407 #endif
3408                     {
3409                         /* Leave the pairs with i > j.
3410                          * x is the major index, so skip half of it.
3411                          */
3412                         cxf = ci_x;
3413                     }
3414
3415                     if (nbl->bSimple)
3416                     {
3417                         set_icell_bb_simple(bb_i, ci, shx, shy, shz,
3418                                             nbl->work->bb_ci);
3419                     }
3420                     else
3421                     {
3422 #ifdef NBNXN_BBXXXX
3423                         set_icell_bbxxxx_supersub(pbb_i, ci, shx, shy, shz,
3424                                                   nbl->work->pbb_ci);
3425 #else
3426                         set_icell_bb_supersub(bb_i, ci, shx, shy, shz,
3427                                               nbl->work->bb_ci);
3428 #endif
3429                     }
3430
3431                     nbs->icell_set_x(cell0_i+ci, shx, shy, shz,
3432                                      nbat->xstride, nbat->x,
3433                                      nbl->work);
3434
3435                     for (int cx = cxf; cx <= cxl; cx++)
3436                     {
3437                         d2zx = d2z;
3438                         if (gridj->c0[XX] + cx*gridj->sx > bx1)
3439                         {
3440                             d2zx += gmx::square(gridj->c0[XX] + cx*gridj->sx - bx1);
3441                         }
3442                         else if (gridj->c0[XX] + (cx+1)*gridj->sx < bx0)
3443                         {
3444                             d2zx += gmx::square(gridj->c0[XX] + (cx+1)*gridj->sx - bx0);
3445                         }
3446
3447 #ifndef NBNXN_SHIFT_BACKWARD
3448                         if (gridi == gridj &&
3449                             cx == 0 && cyf < ci_y)
3450 #else
3451                         if (gridi == gridj &&
3452                             cx == 0 && shift == CENTRAL && cyf < ci_y)
3453 #endif
3454                         {
3455                             /* Leave the pairs with i > j.
3456                              * Skip half of y when i and j have the same x.
3457                              */
3458                             cyf_x = ci_y;
3459                         }
3460                         else
3461                         {
3462                             cyf_x = cyf;
3463                         }
3464
3465                         for (int cy = cyf_x; cy <= cyl; cy++)
3466                         {
3467                             c0 = gridj->cxy_ind[cx*gridj->ncy+cy];
3468                             c1 = gridj->cxy_ind[cx*gridj->ncy+cy+1];
3469 #ifdef NBNXN_SHIFT_BACKWARD
3470                             if (gridi == gridj &&
3471                                 shift == CENTRAL && c0 < ci)
3472                             {
3473                                 c0 = ci;
3474                             }
3475 #endif
3476
3477                             d2zxy = d2zx;
3478                             if (gridj->c0[YY] + cy*gridj->sy > by1)
3479                             {
3480                                 d2zxy += gmx::square(gridj->c0[YY] + cy*gridj->sy - by1);
3481                             }
3482                             else if (gridj->c0[YY] + (cy+1)*gridj->sy < by0)
3483                             {
3484                                 d2zxy += gmx::square(gridj->c0[YY] + (cy+1)*gridj->sy - by0);
3485                             }
3486                             if (c1 > c0 && d2zxy < rl2)
3487                             {
3488                                 cs = c0 + static_cast<int>(bz1_frac*(c1 - c0));
3489                                 if (cs >= c1)
3490                                 {
3491                                     cs = c1 - 1;
3492                                 }
3493
3494                                 d2xy = d2zxy - d2z;
3495
3496                                 /* Find the lowest cell that can possibly
3497                                  * be within range.
3498                                  */
3499                                 cf = cs;
3500                                 while (cf > c0 &&
3501                                        (bbcz_j[cf*NNBSBB_D+1] >= bz0 ||
3502                                         d2xy + gmx::square(bbcz_j[cf*NNBSBB_D+1] - bz0) < rl2))
3503                                 {
3504                                     cf--;
3505                                 }
3506
3507                                 /* Find the highest cell that can possibly
3508                                  * be within range.
3509                                  */
3510                                 cl = cs;
3511                                 while (cl < c1-1 &&
3512                                        (bbcz_j[cl*NNBSBB_D] <= bz1 ||
3513                                         d2xy + gmx::square(bbcz_j[cl*NNBSBB_D] - bz1) < rl2))
3514                                 {
3515                                     cl++;
3516                                 }
3517
3518 #ifdef NBNXN_REFCODE
3519                                 {
3520                                     /* Simple reference code, for debugging,
3521                                      * overrides the more complex code above.
3522                                      */
3523                                     cf = c1;
3524                                     cl = -1;
3525                                     for (int k = c0; k < c1; k++)
3526                                     {
3527                                         if (box_dist2(bx0, bx1, by0, by1, bz0, bz1, bb+k) < rl2 &&
3528                                             k < cf)
3529                                         {
3530                                             cf = k;
3531                                         }
3532                                         if (box_dist2(bx0, bx1, by0, by1, bz0, bz1, bb+k) < rl2 &&
3533                                             k > cl)
3534                                         {
3535                                             cl = k;
3536                                         }
3537                                     }
3538                                 }
3539 #endif
3540
3541                                 if (gridi == gridj)
3542                                 {
3543                                     /* We want each atom/cell pair only once,
3544                                      * only use cj >= ci.
3545                                      */
3546 #ifndef NBNXN_SHIFT_BACKWARD
3547                                     cf = std::max(cf, ci);
3548 #else
3549                                     if (shift == CENTRAL)
3550                                     {
3551                                         cf = std::max(cf, ci);
3552                                     }
3553 #endif
3554                                 }
3555
3556                                 if (cf <= cl)
3557                                 {
3558                                     /* For f buffer flags with simple lists */
3559                                     ncj_old_j = nbl->ncj;
3560
3561                                     switch (nb_kernel_type)
3562                                     {
3563                                         case nbnxnk4x4_PlainC:
3564                                             check_cell_list_space_simple(nbl, cl-cf+1);
3565
3566                                             make_cluster_list_simple(gridj,
3567                                                                      nbl, ci, cf, cl,
3568                                                                      (gridi == gridj && shift == CENTRAL),
3569                                                                      nbat->x,
3570                                                                      rl2, rbb2,
3571                                                                      &ndistc);
3572                                             break;
3573 #ifdef GMX_NBNXN_SIMD_4XN
3574                                         case nbnxnk4xN_SIMD_4xN:
3575                                             check_cell_list_space_simple(nbl, ci_to_cj(na_cj_2log, cl-cf)+2);
3576                                             make_cluster_list_simd_4xn(gridj,
3577                                                                        nbl, ci, cf, cl,
3578                                                                        (gridi == gridj && shift == CENTRAL),
3579                                                                        nbat->x,
3580                                                                        rl2, rbb2,
3581                                                                        &ndistc);
3582                                             break;
3583 #endif
3584 #ifdef GMX_NBNXN_SIMD_2XNN
3585                                         case nbnxnk4xN_SIMD_2xNN:
3586                                             check_cell_list_space_simple(nbl, ci_to_cj(na_cj_2log, cl-cf)+2);
3587                                             make_cluster_list_simd_2xnn(gridj,
3588                                                                         nbl, ci, cf, cl,
3589                                                                         (gridi == gridj && shift == CENTRAL),
3590                                                                         nbat->x,
3591                                                                         rl2, rbb2,
3592                                                                         &ndistc);
3593                                             break;
3594 #endif
3595                                         case nbnxnk8x8x8_PlainC:
3596                                         case nbnxnk8x8x8_GPU:
3597                                             check_cell_list_space_supersub(nbl, cl-cf+1);
3598                                             for (cj = cf; cj <= cl; cj++)
3599                                             {
3600                                                 make_cluster_list_supersub(gridi, gridj,
3601                                                                            nbl, ci, cj,
3602                                                                            (gridi == gridj && shift == CENTRAL && ci == cj),
3603                                                                            nbat->xstride, nbat->x,
3604                                                                            rl2, rbb2,
3605                                                                            &ndistc);
3606                                             }
3607                                             break;
3608                                     }
3609                                     ncpcheck += cl - cf + 1;
3610
3611                                     if (bFBufferFlag && nbl->ncj > ncj_old_j)
3612                                     {
3613                                         int cbf = nbl->cj[ncj_old_j].cj >> gridj_flag_shift;
3614                                         int cbl = nbl->cj[nbl->ncj-1].cj >> gridj_flag_shift;
3615                                         for (int cb = cbf; cb <= cbl; cb++)
3616                                         {
3617                                             bitmask_init_bit(&gridj_flag[cb], th);
3618                                         }
3619                                     }
3620                                 }
3621                             }
3622                         }
3623                     }
3624
3625                     /* Set the exclusions for this ci list */
3626                     if (nbl->bSimple)
3627                     {
3628                         set_ci_top_excls(nbs,
3629                                          nbl,
3630                                          shift == CENTRAL && gridi == gridj,
3631                                          gridj->na_c_2log,
3632                                          na_cj_2log,
3633                                          &(nbl->ci[nbl->nci]),
3634                                          excl);
3635
3636                         if (nbs->bFEP)
3637                         {
3638                             make_fep_list(nbs, nbat, nbl,
3639                                           shift == CENTRAL && gridi == gridj,
3640                                           &(nbl->ci[nbl->nci]),
3641                                           gridi, gridj, nbl_fep);
3642                         }
3643                     }
3644                     else
3645                     {
3646                         set_sci_top_excls(nbs,
3647                                           nbl,
3648                                           shift == CENTRAL && gridi == gridj,
3649                                           gridj->na_c_2log,
3650                                           &(nbl->sci[nbl->nsci]),
3651                                           excl);
3652
3653                         if (nbs->bFEP)
3654                         {
3655                             make_fep_list_supersub(nbs, nbat, nbl,
3656                                                    shift == CENTRAL && gridi == gridj,
3657                                                    &(nbl->sci[nbl->nsci]),
3658                                                    shx, shy, shz,
3659                                                    rl_fep2,
3660                                                    gridi, gridj, nbl_fep);
3661                         }
3662                     }
3663
3664                     /* Close this ci list */
3665                     if (nbl->bSimple)
3666                     {
3667                         close_ci_entry_simple(nbl);
3668                     }
3669                     else
3670                     {
3671                         close_ci_entry_supersub(nbl,
3672                                                 nsubpair_max,
3673                                                 progBal, nsubpair_tot_est,
3674                                                 th, nth);
3675                     }
3676                 }
3677             }
3678         }
3679
3680         if (bFBufferFlag && nbl->ncj > ncj_old_i)
3681         {
3682             bitmask_init_bit(&(work->buffer_flags.flag[(gridi->cell0+ci)>>gridi_flag_shift]), th);
3683         }
3684     }
3685
3686     work->ndistc = ndistc;
3687
3688     nbs_cycle_stop(&work->cc[enbsCCsearch]);
3689
3690     if (debug)
3691     {
3692         fprintf(debug, "number of distance checks %d\n", ndistc);
3693         fprintf(debug, "ncpcheck %s %d\n", gridi == gridj ? "local" : "non-local",
3694                 ncpcheck);
3695
3696         if (nbl->bSimple)
3697         {
3698             print_nblist_statistics_simple(debug, nbl, nbs, rlist);
3699         }
3700         else
3701         {
3702             print_nblist_statistics_supersub(debug, nbl, nbs, rlist);
3703         }
3704
3705         if (nbs->bFEP)
3706         {
3707             fprintf(debug, "nbl FEP list pairs: %d\n", nbl_fep->nrj);
3708         }
3709     }
3710 }
3711
3712 static void reduce_buffer_flags(const nbnxn_search_t        nbs,
3713                                 int                         nsrc,
3714                                 const nbnxn_buffer_flags_t *dest)
3715 {
3716     for (int s = 0; s < nsrc; s++)
3717     {
3718         gmx_bitmask_t * flag = nbs->work[s].buffer_flags.flag;
3719
3720         for (int b = 0; b < dest->nflag; b++)
3721         {
3722             bitmask_union(&(dest->flag[b]), flag[b]);
3723         }
3724     }
3725 }
3726
3727 static void print_reduction_cost(const nbnxn_buffer_flags_t *flags, int nout)
3728 {
3729     int           nelem, nkeep, ncopy, nred, out;
3730     gmx_bitmask_t mask_0;
3731
3732     nelem = 0;
3733     nkeep = 0;
3734     ncopy = 0;
3735     nred  = 0;
3736     bitmask_init_bit(&mask_0, 0);
3737     for (int b = 0; b < flags->nflag; b++)
3738     {
3739         if (bitmask_is_equal(flags->flag[b], mask_0))
3740         {
3741             /* Only flag 0 is set, no copy of reduction required */
3742             nelem++;
3743             nkeep++;
3744         }
3745         else if (!bitmask_is_zero(flags->flag[b]))
3746         {
3747             int c = 0;
3748             for (out = 0; out < nout; out++)
3749             {
3750                 if (bitmask_is_set(flags->flag[b], out))
3751                 {
3752                     c++;
3753                 }
3754             }
3755             nelem += c;
3756             if (c == 1)
3757             {
3758                 ncopy++;
3759             }
3760             else
3761             {
3762                 nred += c;
3763             }
3764         }
3765     }
3766
3767     fprintf(debug, "nbnxn reduction: #flag %d #list %d elem %4.2f, keep %4.2f copy %4.2f red %4.2f\n",
3768             flags->nflag, nout,
3769             nelem/(double)(flags->nflag),
3770             nkeep/(double)(flags->nflag),
3771             ncopy/(double)(flags->nflag),
3772             nred/(double)(flags->nflag));
3773 }
3774
3775 /* Perform a count (linear) sort to sort the smaller lists to the end.
3776  * This avoids load imbalance on the GPU, as large lists will be
3777  * scheduled and executed first and the smaller lists later.
3778  * Load balancing between multi-processors only happens at the end
3779  * and there smaller lists lead to more effective load balancing.
3780  * The sorting is done on the cj4 count, not on the actual pair counts.
3781  * Not only does this make the sort faster, but it also results in
3782  * better load balancing than using a list sorted on exact load.
3783  * This function swaps the pointer in the pair list to avoid a copy operation.
3784  */
3785 static void sort_sci(nbnxn_pairlist_t *nbl)
3786 {
3787     nbnxn_list_work_t *work;
3788     int                m, s0, s1;
3789     nbnxn_sci_t       *sci_sort;
3790
3791     if (nbl->ncj4 <= nbl->nsci)
3792     {
3793         /* nsci = 0 or all sci have size 1, sorting won't change the order */
3794         return;
3795     }
3796
3797     work = nbl->work;
3798
3799     /* We will distinguish differences up to double the average */
3800     m = (2*nbl->ncj4)/nbl->nsci;
3801
3802     if (m + 1 > work->sort_nalloc)
3803     {
3804         work->sort_nalloc = over_alloc_large(m + 1);
3805         srenew(work->sort, work->sort_nalloc);
3806     }
3807
3808     if (work->sci_sort_nalloc != nbl->sci_nalloc)
3809     {
3810         work->sci_sort_nalloc = nbl->sci_nalloc;
3811         nbnxn_realloc_void((void **)&work->sci_sort,
3812                            0,
3813                            work->sci_sort_nalloc*sizeof(*work->sci_sort),
3814                            nbl->alloc, nbl->free);
3815     }
3816
3817     /* Count the entries of each size */
3818     for (int i = 0; i <= m; i++)
3819     {
3820         work->sort[i] = 0;
3821     }
3822     for (int s = 0; s < nbl->nsci; s++)
3823     {
3824         int i = std::min(m, nbl->sci[s].cj4_ind_end - nbl->sci[s].cj4_ind_start);
3825         work->sort[i]++;
3826     }
3827     /* Calculate the offset for each count */
3828     s0            = work->sort[m];
3829     work->sort[m] = 0;
3830     for (int i = m - 1; i >= 0; i--)
3831     {
3832         s1            = work->sort[i];
3833         work->sort[i] = work->sort[i + 1] + s0;
3834         s0            = s1;
3835     }
3836
3837     /* Sort entries directly into place */
3838     sci_sort = work->sci_sort;
3839     for (int s = 0; s < nbl->nsci; s++)
3840     {
3841         int i = std::min(m, nbl->sci[s].cj4_ind_end - nbl->sci[s].cj4_ind_start);
3842         sci_sort[work->sort[i]++] = nbl->sci[s];
3843     }
3844
3845     /* Swap the sci pointers so we use the new, sorted list */
3846     work->sci_sort = nbl->sci;
3847     nbl->sci       = sci_sort;
3848 }
3849
3850 /* Make a local or non-local pair-list, depending on iloc */
3851 void nbnxn_make_pairlist(const nbnxn_search_t  nbs,
3852                          nbnxn_atomdata_t     *nbat,
3853                          const t_blocka       *excl,
3854                          real                  rlist,
3855                          int                   min_ci_balanced,
3856                          nbnxn_pairlist_set_t *nbl_list,
3857                          int                   iloc,
3858                          int                   nb_kernel_type,
3859                          t_nrnb               *nrnb)
3860 {
3861     nbnxn_grid_t      *gridi, *gridj;
3862     gmx_bool           bGPUCPU;
3863     int                nzi, zj0, zj1;
3864     int                nsubpair_target, nsubpair_tot_est;
3865     int                nnbl;
3866     nbnxn_pairlist_t **nbl;
3867     int                ci_block;
3868     gmx_bool           CombineNBLists;
3869     gmx_bool           progBal;
3870     int                np_tot, np_noq, np_hlj, nap;
3871
3872     /* Check if we are running hybrid GPU + CPU nbnxn mode */
3873     bGPUCPU = (!nbs->grid[0].bSimple && nbl_list->bSimple);
3874
3875     nnbl            = nbl_list->nnbl;
3876     nbl             = nbl_list->nbl;
3877     CombineNBLists  = nbl_list->bCombined;
3878
3879     if (debug)
3880     {
3881         fprintf(debug, "ns making %d nblists\n", nnbl);
3882     }
3883
3884     nbat->bUseBufferFlags = (nbat->nout > 1);
3885     /* We should re-init the flags before making the first list */
3886     if (nbat->bUseBufferFlags && (LOCAL_I(iloc) || bGPUCPU))
3887     {
3888         init_buffer_flags(&nbat->buffer_flags, nbat->natoms);
3889     }
3890
3891     if (nbl_list->bSimple)
3892     {
3893 #if GMX_SIMD
3894         switch (nb_kernel_type)
3895         {
3896 #ifdef GMX_NBNXN_SIMD_4XN
3897             case nbnxnk4xN_SIMD_4xN:
3898                 nbs->icell_set_x = icell_set_x_simd_4xn;
3899                 break;
3900 #endif
3901 #ifdef GMX_NBNXN_SIMD_2XNN
3902             case nbnxnk4xN_SIMD_2xNN:
3903                 nbs->icell_set_x = icell_set_x_simd_2xnn;
3904                 break;
3905 #endif
3906             default:
3907                 nbs->icell_set_x = icell_set_x_simple;
3908                 break;
3909         }
3910 #else   // GMX_SIMD
3911         /* MSVC 2013 complains about switch statements without case */
3912         nbs->icell_set_x = icell_set_x_simple;
3913 #endif  // GMX_SIMD
3914     }
3915     else
3916     {
3917         nbs->icell_set_x = icell_set_x_supersub;
3918     }
3919
3920     if (LOCAL_I(iloc))
3921     {
3922         /* Only zone (grid) 0 vs 0 */
3923         nzi = 1;
3924         zj0 = 0;
3925         zj1 = 1;
3926     }
3927     else
3928     {
3929         nzi = nbs->zones->nizone;
3930     }
3931
3932     if (!nbl_list->bSimple && min_ci_balanced > 0)
3933     {
3934         get_nsubpair_target(nbs, iloc, rlist, min_ci_balanced,
3935                             &nsubpair_target, &nsubpair_tot_est);
3936     }
3937     else
3938     {
3939         nsubpair_target  = 0;
3940         nsubpair_tot_est = 0;
3941     }
3942
3943     /* Clear all pair-lists */
3944     for (int th = 0; th < nnbl; th++)
3945     {
3946         clear_pairlist(nbl[th]);
3947
3948         if (nbs->bFEP)
3949         {
3950             clear_pairlist_fep(nbl_list->nbl_fep[th]);
3951         }
3952     }
3953
3954     for (int zi = 0; zi < nzi; zi++)
3955     {
3956         gridi = &nbs->grid[zi];
3957
3958         if (NONLOCAL_I(iloc))
3959         {
3960             zj0 = nbs->zones->izone[zi].j0;
3961             zj1 = nbs->zones->izone[zi].j1;
3962             if (zi == 0)
3963             {
3964                 zj0++;
3965             }
3966         }
3967         for (int zj = zj0; zj < zj1; zj++)
3968         {
3969             gridj = &nbs->grid[zj];
3970
3971             if (debug)
3972             {
3973                 fprintf(debug, "ns search grid %d vs %d\n", zi, zj);
3974             }
3975
3976             nbs_cycle_start(&nbs->cc[enbsCCsearch]);
3977
3978             if (nbl[0]->bSimple && !gridi->bSimple)
3979             {
3980                 /* Hybrid list, determine blocking later */
3981                 ci_block = 0;
3982             }
3983             else
3984             {
3985                 ci_block = get_ci_block_size(gridi, nbs->DomDec, nnbl);
3986             }
3987
3988             /* With GPU: generate progressively smaller lists for
3989              * load balancing for local only or non-local with 2 zones.
3990              */
3991             progBal = (LOCAL_I(iloc) || nbs->zones->n <= 2);
3992
3993 #pragma omp parallel for num_threads(nnbl) schedule(static)
3994             for (int th = 0; th < nnbl; th++)
3995             {
3996                 try
3997                 {
3998                     /* Re-init the thread-local work flag data before making
3999                      * the first list (not an elegant conditional).
4000                      */
4001                     if (nbat->bUseBufferFlags && ((zi == 0 && zj == 0) ||
4002                                                   (bGPUCPU && zi == 0 && zj == 1)))
4003                     {
4004                         init_buffer_flags(&nbs->work[th].buffer_flags, nbat->natoms);
4005                     }
4006
4007                     if (CombineNBLists && th > 0)
4008                     {
4009                         clear_pairlist(nbl[th]);
4010                     }
4011
4012                     /* Divide the i super cell equally over the nblists */
4013                     nbnxn_make_pairlist_part(nbs, gridi, gridj,
4014                                              &nbs->work[th], nbat, excl,
4015                                              rlist,
4016                                              nb_kernel_type,
4017                                              ci_block,
4018                                              nbat->bUseBufferFlags,
4019                                              nsubpair_target,
4020                                              progBal, nsubpair_tot_est,
4021                                              th, nnbl,
4022                                              nbl[th],
4023                                              nbl_list->nbl_fep[th]);
4024                 }
4025                 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4026             }
4027             nbs_cycle_stop(&nbs->cc[enbsCCsearch]);
4028
4029             np_tot = 0;
4030             np_noq = 0;
4031             np_hlj = 0;
4032             for (int th = 0; th < nnbl; th++)
4033             {
4034                 inc_nrnb(nrnb, eNR_NBNXN_DIST2, nbs->work[th].ndistc);
4035
4036                 if (nbl_list->bSimple)
4037                 {
4038                     np_tot += nbl[th]->ncj;
4039                     np_noq += nbl[th]->work->ncj_noq;
4040                     np_hlj += nbl[th]->work->ncj_hlj;
4041                 }
4042                 else
4043                 {
4044                     /* This count ignores potential subsequent pair pruning */
4045                     np_tot += nbl[th]->nci_tot;
4046                 }
4047             }
4048             nap                   = nbl[0]->na_ci*nbl[0]->na_cj;
4049             nbl_list->natpair_ljq = (np_tot - np_noq)*nap - np_hlj*nap/2;
4050             nbl_list->natpair_lj  = np_noq*nap;
4051             nbl_list->natpair_q   = np_hlj*nap/2;
4052
4053             if (CombineNBLists && nnbl > 1)
4054             {
4055                 nbs_cycle_start(&nbs->cc[enbsCCcombine]);
4056
4057                 combine_nblists(nnbl-1, nbl+1, nbl[0]);
4058
4059                 nbs_cycle_stop(&nbs->cc[enbsCCcombine]);
4060             }
4061         }
4062     }
4063
4064     if (!nbl_list->bSimple)
4065     {
4066         /* Sort the entries on size, large ones first */
4067         if (CombineNBLists || nnbl == 1)
4068         {
4069             sort_sci(nbl[0]);
4070         }
4071         else
4072         {
4073 #pragma omp parallel for num_threads(nnbl) schedule(static)
4074             for (int th = 0; th < nnbl; th++)
4075             {
4076                 try
4077                 {
4078                     sort_sci(nbl[th]);
4079                 }
4080                 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
4081             }
4082         }
4083     }
4084
4085     if (nbat->bUseBufferFlags)
4086     {
4087         reduce_buffer_flags(nbs, nnbl, &nbat->buffer_flags);
4088     }
4089
4090     if (nbs->bFEP)
4091     {
4092         /* Balance the free-energy lists over all the threads */
4093         balance_fep_lists(nbs, nbl_list);
4094     }
4095
4096     /* Special performance logging stuff (env.var. GMX_NBNXN_CYCLE) */
4097     if (LOCAL_I(iloc))
4098     {
4099         nbs->search_count++;
4100     }
4101     if (nbs->print_cycles &&
4102         (!nbs->DomDec || (nbs->DomDec && !LOCAL_I(iloc))) &&
4103         nbs->search_count % 100 == 0)
4104     {
4105         nbs_cycle_print(stderr, nbs);
4106     }
4107
4108     if (debug && (CombineNBLists && nnbl > 1))
4109     {
4110         if (nbl[0]->bSimple)
4111         {
4112             print_nblist_statistics_simple(debug, nbl[0], nbs, rlist);
4113         }
4114         else
4115         {
4116             print_nblist_statistics_supersub(debug, nbl[0], nbs, rlist);
4117         }
4118     }
4119
4120     if (debug)
4121     {
4122         if (gmx_debug_at)
4123         {
4124             if (nbl[0]->bSimple)
4125             {
4126                 print_nblist_ci_cj(debug, nbl[0]);
4127             }
4128             else
4129             {
4130                 print_nblist_sci_cj(debug, nbl[0]);
4131             }
4132         }
4133
4134         if (nbat->bUseBufferFlags)
4135         {
4136             print_reduction_cost(&nbat->buffer_flags, nnbl);
4137         }
4138     }
4139 }