eebec5963ff21b506894c643cdc959fc2a506e72
[alexxy/gromacs.git] / src / mdlib / nbnxn_atomdata.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2012, The GROMACS development team,
6  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012,2013, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42
43 #include <math.h>
44 #include <string.h>
45 #include "smalloc.h"
46 #include "macros.h"
47 #include "vec.h"
48 #include "nbnxn_consts.h"
49 #include "nbnxn_internal.h"
50 #include "nbnxn_search.h"
51 #include "nbnxn_atomdata.h"
52 #include "gmx_omp_nthreads.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 /* Determines the combination rule (or none) to be used, stores it,
365  * and sets the LJ parameters required with the rule.
366  */
367 static void set_combination_rule_data(nbnxn_atomdata_t *nbat)
368 {
369     int  nt, i, j;
370     real c6, c12;
371
372     nt = nbat->ntype;
373
374     switch (nbat->comb_rule)
375     {
376         case  ljcrGEOM:
377             nbat->comb_rule = ljcrGEOM;
378
379             for (i = 0; i < nt; i++)
380             {
381                 /* Copy the diagonal from the nbfp matrix */
382                 nbat->nbfp_comb[i*2  ] = sqrt(nbat->nbfp[(i*nt+i)*2  ]);
383                 nbat->nbfp_comb[i*2+1] = sqrt(nbat->nbfp[(i*nt+i)*2+1]);
384             }
385             break;
386         case ljcrLB:
387             for (i = 0; i < nt; i++)
388             {
389                 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
390                 c6  = nbat->nbfp[(i*nt+i)*2  ];
391                 c12 = nbat->nbfp[(i*nt+i)*2+1];
392                 if (c6 > 0 && c12 > 0)
393                 {
394                     /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
395                      * so we get 6*C6 and 12*C12 after combining.
396                      */
397                     nbat->nbfp_comb[i*2  ] = 0.5*pow(c12/c6, 1.0/6.0);
398                     nbat->nbfp_comb[i*2+1] = sqrt(c6*c6/c12);
399                 }
400                 else
401                 {
402                     nbat->nbfp_comb[i*2  ] = 0;
403                     nbat->nbfp_comb[i*2+1] = 0;
404                 }
405             }
406             break;
407         case ljcrNONE:
408             /* nbfp_s4 stores two parameters using a stride of 4,
409              * because this would suit x86 SIMD single-precision
410              * quad-load intrinsics. There's a slight inefficiency in
411              * allocating and initializing nbfp_s4 when it might not
412              * be used, but introducing the conditional code is not
413              * really worth it. */
414             nbat->alloc((void **)&nbat->nbfp_s4, nt*nt*4*sizeof(*nbat->nbfp_s4));
415             for (i = 0; i < nt; i++)
416             {
417                 for (j = 0; j < nt; j++)
418                 {
419                     nbat->nbfp_s4[(i*nt+j)*4+0] = nbat->nbfp[(i*nt+j)*2+0];
420                     nbat->nbfp_s4[(i*nt+j)*4+1] = nbat->nbfp[(i*nt+j)*2+1];
421                     nbat->nbfp_s4[(i*nt+j)*4+2] = 0;
422                     nbat->nbfp_s4[(i*nt+j)*4+3] = 0;
423                 }
424             }
425             break;
426         default:
427             gmx_incons("Unknown combination rule");
428             break;
429     }
430 }
431
432 /* Initializes an nbnxn_atomdata_t data structure */
433 void nbnxn_atomdata_init(FILE *fp,
434                          nbnxn_atomdata_t *nbat,
435                          int nb_kernel_type,
436                          int ntype, const real *nbfp,
437                          int n_energygroups,
438                          int nout,
439                          nbnxn_alloc_t *alloc,
440                          nbnxn_free_t  *free)
441 {
442     int      i, j;
443     real     c6, c12, tol;
444     char    *ptr;
445     gmx_bool simple, bCombGeom, bCombLB;
446
447     if (alloc == NULL)
448     {
449         nbat->alloc = nbnxn_alloc_aligned;
450     }
451     else
452     {
453         nbat->alloc = alloc;
454     }
455     if (free == NULL)
456     {
457         nbat->free = nbnxn_free_aligned;
458     }
459     else
460     {
461         nbat->free = free;
462     }
463
464     if (debug)
465     {
466         fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
467     }
468     nbat->ntype = ntype + 1;
469     nbat->alloc((void **)&nbat->nbfp,
470                 nbat->ntype*nbat->ntype*2*sizeof(*nbat->nbfp));
471     nbat->alloc((void **)&nbat->nbfp_comb, nbat->ntype*2*sizeof(*nbat->nbfp_comb));
472
473     /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
474      * force-field floating point parameters.
475      */
476     tol = 1e-5;
477     ptr = getenv("GMX_LJCOMB_TOL");
478     if (ptr != NULL)
479     {
480         double dbl;
481
482         sscanf(ptr, "%lf", &dbl);
483         tol = dbl;
484     }
485     bCombGeom = TRUE;
486     bCombLB   = TRUE;
487
488     /* Temporarily fill nbat->nbfp_comb with sigma and epsilon
489      * to check for the LB rule.
490      */
491     for (i = 0; i < ntype; i++)
492     {
493         c6  = nbfp[(i*ntype+i)*2  ]/6.0;
494         c12 = nbfp[(i*ntype+i)*2+1]/12.0;
495         if (c6 > 0 && c12 > 0)
496         {
497             nbat->nbfp_comb[i*2  ] = pow(c12/c6, 1.0/6.0);
498             nbat->nbfp_comb[i*2+1] = 0.25*c6*c6/c12;
499         }
500         else if (c6 == 0 && c12 == 0)
501         {
502             nbat->nbfp_comb[i*2  ] = 0;
503             nbat->nbfp_comb[i*2+1] = 0;
504         }
505         else
506         {
507             /* Can not use LB rule with only dispersion or repulsion */
508             bCombLB = FALSE;
509         }
510     }
511
512     for (i = 0; i < nbat->ntype; i++)
513     {
514         for (j = 0; j < nbat->ntype; j++)
515         {
516             if (i < ntype && j < ntype)
517             {
518                 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
519                  * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
520                  */
521                 c6  = nbfp[(i*ntype+j)*2  ];
522                 c12 = nbfp[(i*ntype+j)*2+1];
523                 nbat->nbfp[(i*nbat->ntype+j)*2  ] = c6;
524                 nbat->nbfp[(i*nbat->ntype+j)*2+1] = c12;
525
526                 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
527                 bCombGeom = bCombGeom &&
528                     gmx_within_tol(c6*c6, nbfp[(i*ntype+i)*2  ]*nbfp[(j*ntype+j)*2  ], tol) &&
529                     gmx_within_tol(c12*c12, nbfp[(i*ntype+i)*2+1]*nbfp[(j*ntype+j)*2+1], tol);
530
531                 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
532                 c6     /= 6.0;
533                 c12    /= 12.0;
534                 bCombLB = bCombLB &&
535                     ((c6 == 0 && c12 == 0 &&
536                       (nbat->nbfp_comb[i*2+1] == 0 || nbat->nbfp_comb[j*2+1] == 0)) ||
537                      (c6 > 0 && c12 > 0 &&
538                       gmx_within_tol(pow(c12/c6, 1.0/6.0), 0.5*(nbat->nbfp_comb[i*2]+nbat->nbfp_comb[j*2]), tol) &&
539                       gmx_within_tol(0.25*c6*c6/c12, sqrt(nbat->nbfp_comb[i*2+1]*nbat->nbfp_comb[j*2+1]), tol)));
540             }
541             else
542             {
543                 /* Add zero parameters for the additional dummy atom type */
544                 nbat->nbfp[(i*nbat->ntype+j)*2  ] = 0;
545                 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 0;
546             }
547         }
548     }
549     if (debug)
550     {
551         fprintf(debug, "Combination rules: geometric %d Lorentz-Berthelot %d\n",
552                 bCombGeom, bCombLB);
553     }
554
555     simple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
556
557     if (simple)
558     {
559         /* We prefer the geometic combination rule,
560          * as that gives a slightly faster kernel than the LB rule.
561          */
562         if (bCombGeom)
563         {
564             nbat->comb_rule = ljcrGEOM;
565         }
566         else if (bCombLB)
567         {
568             nbat->comb_rule = ljcrLB;
569         }
570         else
571         {
572             nbat->comb_rule = ljcrNONE;
573
574             nbat->free(nbat->nbfp_comb);
575         }
576
577         if (fp)
578         {
579             if (nbat->comb_rule == ljcrNONE)
580             {
581                 fprintf(fp, "Using full Lennard-Jones parameter combination matrix\n\n");
582             }
583             else
584             {
585                 fprintf(fp, "Using %s Lennard-Jones combination rule\n\n",
586                         nbat->comb_rule == ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
587             }
588         }
589
590         set_combination_rule_data(nbat);
591     }
592     else
593     {
594         nbat->comb_rule = ljcrNONE;
595
596         nbat->free(nbat->nbfp_comb);
597     }
598
599     nbat->natoms  = 0;
600     nbat->type    = NULL;
601     nbat->lj_comb = NULL;
602     if (simple)
603     {
604         int pack_x;
605
606         switch (nb_kernel_type)
607         {
608             case nbnxnk4xN_SIMD_4xN:
609             case nbnxnk4xN_SIMD_2xNN:
610                 pack_x = max(NBNXN_CPU_CLUSTER_I_SIZE,
611                              nbnxn_kernel_to_cj_size(nb_kernel_type));
612                 switch (pack_x)
613                 {
614                     case 4:
615                         nbat->XFormat = nbatX4;
616                         break;
617                     case 8:
618                         nbat->XFormat = nbatX8;
619                         break;
620                     default:
621                         gmx_incons("Unsupported packing width");
622                 }
623                 break;
624             default:
625                 nbat->XFormat = nbatXYZ;
626                 break;
627         }
628
629         nbat->FFormat = nbat->XFormat;
630     }
631     else
632     {
633         nbat->XFormat = nbatXYZQ;
634         nbat->FFormat = nbatXYZ;
635     }
636     nbat->q        = NULL;
637     nbat->nenergrp = n_energygroups;
638     if (!simple)
639     {
640         /* Energy groups not supported yet for super-sub lists */
641         if (n_energygroups > 1 && fp != NULL)
642         {
643             fprintf(fp, "\nNOTE: With GPUs, reporting energy group contributions is not supported\n\n");
644         }
645         nbat->nenergrp = 1;
646     }
647     /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
648     if (nbat->nenergrp > 64)
649     {
650         gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
651     }
652     nbat->neg_2log = 1;
653     while (nbat->nenergrp > (1<<nbat->neg_2log))
654     {
655         nbat->neg_2log++;
656     }
657     nbat->energrp = NULL;
658     nbat->alloc((void **)&nbat->shift_vec, SHIFTS*sizeof(*nbat->shift_vec));
659     nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
660     nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
661     nbat->x       = NULL;
662
663 #ifdef GMX_NBNXN_SIMD
664     if (simple)
665     {
666         /* Set the diagonal cluster pair exclusion mask setup data.
667          * In the kernel we check 0 < j - i to generate the masks.
668          * Here we store j - i for generating the mask for the first i,
669          * we substract 0.5 to avoid rounding issues.
670          * In the kernel we can subtract 1 to generate the subsequent mask.
671          */
672         const int simd_width = GMX_NBNXN_SIMD_BITWIDTH/(sizeof(real)*8);
673         int       simd_4xn_diag_size, real_excl, simd_excl_size, j, s;
674
675         simd_4xn_diag_size = max(NBNXN_CPU_CLUSTER_I_SIZE, simd_width);
676         snew_aligned(nbat->simd_4xn_diag, simd_4xn_diag_size, NBNXN_MEM_ALIGN);
677         for (j = 0; j < simd_4xn_diag_size; j++)
678         {
679             nbat->simd_4xn_diag[j] = j - 0.5;
680         }
681
682         snew_aligned(nbat->simd_2xnn_diag, simd_width, NBNXN_MEM_ALIGN);
683         for (j = 0; j < simd_width/2; j++)
684         {
685             /* The j-cluster size is half the SIMD width */
686             nbat->simd_2xnn_diag[j]              = j - 0.5;
687             /* The next half of the SIMD width is for i + 1 */
688             nbat->simd_2xnn_diag[simd_width/2+j] = j - 1 - 0.5;
689         }
690
691         /* We always use 32-bit integer exclusion masks. When we use
692          * double precision, we fit two integers in a double SIMD register.
693          */
694         real_excl = sizeof(real)/sizeof(*nbat->simd_excl_mask);
695         /* Set bits for use with both 4xN and 2x(N+N) kernels */
696         simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*simd_width*real_excl;
697         snew_aligned(nbat->simd_excl_mask, simd_excl_size*real_excl, NBNXN_MEM_ALIGN);
698         for (j = 0; j < simd_excl_size; j++)
699         {
700             /* Set the consecutive bits for masking pair exclusions.
701              * For double a single-bit mask would be enough.
702              * But using two bits avoids endianness issues.
703              */
704             for (s = 0; s < real_excl; s++)
705             {
706                 /* Set the consecutive bits for masking pair exclusions */
707                 nbat->simd_excl_mask[j*real_excl + s] = (1U << j);
708             }
709         }
710     }
711 #endif
712
713     /* Initialize the output data structures */
714     nbat->nout    = nout;
715     snew(nbat->out, nbat->nout);
716     nbat->nalloc  = 0;
717     for (i = 0; i < nbat->nout; i++)
718     {
719         nbnxn_atomdata_output_init(&nbat->out[i],
720                                    nb_kernel_type,
721                                    nbat->nenergrp, 1<<nbat->neg_2log,
722                                    nbat->alloc);
723     }
724     nbat->buffer_flags.flag        = NULL;
725     nbat->buffer_flags.flag_nalloc = 0;
726 }
727
728 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
729                                        const int *type, int na,
730                                        real *ljparam_at)
731 {
732     int is, k, i;
733
734     /* The LJ params follow the combination rule:
735      * copy the params for the type array to the atom array.
736      */
737     for (is = 0; is < na; is += PACK_X4)
738     {
739         for (k = 0; k < PACK_X4; k++)
740         {
741             i = is + k;
742             ljparam_at[is*2        +k] = ljparam_type[type[i]*2  ];
743             ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
744         }
745     }
746 }
747
748 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
749                                        const int *type, int na,
750                                        real *ljparam_at)
751 {
752     int is, k, i;
753
754     /* The LJ params follow the combination rule:
755      * copy the params for the type array to the atom array.
756      */
757     for (is = 0; is < na; is += PACK_X8)
758     {
759         for (k = 0; k < PACK_X8; k++)
760         {
761             i = is + k;
762             ljparam_at[is*2        +k] = ljparam_type[type[i]*2  ];
763             ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
764         }
765     }
766 }
767
768 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
769 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t    *nbat,
770                                          int                  ngrid,
771                                          const nbnxn_search_t nbs,
772                                          const int           *type)
773 {
774     int                 g, i, ncz, ash;
775     const nbnxn_grid_t *grid;
776
777     for (g = 0; g < ngrid; g++)
778     {
779         grid = &nbs->grid[g];
780
781         /* Loop over all columns and copy and fill */
782         for (i = 0; i < grid->ncx*grid->ncy; i++)
783         {
784             ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
785             ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
786
787             copy_int_to_nbat_int(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
788                                  type, nbat->ntype-1, nbat->type+ash);
789
790             if (nbat->comb_rule != ljcrNONE)
791             {
792                 if (nbat->XFormat == nbatX4)
793                 {
794                     copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
795                                                nbat->type+ash, ncz*grid->na_sc,
796                                                nbat->lj_comb+ash*2);
797                 }
798                 else if (nbat->XFormat == nbatX8)
799                 {
800                     copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
801                                                nbat->type+ash, ncz*grid->na_sc,
802                                                nbat->lj_comb+ash*2);
803                 }
804             }
805         }
806     }
807 }
808
809 /* Sets the charges in nbnxn_atomdata_t *nbat */
810 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t    *nbat,
811                                        int                  ngrid,
812                                        const nbnxn_search_t nbs,
813                                        const real          *charge)
814 {
815     int                 g, cxy, ncz, ash, na, na_round, i, j;
816     real               *q;
817     const nbnxn_grid_t *grid;
818
819     for (g = 0; g < ngrid; g++)
820     {
821         grid = &nbs->grid[g];
822
823         /* Loop over all columns and copy and fill */
824         for (cxy = 0; cxy < grid->ncx*grid->ncy; cxy++)
825         {
826             ash      = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
827             na       = grid->cxy_na[cxy];
828             na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
829
830             if (nbat->XFormat == nbatXYZQ)
831             {
832                 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
833                 for (i = 0; i < na; i++)
834                 {
835                     *q = charge[nbs->a[ash+i]];
836                     q += STRIDE_XYZQ;
837                 }
838                 /* Complete the partially filled last cell with zeros */
839                 for (; i < na_round; i++)
840                 {
841                     *q = 0;
842                     q += STRIDE_XYZQ;
843                 }
844             }
845             else
846             {
847                 q = nbat->q + ash;
848                 for (i = 0; i < na; i++)
849                 {
850                     *q = charge[nbs->a[ash+i]];
851                     q++;
852                 }
853                 /* Complete the partially filled last cell with zeros */
854                 for (; i < na_round; i++)
855                 {
856                     *q = 0;
857                     q++;
858                 }
859             }
860         }
861     }
862 }
863
864 /* Copies the energy group indices to a reordered and packed array */
865 static void copy_egp_to_nbat_egps(const int *a, int na, int na_round,
866                                   int na_c, int bit_shift,
867                                   const int *in, int *innb)
868 {
869     int i, j, sa, at;
870     int comb;
871
872     j = 0;
873     for (i = 0; i < na; i += na_c)
874     {
875         /* Store na_c energy group numbers into one int */
876         comb = 0;
877         for (sa = 0; sa < na_c; sa++)
878         {
879             at = a[i+sa];
880             if (at >= 0)
881             {
882                 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
883             }
884         }
885         innb[j++] = comb;
886     }
887     /* Complete the partially filled last cell with fill */
888     for (; i < na_round; i += na_c)
889     {
890         innb[j++] = 0;
891     }
892 }
893
894 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
895 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t    *nbat,
896                                             int                  ngrid,
897                                             const nbnxn_search_t nbs,
898                                             const int           *atinfo)
899 {
900     int                 g, i, ncz, ash;
901     const nbnxn_grid_t *grid;
902
903     for (g = 0; g < ngrid; g++)
904     {
905         grid = &nbs->grid[g];
906
907         /* Loop over all columns and copy and fill */
908         for (i = 0; i < grid->ncx*grid->ncy; i++)
909         {
910             ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
911             ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
912
913             copy_egp_to_nbat_egps(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
914                                   nbat->na_c, nbat->neg_2log,
915                                   atinfo, nbat->energrp+(ash>>grid->na_c_2log));
916         }
917     }
918 }
919
920 /* Sets all required atom parameter data in nbnxn_atomdata_t */
921 void nbnxn_atomdata_set(nbnxn_atomdata_t    *nbat,
922                         int                  locality,
923                         const nbnxn_search_t nbs,
924                         const t_mdatoms     *mdatoms,
925                         const int           *atinfo)
926 {
927     int ngrid;
928
929     if (locality == eatLocal)
930     {
931         ngrid = 1;
932     }
933     else
934     {
935         ngrid = nbs->ngrid;
936     }
937
938     nbnxn_atomdata_set_atomtypes(nbat, ngrid, nbs, mdatoms->typeA);
939
940     nbnxn_atomdata_set_charges(nbat, ngrid, nbs, mdatoms->chargeA);
941
942     if (nbat->nenergrp > 1)
943     {
944         nbnxn_atomdata_set_energygroups(nbat, ngrid, nbs, atinfo);
945     }
946 }
947
948 /* Copies the shift vector array to nbnxn_atomdata_t */
949 void nbnxn_atomdata_copy_shiftvec(gmx_bool          bDynamicBox,
950                                   rvec             *shift_vec,
951                                   nbnxn_atomdata_t *nbat)
952 {
953     int i;
954
955     nbat->bDynamicBox = bDynamicBox;
956     for (i = 0; i < SHIFTS; i++)
957     {
958         copy_rvec(shift_vec[i], nbat->shift_vec[i]);
959     }
960 }
961
962 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
963 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
964                                      int                  locality,
965                                      gmx_bool             FillLocal,
966                                      rvec                *x,
967                                      nbnxn_atomdata_t    *nbat)
968 {
969     int g0 = 0, g1 = 0;
970     int nth, th;
971
972     switch (locality)
973     {
974         case eatAll:
975             g0 = 0;
976             g1 = nbs->ngrid;
977             break;
978         case eatLocal:
979             g0 = 0;
980             g1 = 1;
981             break;
982         case eatNonlocal:
983             g0 = 1;
984             g1 = nbs->ngrid;
985             break;
986     }
987
988     if (FillLocal)
989     {
990         nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
991     }
992
993     nth = gmx_omp_nthreads_get(emntPairsearch);
994
995 #pragma omp parallel for num_threads(nth) schedule(static)
996     for (th = 0; th < nth; th++)
997     {
998         int g;
999
1000         for (g = g0; g < g1; g++)
1001         {
1002             const nbnxn_grid_t *grid;
1003             int                 cxy0, cxy1, cxy;
1004
1005             grid = &nbs->grid[g];
1006
1007             cxy0 = (grid->ncx*grid->ncy* th   +nth-1)/nth;
1008             cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
1009
1010             for (cxy = cxy0; cxy < cxy1; cxy++)
1011             {
1012                 int na, ash, na_fill;
1013
1014                 na  = grid->cxy_na[cxy];
1015                 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1016
1017                 if (g == 0 && FillLocal)
1018                 {
1019                     na_fill =
1020                         (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
1021                 }
1022                 else
1023                 {
1024                     /* We fill only the real particle locations.
1025                      * We assume the filling entries at the end have been
1026                      * properly set before during ns.
1027                      */
1028                     na_fill = na;
1029                 }
1030                 copy_rvec_to_nbat_real(nbs->a+ash, na, na_fill, x,
1031                                        nbat->XFormat, nbat->x, ash,
1032                                        0, 0, 0);
1033             }
1034         }
1035     }
1036 }
1037
1038 static void
1039 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
1040                            int i0, int i1)
1041 {
1042     int i;
1043
1044     for (i = i0; i < i1; i++)
1045     {
1046         dest[i] = 0;
1047     }
1048 }
1049
1050 static void
1051 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1052                             gmx_bool bDestSet,
1053                             real ** gmx_restrict src,
1054                             int nsrc,
1055                             int i0, int i1)
1056 {
1057     int i, s;
1058
1059     if (bDestSet)
1060     {
1061         /* The destination buffer contains data, add to it */
1062         for (i = i0; i < i1; i++)
1063         {
1064             for (s = 0; s < nsrc; s++)
1065             {
1066                 dest[i] += src[s][i];
1067             }
1068         }
1069     }
1070     else
1071     {
1072         /* The destination buffer is unitialized, set it first */
1073         for (i = i0; i < i1; i++)
1074         {
1075             dest[i] = src[0][i];
1076             for (s = 1; s < nsrc; s++)
1077             {
1078                 dest[i] += src[s][i];
1079             }
1080         }
1081     }
1082 }
1083
1084 static void
1085 nbnxn_atomdata_reduce_reals_simd(real * gmx_restrict dest,
1086                                  gmx_bool bDestSet,
1087                                  real ** gmx_restrict src,
1088                                  int nsrc,
1089                                  int i0, int i1)
1090 {
1091 #ifdef GMX_NBNXN_SIMD
1092 /* The SIMD width here is actually independent of that in the kernels,
1093  * but we use the same width for simplicity (usually optimal anyhow).
1094  */
1095 #ifdef GMX_NBNXN_HALF_WIDTH_SIMD
1096 #define GMX_USE_HALF_WIDTH_SIMD_HERE
1097 #endif
1098 #include "gmx_simd_macros.h"
1099
1100     int       i, s;
1101     gmx_mm_pr dest_SSE, src_SSE;
1102
1103     if (bDestSet)
1104     {
1105         for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1106         {
1107             dest_SSE = gmx_load_pr(dest+i);
1108             for (s = 0; s < nsrc; s++)
1109             {
1110                 src_SSE  = gmx_load_pr(src[s]+i);
1111                 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1112             }
1113             gmx_store_pr(dest+i, dest_SSE);
1114         }
1115     }
1116     else
1117     {
1118         for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1119         {
1120             dest_SSE = gmx_load_pr(src[0]+i);
1121             for (s = 1; s < nsrc; s++)
1122             {
1123                 src_SSE  = gmx_load_pr(src[s]+i);
1124                 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1125             }
1126             gmx_store_pr(dest+i, dest_SSE);
1127         }
1128     }
1129 #endif
1130 }
1131
1132 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1133 static void
1134 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1135                                     const nbnxn_atomdata_t *nbat,
1136                                     nbnxn_atomdata_output_t *out,
1137                                     int nfa,
1138                                     int a0, int a1,
1139                                     rvec *f)
1140 {
1141     int         a, i, fa;
1142     const int  *cell;
1143     const real *fnb;
1144
1145     cell = nbs->cell;
1146
1147     /* Loop over all columns and copy and fill */
1148     switch (nbat->FFormat)
1149     {
1150         case nbatXYZ:
1151         case nbatXYZQ:
1152             if (nfa == 1)
1153             {
1154                 fnb = out[0].f;
1155
1156                 for (a = a0; a < a1; a++)
1157                 {
1158                     i = cell[a]*nbat->fstride;
1159
1160                     f[a][XX] += fnb[i];
1161                     f[a][YY] += fnb[i+1];
1162                     f[a][ZZ] += fnb[i+2];
1163                 }
1164             }
1165             else
1166             {
1167                 for (a = a0; a < a1; a++)
1168                 {
1169                     i = cell[a]*nbat->fstride;
1170
1171                     for (fa = 0; fa < nfa; fa++)
1172                     {
1173                         f[a][XX] += out[fa].f[i];
1174                         f[a][YY] += out[fa].f[i+1];
1175                         f[a][ZZ] += out[fa].f[i+2];
1176                     }
1177                 }
1178             }
1179             break;
1180         case nbatX4:
1181             if (nfa == 1)
1182             {
1183                 fnb = out[0].f;
1184
1185                 for (a = a0; a < a1; a++)
1186                 {
1187                     i = X4_IND_A(cell[a]);
1188
1189                     f[a][XX] += fnb[i+XX*PACK_X4];
1190                     f[a][YY] += fnb[i+YY*PACK_X4];
1191                     f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1192                 }
1193             }
1194             else
1195             {
1196                 for (a = a0; a < a1; a++)
1197                 {
1198                     i = X4_IND_A(cell[a]);
1199
1200                     for (fa = 0; fa < nfa; fa++)
1201                     {
1202                         f[a][XX] += out[fa].f[i+XX*PACK_X4];
1203                         f[a][YY] += out[fa].f[i+YY*PACK_X4];
1204                         f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1205                     }
1206                 }
1207             }
1208             break;
1209         case nbatX8:
1210             if (nfa == 1)
1211             {
1212                 fnb = out[0].f;
1213
1214                 for (a = a0; a < a1; a++)
1215                 {
1216                     i = X8_IND_A(cell[a]);
1217
1218                     f[a][XX] += fnb[i+XX*PACK_X8];
1219                     f[a][YY] += fnb[i+YY*PACK_X8];
1220                     f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1221                 }
1222             }
1223             else
1224             {
1225                 for (a = a0; a < a1; a++)
1226                 {
1227                     i = X8_IND_A(cell[a]);
1228
1229                     for (fa = 0; fa < nfa; fa++)
1230                     {
1231                         f[a][XX] += out[fa].f[i+XX*PACK_X8];
1232                         f[a][YY] += out[fa].f[i+YY*PACK_X8];
1233                         f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1234                     }
1235                 }
1236             }
1237             break;
1238         default:
1239             gmx_incons("Unsupported nbnxn_atomdata_t format");
1240     }
1241 }
1242
1243 /* Add the force array(s) from nbnxn_atomdata_t to f */
1244 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t    nbs,
1245                                     int                     locality,
1246                                     const nbnxn_atomdata_t *nbat,
1247                                     rvec                   *f)
1248 {
1249     int a0 = 0, na = 0;
1250     int nth, th;
1251
1252     nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1253
1254     switch (locality)
1255     {
1256         case eatAll:
1257             a0 = 0;
1258             na = nbs->natoms_nonlocal;
1259             break;
1260         case eatLocal:
1261             a0 = 0;
1262             na = nbs->natoms_local;
1263             break;
1264         case eatNonlocal:
1265             a0 = nbs->natoms_local;
1266             na = nbs->natoms_nonlocal - nbs->natoms_local;
1267             break;
1268     }
1269
1270     nth = gmx_omp_nthreads_get(emntNonbonded);
1271
1272     if (nbat->nout > 1)
1273     {
1274         if (locality != eatAll)
1275         {
1276             gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1277         }
1278
1279         /* Reduce the force thread output buffers into buffer 0, before adding
1280          * them to the, differently ordered, "real" force buffer.
1281          */
1282 #pragma omp parallel for num_threads(nth) schedule(static)
1283         for (th = 0; th < nth; th++)
1284         {
1285             const nbnxn_buffer_flags_t *flags;
1286             int   b0, b1, b;
1287             int   i0, i1;
1288             int   nfptr;
1289             real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1290             int   out;
1291
1292             flags = &nbat->buffer_flags;
1293
1294             /* Calculate the cell-block range for our thread */
1295             b0 = (flags->nflag* th   )/nth;
1296             b1 = (flags->nflag*(th+1))/nth;
1297
1298             for (b = b0; b < b1; b++)
1299             {
1300                 i0 =  b   *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1301                 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1302
1303                 nfptr = 0;
1304                 for (out = 1; out < nbat->nout; out++)
1305                 {
1306                     if (flags->flag[b] & (1U<<out))
1307                     {
1308                         fptr[nfptr++] = nbat->out[out].f;
1309                     }
1310                 }
1311                 if (nfptr > 0)
1312                 {
1313 #ifdef GMX_NBNXN_SIMD
1314                     nbnxn_atomdata_reduce_reals_simd
1315 #else
1316                     nbnxn_atomdata_reduce_reals
1317 #endif
1318                         (nbat->out[0].f,
1319                         flags->flag[b] & (1U<<0),
1320                         fptr, nfptr,
1321                         i0, i1);
1322                 }
1323                 else if (!(flags->flag[b] & (1U<<0)))
1324                 {
1325                     nbnxn_atomdata_clear_reals(nbat->out[0].f,
1326                                                i0, i1);
1327                 }
1328             }
1329         }
1330     }
1331
1332 #pragma omp parallel for num_threads(nth) schedule(static)
1333     for (th = 0; th < nth; th++)
1334     {
1335         nbnxn_atomdata_add_nbat_f_to_f_part(nbs, nbat,
1336                                             nbat->out,
1337                                             1,
1338                                             a0+((th+0)*na)/nth,
1339                                             a0+((th+1)*na)/nth,
1340                                             f);
1341     }
1342
1343     nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1344 }
1345
1346 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1347 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1348                                               rvec                   *fshift)
1349 {
1350     const nbnxn_atomdata_output_t *out;
1351     int  th;
1352     int  s;
1353     rvec sum;
1354
1355     out = nbat->out;
1356
1357     for (s = 0; s < SHIFTS; s++)
1358     {
1359         clear_rvec(sum);
1360         for (th = 0; th < nbat->nout; th++)
1361         {
1362             sum[XX] += out[th].fshift[s*DIM+XX];
1363             sum[YY] += out[th].fshift[s*DIM+YY];
1364             sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1365         }
1366         rvec_inc(fshift[s], sum);
1367     }
1368 }