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