Remove t_mdatoms from Nbnxm setAtomProperties()
[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.flag        = nullptr;
679     nbat->buffer_flags.flag_nalloc = 0;
680
681     const int nth = gmx_omp_nthreads_get(emntNonbonded);
682
683     const char* ptr = getenv("GMX_USE_TREEREDUCE");
684     if (ptr != nullptr)
685     {
686         nbat->bUseTreeReduce = (strtol(ptr, nullptr, 10) != 0);
687     }
688 #if defined __MIC__
689     else if (nth > 8) /*on the CPU we currently don't benefit even at 32*/
690     {
691         nbat->bUseTreeReduce = 1;
692     }
693 #endif
694     else
695     {
696         nbat->bUseTreeReduce = false;
697     }
698     if (nbat->bUseTreeReduce)
699     {
700         GMX_LOG(mdlog.info).asParagraph().appendText("Using tree force reduction");
701
702         nbat->syncStep = new tMPI_Atomic[nth];
703     }
704 }
705
706 template<int packSize>
707 static void copy_lj_to_nbat_lj_comb(gmx::ArrayRef<const real> ljparam_type, const int* type, int na, real* ljparam_at)
708 {
709     /* The LJ params follow the combination rule:
710      * copy the params for the type array to the atom array.
711      */
712     for (int is = 0; is < na; is += packSize)
713     {
714         for (int k = 0; k < packSize; k++)
715         {
716             int i                             = is + k;
717             ljparam_at[is * 2 + k]            = ljparam_type[type[i] * 2];
718             ljparam_at[is * 2 + packSize + k] = ljparam_type[type[i] * 2 + 1];
719         }
720     }
721 }
722
723 /* Sets the atom type in nbnxn_atomdata_t */
724 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t::Params* params,
725                                          const Nbnxm::GridSet&     gridSet,
726                                          ArrayRef<const int>       atomTypes)
727 {
728     params->type.resize(gridSet.numGridAtomsTotal());
729
730     for (const Nbnxm::Grid& grid : gridSet.grids())
731     {
732         /* Loop over all columns and copy and fill */
733         for (int i = 0; i < grid.numColumns(); i++)
734         {
735             const int numAtoms   = grid.paddedNumAtomsInColumn(i);
736             const int atomOffset = grid.firstAtomInColumn(i);
737
738             copy_int_to_nbat_int(gridSet.atomIndices().data() + atomOffset,
739                                  grid.numAtomsInColumn(i), numAtoms, atomTypes.data(),
740                                  params->numTypes - 1, params->type.data() + atomOffset);
741         }
742     }
743 }
744
745 /* Sets the LJ combination rule parameters in nbnxn_atomdata_t */
746 static void nbnxn_atomdata_set_ljcombparams(nbnxn_atomdata_t::Params* params,
747                                             const int                 XFormat,
748                                             const Nbnxm::GridSet&     gridSet)
749 {
750     params->lj_comb.resize(gridSet.numGridAtomsTotal() * 2);
751
752     if (params->comb_rule != ljcrNONE)
753     {
754         for (const Nbnxm::Grid& grid : gridSet.grids())
755         {
756             /* Loop over all columns and copy and fill */
757             for (int i = 0; i < grid.numColumns(); i++)
758             {
759                 const int numAtoms   = grid.paddedNumAtomsInColumn(i);
760                 const int atomOffset = grid.firstAtomInColumn(i);
761
762                 if (XFormat == nbatX4)
763                 {
764                     copy_lj_to_nbat_lj_comb<c_packX4>(params->nbfp_comb,
765                                                       params->type.data() + atomOffset, numAtoms,
766                                                       params->lj_comb.data() + atomOffset * 2);
767                 }
768                 else if (XFormat == nbatX8)
769                 {
770                     copy_lj_to_nbat_lj_comb<c_packX8>(params->nbfp_comb,
771                                                       params->type.data() + atomOffset, numAtoms,
772                                                       params->lj_comb.data() + atomOffset * 2);
773                 }
774                 else if (XFormat == nbatXYZQ)
775                 {
776                     copy_lj_to_nbat_lj_comb<1>(params->nbfp_comb, params->type.data() + atomOffset,
777                                                numAtoms, params->lj_comb.data() + atomOffset * 2);
778                 }
779             }
780         }
781     }
782 }
783
784 /* Sets the charges in nbnxn_atomdata_t *nbat */
785 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t*     nbat,
786                                        const Nbnxm::GridSet& gridSet,
787                                        ArrayRef<const real>  charges)
788 {
789     if (nbat->XFormat != nbatXYZQ)
790     {
791         nbat->paramsDeprecated().q.resize(nbat->numAtoms());
792     }
793
794     for (const Nbnxm::Grid& grid : gridSet.grids())
795     {
796         /* Loop over all columns and copy and fill */
797         for (int cxy = 0; cxy < grid.numColumns(); cxy++)
798         {
799             const int atomOffset     = grid.firstAtomInColumn(cxy);
800             const int numAtoms       = grid.numAtomsInColumn(cxy);
801             const int paddedNumAtoms = grid.paddedNumAtomsInColumn(cxy);
802
803             if (nbat->XFormat == nbatXYZQ)
804             {
805                 real* q = nbat->x().data() + atomOffset * STRIDE_XYZQ + ZZ + 1;
806                 int   i;
807                 for (i = 0; i < numAtoms; i++)
808                 {
809                     *q = charges[gridSet.atomIndices()[atomOffset + i]];
810                     q += STRIDE_XYZQ;
811                 }
812                 /* Complete the partially filled last cell with zeros */
813                 for (; i < paddedNumAtoms; i++)
814                 {
815                     *q = 0;
816                     q += STRIDE_XYZQ;
817                 }
818             }
819             else
820             {
821                 real* q = nbat->paramsDeprecated().q.data() + atomOffset;
822                 int   i;
823                 for (i = 0; i < numAtoms; i++)
824                 {
825                     *q = charges[gridSet.atomIndices()[atomOffset + i]];
826                     q++;
827                 }
828                 /* Complete the partially filled last cell with zeros */
829                 for (; i < paddedNumAtoms; i++)
830                 {
831                     *q = 0;
832                     q++;
833                 }
834             }
835         }
836     }
837 }
838
839 /* Set the charges of perturbed atoms in nbnxn_atomdata_t to 0.
840  * This is to automatically remove the RF/PME self term in the nbnxn kernels.
841  * Part of the zero interactions are still calculated in the normal kernels.
842  * All perturbed interactions are calculated in the free energy kernel,
843  * using the original charge and LJ data, not nbnxn_atomdata_t.
844  */
845 static void nbnxn_atomdata_mask_fep(nbnxn_atomdata_t* nbat, const Nbnxm::GridSet& gridSet)
846 {
847     nbnxn_atomdata_t::Params& params = nbat->paramsDeprecated();
848     real*                     q;
849     int                       stride_q;
850
851     if (nbat->XFormat == nbatXYZQ)
852     {
853         q        = nbat->x().data() + ZZ + 1;
854         stride_q = STRIDE_XYZQ;
855     }
856     else
857     {
858         q        = params.q.data();
859         stride_q = 1;
860     }
861
862     for (const Nbnxm::Grid& grid : gridSet.grids())
863     {
864         int nsubc;
865         if (grid.geometry().isSimple)
866         {
867             nsubc = 1;
868         }
869         else
870         {
871             nsubc = c_gpuNumClusterPerCell;
872         }
873
874         int c_offset = grid.firstAtomInColumn(0);
875
876         /* Loop over all columns and copy and fill */
877         for (int c = 0; c < grid.numCells() * nsubc; c++)
878         {
879             /* Does this cluster contain perturbed particles? */
880             if (grid.clusterIsPerturbed(c))
881             {
882                 const int numAtomsPerCluster = grid.geometry().numAtomsICluster;
883                 for (int i = 0; i < numAtomsPerCluster; i++)
884                 {
885                     /* Is this a perturbed particle? */
886                     if (grid.atomIsPerturbed(c, i))
887                     {
888                         int ind = c_offset + c * numAtomsPerCluster + i;
889                         /* Set atom type and charge to non-interacting */
890                         params.type[ind]  = params.numTypes - 1;
891                         q[ind * stride_q] = 0;
892                     }
893                 }
894             }
895         }
896     }
897 }
898
899 /* Copies the energy group indices to a reordered and packed array */
900 static void
901 copy_egp_to_nbat_egps(const int* a, int na, int na_round, int na_c, int bit_shift, const int* in, int* innb)
902 {
903     int i;
904     int comb;
905
906     int j = 0;
907     for (i = 0; i < na; i += na_c)
908     {
909         /* Store na_c energy group numbers into one int */
910         comb = 0;
911         for (int sa = 0; sa < na_c; sa++)
912         {
913             int at = a[i + sa];
914             if (at >= 0)
915             {
916                 comb |= (GET_CGINFO_GID(in[at]) << (sa * bit_shift));
917             }
918         }
919         innb[j++] = comb;
920     }
921     /* Complete the partially filled last cell with fill */
922     for (; i < na_round; i += na_c)
923     {
924         innb[j++] = 0;
925     }
926 }
927
928 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
929 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t::Params* params,
930                                             const Nbnxm::GridSet&     gridSet,
931                                             ArrayRef<const int>       atomInfo)
932 {
933     if (params->nenergrp == 1)
934     {
935         return;
936     }
937
938     params->energrp.resize(gridSet.numGridAtomsTotal());
939
940     for (const Nbnxm::Grid& grid : gridSet.grids())
941     {
942         /* Loop over all columns and copy and fill */
943         for (int i = 0; i < grid.numColumns(); i++)
944         {
945             const int numAtoms   = grid.paddedNumAtomsInColumn(i);
946             const int atomOffset = grid.firstAtomInColumn(i);
947
948             copy_egp_to_nbat_egps(gridSet.atomIndices().data() + atomOffset, grid.numAtomsInColumn(i),
949                                   numAtoms, c_nbnxnCpuIClusterSize, params->neg_2log, atomInfo.data(),
950                                   params->energrp.data() + grid.atomToCluster(atomOffset));
951         }
952     }
953 }
954
955 /* Sets all required atom parameter data in nbnxn_atomdata_t */
956 void nbnxn_atomdata_set(nbnxn_atomdata_t*     nbat,
957                         const Nbnxm::GridSet& gridSet,
958                         ArrayRef<const int>   atomTypes,
959                         ArrayRef<const real>  atomCharges,
960                         ArrayRef<const int>   atomInfo)
961 {
962     nbnxn_atomdata_t::Params& params = nbat->paramsDeprecated();
963
964     nbnxn_atomdata_set_atomtypes(&params, gridSet, atomTypes);
965
966     nbnxn_atomdata_set_charges(nbat, gridSet, atomCharges);
967
968     if (gridSet.haveFep())
969     {
970         nbnxn_atomdata_mask_fep(nbat, gridSet);
971     }
972
973     /* This must be done after masking types for FEP */
974     nbnxn_atomdata_set_ljcombparams(&params, nbat->XFormat, gridSet);
975
976     nbnxn_atomdata_set_energygroups(&params, gridSet, atomInfo);
977 }
978
979 /* Copies the shift vector array to nbnxn_atomdata_t */
980 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox, rvec* shift_vec, nbnxn_atomdata_t* nbat)
981 {
982     int i;
983
984     nbat->bDynamicBox = bDynamicBox;
985     for (i = 0; i < SHIFTS; i++)
986     {
987         copy_rvec(shift_vec[i], nbat->shift_vec[i]);
988     }
989 }
990
991 // This is slightly different from nbnxn_get_atom_range(...) at the end of the file
992 // TODO: Combine if possible
993 static void getAtomRanges(const Nbnxm::GridSet&   gridSet,
994                           const gmx::AtomLocality locality,
995                           int*                    gridBegin,
996                           int*                    gridEnd)
997 {
998     switch (locality)
999     {
1000         case gmx::AtomLocality::All:
1001             *gridBegin = 0;
1002             *gridEnd   = gridSet.grids().size();
1003             break;
1004         case gmx::AtomLocality::Local:
1005             *gridBegin = 0;
1006             *gridEnd   = 1;
1007             break;
1008         case gmx::AtomLocality::NonLocal:
1009             *gridBegin = 1;
1010             *gridEnd   = gridSet.grids().size();
1011             break;
1012         case gmx::AtomLocality::Count:
1013             GMX_ASSERT(false, "Count is invalid locality specifier");
1014             break;
1015     }
1016 }
1017
1018 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1019 void nbnxn_atomdata_copy_x_to_nbat_x(const Nbnxm::GridSet&   gridSet,
1020                                      const gmx::AtomLocality locality,
1021                                      bool                    fillLocal,
1022                                      const rvec*             coordinates,
1023                                      nbnxn_atomdata_t*       nbat)
1024 {
1025
1026     int gridBegin = 0;
1027     int gridEnd   = 0;
1028     getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
1029
1030     if (fillLocal)
1031     {
1032         nbat->natoms_local = gridSet.grids()[0].atomIndexEnd();
1033     }
1034
1035     const int nth = gmx_omp_nthreads_get(emntPairsearch);
1036 #pragma omp parallel for num_threads(nth) schedule(static)
1037     for (int th = 0; th < nth; th++)
1038     {
1039         try
1040         {
1041             for (int g = gridBegin; g < gridEnd; g++)
1042             {
1043                 const Nbnxm::Grid& grid       = gridSet.grids()[g];
1044                 const int          numCellsXY = grid.numColumns();
1045
1046                 const int cxy0 = (numCellsXY * th + nth - 1) / nth;
1047                 const int cxy1 = (numCellsXY * (th + 1) + nth - 1) / nth;
1048
1049                 for (int cxy = cxy0; cxy < cxy1; cxy++)
1050                 {
1051                     const int na  = grid.numAtomsInColumn(cxy);
1052                     const int ash = grid.firstAtomInColumn(cxy);
1053
1054                     int na_fill;
1055                     if (g == 0 && fillLocal)
1056                     {
1057                         na_fill = grid.paddedNumAtomsInColumn(cxy);
1058                     }
1059                     else
1060                     {
1061                         /* We fill only the real particle locations.
1062                          * We assume the filling entries at the end have been
1063                          * properly set before during pair-list generation.
1064                          */
1065                         na_fill = na;
1066                     }
1067                     copy_rvec_to_nbat_real(gridSet.atomIndices().data() + ash, na, na_fill,
1068                                            coordinates, nbat->XFormat, nbat->x().data(), ash);
1069                 }
1070             }
1071         }
1072         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1073     }
1074 }
1075
1076 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t on the GPU*/
1077 void nbnxn_atomdata_x_to_nbat_x_gpu(const Nbnxm::GridSet&   gridSet,
1078                                     const gmx::AtomLocality locality,
1079                                     bool                    fillLocal,
1080                                     NbnxmGpu*               gpu_nbv,
1081                                     DeviceBuffer<RVec>      d_x,
1082                                     GpuEventSynchronizer*   xReadyOnDevice)
1083 {
1084
1085     int gridBegin = 0;
1086     int gridEnd   = 0;
1087     getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
1088
1089     for (int g = gridBegin; g < gridEnd; g++)
1090     {
1091         nbnxn_gpu_x_to_nbat_x(gridSet.grids()[g], fillLocal && g == 0, gpu_nbv, d_x, xReadyOnDevice,
1092                               locality, g, gridSet.numColumnsMax());
1093     }
1094 }
1095
1096 static void nbnxn_atomdata_clear_reals(gmx::ArrayRef<real> dest, int i0, int i1)
1097 {
1098     for (int i = i0; i < i1; i++)
1099     {
1100         dest[i] = 0;
1101     }
1102 }
1103
1104 gmx_unused static void nbnxn_atomdata_reduce_reals(real* gmx_restrict dest,
1105                                                    gmx_bool           bDestSet,
1106                                                    const real** gmx_restrict src,
1107                                                    int                       nsrc,
1108                                                    int                       i0,
1109                                                    int                       i1)
1110 {
1111     if (bDestSet)
1112     {
1113         /* The destination buffer contains data, add to it */
1114         for (int i = i0; i < i1; i++)
1115         {
1116             for (int s = 0; s < nsrc; s++)
1117             {
1118                 dest[i] += src[s][i];
1119             }
1120         }
1121     }
1122     else
1123     {
1124         /* The destination buffer is unitialized, set it first */
1125         for (int i = i0; i < i1; i++)
1126         {
1127             dest[i] = src[0][i];
1128             for (int s = 1; s < nsrc; s++)
1129             {
1130                 dest[i] += src[s][i];
1131             }
1132         }
1133     }
1134 }
1135
1136 gmx_unused static void nbnxn_atomdata_reduce_reals_simd(real gmx_unused* gmx_restrict dest,
1137                                                         gmx_bool gmx_unused bDestSet,
1138                                                         const gmx_unused real** gmx_restrict src,
1139                                                         int gmx_unused nsrc,
1140                                                         int gmx_unused i0,
1141                                                         int gmx_unused i1)
1142 {
1143 #if GMX_SIMD
1144     /* The SIMD width here is actually independent of that in the kernels,
1145      * but we use the same width for simplicity (usually optimal anyhow).
1146      */
1147     SimdReal dest_SSE, src_SSE;
1148
1149     if (bDestSet)
1150     {
1151         for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1152         {
1153             dest_SSE = load<SimdReal>(dest + i);
1154             for (int s = 0; s < nsrc; s++)
1155             {
1156                 src_SSE  = load<SimdReal>(src[s] + i);
1157                 dest_SSE = dest_SSE + src_SSE;
1158             }
1159             store(dest + i, dest_SSE);
1160         }
1161     }
1162     else
1163     {
1164         for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1165         {
1166             dest_SSE = load<SimdReal>(src[0] + i);
1167             for (int s = 1; s < nsrc; s++)
1168             {
1169                 src_SSE  = load<SimdReal>(src[s] + i);
1170                 dest_SSE = dest_SSE + src_SSE;
1171             }
1172             store(dest + i, dest_SSE);
1173         }
1174     }
1175 #endif
1176 }
1177
1178 /* Add part of the force array(s) from nbnxn_atomdata_t to f
1179  *
1180  * Note: Adding restrict to f makes this function 50% slower with gcc 7.3
1181  */
1182 static void nbnxn_atomdata_add_nbat_f_to_f_part(const Nbnxm::GridSet&          gridSet,
1183                                                 const nbnxn_atomdata_t&        nbat,
1184                                                 const nbnxn_atomdata_output_t& out,
1185                                                 const int                      a0,
1186                                                 const int                      a1,
1187                                                 rvec*                          f)
1188 {
1189     gmx::ArrayRef<const int> cell = gridSet.cells();
1190     // Note: Using ArrayRef instead makes this code 25% slower with gcc 7.3
1191     const real* fnb = out.f.data();
1192
1193     /* Loop over all columns and copy and fill */
1194     switch (nbat.FFormat)
1195     {
1196         case nbatXYZ:
1197         case nbatXYZQ:
1198             for (int a = a0; a < a1; a++)
1199             {
1200                 int i = cell[a] * nbat.fstride;
1201
1202                 f[a][XX] += fnb[i];
1203                 f[a][YY] += fnb[i + 1];
1204                 f[a][ZZ] += fnb[i + 2];
1205             }
1206             break;
1207         case nbatX4:
1208             for (int a = a0; a < a1; a++)
1209             {
1210                 int i = atom_to_x_index<c_packX4>(cell[a]);
1211
1212                 f[a][XX] += fnb[i + XX * c_packX4];
1213                 f[a][YY] += fnb[i + YY * c_packX4];
1214                 f[a][ZZ] += fnb[i + ZZ * c_packX4];
1215             }
1216             break;
1217         case nbatX8:
1218             for (int a = a0; a < a1; a++)
1219             {
1220                 int i = atom_to_x_index<c_packX8>(cell[a]);
1221
1222                 f[a][XX] += fnb[i + XX * c_packX8];
1223                 f[a][YY] += fnb[i + YY * c_packX8];
1224                 f[a][ZZ] += fnb[i + ZZ * c_packX8];
1225             }
1226             break;
1227         default: gmx_incons("Unsupported nbnxn_atomdata_t format");
1228     }
1229 }
1230
1231 static inline unsigned char reverse_bits(unsigned char b)
1232 {
1233     /* http://graphics.stanford.edu/~seander/bithacks.html#ReverseByteWith64BitsDiv */
1234     return (b * 0x0202020202ULL & 0x010884422010ULL) % 1023;
1235 }
1236
1237 static void nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbnxn_atomdata_t* nbat, int nth)
1238 {
1239     const nbnxn_buffer_flags_t* flags = &nbat->buffer_flags;
1240
1241     int next_pow2 = 1 << (gmx::log2I(nth - 1) + 1);
1242
1243     const int numOutputBuffers = nbat->out.size();
1244     GMX_ASSERT(numOutputBuffers == nth,
1245                "tree-reduce currently only works for numOutputBuffers==nth");
1246
1247     memset(nbat->syncStep, 0, sizeof(*(nbat->syncStep)) * nth);
1248
1249 #pragma omp parallel num_threads(nth)
1250     {
1251         try
1252         {
1253             int b0, b1, b;
1254             int i0, i1;
1255             int group_size, th;
1256
1257             th = gmx_omp_get_thread_num();
1258
1259             for (group_size = 2; group_size < 2 * next_pow2; group_size *= 2)
1260             {
1261                 int index[2], group_pos, partner_pos, wu;
1262                 int partner_th = th ^ (group_size / 2);
1263
1264                 if (group_size > 2)
1265                 {
1266 #ifdef TMPI_ATOMICS
1267                     /* wait on partner thread - replaces full barrier */
1268                     int sync_th, sync_group_size;
1269
1270 #    if defined(__clang__) && __clang_major__ >= 8
1271                     // Suppress warnings that the use of memory_barrier may be excessive
1272                     // Only exists beginning with clang-8
1273 #        pragma clang diagnostic push
1274 #        pragma clang diagnostic ignored "-Watomic-implicit-seq-cst"
1275 #    endif
1276
1277                     tMPI_Atomic_memory_barrier(); /* guarantee data is saved before marking work as done */
1278                     tMPI_Atomic_set(&(nbat->syncStep[th]), group_size / 2); /* mark previous step as completed */
1279
1280                     /* find thread to sync with. Equal to partner_th unless nth is not a power of two. */
1281                     for (sync_th = partner_th, sync_group_size = group_size;
1282                          sync_th >= nth && sync_group_size > 2; sync_group_size /= 2)
1283                     {
1284                         sync_th &= ~(sync_group_size / 4);
1285                     }
1286                     if (sync_th < nth) /* otherwise nothing to sync index[1] will be >=nout */
1287                     {
1288                         /* wait on the thread which computed input data in previous step */
1289                         while (tMPI_Atomic_get(static_cast<volatile tMPI_Atomic_t*>(&(nbat->syncStep[sync_th])))
1290                                < group_size / 2)
1291                         {
1292                             gmx_pause();
1293                         }
1294                         /* guarantee that no later load happens before wait loop is finisehd */
1295                         tMPI_Atomic_memory_barrier();
1296                     }
1297 #    if defined(__clang__) && __clang_major__ >= 8
1298 #        pragma clang diagnostic pop
1299 #    endif
1300 #else /* TMPI_ATOMICS */
1301 #    pragma omp barrier
1302 #endif
1303                 }
1304
1305                 /* Calculate buffers to sum (result goes into first buffer) */
1306                 group_pos = th % group_size;
1307                 index[0]  = th - group_pos;
1308                 index[1]  = index[0] + group_size / 2;
1309
1310                 /* If no second buffer, nothing to do */
1311                 if (index[1] >= numOutputBuffers && group_size > 2)
1312                 {
1313                     continue;
1314                 }
1315
1316 #if NBNXN_BUFFERFLAG_MAX_THREADS > 256
1317 #    error reverse_bits assumes max 256 threads
1318 #endif
1319                 /* Position is permuted so that one of the 2 vectors being added was computed on the same thread in the previous step.
1320                    This improves locality and enables to sync with just a single thread between steps (=the levels in the btree).
1321                    The permutation which allows this corresponds to reversing the bits of the group position.
1322                  */
1323                 group_pos = reverse_bits(group_pos) / (256 / group_size);
1324
1325                 partner_pos = group_pos ^ 1;
1326
1327                 /* loop over two work-units (own and partner) */
1328                 for (wu = 0; wu < 2; wu++)
1329                 {
1330                     if (wu == 1)
1331                     {
1332                         if (partner_th < nth)
1333                         {
1334                             break; /* partner exists we don't have to do his work */
1335                         }
1336                         else
1337                         {
1338                             group_pos = partner_pos;
1339                         }
1340                     }
1341
1342                     /* Calculate the cell-block range for our thread */
1343                     b0 = (flags->nflag * group_pos) / group_size;
1344                     b1 = (flags->nflag * (group_pos + 1)) / group_size;
1345
1346                     for (b = b0; b < b1; b++)
1347                     {
1348                         i0 = b * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1349                         i1 = (b + 1) * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1350
1351                         if (bitmask_is_set(flags->flag[b], index[1]) || group_size > 2)
1352                         {
1353                             const real* fIndex1 = nbat->out[index[1]].f.data();
1354 #if GMX_SIMD
1355                             nbnxn_atomdata_reduce_reals_simd
1356 #else
1357                             nbnxn_atomdata_reduce_reals
1358 #endif
1359                                     (nbat->out[index[0]].f.data(),
1360                                      bitmask_is_set(flags->flag[b], index[0]) || group_size > 2,
1361                                      &fIndex1, 1, i0, i1);
1362                         }
1363                         else if (!bitmask_is_set(flags->flag[b], index[0]))
1364                         {
1365                             nbnxn_atomdata_clear_reals(nbat->out[index[0]].f, i0, i1);
1366                         }
1367                     }
1368                 }
1369             }
1370         }
1371         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1372     }
1373 }
1374
1375
1376 static void nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbnxn_atomdata_t* nbat, int nth)
1377 {
1378 #pragma omp parallel for num_threads(nth) schedule(static)
1379     for (int th = 0; th < nth; th++)
1380     {
1381         try
1382         {
1383             const nbnxn_buffer_flags_t* flags;
1384             int                         nfptr;
1385             const real*                 fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1386
1387             flags = &nbat->buffer_flags;
1388
1389             /* Calculate the cell-block range for our thread */
1390             int b0 = (flags->nflag * th) / nth;
1391             int b1 = (flags->nflag * (th + 1)) / nth;
1392
1393             for (int b = b0; b < b1; b++)
1394             {
1395                 int i0 = b * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1396                 int i1 = (b + 1) * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1397
1398                 nfptr = 0;
1399                 for (gmx::index out = 1; out < gmx::ssize(nbat->out); out++)
1400                 {
1401                     if (bitmask_is_set(flags->flag[b], out))
1402                     {
1403                         fptr[nfptr++] = nbat->out[out].f.data();
1404                     }
1405                 }
1406                 if (nfptr > 0)
1407                 {
1408 #if GMX_SIMD
1409                     nbnxn_atomdata_reduce_reals_simd
1410 #else
1411                     nbnxn_atomdata_reduce_reals
1412 #endif
1413                             (nbat->out[0].f.data(), bitmask_is_set(flags->flag[b], 0), fptr, nfptr, i0, i1);
1414                 }
1415                 else if (!bitmask_is_set(flags->flag[b], 0))
1416                 {
1417                     nbnxn_atomdata_clear_reals(nbat->out[0].f, i0, i1);
1418                 }
1419             }
1420         }
1421         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1422     }
1423 }
1424
1425
1426 /* Add the force array(s) from nbnxn_atomdata_t to f */
1427 void reduceForces(nbnxn_atomdata_t* nbat, const gmx::AtomLocality locality, const Nbnxm::GridSet& gridSet, rvec* f)
1428 {
1429     int a0 = 0;
1430     int na = 0;
1431
1432     nbnxn_get_atom_range(locality, gridSet, &a0, &na);
1433
1434     if (na == 0)
1435     {
1436         /* The are no atoms for this reduction, avoid some overhead */
1437         return;
1438     }
1439
1440     int nth = gmx_omp_nthreads_get(emntNonbonded);
1441
1442     if (nbat->out.size() > 1)
1443     {
1444         if (locality != gmx::AtomLocality::All)
1445         {
1446             gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1447         }
1448
1449         /* Reduce the force thread output buffers into buffer 0, before adding
1450          * them to the, differently ordered, "real" force buffer.
1451          */
1452         if (nbat->bUseTreeReduce)
1453         {
1454             nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbat, nth);
1455         }
1456         else
1457         {
1458             nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbat, nth);
1459         }
1460     }
1461 #pragma omp parallel for num_threads(nth) schedule(static)
1462     for (int th = 0; th < nth; th++)
1463     {
1464         try
1465         {
1466             nbnxn_atomdata_add_nbat_f_to_f_part(gridSet, *nbat, nbat->out[0], a0 + ((th + 0) * na) / nth,
1467                                                 a0 + ((th + 1) * na) / nth, f);
1468         }
1469         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1470     }
1471 }
1472
1473 /* Add the force array(s) from nbnxn_atomdata_t to f */
1474 void reduceForcesGpu(const gmx::AtomLocality                    locality,
1475                      DeviceBuffer<RVec>                         totalForcesDevice,
1476                      const Nbnxm::GridSet&                      gridSet,
1477                      void*                                      pmeForcesDevice,
1478                      gmx::ArrayRef<GpuEventSynchronizer* const> dependencyList,
1479                      NbnxmGpu*                                  gpu_nbv,
1480                      bool                                       useGpuFPmeReduction,
1481                      bool                                       accumulateForce)
1482 {
1483     int atomsStart = 0;
1484     int numAtoms   = 0;
1485
1486     nbnxn_get_atom_range(locality, gridSet, &atomsStart, &numAtoms);
1487
1488     if (numAtoms == 0)
1489     {
1490         /* The are no atoms for this reduction, avoid some overhead */
1491         return;
1492     }
1493
1494     Nbnxm::nbnxn_gpu_add_nbat_f_to_f(locality, totalForcesDevice, gpu_nbv, pmeForcesDevice, dependencyList,
1495                                      atomsStart, numAtoms, useGpuFPmeReduction, accumulateForce);
1496 }
1497
1498 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t& nbat, gmx::ArrayRef<gmx::RVec> fshift)
1499 {
1500     gmx::ArrayRef<const nbnxn_atomdata_output_t> outputBuffers = nbat.out;
1501
1502     for (int s = 0; s < SHIFTS; s++)
1503     {
1504         rvec sum;
1505         clear_rvec(sum);
1506         for (const nbnxn_atomdata_output_t& out : outputBuffers)
1507         {
1508             sum[XX] += out.fshift[s * DIM + XX];
1509             sum[YY] += out.fshift[s * DIM + YY];
1510             sum[ZZ] += out.fshift[s * DIM + ZZ];
1511         }
1512         fshift[s] += sum;
1513     }
1514 }
1515
1516 void nbnxn_get_atom_range(const gmx::AtomLocality atomLocality,
1517                           const Nbnxm::GridSet&   gridSet,
1518                           int*                    atomStart,
1519                           int*                    nAtoms)
1520 {
1521
1522     switch (atomLocality)
1523     {
1524         case gmx::AtomLocality::All:
1525             *atomStart = 0;
1526             *nAtoms    = gridSet.numRealAtomsTotal();
1527             break;
1528         case gmx::AtomLocality::Local:
1529             *atomStart = 0;
1530             *nAtoms    = gridSet.numRealAtomsLocal();
1531             break;
1532         case gmx::AtomLocality::NonLocal:
1533             *atomStart = gridSet.numRealAtomsLocal();
1534             *nAtoms    = gridSet.numRealAtomsTotal() - gridSet.numRealAtomsLocal();
1535             break;
1536         case gmx::AtomLocality::Count:
1537             GMX_ASSERT(false, "Count is invalid locality specifier");
1538             break;
1539     }
1540 }