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