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