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