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