F buffer operations in CUDA
[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 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1002 template <bool useGpu>
1003 void nbnxn_atomdata_copy_x_to_nbat_x(const Nbnxm::GridSet     &gridSet,
1004                                      const Nbnxm::AtomLocality locality,
1005                                      gmx_bool                  FillLocal,
1006                                      const rvec               *x,
1007                                      nbnxn_atomdata_t         *nbat,
1008                                      gmx_nbnxn_gpu_t          *gpu_nbv,
1009                                      void                     *xPmeDevicePtr)
1010 {
1011     int gridBegin = 0;
1012     int gridEnd   = 0;
1013
1014     switch (locality)
1015     {
1016         case Nbnxm::AtomLocality::All:
1017             gridBegin = 0;
1018             gridEnd   = gridSet.grids().size();
1019             break;
1020         case Nbnxm::AtomLocality::Local:
1021             gridBegin = 0;
1022             gridEnd   = 1;
1023             break;
1024         case Nbnxm::AtomLocality::NonLocal:
1025             gridBegin = 1;
1026             gridEnd   = gridSet.grids().size();
1027             break;
1028         case Nbnxm::AtomLocality::Count:
1029             GMX_ASSERT(false, "Count is invalid locality specifier");
1030             break;
1031     }
1032
1033     if (FillLocal)
1034     {
1035         nbat->natoms_local = gridSet.grids()[0].atomIndexEnd();
1036     }
1037
1038     if (useGpu)
1039     {
1040         for (int g = gridBegin; g < gridEnd; g++)
1041         {
1042             nbnxn_gpu_x_to_nbat_x(gridSet.grids()[g],
1043                                   FillLocal && g == 0,
1044                                   gpu_nbv,
1045                                   xPmeDevicePtr,
1046                                   locality,
1047                                   x, g, gridSet.numColumnsMax());
1048         }
1049     }
1050     else
1051     {
1052         const int nth = gmx_omp_nthreads_get(emntPairsearch);
1053 #pragma omp parallel for num_threads(nth) schedule(static)
1054         for (int th = 0; th < nth; th++)
1055         {
1056             try
1057             {
1058                 for (int g = gridBegin; g < gridEnd; g++)
1059                 {
1060                     const Nbnxm::Grid  &grid       = gridSet.grids()[g];
1061                     const int           numCellsXY = grid.numColumns();
1062
1063                     const int           cxy0 = (numCellsXY* th      + nth - 1)/nth;
1064                     const int           cxy1 = (numCellsXY*(th + 1) + nth - 1)/nth;
1065
1066                     for (int cxy = cxy0; cxy < cxy1; cxy++)
1067                     {
1068                         const int na  = grid.numAtomsInColumn(cxy);
1069                         const int ash = grid.firstAtomInColumn(cxy);
1070
1071                         int       na_fill;
1072                         if (g == 0 && FillLocal)
1073                         {
1074                             na_fill = grid.paddedNumAtomsInColumn(cxy);
1075                         }
1076                         else
1077                         {
1078                             /* We fill only the real particle locations.
1079                              * We assume the filling entries at the end have been
1080                              * properly set before during pair-list generation.
1081                              */
1082                             na_fill = na;
1083                         }
1084                         copy_rvec_to_nbat_real(gridSet.atomIndices().data() + ash,
1085                                                na, na_fill, x,
1086                                                nbat->XFormat, nbat->x().data(), ash);
1087                     }
1088                 }
1089             }
1090             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1091         }
1092     }
1093 }
1094
1095 template
1096 void nbnxn_atomdata_copy_x_to_nbat_x<true>(const Nbnxm::GridSet &,
1097                                            const Nbnxm::AtomLocality,
1098                                            gmx_bool,
1099                                            const rvec*,
1100                                            nbnxn_atomdata_t *,
1101                                            gmx_nbnxn_gpu_t*,
1102                                            void *);
1103 template
1104 void nbnxn_atomdata_copy_x_to_nbat_x<false>(const Nbnxm::GridSet &,
1105                                             const Nbnxm::AtomLocality,
1106                                             gmx_bool,
1107                                             const rvec*,
1108                                             nbnxn_atomdata_t *,
1109                                             gmx_nbnxn_gpu_t*,
1110                                             void *);
1111
1112 static void
1113 nbnxn_atomdata_clear_reals(gmx::ArrayRef<real> dest,
1114                            int i0, int i1)
1115 {
1116     for (int i = i0; i < i1; i++)
1117     {
1118         dest[i] = 0;
1119     }
1120 }
1121
1122 gmx_unused static void
1123 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1124                             gmx_bool bDestSet,
1125                             const real ** gmx_restrict src,
1126                             int nsrc,
1127                             int i0, int i1)
1128 {
1129     if (bDestSet)
1130     {
1131         /* The destination buffer contains data, add to it */
1132         for (int i = i0; i < i1; i++)
1133         {
1134             for (int s = 0; s < nsrc; s++)
1135             {
1136                 dest[i] += src[s][i];
1137             }
1138         }
1139     }
1140     else
1141     {
1142         /* The destination buffer is unitialized, set it first */
1143         for (int i = i0; i < i1; i++)
1144         {
1145             dest[i] = src[0][i];
1146             for (int s = 1; s < nsrc; s++)
1147             {
1148                 dest[i] += src[s][i];
1149             }
1150         }
1151     }
1152 }
1153
1154 gmx_unused static void
1155 nbnxn_atomdata_reduce_reals_simd(real gmx_unused * gmx_restrict dest,
1156                                  gmx_bool gmx_unused bDestSet,
1157                                  const gmx_unused real ** gmx_restrict src,
1158                                  int gmx_unused nsrc,
1159                                  int gmx_unused i0, int gmx_unused i1)
1160 {
1161 #if GMX_SIMD
1162 /* The SIMD width here is actually independent of that in the kernels,
1163  * but we use the same width for simplicity (usually optimal anyhow).
1164  */
1165     SimdReal dest_SSE, src_SSE;
1166
1167     if (bDestSet)
1168     {
1169         for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1170         {
1171             dest_SSE = load<SimdReal>(dest+i);
1172             for (int s = 0; s < nsrc; s++)
1173             {
1174                 src_SSE  = load<SimdReal>(src[s]+i);
1175                 dest_SSE = dest_SSE + src_SSE;
1176             }
1177             store(dest+i, dest_SSE);
1178         }
1179     }
1180     else
1181     {
1182         for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1183         {
1184             dest_SSE = load<SimdReal>(src[0]+i);
1185             for (int s = 1; s < nsrc; s++)
1186             {
1187                 src_SSE  = load<SimdReal>(src[s]+i);
1188                 dest_SSE = dest_SSE + src_SSE;
1189             }
1190             store(dest+i, dest_SSE);
1191         }
1192     }
1193 #endif
1194 }
1195
1196 /* Add part of the force array(s) from nbnxn_atomdata_t to f
1197  *
1198  * Note: Adding restrict to f makes this function 50% slower with gcc 7.3
1199  */
1200 static void
1201 nbnxn_atomdata_add_nbat_f_to_f_part(const Nbnxm::GridSet          &gridSet,
1202                                     const nbnxn_atomdata_t        &nbat,
1203                                     const nbnxn_atomdata_output_t &out,
1204                                     const int                      a0,
1205                                     const int                      a1,
1206                                     rvec *                         f)
1207 {
1208     gmx::ArrayRef<const int>  cell = gridSet.cells();
1209     // Note: Using ArrayRef instead makes this code 25% slower with gcc 7.3
1210     const real               *fnb  = out.f.data();
1211
1212     /* Loop over all columns and copy and fill */
1213     switch (nbat.FFormat)
1214     {
1215         case nbatXYZ:
1216         case nbatXYZQ:
1217             for (int a = a0; a < a1; a++)
1218             {
1219                 int i = cell[a]*nbat.fstride;
1220
1221                 f[a][XX] += fnb[i];
1222                 f[a][YY] += fnb[i + 1];
1223                 f[a][ZZ] += fnb[i + 2];
1224             }
1225             break;
1226         case nbatX4:
1227             for (int a = a0; a < a1; a++)
1228             {
1229                 int i = atom_to_x_index<c_packX4>(cell[a]);
1230
1231                 f[a][XX] += fnb[i + XX*c_packX4];
1232                 f[a][YY] += fnb[i + YY*c_packX4];
1233                 f[a][ZZ] += fnb[i + ZZ*c_packX4];
1234             }
1235             break;
1236         case nbatX8:
1237             for (int a = a0; a < a1; a++)
1238             {
1239                 int i = atom_to_x_index<c_packX8>(cell[a]);
1240
1241                 f[a][XX] += fnb[i + XX*c_packX8];
1242                 f[a][YY] += fnb[i + YY*c_packX8];
1243                 f[a][ZZ] += fnb[i + ZZ*c_packX8];
1244             }
1245             break;
1246         default:
1247             gmx_incons("Unsupported nbnxn_atomdata_t format");
1248     }
1249 }
1250
1251 static inline unsigned char reverse_bits(unsigned char b)
1252 {
1253     /* http://graphics.stanford.edu/~seander/bithacks.html#ReverseByteWith64BitsDiv */
1254     return (b * 0x0202020202ULL & 0x010884422010ULL) % 1023;
1255 }
1256
1257 static void nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbnxn_atomdata_t *nbat,
1258                                                       int               nth)
1259 {
1260     const nbnxn_buffer_flags_t *flags = &nbat->buffer_flags;
1261
1262     int                         next_pow2 = 1<<(gmx::log2I(nth-1)+1);
1263
1264     const int                   numOutputBuffers = nbat->out.size();
1265     GMX_ASSERT(numOutputBuffers == nth,
1266                "tree-reduce currently only works for numOutputBuffers==nth");
1267
1268     memset(nbat->syncStep, 0, sizeof(*(nbat->syncStep))*nth);
1269
1270 #pragma omp parallel num_threads(nth)
1271     {
1272         try
1273         {
1274             int   b0, b1, b;
1275             int   i0, i1;
1276             int   group_size, th;
1277
1278             th = gmx_omp_get_thread_num();
1279
1280             for (group_size = 2; group_size < 2*next_pow2; group_size *= 2)
1281             {
1282                 int index[2], group_pos, partner_pos, wu;
1283                 int partner_th = th ^ (group_size/2);
1284
1285                 if (group_size > 2)
1286                 {
1287 #ifdef TMPI_ATOMICS
1288                     /* wait on partner thread - replaces full barrier */
1289                     int sync_th, sync_group_size;
1290
1291                     tMPI_Atomic_memory_barrier();                         /* gurantee data is saved before marking work as done */
1292                     tMPI_Atomic_set(&(nbat->syncStep[th]), group_size/2); /* mark previous step as completed */
1293
1294                     /* find thread to sync with. Equal to partner_th unless nth is not a power of two. */
1295                     for (sync_th = partner_th, sync_group_size = group_size; sync_th >= nth && sync_group_size > 2; sync_group_size /= 2)
1296                     {
1297                         sync_th &= ~(sync_group_size/4);
1298                     }
1299                     if (sync_th < nth) /* otherwise nothing to sync index[1] will be >=nout */
1300                     {
1301                         /* wait on the thread which computed input data in previous step */
1302                         while (tMPI_Atomic_get(static_cast<volatile tMPI_Atomic_t*>(&(nbat->syncStep[sync_th]))) < group_size/2)
1303                         {
1304                             gmx_pause();
1305                         }
1306                         /* guarantee that no later load happens before wait loop is finisehd */
1307                         tMPI_Atomic_memory_barrier();
1308                     }
1309 #else               /* TMPI_ATOMICS */
1310 #pragma omp barrier
1311 #endif
1312                 }
1313
1314                 /* Calculate buffers to sum (result goes into first buffer) */
1315                 group_pos = th % group_size;
1316                 index[0]  = th - group_pos;
1317                 index[1]  = index[0] + group_size/2;
1318
1319                 /* If no second buffer, nothing to do */
1320                 if (index[1] >= numOutputBuffers && group_size > 2)
1321                 {
1322                     continue;
1323                 }
1324
1325 #if NBNXN_BUFFERFLAG_MAX_THREADS > 256
1326 #error reverse_bits assumes max 256 threads
1327 #endif
1328                 /* Position is permuted so that one of the 2 vectors being added was computed on the same thread in the previous step.
1329                    This improves locality and enables to sync with just a single thread between steps (=the levels in the btree).
1330                    The permutation which allows this corresponds to reversing the bits of the group position.
1331                  */
1332                 group_pos = reverse_bits(group_pos)/(256/group_size);
1333
1334                 partner_pos = group_pos ^ 1;
1335
1336                 /* loop over two work-units (own and partner) */
1337                 for (wu = 0; wu < 2; wu++)
1338                 {
1339                     if (wu == 1)
1340                     {
1341                         if (partner_th < nth)
1342                         {
1343                             break; /* partner exists we don't have to do his work */
1344                         }
1345                         else
1346                         {
1347                             group_pos = partner_pos;
1348                         }
1349                     }
1350
1351                     /* Calculate the cell-block range for our thread */
1352                     b0 = (flags->nflag* group_pos   )/group_size;
1353                     b1 = (flags->nflag*(group_pos+1))/group_size;
1354
1355                     for (b = b0; b < b1; b++)
1356                     {
1357                         i0 =  b   *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1358                         i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1359
1360                         if (bitmask_is_set(flags->flag[b], index[1]) || group_size > 2)
1361                         {
1362                             const real *fIndex1 = nbat->out[index[1]].f.data();
1363 #if GMX_SIMD
1364                             nbnxn_atomdata_reduce_reals_simd
1365 #else
1366                             nbnxn_atomdata_reduce_reals
1367 #endif
1368                                 (nbat->out[index[0]].f.data(),
1369                                 bitmask_is_set(flags->flag[b], index[0]) || group_size > 2,
1370                                 &fIndex1, 1, i0, i1);
1371
1372                         }
1373                         else if (!bitmask_is_set(flags->flag[b], index[0]))
1374                         {
1375                             nbnxn_atomdata_clear_reals(nbat->out[index[0]].f,
1376                                                        i0, i1);
1377                         }
1378                     }
1379                 }
1380             }
1381         }
1382         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1383     }
1384 }
1385
1386
1387 static void nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbnxn_atomdata_t *nbat,
1388                                                      int               nth)
1389 {
1390 #pragma omp parallel for num_threads(nth) schedule(static)
1391     for (int th = 0; th < nth; th++)
1392     {
1393         try
1394         {
1395             const nbnxn_buffer_flags_t *flags;
1396             int                         nfptr;
1397             const real                 *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1398
1399             flags = &nbat->buffer_flags;
1400
1401             /* Calculate the cell-block range for our thread */
1402             int b0 = (flags->nflag* th   )/nth;
1403             int b1 = (flags->nflag*(th+1))/nth;
1404
1405             for (int b = b0; b < b1; b++)
1406             {
1407                 int i0 =  b   *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1408                 int i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1409
1410                 nfptr = 0;
1411                 for (int out = 1; out < gmx::ssize(nbat->out); out++)
1412                 {
1413                     if (bitmask_is_set(flags->flag[b], out))
1414                     {
1415                         fptr[nfptr++] = nbat->out[out].f.data();
1416                     }
1417                 }
1418                 if (nfptr > 0)
1419                 {
1420 #if GMX_SIMD
1421                     nbnxn_atomdata_reduce_reals_simd
1422 #else
1423                     nbnxn_atomdata_reduce_reals
1424 #endif
1425                         (nbat->out[0].f.data(),
1426                         bitmask_is_set(flags->flag[b], 0),
1427                         fptr, nfptr,
1428                         i0, i1);
1429                 }
1430                 else if (!bitmask_is_set(flags->flag[b], 0))
1431                 {
1432                     nbnxn_atomdata_clear_reals(nbat->out[0].f,
1433                                                i0, i1);
1434                 }
1435             }
1436         }
1437         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1438     }
1439 }
1440
1441
1442 /* Add the force array(s) from nbnxn_atomdata_t to f */
1443 template <bool  useGpu>
1444 void reduceForces(nbnxn_atomdata_t                   *nbat,
1445                   const Nbnxm::AtomLocality           locality,
1446                   const Nbnxm::GridSet               &gridSet,
1447                   rvec                               *f,
1448                   gmx_nbnxn_gpu_t                    *gpu_nbv,
1449                   GpuBufferOpsAccumulateForce         accumulateForce)
1450 {
1451     int a0 = 0;
1452     int na = 0;
1453
1454     nbnxn_get_atom_range(locality, gridSet, &a0, &na);
1455
1456     if (na == 0)
1457     {
1458         /* The are no atoms for this reduction, avoid some overhead */
1459         return;
1460     }
1461
1462     if (useGpu)
1463     {
1464         Nbnxm::nbnxn_gpu_add_nbat_f_to_f(locality,
1465                                          gpu_nbv,
1466                                          a0, na,
1467                                          accumulateForce);
1468
1469     }
1470     else
1471     {
1472         int nth = gmx_omp_nthreads_get(emntNonbonded);
1473
1474         if (nbat->out.size() > 1)
1475         {
1476             if (locality != Nbnxm::AtomLocality::All)
1477             {
1478                 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1479             }
1480
1481             /* Reduce the force thread output buffers into buffer 0, before adding
1482              * them to the, differently ordered, "real" force buffer.
1483              */
1484             if (nbat->bUseTreeReduce)
1485             {
1486                 nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbat, nth);
1487             }
1488             else
1489             {
1490                 nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbat, nth);
1491             }
1492         }
1493 #pragma omp parallel for num_threads(nth) schedule(static)
1494         for (int th = 0; th < nth; th++)
1495         {
1496             try
1497             {
1498                 nbnxn_atomdata_add_nbat_f_to_f_part(gridSet, *nbat,
1499                                                     nbat->out[0],
1500                                                     a0 + ((th + 0)*na)/nth,
1501                                                     a0 + ((th + 1)*na)/nth,
1502                                                     f);
1503             }
1504             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1505         }
1506     }
1507 }
1508 template
1509 void reduceForces<true>(nbnxn_atomdata_t             *nbat,
1510                         const Nbnxm::AtomLocality     locality,
1511                         const Nbnxm::GridSet         &gridSet,
1512                         rvec                         *f,
1513                         gmx_nbnxn_gpu_t              *gpu_nbv,
1514                         GpuBufferOpsAccumulateForce   accumulateForce);
1515
1516 template
1517 void reduceForces<false>(nbnxn_atomdata_t             *nbat,
1518                          const Nbnxm::AtomLocality     locality,
1519                          const Nbnxm::GridSet         &gridSet,
1520                          rvec                         *f,
1521                          gmx_nbnxn_gpu_t              *gpu_nbv,
1522                          GpuBufferOpsAccumulateForce   accumulateForce);
1523
1524 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1525                                               rvec                   *fshift)
1526 {
1527     gmx::ArrayRef<const nbnxn_atomdata_output_t> outputBuffers = nbat->out;
1528
1529     for (int s = 0; s < SHIFTS; s++)
1530     {
1531         rvec sum;
1532         clear_rvec(sum);
1533         for (const nbnxn_atomdata_output_t &out : outputBuffers)
1534         {
1535             sum[XX] += out.fshift[s*DIM+XX];
1536             sum[YY] += out.fshift[s*DIM+YY];
1537             sum[ZZ] += out.fshift[s*DIM+ZZ];
1538         }
1539         rvec_inc(fshift[s], sum);
1540     }
1541 }
1542
1543 void nbnxn_get_atom_range(const Nbnxm::AtomLocality        atomLocality,
1544                           const Nbnxm::GridSet            &gridSet,
1545                           int                             *atomStart,
1546                           int                             *nAtoms)
1547 {
1548
1549     switch (atomLocality)
1550     {
1551         case Nbnxm::AtomLocality::All:
1552             *atomStart = 0;
1553             *nAtoms    = gridSet.numRealAtomsTotal();
1554             break;
1555         case Nbnxm::AtomLocality::Local:
1556             *atomStart = 0;
1557             *nAtoms    = gridSet.numRealAtomsLocal();
1558             break;
1559         case Nbnxm::AtomLocality::NonLocal:
1560             *atomStart = gridSet.numRealAtomsLocal();
1561             *nAtoms    = gridSet.numRealAtomsTotal() - gridSet.numRealAtomsLocal();
1562             break;
1563         case Nbnxm::AtomLocality::Count:
1564             GMX_ASSERT(false, "Count is invalid locality specifier");
1565             break;
1566     }
1567 }