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