Decouple GPU force buffer management from buffer ops in NBNXM
[alexxy/gromacs.git] / src / gromacs / nbnxm / atomdata.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
36 #include "gmxpre.h"
37
38 #include "atomdata.h"
39
40 #include <cassert>
41 #include <cmath>
42 #include <cstdlib>
43 #include <cstring>
44
45 #include <algorithm>
46
47 #include "thread_mpi/atomic.h"
48
49 #include "gromacs/math/functions.h"
50 #include "gromacs/math/utilities.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/mdlib/gmx_omp_nthreads.h"
53 #include "gromacs/mdtypes/forcerec.h" // only for GET_CGINFO_*
54 #include "gromacs/mdtypes/mdatom.h"
55 #include "gromacs/nbnxm/nbnxm.h"
56 #include "gromacs/nbnxm/nbnxm_geometry.h"
57 #include "gromacs/nbnxm/pairlist.h"
58 #include "gromacs/pbcutil/ishift.h"
59 #include "gromacs/simd/simd.h"
60 #include "gromacs/utility/exceptions.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/gmxomp.h"
63 #include "gromacs/utility/logger.h"
64 #include "gromacs/utility/smalloc.h"
65 #include "gromacs/utility/strconvert.h"
66 #include "gromacs/utility/stringutil.h"
67
68 #include "grid.h"
69 #include "gridset.h"
70
71 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
72
73 void nbnxn_atomdata_t::resizeCoordinateBuffer(int numAtoms)
74 {
75     numAtoms_ = numAtoms;
76
77     x_.resize(numAtoms*xstride);
78 }
79
80 void nbnxn_atomdata_t::resizeForceBuffers()
81 {
82     /* Force buffers need padding up to a multiple of the buffer flag size */
83     const int paddedSize = (numAtoms() + NBNXN_BUFFERFLAG_SIZE - 1)/NBNXN_BUFFERFLAG_SIZE*NBNXN_BUFFERFLAG_SIZE;
84
85     /* Should we let each thread allocate it's own data instead? */
86     for (nbnxn_atomdata_output_t &outBuffer : out)
87     {
88         outBuffer.f.resize(paddedSize*fstride);
89     }
90 }
91
92 /* Initializes an nbnxn_atomdata_output_t data structure */
93 nbnxn_atomdata_output_t::nbnxn_atomdata_output_t(Nbnxm::KernelType  kernelType,
94                                                  int                numEnergyGroups,
95                                                  int                simdEnergyBufferStride,
96                                                  gmx::PinningPolicy pinningPolicy) :
97     f({}, {pinningPolicy}),
98     fshift({}, {pinningPolicy}),
99     Vvdw({}, {pinningPolicy}),
100     Vc({}, {pinningPolicy})
101 {
102     fshift.resize(SHIFTS*DIM);
103     Vvdw.resize(numEnergyGroups*numEnergyGroups);
104     Vc.resize(numEnergyGroups*numEnergyGroups);
105
106     if (Nbnxm::kernelTypeIsSimd(kernelType))
107     {
108         int cj_size     = Nbnxm::JClusterSizePerKernelType[kernelType];
109         int numElements = numEnergyGroups*numEnergyGroups*simdEnergyBufferStride*(cj_size/2)*cj_size;
110         VSvdw.resize(numElements);
111         VSc.resize(numElements);
112     }
113 }
114
115 static void copy_int_to_nbat_int(const int *a, int na, int na_round,
116                                  const int *in, int fill, int *innb)
117 {
118     int i, j;
119
120     j = 0;
121     for (i = 0; i < na; i++)
122     {
123         innb[j++] = in[a[i]];
124     }
125     /* Complete the partially filled last cell with fill */
126     for (; i < na_round; i++)
127     {
128         innb[j++] = fill;
129     }
130 }
131
132 void copy_rvec_to_nbat_real(const int *a, int na, int na_round,
133                             const rvec *x, int nbatFormat,
134                             real *xnb, int a0)
135 {
136     /* We complete partially filled cells, can only be the last one in each
137      * column, with coordinates farAway. The actual coordinate value does
138      * not influence the results, since these filler particles do not interact.
139      * Clusters with normal atoms + fillers have a bounding box based only
140      * on the coordinates of the atoms. Clusters with only fillers have as
141      * the bounding box the coordinates of the first filler. Such clusters
142      * are not considered as i-entries, but they are considered as j-entries.
143      * So for performance it is better to have their bounding boxes far away,
144      * such that filler only clusters don't end up in the pair list.
145      */
146     const real farAway = -1000000;
147
148     int        i, j, c;
149
150     switch (nbatFormat)
151     {
152         case nbatXYZ:
153             j = a0*STRIDE_XYZ;
154             for (i = 0; i < na; i++)
155             {
156                 xnb[j++] = x[a[i]][XX];
157                 xnb[j++] = x[a[i]][YY];
158                 xnb[j++] = x[a[i]][ZZ];
159             }
160             /* Complete the partially filled last cell with farAway elements */
161             for (; i < na_round; i++)
162             {
163                 xnb[j++] = farAway;
164                 xnb[j++] = farAway;
165                 xnb[j++] = farAway;
166             }
167             break;
168         case nbatXYZQ:
169             j = a0*STRIDE_XYZQ;
170             for (i = 0; i < na; i++)
171             {
172                 xnb[j++] = x[a[i]][XX];
173                 xnb[j++] = x[a[i]][YY];
174                 xnb[j++] = x[a[i]][ZZ];
175                 j++;
176             }
177             /* Complete the partially filled last cell with zeros */
178             for (; i < na_round; i++)
179             {
180                 xnb[j++] = farAway;
181                 xnb[j++] = farAway;
182                 xnb[j++] = farAway;
183                 j++;
184             }
185             break;
186         case nbatX4:
187             j = atom_to_x_index<c_packX4>(a0);
188             c = a0 & (c_packX4-1);
189             for (i = 0; i < na; i++)
190             {
191                 xnb[j+XX*c_packX4] = x[a[i]][XX];
192                 xnb[j+YY*c_packX4] = x[a[i]][YY];
193                 xnb[j+ZZ*c_packX4] = x[a[i]][ZZ];
194                 j++;
195                 c++;
196                 if (c == c_packX4)
197                 {
198                     j += (DIM-1)*c_packX4;
199                     c  = 0;
200                 }
201             }
202             /* Complete the partially filled last cell with zeros */
203             for (; i < na_round; i++)
204             {
205                 xnb[j+XX*c_packX4] = farAway;
206                 xnb[j+YY*c_packX4] = farAway;
207                 xnb[j+ZZ*c_packX4] = farAway;
208                 j++;
209                 c++;
210                 if (c == c_packX4)
211                 {
212                     j += (DIM-1)*c_packX4;
213                     c  = 0;
214                 }
215             }
216             break;
217         case nbatX8:
218             j = atom_to_x_index<c_packX8>(a0);
219             c = a0 & (c_packX8 - 1);
220             for (i = 0; i < na; i++)
221             {
222                 xnb[j+XX*c_packX8] = x[a[i]][XX];
223                 xnb[j+YY*c_packX8] = x[a[i]][YY];
224                 xnb[j+ZZ*c_packX8] = x[a[i]][ZZ];
225                 j++;
226                 c++;
227                 if (c == c_packX8)
228                 {
229                     j += (DIM-1)*c_packX8;
230                     c  = 0;
231                 }
232             }
233             /* Complete the partially filled last cell with zeros */
234             for (; i < na_round; i++)
235             {
236                 xnb[j+XX*c_packX8] = farAway;
237                 xnb[j+YY*c_packX8] = farAway;
238                 xnb[j+ZZ*c_packX8] = farAway;
239                 j++;
240                 c++;
241                 if (c == c_packX8)
242                 {
243                     j += (DIM-1)*c_packX8;
244                     c  = 0;
245                 }
246             }
247             break;
248         default:
249             gmx_incons("Unsupported nbnxn_atomdata_t format");
250     }
251 }
252
253 /* Stores the LJ parameter data in a format convenient for different kernels */
254 static void set_lj_parameter_data(nbnxn_atomdata_t::Params *params, gmx_bool bSIMD)
255 {
256     int  nt = params->numTypes;
257
258     if (bSIMD)
259     {
260 #if GMX_SIMD
261         /* nbfp_aligned stores two parameters using the stride most suitable
262          * for the present SIMD architecture, as specified by the constant
263          * c_simdBestPairAlignment from the SIMD header.
264          * There's a slight inefficiency in allocating and initializing nbfp_aligned
265          * when it might not be used, but introducing the conditional code is not
266          * really worth it.
267          */
268         params->nbfp_aligned.resize(nt*nt*c_simdBestPairAlignment);
269
270         for (int i = 0; i < nt; i++)
271         {
272             for (int j = 0; j < nt; j++)
273             {
274                 params->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+0] = params->nbfp[(i*nt+j)*2+0];
275                 params->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+1] = params->nbfp[(i*nt+j)*2+1];
276                 if (c_simdBestPairAlignment > 2)
277                 {
278                     params->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+2] = 0;
279                     params->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+3] = 0;
280                 }
281             }
282         }
283 #endif
284     }
285
286     /* We use combination rule data for SIMD combination rule kernels
287      * and with LJ-PME kernels. We then only need parameters per atom type,
288      * not per pair of atom types.
289      */
290     params->nbfp_comb.resize(nt*2);
291     switch (params->comb_rule)
292     {
293         case ljcrGEOM:
294             params->comb_rule = ljcrGEOM;
295
296             for (int i = 0; i < nt; i++)
297             {
298                 /* Store the sqrt of the diagonal from the nbfp matrix */
299                 params->nbfp_comb[i*2  ] = std::sqrt(params->nbfp[(i*nt+i)*2  ]);
300                 params->nbfp_comb[i*2+1] = std::sqrt(params->nbfp[(i*nt+i)*2+1]);
301             }
302             break;
303         case ljcrLB:
304             for (int i = 0; i < nt; i++)
305             {
306                 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
307                 const real c6  = params->nbfp[(i*nt+i)*2  ];
308                 const real c12 = params->nbfp[(i*nt+i)*2+1];
309                 if (c6 > 0 && c12 > 0)
310                 {
311                     /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
312                      * so we get 6*C6 and 12*C12 after combining.
313                      */
314                     params->nbfp_comb[i*2  ] = 0.5*gmx::sixthroot(c12/c6);
315                     params->nbfp_comb[i*2+1] = std::sqrt(c6*c6/c12);
316                 }
317                 else
318                 {
319                     params->nbfp_comb[i*2  ] = 0;
320                     params->nbfp_comb[i*2+1] = 0;
321                 }
322             }
323             break;
324         case ljcrNONE:
325             /* We always store the full matrix (see code above) */
326             break;
327         default:
328             gmx_incons("Unknown combination rule");
329     }
330 }
331
332 nbnxn_atomdata_t::SimdMasks::SimdMasks()
333 {
334 #if GMX_SIMD
335     constexpr int simd_width = GMX_SIMD_REAL_WIDTH;
336     /* Set the diagonal cluster pair exclusion mask setup data.
337      * In the kernel we check 0 < j - i to generate the masks.
338      * Here we store j - i for generating the mask for the first i,
339      * we substract 0.5 to avoid rounding issues.
340      * In the kernel we can subtract 1 to generate the subsequent mask.
341      */
342     const int simd_4xn_diag_size = std::max(c_nbnxnCpuIClusterSize, simd_width);
343     diagonal_4xn_j_minus_i.resize(simd_4xn_diag_size);
344     for (int j = 0; j < simd_4xn_diag_size; j++)
345     {
346         diagonal_4xn_j_minus_i[j] = j - 0.5;
347     }
348
349     diagonal_2xnn_j_minus_i.resize(simd_width);
350     for (int j = 0; j < simd_width/2; j++)
351     {
352         /* The j-cluster size is half the SIMD width */
353         diagonal_2xnn_j_minus_i[j]                = j - 0.5;
354         /* The next half of the SIMD width is for i + 1 */
355         diagonal_2xnn_j_minus_i[simd_width/2 + j] = j - 1 - 0.5;
356     }
357
358     /* We use up to 32 bits for exclusion masking.
359      * The same masks are used for the 4xN and 2x(N+N) kernels.
360      * The masks are read either into integer SIMD registers or into
361      * real SIMD registers (together with a cast).
362      * In single precision this means the real and integer SIMD registers
363      * are of equal size.
364      */
365     const int simd_excl_size = c_nbnxnCpuIClusterSize*simd_width;
366 #if GMX_DOUBLE && !GMX_SIMD_HAVE_INT32_LOGICAL
367     exclusion_filter64.resize(simd_excl_size);
368 #else
369     exclusion_filter.resize(simd_excl_size);
370 #endif
371
372     for (int j = 0; j < simd_excl_size; j++)
373     {
374         /* Set the consecutive bits for masking pair exclusions */
375 #if GMX_DOUBLE && !GMX_SIMD_HAVE_INT32_LOGICAL
376         exclusion_filter64[j] = (1U << j);
377 #else
378         exclusion_filter[j]   = (1U << j);
379 #endif
380     }
381
382     if (!GMX_SIMD_HAVE_LOGICAL && !GMX_SIMD_HAVE_INT32_LOGICAL) // NOLINT(misc-redundant-expression)
383     {
384         // If the SIMD implementation has no bitwise logical operation support
385         // whatsoever we cannot use the normal masking. Instead,
386         // we generate a vector of all 2^4 possible ways an i atom
387         // interacts with its 4 j atoms. Each array entry contains
388         // GMX_SIMD_REAL_WIDTH values that are read with a single aligned SIMD load.
389         // Since there is no logical value representation in this case, we use
390         // any nonzero value to indicate 'true', while zero mean 'false'.
391         // This can then be converted to a SIMD boolean internally in the SIMD
392         // module by comparing to zero.
393         // Each array entry encodes how this i atom will interact with the 4 j atoms.
394         // Matching code exists in set_ci_top_excls() to generate indices into this array.
395         // Those indices are used in the kernels.
396
397         const int  simd_excl_size = c_nbnxnCpuIClusterSize*c_nbnxnCpuIClusterSize;
398         const real simdFalse      = 0.0;
399         const real simdTrue       = 1.0;
400
401         interaction_array.resize(simd_excl_size * GMX_SIMD_REAL_WIDTH);
402         for (int j = 0; j < simd_excl_size; j++)
403         {
404             const int index = j * GMX_SIMD_REAL_WIDTH;
405             for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
406             {
407                 interaction_array[index + i] = (j & (1 << i)) ? simdTrue : simdFalse;
408             }
409         }
410     }
411 #endif
412 }
413
414 nbnxn_atomdata_t::Params::Params(gmx::PinningPolicy pinningPolicy) :
415     numTypes(0),
416     nbfp({}, {pinningPolicy}),
417     nbfp_comb({}, {pinningPolicy}),
418     type({}, {pinningPolicy}),
419     lj_comb({}, {pinningPolicy}),
420     q({}, {pinningPolicy}),
421     nenergrp(0),
422     neg_2log(0),
423     energrp({}, {pinningPolicy})
424 {
425 }
426
427 nbnxn_atomdata_t::nbnxn_atomdata_t(gmx::PinningPolicy pinningPolicy) :
428     params_(pinningPolicy),
429     numAtoms_(0),
430     natoms_local(0),
431     shift_vec({}, {pinningPolicy}),
432     x_({}, {pinningPolicy}),
433     simdMasks(),
434     bUseBufferFlags(FALSE),
435     bUseTreeReduce(FALSE)
436 {
437 }
438
439 /* Initializes an nbnxn_atomdata_t::Params data structure */
440 static void nbnxn_atomdata_params_init(const gmx::MDLogger &mdlog,
441                                        nbnxn_atomdata_t::Params *params,
442                                        const Nbnxm::KernelType kernelType,
443                                        int enbnxninitcombrule,
444                                        int ntype, const real *nbfp,
445                                        int n_energygroups)
446 {
447     real     c6, c12, tol;
448     char    *ptr;
449     gmx_bool simple, bCombGeom, bCombLB, bSIMD;
450
451     if (debug)
452     {
453         fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
454     }
455     params->numTypes = ntype + 1;
456     params->nbfp.resize(params->numTypes*params->numTypes*2);
457     params->nbfp_comb.resize(params->numTypes*2);
458
459     /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
460      * force-field floating point parameters.
461      */
462     tol = 1e-5;
463     ptr = getenv("GMX_LJCOMB_TOL");
464     if (ptr != nullptr)
465     {
466         double dbl;
467
468         sscanf(ptr, "%lf", &dbl);
469         tol = dbl;
470     }
471     bCombGeom = TRUE;
472     bCombLB   = TRUE;
473
474     /* Temporarily fill params->nbfp_comb with sigma and epsilon
475      * to check for the LB rule.
476      */
477     for (int i = 0; i < ntype; i++)
478     {
479         c6  = nbfp[(i*ntype+i)*2    ]/6.0;
480         c12 = nbfp[(i*ntype+i)*2 + 1]/12.0;
481         if (c6 > 0 && c12 > 0)
482         {
483             params->nbfp_comb[i*2    ] = gmx::sixthroot(c12/c6);
484             params->nbfp_comb[i*2 + 1] = 0.25*c6*c6/c12;
485         }
486         else if (c6 == 0 && c12 == 0)
487         {
488             params->nbfp_comb[i*2    ] = 0;
489             params->nbfp_comb[i*2 + 1] = 0;
490         }
491         else
492         {
493             /* Can not use LB rule with only dispersion or repulsion */
494             bCombLB = FALSE;
495         }
496     }
497
498     for (int i = 0; i < params->numTypes; i++)
499     {
500         for (int j = 0; j < params->numTypes; j++)
501         {
502             if (i < ntype && j < ntype)
503             {
504                 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
505                  * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
506                  */
507                 c6  = nbfp[(i*ntype+j)*2    ];
508                 c12 = nbfp[(i*ntype+j)*2 + 1];
509                 params->nbfp[(i*params->numTypes+j)*2    ] = c6;
510                 params->nbfp[(i*params->numTypes+j)*2 + 1] = c12;
511
512                 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
513                 bCombGeom = bCombGeom &&
514                     gmx_within_tol(c6*c6, nbfp[(i*ntype+i)*2  ]*nbfp[(j*ntype+j)*2  ], tol) &&
515                     gmx_within_tol(c12*c12, nbfp[(i*ntype+i)*2 + 1]*nbfp[(j*ntype+j)*2 + 1], tol);
516
517                 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
518                 c6     /= 6.0;
519                 c12    /= 12.0;
520                 bCombLB = bCombLB &&
521                     ((c6 == 0 && c12 == 0 &&
522                       (params->nbfp_comb[i*2 + 1] == 0 || params->nbfp_comb[j*2 + 1] == 0)) ||
523                      (c6 > 0 && c12 > 0 &&
524                       gmx_within_tol(gmx::sixthroot(c12/c6),
525                                      0.5*(params->nbfp_comb[i*2]+params->nbfp_comb[j*2]), tol) &&
526                       gmx_within_tol(0.25*c6*c6/c12, std::sqrt(params->nbfp_comb[i*2 + 1]*params->nbfp_comb[j*2 + 1]), tol)));
527             }
528             else
529             {
530                 /* Add zero parameters for the additional dummy atom type */
531                 params->nbfp[(i*params->numTypes + j)*2  ] = 0;
532                 params->nbfp[(i*params->numTypes + j)*2+1] = 0;
533             }
534         }
535     }
536     if (debug)
537     {
538         fprintf(debug, "Combination rules: geometric %s Lorentz-Berthelot %s\n",
539                 gmx::boolToString(bCombGeom), gmx::boolToString(bCombLB));
540     }
541
542     simple = Nbnxm::kernelTypeUsesSimplePairlist(kernelType);
543
544     switch (enbnxninitcombrule)
545     {
546         case enbnxninitcombruleDETECT:
547             /* We prefer the geometic combination rule,
548              * as that gives a slightly faster kernel than the LB rule.
549              */
550             if (bCombGeom)
551             {
552                 params->comb_rule = ljcrGEOM;
553             }
554             else if (bCombLB)
555             {
556                 params->comb_rule = ljcrLB;
557             }
558             else
559             {
560                 params->comb_rule = ljcrNONE;
561
562                 params->nbfp_comb.clear();
563             }
564
565             {
566                 std::string mesg;
567                 if (params->comb_rule == ljcrNONE)
568                 {
569                     mesg = "Using full Lennard-Jones parameter combination matrix";
570                 }
571                 else
572                 {
573                     mesg = gmx::formatString("Using %s Lennard-Jones combination rule",
574                                              params->comb_rule == ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
575                 }
576                 GMX_LOG(mdlog.info).asParagraph().appendText(mesg);
577             }
578             break;
579         case enbnxninitcombruleGEOM:
580             params->comb_rule = ljcrGEOM;
581             break;
582         case enbnxninitcombruleLB:
583             params->comb_rule = ljcrLB;
584             break;
585         case enbnxninitcombruleNONE:
586             params->comb_rule = ljcrNONE;
587
588             params->nbfp_comb.clear();
589             break;
590         default:
591             gmx_incons("Unknown enbnxninitcombrule");
592     }
593
594     bSIMD = Nbnxm::kernelTypeIsSimd(kernelType);
595
596     set_lj_parameter_data(params, bSIMD);
597
598     params->nenergrp = n_energygroups;
599     if (!simple)
600     {
601         // We now check for energy groups already when starting mdrun
602         GMX_RELEASE_ASSERT(n_energygroups == 1, "GPU kernels do not support energy groups");
603     }
604     /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
605     if (params->nenergrp > 64)
606     {
607         gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
608     }
609     params->neg_2log = 1;
610     while (params->nenergrp > (1<<params->neg_2log))
611     {
612         params->neg_2log++;
613     }
614 }
615
616 /* Initializes an nbnxn_atomdata_t data structure */
617 void nbnxn_atomdata_init(const gmx::MDLogger &mdlog,
618                          nbnxn_atomdata_t *nbat,
619                          const Nbnxm::KernelType kernelType,
620                          int enbnxninitcombrule,
621                          int ntype, const real *nbfp,
622                          int n_energygroups,
623                          int nout)
624 {
625     nbnxn_atomdata_params_init(mdlog, &nbat->paramsDeprecated(), kernelType,
626                                enbnxninitcombrule, ntype, nbfp, n_energygroups);
627
628     const bool simple = Nbnxm::kernelTypeUsesSimplePairlist(kernelType);
629     const bool bSIMD  = Nbnxm::kernelTypeIsSimd(kernelType);
630
631     if (simple)
632     {
633         int pack_x;
634
635         if (bSIMD)
636         {
637             pack_x = std::max(c_nbnxnCpuIClusterSize,
638                               Nbnxm::JClusterSizePerKernelType[kernelType]);
639             switch (pack_x)
640             {
641                 case 4:
642                     nbat->XFormat = nbatX4;
643                     break;
644                 case 8:
645                     nbat->XFormat = nbatX8;
646                     break;
647                 default:
648                     gmx_incons("Unsupported packing width");
649             }
650         }
651         else
652         {
653             nbat->XFormat = nbatXYZ;
654         }
655
656         nbat->FFormat = nbat->XFormat;
657     }
658     else
659     {
660         nbat->XFormat = nbatXYZQ;
661         nbat->FFormat = nbatXYZ;
662     }
663
664     nbat->shift_vec.resize(SHIFTS);
665
666     nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
667     nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
668
669     /* Initialize the output data structures */
670     for (int i = 0; i < nout; i++)
671     {
672         const auto &pinningPolicy = nbat->params().type.get_allocator().pinningPolicy();
673         nbat->out.emplace_back(kernelType, nbat->params().nenergrp, 1 << nbat->params().neg_2log,
674                                pinningPolicy);
675     }
676
677     nbat->buffer_flags.flag        = nullptr;
678     nbat->buffer_flags.flag_nalloc = 0;
679
680     const int   nth = gmx_omp_nthreads_get(emntNonbonded);
681
682     const char *ptr = getenv("GMX_USE_TREEREDUCE");
683     if (ptr != nullptr)
684     {
685         nbat->bUseTreeReduce = (strtol(ptr, nullptr, 10) != 0);
686     }
687 #if defined __MIC__
688     else if (nth > 8) /*on the CPU we currently don't benefit even at 32*/
689     {
690         nbat->bUseTreeReduce = 1;
691     }
692 #endif
693     else
694     {
695         nbat->bUseTreeReduce = false;
696     }
697     if (nbat->bUseTreeReduce)
698     {
699         GMX_LOG(mdlog.info).asParagraph().appendText("Using tree force reduction");
700
701         nbat->syncStep = new tMPI_Atomic[nth];
702     }
703 }
704
705 template<int packSize>
706 static void copy_lj_to_nbat_lj_comb(gmx::ArrayRef<const real> ljparam_type,
707                                     const int *type, int na,
708                                     real *ljparam_at)
709 {
710     /* The LJ params follow the combination rule:
711      * copy the params for the type array to the atom array.
712      */
713     for (int is = 0; is < na; is += packSize)
714     {
715         for (int k = 0; k < packSize; k++)
716         {
717             int i = is + k;
718             ljparam_at[is*2            + k] = ljparam_type[type[i]*2    ];
719             ljparam_at[is*2 + packSize + k] = ljparam_type[type[i]*2 + 1];
720         }
721     }
722 }
723
724 /* Sets the atom type in nbnxn_atomdata_t */
725 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t::Params *params,
726                                          const Nbnxm::GridSet     &gridSet,
727                                          const int                *type)
728 {
729     params->type.resize(gridSet.numGridAtomsTotal());
730
731     for (const Nbnxm::Grid &grid : gridSet.grids())
732     {
733         /* Loop over all columns and copy and fill */
734         for (int i = 0; i < grid.numColumns(); i++)
735         {
736             const int numAtoms   = grid.paddedNumAtomsInColumn(i);
737             const int atomOffset = grid.firstAtomInColumn(i);
738
739             copy_int_to_nbat_int(gridSet.atomIndices().data() + atomOffset,
740                                  grid.numAtomsInColumn(i), numAtoms,
741                                  type, params->numTypes - 1, params->type.data() + atomOffset);
742         }
743     }
744 }
745
746 /* Sets the LJ combination rule parameters in nbnxn_atomdata_t */
747 static void nbnxn_atomdata_set_ljcombparams(nbnxn_atomdata_t::Params *params,
748                                             const int                 XFormat,
749                                             const Nbnxm::GridSet     &gridSet)
750 {
751     params->lj_comb.resize(gridSet.numGridAtomsTotal()*2);
752
753     if (params->comb_rule != ljcrNONE)
754     {
755         for (const Nbnxm::Grid &grid : gridSet.grids())
756         {
757             /* Loop over all columns and copy and fill */
758             for (int i = 0; i < grid.numColumns(); i++)
759             {
760                 const int numAtoms   = grid.paddedNumAtomsInColumn(i);
761                 const int atomOffset = grid.firstAtomInColumn(i);
762
763                 if (XFormat == nbatX4)
764                 {
765                     copy_lj_to_nbat_lj_comb<c_packX4>(params->nbfp_comb,
766                                                       params->type.data() + atomOffset,
767                                                       numAtoms,
768                                                       params->lj_comb.data() + atomOffset*2);
769                 }
770                 else if (XFormat == nbatX8)
771                 {
772                     copy_lj_to_nbat_lj_comb<c_packX8>(params->nbfp_comb,
773                                                       params->type.data() + atomOffset,
774                                                       numAtoms,
775                                                       params->lj_comb.data() + atomOffset*2);
776                 }
777                 else if (XFormat == nbatXYZQ)
778                 {
779                     copy_lj_to_nbat_lj_comb<1>(params->nbfp_comb,
780                                                params->type.data() + atomOffset,
781                                                numAtoms,
782                                                params->lj_comb.data() + atomOffset*2);
783                 }
784             }
785         }
786     }
787 }
788
789 /* Sets the charges in nbnxn_atomdata_t *nbat */
790 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t     *nbat,
791                                        const Nbnxm::GridSet &gridSet,
792                                        const real           *charge)
793 {
794     if (nbat->XFormat != nbatXYZQ)
795     {
796         nbat->paramsDeprecated().q.resize(nbat->numAtoms());
797     }
798
799     for (const Nbnxm::Grid &grid : gridSet.grids())
800     {
801         /* Loop over all columns and copy and fill */
802         for (int cxy = 0; cxy < grid.numColumns(); cxy++)
803         {
804             const int atomOffset     = grid.firstAtomInColumn(cxy);
805             const int numAtoms       = grid.numAtomsInColumn(cxy);
806             const int paddedNumAtoms = grid.paddedNumAtomsInColumn(cxy);
807
808             if (nbat->XFormat == nbatXYZQ)
809             {
810                 real *q = nbat->x().data() + atomOffset*STRIDE_XYZQ + ZZ + 1;
811                 int   i;
812                 for (i = 0; i < numAtoms; i++)
813                 {
814                     *q = charge[gridSet.atomIndices()[atomOffset + i]];
815                     q += STRIDE_XYZQ;
816                 }
817                 /* Complete the partially filled last cell with zeros */
818                 for (; i < paddedNumAtoms; i++)
819                 {
820                     *q = 0;
821                     q += STRIDE_XYZQ;
822                 }
823             }
824             else
825             {
826                 real *q = nbat->paramsDeprecated().q.data() + atomOffset;
827                 int   i;
828                 for (i = 0; i < numAtoms; i++)
829                 {
830                     *q = charge[gridSet.atomIndices()[atomOffset + i]];
831                     q++;
832                 }
833                 /* Complete the partially filled last cell with zeros */
834                 for (; i < paddedNumAtoms; i++)
835                 {
836                     *q = 0;
837                     q++;
838                 }
839             }
840         }
841     }
842 }
843
844 /* Set the charges of perturbed atoms in nbnxn_atomdata_t to 0.
845  * This is to automatically remove the RF/PME self term in the nbnxn kernels.
846  * Part of the zero interactions are still calculated in the normal kernels.
847  * All perturbed interactions are calculated in the free energy kernel,
848  * using the original charge and LJ data, not nbnxn_atomdata_t.
849  */
850 static void nbnxn_atomdata_mask_fep(nbnxn_atomdata_t     *nbat,
851                                     const Nbnxm::GridSet &gridSet)
852 {
853     nbnxn_atomdata_t::Params &params = nbat->paramsDeprecated();
854     real                     *q;
855     int                       stride_q;
856
857     if (nbat->XFormat == nbatXYZQ)
858     {
859         q        = nbat->x().data() + ZZ + 1;
860         stride_q = STRIDE_XYZQ;
861     }
862     else
863     {
864         q        = params.q.data();
865         stride_q = 1;
866     }
867
868     for (const Nbnxm::Grid &grid : gridSet.grids())
869     {
870         int nsubc;
871         if (grid.geometry().isSimple)
872         {
873             nsubc = 1;
874         }
875         else
876         {
877             nsubc = c_gpuNumClusterPerCell;
878         }
879
880         int c_offset = grid.firstAtomInColumn(0);
881
882         /* Loop over all columns and copy and fill */
883         for (int c = 0; c < grid.numCells()*nsubc; c++)
884         {
885             /* Does this cluster contain perturbed particles? */
886             if (grid.clusterIsPerturbed(c))
887             {
888                 const int numAtomsPerCluster = grid.geometry().numAtomsICluster;
889                 for (int i = 0; i < numAtomsPerCluster; i++)
890                 {
891                     /* Is this a perturbed particle? */
892                     if (grid.atomIsPerturbed(c, i))
893                     {
894                         int ind = c_offset + c*numAtomsPerCluster + i;
895                         /* Set atom type and charge to non-interacting */
896                         params.type[ind] = params.numTypes - 1;
897                         q[ind*stride_q]  = 0;
898                     }
899                 }
900             }
901         }
902     }
903 }
904
905 /* Copies the energy group indices to a reordered and packed array */
906 static void copy_egp_to_nbat_egps(const int *a, int na, int na_round,
907                                   int na_c, int bit_shift,
908                                   const int *in, int *innb)
909 {
910     int i;
911     int comb;
912
913     int j = 0;
914     for (i = 0; i < na; i += na_c)
915     {
916         /* Store na_c energy group numbers into one int */
917         comb = 0;
918         for (int sa = 0; sa < na_c; sa++)
919         {
920             int at = a[i+sa];
921             if (at >= 0)
922             {
923                 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
924             }
925         }
926         innb[j++] = comb;
927     }
928     /* Complete the partially filled last cell with fill */
929     for (; i < na_round; i += na_c)
930     {
931         innb[j++] = 0;
932     }
933 }
934
935 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
936 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t::Params *params,
937                                             const Nbnxm::GridSet     &gridSet,
938                                             const int                *atinfo)
939 {
940     if (params->nenergrp == 1)
941     {
942         return;
943     }
944
945     params->energrp.resize(gridSet.numGridAtomsTotal());
946
947     for (const Nbnxm::Grid &grid : gridSet.grids())
948     {
949         /* Loop over all columns and copy and fill */
950         for (int i = 0; i < grid.numColumns(); i++)
951         {
952             const int numAtoms   = grid.paddedNumAtomsInColumn(i);
953             const int atomOffset = grid.firstAtomInColumn(i);
954
955             copy_egp_to_nbat_egps(gridSet.atomIndices().data() + atomOffset,
956                                   grid.numAtomsInColumn(i), numAtoms,
957                                   c_nbnxnCpuIClusterSize, params->neg_2log,
958                                   atinfo,
959                                   params->energrp.data() + grid.atomToCluster(atomOffset));
960         }
961     }
962 }
963
964 /* Sets all required atom parameter data in nbnxn_atomdata_t */
965 void nbnxn_atomdata_set(nbnxn_atomdata_t     *nbat,
966                         const Nbnxm::GridSet &gridSet,
967                         const t_mdatoms      *mdatoms,
968                         const int            *atinfo)
969 {
970     nbnxn_atomdata_t::Params &params = nbat->paramsDeprecated();
971
972     nbnxn_atomdata_set_atomtypes(&params, gridSet, mdatoms->typeA);
973
974     nbnxn_atomdata_set_charges(nbat, gridSet, mdatoms->chargeA);
975
976     if (gridSet.haveFep())
977     {
978         nbnxn_atomdata_mask_fep(nbat, gridSet);
979     }
980
981     /* This must be done after masking types for FEP */
982     nbnxn_atomdata_set_ljcombparams(&params, nbat->XFormat, gridSet);
983
984     nbnxn_atomdata_set_energygroups(&params, gridSet, atinfo);
985 }
986
987 /* Copies the shift vector array to nbnxn_atomdata_t */
988 void nbnxn_atomdata_copy_shiftvec(gmx_bool          bDynamicBox,
989                                   rvec             *shift_vec,
990                                   nbnxn_atomdata_t *nbat)
991 {
992     int i;
993
994     nbat->bDynamicBox = bDynamicBox;
995     for (i = 0; i < SHIFTS; i++)
996     {
997         copy_rvec(shift_vec[i], nbat->shift_vec[i]);
998     }
999 }
1000
1001 // This is slightly different from nbnxn_get_atom_range(...) at the end of the file
1002 // TODO: Combine if possible
1003 static void getAtomRanges(const Nbnxm::GridSet      &gridSet,
1004                           const Nbnxm::AtomLocality  locality,
1005                           int                       *gridBegin,
1006                           int                       *gridEnd)
1007 {
1008     switch (locality)
1009     {
1010         case Nbnxm::AtomLocality::All:
1011             *gridBegin = 0;
1012             *gridEnd   = gridSet.grids().size();
1013             break;
1014         case Nbnxm::AtomLocality::Local:
1015             *gridBegin = 0;
1016             *gridEnd   = 1;
1017             break;
1018         case Nbnxm::AtomLocality::NonLocal:
1019             *gridBegin = 1;
1020             *gridEnd   = gridSet.grids().size();
1021             break;
1022         case Nbnxm::AtomLocality::Count:
1023             GMX_ASSERT(false, "Count is invalid locality specifier");
1024             break;
1025     }
1026 }
1027
1028 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1029 void nbnxn_atomdata_copy_x_to_nbat_x(const Nbnxm::GridSet     &gridSet,
1030                                      const Nbnxm::AtomLocality locality,
1031                                      bool                      fillLocal,
1032                                      const rvec               *coordinates,
1033                                      nbnxn_atomdata_t         *nbat)
1034 {
1035
1036     int gridBegin = 0;
1037     int gridEnd   = 0;
1038     getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
1039
1040     if (fillLocal)
1041     {
1042         nbat->natoms_local = gridSet.grids()[0].atomIndexEnd();
1043     }
1044
1045     const int nth = gmx_omp_nthreads_get(emntPairsearch);
1046 #pragma omp parallel for num_threads(nth) schedule(static)
1047     for (int th = 0; th < nth; th++)
1048     {
1049         try
1050         {
1051             for (int g = gridBegin; g < gridEnd; g++)
1052             {
1053                 const Nbnxm::Grid  &grid       = gridSet.grids()[g];
1054                 const int           numCellsXY = grid.numColumns();
1055
1056                 const int           cxy0 = (numCellsXY* th      + nth - 1)/nth;
1057                 const int           cxy1 = (numCellsXY*(th + 1) + nth - 1)/nth;
1058
1059                 for (int cxy = cxy0; cxy < cxy1; cxy++)
1060                 {
1061                     const int na  = grid.numAtomsInColumn(cxy);
1062                     const int ash = grid.firstAtomInColumn(cxy);
1063
1064                     int       na_fill;
1065                     if (g == 0 && fillLocal)
1066                     {
1067                         na_fill = grid.paddedNumAtomsInColumn(cxy);
1068                     }
1069                     else
1070                     {
1071                         /* We fill only the real particle locations.
1072                          * We assume the filling entries at the end have been
1073                          * properly set before during pair-list generation.
1074                          */
1075                         na_fill = na;
1076                     }
1077                     copy_rvec_to_nbat_real(gridSet.atomIndices().data() + ash,
1078                                            na, na_fill, coordinates,
1079                                            nbat->XFormat, nbat->x().data(), ash);
1080                 }
1081             }
1082         }
1083         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1084     }
1085 }
1086
1087 void nbnxn_atomdata_copy_x_to_gpu(const Nbnxm::GridSet     &gridSet,
1088                                   const Nbnxm::AtomLocality locality,
1089                                   bool                      fillLocal,
1090                                   nbnxn_atomdata_t         *nbat,
1091                                   gmx_nbnxn_gpu_t          *gpu_nbv,
1092                                   const rvec               *coordinatesHost)
1093 {
1094     int gridBegin = 0;
1095     int gridEnd   = 0;
1096     getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
1097
1098     if (fillLocal)
1099     {
1100         nbat->natoms_local = gridSet.grids()[0].atomIndexEnd();
1101     }
1102
1103     for (int g = gridBegin; g < gridEnd; g++)
1104     {
1105         nbnxn_gpu_copy_x_to_gpu(gridSet.grids()[g],
1106                                 fillLocal && g == 0,
1107                                 gpu_nbv,
1108                                 locality,
1109                                 coordinatesHost,
1110                                 g,
1111                                 gridSet.numColumnsMax());
1112     }
1113 }
1114
1115 DeviceBuffer<float> nbnxn_atomdata_get_x_gpu(gmx_nbnxn_gpu_t *gpu_nbv)
1116 {
1117     return Nbnxm::nbnxn_gpu_get_x_gpu(gpu_nbv);
1118 }
1119
1120 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t on the GPU*/
1121 void nbnxn_atomdata_x_to_nbat_x_gpu(const Nbnxm::GridSet     &gridSet,
1122                                     const Nbnxm::AtomLocality locality,
1123                                     bool                      fillLocal,
1124                                     gmx_nbnxn_gpu_t          *gpu_nbv,
1125                                     DeviceBuffer<float>       coordinatesDevice)
1126 {
1127
1128     int gridBegin = 0;
1129     int gridEnd   = 0;
1130     getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
1131
1132     for (int g = gridBegin; g < gridEnd; g++)
1133     {
1134         nbnxn_gpu_x_to_nbat_x(gridSet.grids()[g],
1135                               fillLocal && g == 0,
1136                               gpu_nbv,
1137                               coordinatesDevice,
1138                               locality,
1139                               g,
1140                               gridSet.numColumnsMax());
1141     }
1142 }
1143
1144 static void
1145 nbnxn_atomdata_clear_reals(gmx::ArrayRef<real> dest,
1146                            int i0, int i1)
1147 {
1148     for (int i = i0; i < i1; i++)
1149     {
1150         dest[i] = 0;
1151     }
1152 }
1153
1154 gmx_unused static void
1155 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1156                             gmx_bool bDestSet,
1157                             const real ** gmx_restrict src,
1158                             int nsrc,
1159                             int i0, int i1)
1160 {
1161     if (bDestSet)
1162     {
1163         /* The destination buffer contains data, add to it */
1164         for (int i = i0; i < i1; i++)
1165         {
1166             for (int s = 0; s < nsrc; s++)
1167             {
1168                 dest[i] += src[s][i];
1169             }
1170         }
1171     }
1172     else
1173     {
1174         /* The destination buffer is unitialized, set it first */
1175         for (int i = i0; i < i1; i++)
1176         {
1177             dest[i] = src[0][i];
1178             for (int s = 1; s < nsrc; s++)
1179             {
1180                 dest[i] += src[s][i];
1181             }
1182         }
1183     }
1184 }
1185
1186 gmx_unused static void
1187 nbnxn_atomdata_reduce_reals_simd(real gmx_unused * gmx_restrict dest,
1188                                  gmx_bool gmx_unused bDestSet,
1189                                  const gmx_unused real ** gmx_restrict src,
1190                                  int gmx_unused nsrc,
1191                                  int gmx_unused i0, int gmx_unused i1)
1192 {
1193 #if GMX_SIMD
1194 /* The SIMD width here is actually independent of that in the kernels,
1195  * but we use the same width for simplicity (usually optimal anyhow).
1196  */
1197     SimdReal dest_SSE, src_SSE;
1198
1199     if (bDestSet)
1200     {
1201         for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1202         {
1203             dest_SSE = load<SimdReal>(dest+i);
1204             for (int s = 0; s < nsrc; s++)
1205             {
1206                 src_SSE  = load<SimdReal>(src[s]+i);
1207                 dest_SSE = dest_SSE + src_SSE;
1208             }
1209             store(dest+i, dest_SSE);
1210         }
1211     }
1212     else
1213     {
1214         for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1215         {
1216             dest_SSE = load<SimdReal>(src[0]+i);
1217             for (int s = 1; s < nsrc; s++)
1218             {
1219                 src_SSE  = load<SimdReal>(src[s]+i);
1220                 dest_SSE = dest_SSE + src_SSE;
1221             }
1222             store(dest+i, dest_SSE);
1223         }
1224     }
1225 #endif
1226 }
1227
1228 /* Add part of the force array(s) from nbnxn_atomdata_t to f
1229  *
1230  * Note: Adding restrict to f makes this function 50% slower with gcc 7.3
1231  */
1232 static void
1233 nbnxn_atomdata_add_nbat_f_to_f_part(const Nbnxm::GridSet          &gridSet,
1234                                     const nbnxn_atomdata_t        &nbat,
1235                                     const nbnxn_atomdata_output_t &out,
1236                                     const int                      a0,
1237                                     const int                      a1,
1238                                     rvec *                         f)
1239 {
1240     gmx::ArrayRef<const int>  cell = gridSet.cells();
1241     // Note: Using ArrayRef instead makes this code 25% slower with gcc 7.3
1242     const real               *fnb  = out.f.data();
1243
1244     /* Loop over all columns and copy and fill */
1245     switch (nbat.FFormat)
1246     {
1247         case nbatXYZ:
1248         case nbatXYZQ:
1249             for (int a = a0; a < a1; a++)
1250             {
1251                 int i = cell[a]*nbat.fstride;
1252
1253                 f[a][XX] += fnb[i];
1254                 f[a][YY] += fnb[i + 1];
1255                 f[a][ZZ] += fnb[i + 2];
1256             }
1257             break;
1258         case nbatX4:
1259             for (int a = a0; a < a1; a++)
1260             {
1261                 int i = atom_to_x_index<c_packX4>(cell[a]);
1262
1263                 f[a][XX] += fnb[i + XX*c_packX4];
1264                 f[a][YY] += fnb[i + YY*c_packX4];
1265                 f[a][ZZ] += fnb[i + ZZ*c_packX4];
1266             }
1267             break;
1268         case nbatX8:
1269             for (int a = a0; a < a1; a++)
1270             {
1271                 int i = atom_to_x_index<c_packX8>(cell[a]);
1272
1273                 f[a][XX] += fnb[i + XX*c_packX8];
1274                 f[a][YY] += fnb[i + YY*c_packX8];
1275                 f[a][ZZ] += fnb[i + ZZ*c_packX8];
1276             }
1277             break;
1278         default:
1279             gmx_incons("Unsupported nbnxn_atomdata_t format");
1280     }
1281 }
1282
1283 static inline unsigned char reverse_bits(unsigned char b)
1284 {
1285     /* http://graphics.stanford.edu/~seander/bithacks.html#ReverseByteWith64BitsDiv */
1286     return (b * 0x0202020202ULL & 0x010884422010ULL) % 1023;
1287 }
1288
1289 static void nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbnxn_atomdata_t *nbat,
1290                                                       int               nth)
1291 {
1292     const nbnxn_buffer_flags_t *flags = &nbat->buffer_flags;
1293
1294     int                         next_pow2 = 1<<(gmx::log2I(nth-1)+1);
1295
1296     const int                   numOutputBuffers = nbat->out.size();
1297     GMX_ASSERT(numOutputBuffers == nth,
1298                "tree-reduce currently only works for numOutputBuffers==nth");
1299
1300     memset(nbat->syncStep, 0, sizeof(*(nbat->syncStep))*nth);
1301
1302 #pragma omp parallel num_threads(nth)
1303     {
1304         try
1305         {
1306             int   b0, b1, b;
1307             int   i0, i1;
1308             int   group_size, th;
1309
1310             th = gmx_omp_get_thread_num();
1311
1312             for (group_size = 2; group_size < 2*next_pow2; group_size *= 2)
1313             {
1314                 int index[2], group_pos, partner_pos, wu;
1315                 int partner_th = th ^ (group_size/2);
1316
1317                 if (group_size > 2)
1318                 {
1319 #ifdef TMPI_ATOMICS
1320                     /* wait on partner thread - replaces full barrier */
1321                     int sync_th, sync_group_size;
1322
1323                     tMPI_Atomic_memory_barrier();                         /* gurantee data is saved before marking work as done */
1324                     tMPI_Atomic_set(&(nbat->syncStep[th]), group_size/2); /* mark previous step as completed */
1325
1326                     /* find thread to sync with. Equal to partner_th unless nth is not a power of two. */
1327                     for (sync_th = partner_th, sync_group_size = group_size; sync_th >= nth && sync_group_size > 2; sync_group_size /= 2)
1328                     {
1329                         sync_th &= ~(sync_group_size/4);
1330                     }
1331                     if (sync_th < nth) /* otherwise nothing to sync index[1] will be >=nout */
1332                     {
1333                         /* wait on the thread which computed input data in previous step */
1334                         while (tMPI_Atomic_get(static_cast<volatile tMPI_Atomic_t*>(&(nbat->syncStep[sync_th]))) < group_size/2)
1335                         {
1336                             gmx_pause();
1337                         }
1338                         /* guarantee that no later load happens before wait loop is finisehd */
1339                         tMPI_Atomic_memory_barrier();
1340                     }
1341 #else               /* TMPI_ATOMICS */
1342 #pragma omp barrier
1343 #endif
1344                 }
1345
1346                 /* Calculate buffers to sum (result goes into first buffer) */
1347                 group_pos = th % group_size;
1348                 index[0]  = th - group_pos;
1349                 index[1]  = index[0] + group_size/2;
1350
1351                 /* If no second buffer, nothing to do */
1352                 if (index[1] >= numOutputBuffers && group_size > 2)
1353                 {
1354                     continue;
1355                 }
1356
1357 #if NBNXN_BUFFERFLAG_MAX_THREADS > 256
1358 #error reverse_bits assumes max 256 threads
1359 #endif
1360                 /* Position is permuted so that one of the 2 vectors being added was computed on the same thread in the previous step.
1361                    This improves locality and enables to sync with just a single thread between steps (=the levels in the btree).
1362                    The permutation which allows this corresponds to reversing the bits of the group position.
1363                  */
1364                 group_pos = reverse_bits(group_pos)/(256/group_size);
1365
1366                 partner_pos = group_pos ^ 1;
1367
1368                 /* loop over two work-units (own and partner) */
1369                 for (wu = 0; wu < 2; wu++)
1370                 {
1371                     if (wu == 1)
1372                     {
1373                         if (partner_th < nth)
1374                         {
1375                             break; /* partner exists we don't have to do his work */
1376                         }
1377                         else
1378                         {
1379                             group_pos = partner_pos;
1380                         }
1381                     }
1382
1383                     /* Calculate the cell-block range for our thread */
1384                     b0 = (flags->nflag* group_pos   )/group_size;
1385                     b1 = (flags->nflag*(group_pos+1))/group_size;
1386
1387                     for (b = b0; b < b1; b++)
1388                     {
1389                         i0 =  b   *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1390                         i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1391
1392                         if (bitmask_is_set(flags->flag[b], index[1]) || group_size > 2)
1393                         {
1394                             const real *fIndex1 = nbat->out[index[1]].f.data();
1395 #if GMX_SIMD
1396                             nbnxn_atomdata_reduce_reals_simd
1397 #else
1398                             nbnxn_atomdata_reduce_reals
1399 #endif
1400                                 (nbat->out[index[0]].f.data(),
1401                                 bitmask_is_set(flags->flag[b], index[0]) || group_size > 2,
1402                                 &fIndex1, 1, i0, i1);
1403
1404                         }
1405                         else if (!bitmask_is_set(flags->flag[b], index[0]))
1406                         {
1407                             nbnxn_atomdata_clear_reals(nbat->out[index[0]].f,
1408                                                        i0, i1);
1409                         }
1410                     }
1411                 }
1412             }
1413         }
1414         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1415     }
1416 }
1417
1418
1419 static void nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbnxn_atomdata_t *nbat,
1420                                                      int               nth)
1421 {
1422 #pragma omp parallel for num_threads(nth) schedule(static)
1423     for (int th = 0; th < nth; th++)
1424     {
1425         try
1426         {
1427             const nbnxn_buffer_flags_t *flags;
1428             int                         nfptr;
1429             const real                 *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1430
1431             flags = &nbat->buffer_flags;
1432
1433             /* Calculate the cell-block range for our thread */
1434             int b0 = (flags->nflag* th   )/nth;
1435             int b1 = (flags->nflag*(th+1))/nth;
1436
1437             for (int b = b0; b < b1; b++)
1438             {
1439                 int i0 =  b   *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1440                 int i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1441
1442                 nfptr = 0;
1443                 for (gmx::index out = 1; out < gmx::ssize(nbat->out); out++)
1444                 {
1445                     if (bitmask_is_set(flags->flag[b], out))
1446                     {
1447                         fptr[nfptr++] = nbat->out[out].f.data();
1448                     }
1449                 }
1450                 if (nfptr > 0)
1451                 {
1452 #if GMX_SIMD
1453                     nbnxn_atomdata_reduce_reals_simd
1454 #else
1455                     nbnxn_atomdata_reduce_reals
1456 #endif
1457                         (nbat->out[0].f.data(),
1458                         bitmask_is_set(flags->flag[b], 0),
1459                         fptr, nfptr,
1460                         i0, i1);
1461                 }
1462                 else if (!bitmask_is_set(flags->flag[b], 0))
1463                 {
1464                     nbnxn_atomdata_clear_reals(nbat->out[0].f,
1465                                                i0, i1);
1466                 }
1467             }
1468         }
1469         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1470     }
1471 }
1472
1473
1474 /* Add the force array(s) from nbnxn_atomdata_t to f */
1475 void reduceForces(nbnxn_atomdata_t                *nbat,
1476                   const Nbnxm::AtomLocality        locality,
1477                   const Nbnxm::GridSet            &gridSet,
1478                   rvec                            *f)
1479 {
1480     int a0 = 0;
1481     int na = 0;
1482
1483     nbnxn_get_atom_range(locality, gridSet, &a0, &na);
1484
1485     if (na == 0)
1486     {
1487         /* The are no atoms for this reduction, avoid some overhead */
1488         return;
1489     }
1490
1491     int nth = gmx_omp_nthreads_get(emntNonbonded);
1492
1493     if (nbat->out.size() > 1)
1494     {
1495         if (locality != Nbnxm::AtomLocality::All)
1496         {
1497             gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1498         }
1499
1500         /* Reduce the force thread output buffers into buffer 0, before adding
1501          * them to the, differently ordered, "real" force buffer.
1502          */
1503         if (nbat->bUseTreeReduce)
1504         {
1505             nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbat, nth);
1506         }
1507         else
1508         {
1509             nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbat, nth);
1510         }
1511     }
1512 #pragma omp parallel for num_threads(nth) schedule(static)
1513     for (int th = 0; th < nth; th++)
1514     {
1515         try
1516         {
1517             nbnxn_atomdata_add_nbat_f_to_f_part(gridSet, *nbat,
1518                                                 nbat->out[0],
1519                                                 a0 + ((th + 0)*na)/nth,
1520                                                 a0 + ((th + 1)*na)/nth,
1521                                                 f);
1522         }
1523         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1524     }
1525 }
1526
1527 /* Add the force array(s) from nbnxn_atomdata_t to f */
1528 void reduceForcesGpu(const Nbnxm::AtomLocality        locality,
1529                      DeviceBuffer<float>              totalForcesDevice,
1530                      const Nbnxm::GridSet            &gridSet,
1531                      void                            *pmeForcesDevice,
1532                      GpuEventSynchronizer            *pmeForcesReady,
1533                      gmx_nbnxn_gpu_t                 *gpu_nbv,
1534                      bool                             useGpuFPmeReduction,
1535                      bool                             accumulateForce)
1536 {
1537     int atomsStart = 0;
1538     int numAtoms   = 0;
1539
1540     nbnxn_get_atom_range(locality, gridSet, &atomsStart, &numAtoms);
1541
1542     if (numAtoms == 0)
1543     {
1544         /* The are no atoms for this reduction, avoid some overhead */
1545         return;
1546     }
1547
1548     Nbnxm::nbnxn_gpu_add_nbat_f_to_f(locality,
1549                                      totalForcesDevice,
1550                                      gpu_nbv,
1551                                      pmeForcesDevice,
1552                                      pmeForcesReady,
1553                                      atomsStart, numAtoms,
1554                                      useGpuFPmeReduction,
1555                                      accumulateForce);
1556 }
1557
1558 DeviceBuffer<float> nbnxn_atomdata_get_f_gpu(gmx_nbnxn_gpu_t *gpu_nbv)
1559 {
1560     return Nbnxm::nbnxn_gpu_get_f_gpu(gpu_nbv);
1561 }
1562
1563 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t   &nbat,
1564                                               gmx::ArrayRef<gmx::RVec>  fshift)
1565 {
1566     gmx::ArrayRef<const nbnxn_atomdata_output_t> outputBuffers = nbat.out;
1567
1568     for (int s = 0; s < SHIFTS; s++)
1569     {
1570         rvec sum;
1571         clear_rvec(sum);
1572         for (const nbnxn_atomdata_output_t &out : outputBuffers)
1573         {
1574             sum[XX] += out.fshift[s*DIM+XX];
1575             sum[YY] += out.fshift[s*DIM+YY];
1576             sum[ZZ] += out.fshift[s*DIM+ZZ];
1577         }
1578         fshift[s] += sum;
1579     }
1580 }
1581
1582 void nbnxn_get_atom_range(const Nbnxm::AtomLocality        atomLocality,
1583                           const Nbnxm::GridSet            &gridSet,
1584                           int                             *atomStart,
1585                           int                             *nAtoms)
1586 {
1587
1588     switch (atomLocality)
1589     {
1590         case Nbnxm::AtomLocality::All:
1591             *atomStart = 0;
1592             *nAtoms    = gridSet.numRealAtomsTotal();
1593             break;
1594         case Nbnxm::AtomLocality::Local:
1595             *atomStart = 0;
1596             *nAtoms    = gridSet.numRealAtomsLocal();
1597             break;
1598         case Nbnxm::AtomLocality::NonLocal:
1599             *atomStart = gridSet.numRealAtomsLocal();
1600             *nAtoms    = gridSet.numRealAtomsTotal() - gridSet.numRealAtomsLocal();
1601             break;
1602         case Nbnxm::AtomLocality::Count:
1603             GMX_ASSERT(false, "Count is invalid locality specifier");
1604             break;
1605     }
1606 }