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