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