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