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