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