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