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