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