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