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