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