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