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