bce2b142aecf82791e308f0aa93bf12a1af7d0c6
[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(c6 * c6, nbfp[(i * ntype + i) * 2] * nbfp[(j * ntype + j) * 2], tol)
508                         && gmx_within_tol(c12 * c12,
509                                           nbfp[(i * ntype + i) * 2 + 1] * nbfp[(j * ntype + j) * 2 + 1],
510                                           tol);
511
512                 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
513                 c6 /= 6.0;
514                 c12 /= 12.0;
515                 bCombLB =
516                         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(gmx::sixthroot(c12 / c6),
521                                                   0.5 * (params->nbfp_comb[i * 2] + params->nbfp_comb[j * 2]),
522                                                   tol)
523                                 && gmx_within_tol(0.25 * c6 * c6 / c12,
524                                                   std::sqrt(params->nbfp_comb[i * 2 + 1]
525                                                             * params->nbfp_comb[j * 2 + 1]),
526                                                   tol)));
527             }
528             else
529             {
530                 /* Add zero parameters for the additional dummy atom type */
531                 params->nbfp[(i * params->numTypes + j) * 2]     = 0;
532                 params->nbfp[(i * params->numTypes + j) * 2 + 1] = 0;
533             }
534         }
535     }
536     if (debug)
537     {
538         fprintf(debug,
539                 "Combination rules: geometric %s Lorentz-Berthelot %s\n",
540                 gmx::boolToString(bCombGeom),
541                 gmx::boolToString(bCombLB));
542     }
543
544     const bool simple = Nbnxm::kernelTypeUsesSimplePairlist(kernelType);
545
546     switch (enbnxninitcombrule)
547     {
548         case enbnxninitcombruleDETECT:
549             /* We prefer the geometric combination rule,
550              * as that gives a slightly faster kernel than the LB rule.
551              */
552             if (bCombGeom)
553             {
554                 params->ljCombinationRule = LJCombinationRule::Geometric;
555             }
556             else if (bCombLB)
557             {
558                 params->ljCombinationRule = LJCombinationRule::LorentzBerthelot;
559             }
560             else
561             {
562                 params->ljCombinationRule = LJCombinationRule::None;
563
564                 params->nbfp_comb.clear();
565             }
566
567             {
568                 std::string mesg;
569                 if (params->ljCombinationRule == LJCombinationRule::None)
570                 {
571                     mesg = "Using full Lennard-Jones parameter combination matrix";
572                 }
573                 else
574                 {
575                     mesg = gmx::formatString("Using %s Lennard-Jones combination rule",
576                                              enumValueToString(params->ljCombinationRule));
577                 }
578                 GMX_LOG(mdlog.info).asParagraph().appendText(mesg);
579             }
580             break;
581         case enbnxninitcombruleGEOM:
582             params->ljCombinationRule = LJCombinationRule::Geometric;
583             break;
584         case enbnxninitcombruleLB:
585             params->ljCombinationRule = LJCombinationRule::LorentzBerthelot;
586             break;
587         case enbnxninitcombruleNONE:
588             params->ljCombinationRule = LJCombinationRule::None;
589
590             params->nbfp_comb.clear();
591             break;
592         default: gmx_incons("Unknown enbnxninitcombrule");
593     }
594
595     const bool bSIMD = Nbnxm::kernelTypeIsSimd(kernelType);
596
597     set_lj_parameter_data(params, bSIMD);
598
599     params->nenergrp = n_energygroups;
600     if (!simple)
601     {
602         // We now check for energy groups already when starting mdrun
603         GMX_RELEASE_ASSERT(n_energygroups == 1, "GPU kernels do not support energy groups");
604     }
605     /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
606     if (params->nenergrp > 64)
607     {
608         gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
609     }
610     params->neg_2log = 1;
611     while (params->nenergrp > (1 << params->neg_2log))
612     {
613         params->neg_2log++;
614     }
615 }
616
617 /* Initializes an nbnxn_atomdata_t data structure */
618 nbnxn_atomdata_t::nbnxn_atomdata_t(gmx::PinningPolicy      pinningPolicy,
619                                    const gmx::MDLogger&    mdlog,
620                                    const Nbnxm::KernelType kernelType,
621                                    int                     enbnxninitcombrule,
622                                    int                     ntype,
623                                    ArrayRef<const real>    nbfp,
624                                    int                     n_energygroups,
625                                    int                     nout) :
626     params_(pinningPolicy),
627     numAtoms_(0),
628     natoms_local(0),
629     shift_vec({}, { pinningPolicy }),
630     x_({}, { pinningPolicy }),
631     simdMasks(),
632     bUseBufferFlags(FALSE)
633 {
634     nbnxn_atomdata_params_init(
635             mdlog, &paramsDeprecated(), kernelType, enbnxninitcombrule, ntype, nbfp, n_energygroups);
636
637     const bool simple = Nbnxm::kernelTypeUsesSimplePairlist(kernelType);
638     const bool bSIMD  = Nbnxm::kernelTypeIsSimd(kernelType);
639
640     if (simple)
641     {
642         if (bSIMD)
643         {
644             int pack_x = std::max(c_nbnxnCpuIClusterSize, Nbnxm::JClusterSizePerKernelType[kernelType]);
645             switch (pack_x)
646             {
647                 case 4: XFormat = nbatX4; break;
648                 case 8: XFormat = nbatX8; break;
649                 default: gmx_incons("Unsupported packing width");
650             }
651         }
652         else
653         {
654             XFormat = nbatXYZ;
655         }
656
657         FFormat = XFormat;
658     }
659     else
660     {
661         XFormat = nbatXYZQ;
662         FFormat = nbatXYZ;
663     }
664
665     shift_vec.resize(gmx::c_numShiftVectors);
666
667     xstride = (XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
668     fstride = (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 = params().type.get_allocator().pinningPolicy();
674         out.emplace_back(kernelType, params().nenergrp, 1 << params().neg_2log, pinningPolicy);
675     }
676
677     buffer_flags.clear();
678 }
679
680 template<int packSize>
681 static void copy_lj_to_nbat_lj_comb(gmx::ArrayRef<const real> ljparam_type, const int* type, int na, real* ljparam_at)
682 {
683     /* The LJ params follow the combination rule:
684      * copy the params for the type array to the atom array.
685      */
686     for (int is = 0; is < na; is += packSize)
687     {
688         for (int k = 0; k < packSize; k++)
689         {
690             int i                             = is + k;
691             ljparam_at[is * 2 + k]            = ljparam_type[type[i] * 2];
692             ljparam_at[is * 2 + packSize + k] = ljparam_type[type[i] * 2 + 1];
693         }
694     }
695 }
696
697 /* Sets the atom type in nbnxn_atomdata_t */
698 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t::Params* params,
699                                          const Nbnxm::GridSet&     gridSet,
700                                          ArrayRef<const int>       atomTypes)
701 {
702     params->type.resize(gridSet.numGridAtomsTotal());
703
704     for (const Nbnxm::Grid& grid : gridSet.grids())
705     {
706         /* Loop over all columns and copy and fill */
707         for (int i = 0; i < grid.numColumns(); i++)
708         {
709             const int numAtoms   = grid.paddedNumAtomsInColumn(i);
710             const int atomOffset = grid.firstAtomInColumn(i);
711
712             copy_int_to_nbat_int(gridSet.atomIndices().data() + atomOffset,
713                                  grid.numAtomsInColumn(i),
714                                  numAtoms,
715                                  atomTypes.data(),
716                                  params->numTypes - 1,
717                                  params->type.data() + atomOffset);
718         }
719     }
720 }
721
722 /* Sets the LJ combination rule parameters in nbnxn_atomdata_t */
723 static void nbnxn_atomdata_set_ljcombparams(nbnxn_atomdata_t::Params* params,
724                                             const int                 XFormat,
725                                             const Nbnxm::GridSet&     gridSet)
726 {
727     params->lj_comb.resize(gridSet.numGridAtomsTotal() * 2);
728
729     if (params->ljCombinationRule != LJCombinationRule::None)
730     {
731         for (const Nbnxm::Grid& grid : gridSet.grids())
732         {
733             /* Loop over all columns and copy and fill */
734             for (int i = 0; i < grid.numColumns(); i++)
735             {
736                 const int numAtoms   = grid.paddedNumAtomsInColumn(i);
737                 const int atomOffset = grid.firstAtomInColumn(i);
738
739                 if (XFormat == nbatX4)
740                 {
741                     copy_lj_to_nbat_lj_comb<c_packX4>(params->nbfp_comb,
742                                                       params->type.data() + atomOffset,
743                                                       numAtoms,
744                                                       params->lj_comb.data() + atomOffset * 2);
745                 }
746                 else if (XFormat == nbatX8)
747                 {
748                     copy_lj_to_nbat_lj_comb<c_packX8>(params->nbfp_comb,
749                                                       params->type.data() + atomOffset,
750                                                       numAtoms,
751                                                       params->lj_comb.data() + atomOffset * 2);
752                 }
753                 else if (XFormat == nbatXYZQ)
754                 {
755                     copy_lj_to_nbat_lj_comb<1>(params->nbfp_comb,
756                                                params->type.data() + atomOffset,
757                                                numAtoms,
758                                                params->lj_comb.data() + atomOffset * 2);
759                 }
760             }
761         }
762     }
763 }
764
765 /* Sets the charges in nbnxn_atomdata_t *nbat */
766 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t*     nbat,
767                                        const Nbnxm::GridSet& gridSet,
768                                        ArrayRef<const real>  charges)
769 {
770     if (nbat->XFormat != nbatXYZQ)
771     {
772         nbat->paramsDeprecated().q.resize(nbat->numAtoms());
773     }
774
775     for (const Nbnxm::Grid& grid : gridSet.grids())
776     {
777         /* Loop over all columns and copy and fill */
778         for (int cxy = 0; cxy < grid.numColumns(); cxy++)
779         {
780             const int atomOffset     = grid.firstAtomInColumn(cxy);
781             const int numAtoms       = grid.numAtomsInColumn(cxy);
782             const int paddedNumAtoms = grid.paddedNumAtomsInColumn(cxy);
783
784             if (nbat->XFormat == nbatXYZQ)
785             {
786                 real* q = nbat->x().data() + atomOffset * STRIDE_XYZQ + ZZ + 1;
787                 for (int i = 0; i < numAtoms; i++)
788                 {
789                     *q = charges[gridSet.atomIndices()[atomOffset + i]];
790                     q += STRIDE_XYZQ;
791                 }
792                 /* Complete the partially filled last cell with zeros */
793                 for (int i = numAtoms; i < paddedNumAtoms; i++)
794                 {
795                     *q = 0;
796                     q += STRIDE_XYZQ;
797                 }
798             }
799             else
800             {
801                 real* q = nbat->paramsDeprecated().q.data() + atomOffset;
802                 for (int i = 0; i < numAtoms; i++)
803                 {
804                     *q = charges[gridSet.atomIndices()[atomOffset + i]];
805                     q++;
806                 }
807                 /* Complete the partially filled last cell with zeros */
808                 for (int i = numAtoms; i < paddedNumAtoms; i++)
809                 {
810                     *q = 0;
811                     q++;
812                 }
813             }
814         }
815     }
816 }
817
818 /* Set the charges of perturbed atoms in nbnxn_atomdata_t to 0.
819  * This is to automatically remove the RF/PME self term in the nbnxn kernels.
820  * Part of the zero interactions are still calculated in the normal kernels.
821  * All perturbed interactions are calculated in the free energy kernel,
822  * using the original charge and LJ data, not nbnxn_atomdata_t.
823  */
824 static void nbnxn_atomdata_mask_fep(nbnxn_atomdata_t* nbat, const Nbnxm::GridSet& gridSet)
825 {
826     nbnxn_atomdata_t::Params& params = nbat->paramsDeprecated();
827
828     const bool formatIsXYZQ = (nbat->XFormat == nbatXYZQ);
829
830     real* q        = formatIsXYZQ ? (nbat->x().data() + ZZ + 1) : params.q.data();
831     int   stride_q = formatIsXYZQ ? STRIDE_XYZQ : 1;
832
833     for (const Nbnxm::Grid& grid : gridSet.grids())
834     {
835         const int nsubc = (grid.geometry().isSimple) ? 1 : c_gpuNumClusterPerCell;
836
837         const int c_offset = grid.firstAtomInColumn(0);
838
839         /* Loop over all columns and copy and fill */
840         for (int c = 0; c < grid.numCells() * nsubc; c++)
841         {
842             /* Does this cluster contain perturbed particles? */
843             if (grid.clusterIsPerturbed(c))
844             {
845                 const int numAtomsPerCluster = grid.geometry().numAtomsICluster;
846                 for (int i = 0; i < numAtomsPerCluster; i++)
847                 {
848                     /* Is this a perturbed particle? */
849                     if (grid.atomIsPerturbed(c, i))
850                     {
851                         int ind = c_offset + c * numAtomsPerCluster + i;
852                         /* Set atom type and charge to non-interacting */
853                         params.type[ind]  = params.numTypes - 1;
854                         q[ind * stride_q] = 0;
855                     }
856                 }
857             }
858         }
859     }
860 }
861
862 /* Copies the energy group indices to a reordered and packed array */
863 static void
864 copy_egp_to_nbat_egps(const int* a, int na, int na_round, int na_c, int bit_shift, const int* in, int* innb)
865 {
866     int i = 0, j = 0;
867     for (; i < na; i += na_c)
868     {
869         /* Store na_c energy group numbers into one int */
870         int comb = 0;
871         for (int sa = 0; sa < na_c; sa++)
872         {
873             int at = a[i + sa];
874             if (at >= 0)
875             {
876                 comb |= (GET_CGINFO_GID(in[at]) << (sa * bit_shift));
877             }
878         }
879         innb[j++] = comb;
880     }
881     /* Complete the partially filled last cell with fill */
882     for (; i < na_round; i += na_c)
883     {
884         innb[j++] = 0;
885     }
886 }
887
888 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
889 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t::Params* params,
890                                             const Nbnxm::GridSet&     gridSet,
891                                             ArrayRef<const int>       atomInfo)
892 {
893     if (params->nenergrp == 1)
894     {
895         return;
896     }
897
898     params->energrp.resize(gridSet.numGridAtomsTotal());
899
900     for (const Nbnxm::Grid& grid : gridSet.grids())
901     {
902         /* Loop over all columns and copy and fill */
903         for (int i = 0; i < grid.numColumns(); i++)
904         {
905             const int numAtoms   = grid.paddedNumAtomsInColumn(i);
906             const int atomOffset = grid.firstAtomInColumn(i);
907
908             copy_egp_to_nbat_egps(gridSet.atomIndices().data() + atomOffset,
909                                   grid.numAtomsInColumn(i),
910                                   numAtoms,
911                                   c_nbnxnCpuIClusterSize,
912                                   params->neg_2log,
913                                   atomInfo.data(),
914                                   params->energrp.data() + grid.atomToCluster(atomOffset));
915         }
916     }
917 }
918
919 /* Sets all required atom parameter data in nbnxn_atomdata_t */
920 void nbnxn_atomdata_set(nbnxn_atomdata_t*     nbat,
921                         const Nbnxm::GridSet& gridSet,
922                         ArrayRef<const int>   atomTypes,
923                         ArrayRef<const real>  atomCharges,
924                         ArrayRef<const int>   atomInfo)
925 {
926     nbnxn_atomdata_t::Params& params = nbat->paramsDeprecated();
927
928     nbnxn_atomdata_set_atomtypes(&params, gridSet, atomTypes);
929
930     nbnxn_atomdata_set_charges(nbat, gridSet, atomCharges);
931
932     if (gridSet.haveFep())
933     {
934         nbnxn_atomdata_mask_fep(nbat, gridSet);
935     }
936
937     /* This must be done after masking types for FEP */
938     nbnxn_atomdata_set_ljcombparams(&params, nbat->XFormat, gridSet);
939
940     nbnxn_atomdata_set_energygroups(&params, gridSet, atomInfo);
941 }
942
943 /* Copies the shift vector array to nbnxn_atomdata_t */
944 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox, gmx::ArrayRef<gmx::RVec> shift_vec, nbnxn_atomdata_t* nbat)
945 {
946     nbat->bDynamicBox = bDynamicBox;
947     std::copy(shift_vec.begin(), shift_vec.end(), nbat->shift_vec.begin());
948 }
949
950 // This is slightly different from nbnxn_get_atom_range(...) at the end of the file
951 // TODO: Combine if possible
952 static void getAtomRanges(const Nbnxm::GridSet&   gridSet,
953                           const gmx::AtomLocality locality,
954                           int*                    gridBegin,
955                           int*                    gridEnd)
956 {
957     switch (locality)
958     {
959         case gmx::AtomLocality::All:
960             *gridBegin = 0;
961             *gridEnd   = gridSet.grids().size();
962             break;
963         case gmx::AtomLocality::Local:
964             *gridBegin = 0;
965             *gridEnd   = 1;
966             break;
967         case gmx::AtomLocality::NonLocal:
968             *gridBegin = 1;
969             *gridEnd   = gridSet.grids().size();
970             break;
971         case gmx::AtomLocality::Count:
972             GMX_ASSERT(false, "Count is invalid locality specifier");
973             break;
974     }
975 }
976
977 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
978 void nbnxn_atomdata_copy_x_to_nbat_x(const Nbnxm::GridSet&   gridSet,
979                                      const gmx::AtomLocality locality,
980                                      const rvec*             coordinates,
981                                      nbnxn_atomdata_t*       nbat)
982 {
983
984     int gridBegin = 0;
985     int gridEnd   = 0;
986     getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
987
988     const int nth = gmx_omp_nthreads_get(ModuleMultiThread::Pairsearch);
989 #pragma omp parallel for num_threads(nth) schedule(static)
990     for (int th = 0; th < nth; th++)
991     {
992         try
993         {
994             for (int g = gridBegin; g < gridEnd; g++)
995             {
996                 const Nbnxm::Grid& grid       = gridSet.grids()[g];
997                 const int          numCellsXY = grid.numColumns();
998
999                 const int cxy0 = (numCellsXY * th + nth - 1) / nth;
1000                 const int cxy1 = (numCellsXY * (th + 1) + nth - 1) / nth;
1001
1002                 for (int cxy = cxy0; cxy < cxy1; cxy++)
1003                 {
1004                     const int na  = grid.numAtomsInColumn(cxy);
1005                     const int ash = grid.firstAtomInColumn(cxy);
1006
1007                     copy_rvec_to_nbat_real(gridSet.atomIndices().data() + ash,
1008                                            na,
1009                                            na,
1010                                            coordinates,
1011                                            nbat->XFormat,
1012                                            nbat->x().data(),
1013                                            ash);
1014                 }
1015             }
1016         }
1017         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1018     }
1019 }
1020
1021 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t on the GPU*/
1022 void nbnxn_atomdata_x_to_nbat_x_gpu(const Nbnxm::GridSet&   gridSet,
1023                                     const gmx::AtomLocality locality,
1024                                     NbnxmGpu*               gpu_nbv,
1025                                     DeviceBuffer<RVec>      d_x,
1026                                     GpuEventSynchronizer*   xReadyOnDevice)
1027 {
1028     GMX_ASSERT(xReadyOnDevice != nullptr, "Need a valid GpuEventSynchronizer object");
1029
1030     int gridBegin = 0;
1031     int gridEnd   = 0;
1032     getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
1033
1034     for (int g = gridBegin; g < gridEnd; g++)
1035     {
1036         nbnxn_gpu_x_to_nbat_x(gridSet.grids()[g],
1037                               gpu_nbv,
1038                               d_x,
1039                               (g == gridBegin) ? xReadyOnDevice : nullptr, // Sync on first iteration only
1040                               locality,
1041                               g,
1042                               gridSet.numColumnsMax(),
1043                               (g == gridEnd - 1));
1044     }
1045 }
1046
1047 static void nbnxn_atomdata_clear_reals(gmx::ArrayRef<real> dest, int i0, int i1)
1048 {
1049     for (int i = i0; i < i1; i++)
1050     {
1051         dest[i] = 0;
1052     }
1053 }
1054
1055 gmx_unused static void nbnxn_atomdata_reduce_reals(real* gmx_restrict dest,
1056                                                    gmx_bool           bDestSet,
1057                                                    const real** gmx_restrict src,
1058                                                    int                       nsrc,
1059                                                    int                       i0,
1060                                                    int                       i1)
1061 {
1062     if (bDestSet)
1063     {
1064         /* The destination buffer contains data, add to it */
1065         for (int i = i0; i < i1; i++)
1066         {
1067             for (int s = 0; s < nsrc; s++)
1068             {
1069                 dest[i] += src[s][i];
1070             }
1071         }
1072     }
1073     else
1074     {
1075         /* The destination buffer is uninitialized, set it first */
1076         for (int i = i0; i < i1; i++)
1077         {
1078             dest[i] = src[0][i];
1079             for (int s = 1; s < nsrc; s++)
1080             {
1081                 dest[i] += src[s][i];
1082             }
1083         }
1084     }
1085 }
1086
1087 gmx_unused static void nbnxn_atomdata_reduce_reals_simd(real gmx_unused* gmx_restrict dest,
1088                                                         gmx_bool gmx_unused bDestSet,
1089                                                         const gmx_unused real** gmx_restrict src,
1090                                                         int gmx_unused nsrc,
1091                                                         int gmx_unused i0,
1092                                                         int gmx_unused i1)
1093 {
1094 #if GMX_SIMD
1095     /* The SIMD width here is actually independent of that in the kernels,
1096      * but we use the same width for simplicity (usually optimal anyhow).
1097      */
1098     SimdReal dest_SSE, src_SSE;
1099
1100     if (bDestSet)
1101     {
1102         for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1103         {
1104             dest_SSE = load<SimdReal>(dest + i);
1105             for (int s = 0; s < nsrc; s++)
1106             {
1107                 src_SSE  = load<SimdReal>(src[s] + i);
1108                 dest_SSE = dest_SSE + src_SSE;
1109             }
1110             store(dest + i, dest_SSE);
1111         }
1112     }
1113     else
1114     {
1115         for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1116         {
1117             dest_SSE = load<SimdReal>(src[0] + i);
1118             for (int s = 1; s < nsrc; s++)
1119             {
1120                 src_SSE  = load<SimdReal>(src[s] + i);
1121                 dest_SSE = dest_SSE + src_SSE;
1122             }
1123             store(dest + i, dest_SSE);
1124         }
1125     }
1126 #endif
1127 }
1128
1129 /* Add part of the force array(s) from nbnxn_atomdata_t to f
1130  *
1131  * Note: Adding restrict to f makes this function 50% slower with gcc 7.3
1132  */
1133 static void nbnxn_atomdata_add_nbat_f_to_f_part(const Nbnxm::GridSet&          gridSet,
1134                                                 const nbnxn_atomdata_t&        nbat,
1135                                                 const nbnxn_atomdata_output_t& out,
1136                                                 const int                      a0,
1137                                                 const int                      a1,
1138                                                 rvec*                          f)
1139 {
1140     gmx::ArrayRef<const int> cell = gridSet.cells();
1141     // Note: Using ArrayRef instead makes this code 25% slower with gcc 7.3
1142     const real* fnb = out.f.data();
1143
1144     /* Loop over all columns and copy and fill */
1145     switch (nbat.FFormat)
1146     {
1147         case nbatXYZ:
1148         case nbatXYZQ:
1149             for (int a = a0; a < a1; a++)
1150             {
1151                 int i = cell[a] * nbat.fstride;
1152
1153                 f[a][XX] += fnb[i];
1154                 f[a][YY] += fnb[i + 1];
1155                 f[a][ZZ] += fnb[i + 2];
1156             }
1157             break;
1158         case nbatX4:
1159             for (int a = a0; a < a1; a++)
1160             {
1161                 int i = atom_to_x_index<c_packX4>(cell[a]);
1162
1163                 f[a][XX] += fnb[i + XX * c_packX4];
1164                 f[a][YY] += fnb[i + YY * c_packX4];
1165                 f[a][ZZ] += fnb[i + ZZ * c_packX4];
1166             }
1167             break;
1168         case nbatX8:
1169             for (int a = a0; a < a1; a++)
1170             {
1171                 int i = atom_to_x_index<c_packX8>(cell[a]);
1172
1173                 f[a][XX] += fnb[i + XX * c_packX8];
1174                 f[a][YY] += fnb[i + YY * c_packX8];
1175                 f[a][ZZ] += fnb[i + ZZ * c_packX8];
1176             }
1177             break;
1178         default: gmx_incons("Unsupported nbnxn_atomdata_t format");
1179     }
1180 }
1181
1182 static void nbnxn_atomdata_add_nbat_f_to_f_reduce(nbnxn_atomdata_t* nbat, int nth)
1183 {
1184 #pragma omp parallel for num_threads(nth) schedule(static)
1185     for (int th = 0; th < nth; th++)
1186     {
1187         try
1188         {
1189             const real* fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1190
1191             gmx::ArrayRef<const gmx_bitmask_t> flags = nbat->buffer_flags;
1192
1193             /* Calculate the cell-block range for our thread */
1194             const int b0 = (flags.size() * th) / nth;
1195             const int b1 = (flags.size() * (th + 1)) / nth;
1196
1197             for (int b = b0; b < b1; b++)
1198             {
1199                 const int i0 = b * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1200                 const int i1 = (b + 1) * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1201
1202                 int nfptr = 0;
1203                 for (gmx::index out = 1; out < gmx::ssize(nbat->out); out++)
1204                 {
1205                     if (bitmask_is_set(flags[b], out))
1206                     {
1207                         fptr[nfptr++] = nbat->out[out].f.data();
1208                     }
1209                 }
1210                 if (nfptr > 0)
1211                 {
1212 #if GMX_SIMD
1213                     nbnxn_atomdata_reduce_reals_simd
1214 #else
1215                     nbnxn_atomdata_reduce_reals
1216 #endif
1217                             (nbat->out[0].f.data(), bitmask_is_set(flags[b], 0), fptr, nfptr, i0, i1);
1218                 }
1219                 else if (!bitmask_is_set(flags[b], 0))
1220                 {
1221                     nbnxn_atomdata_clear_reals(nbat->out[0].f, i0, i1);
1222                 }
1223             }
1224         }
1225         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1226     }
1227 }
1228
1229
1230 /* Add the force array(s) from nbnxn_atomdata_t to f */
1231 void reduceForces(nbnxn_atomdata_t* nbat, const gmx::AtomLocality locality, const Nbnxm::GridSet& gridSet, rvec* f)
1232 {
1233     int a0 = 0;
1234     int na = 0;
1235
1236     nbnxn_get_atom_range(locality, gridSet, &a0, &na);
1237
1238     if (na == 0)
1239     {
1240         /* The are no atoms for this reduction, avoid some overhead */
1241         return;
1242     }
1243
1244     int nth = gmx_omp_nthreads_get(ModuleMultiThread::Nonbonded);
1245
1246     if (nbat->out.size() > 1)
1247     {
1248         if (locality != gmx::AtomLocality::All)
1249         {
1250             gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1251         }
1252
1253         /* Reduce the force thread output buffers into buffer 0, before adding
1254          * them to the, differently ordered, "real" force buffer.
1255          */
1256         nbnxn_atomdata_add_nbat_f_to_f_reduce(nbat, nth);
1257     }
1258 #pragma omp parallel for num_threads(nth) schedule(static)
1259     for (int th = 0; th < nth; th++)
1260     {
1261         try
1262         {
1263             nbnxn_atomdata_add_nbat_f_to_f_part(
1264                     gridSet, *nbat, nbat->out[0], a0 + ((th + 0) * na) / nth, a0 + ((th + 1) * na) / nth, f);
1265         }
1266         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1267     }
1268 }
1269
1270 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t& nbat, gmx::ArrayRef<gmx::RVec> fshift)
1271 {
1272     gmx::ArrayRef<const nbnxn_atomdata_output_t> outputBuffers = nbat.out;
1273
1274     for (int s = 0; s < gmx::c_numShiftVectors; s++)
1275     {
1276         rvec sum;
1277         clear_rvec(sum);
1278         for (const nbnxn_atomdata_output_t& out : outputBuffers)
1279         {
1280             sum[XX] += out.fshift[s * DIM + XX];
1281             sum[YY] += out.fshift[s * DIM + YY];
1282             sum[ZZ] += out.fshift[s * DIM + ZZ];
1283         }
1284         fshift[s] += sum;
1285     }
1286 }
1287
1288 void nbnxn_get_atom_range(const gmx::AtomLocality atomLocality,
1289                           const Nbnxm::GridSet&   gridSet,
1290                           int*                    atomStart,
1291                           int*                    nAtoms)
1292 {
1293
1294     switch (atomLocality)
1295     {
1296         case gmx::AtomLocality::All:
1297             *atomStart = 0;
1298             *nAtoms    = gridSet.numRealAtomsTotal();
1299             break;
1300         case gmx::AtomLocality::Local:
1301             *atomStart = 0;
1302             *nAtoms    = gridSet.numRealAtomsLocal();
1303             break;
1304         case gmx::AtomLocality::NonLocal:
1305             *atomStart = gridSet.numRealAtomsLocal();
1306             *nAtoms    = gridSet.numRealAtomsTotal() - gridSet.numRealAtomsLocal();
1307             break;
1308         case gmx::AtomLocality::Count:
1309             GMX_ASSERT(false, "Count is invalid locality specifier");
1310             break;
1311     }
1312 }