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