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