a82a3b3275a66b23c88cfbbc7986a3e3b0ba0aed
[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 interaction 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 (i=0);
669          * we substract 0.5 to avoid rounding issues.
670          * In the kernel we can subtract 1 to generate the mask for the next i.
671          */
672         const int simd_width = GMX_SIMD_WIDTH_HERE;
673         int       simd_4xn_diag_ind_size, simd_interaction_size, j;
674
675         simd_4xn_diag_ind_size = max(NBNXN_CPU_CLUSTER_I_SIZE, simd_width);
676         snew_aligned(nbat->simd_4xn_diagonal_j_minus_i,
677                      simd_4xn_diag_ind_size, NBNXN_MEM_ALIGN);
678         for (j = 0; j < simd_4xn_diag_ind_size; j++)
679         {
680             nbat->simd_4xn_diagonal_j_minus_i[j] = j - 0.5;
681         }
682
683         snew_aligned(nbat->simd_2xnn_diagonal_j_minus_i,
684                      simd_width, NBNXN_MEM_ALIGN);
685         for (j = 0; j < simd_width/2; j++)
686         {
687             /* The j-cluster size is half the SIMD width */
688             nbat->simd_2xnn_diagonal_j_minus_i[j]              = j - 0.5;
689             /* The next half of the SIMD width is for i + 1 */
690             nbat->simd_2xnn_diagonal_j_minus_i[simd_width/2+j] = j - 1 - 0.5;
691         }
692
693         /* We use up to 32 bits for exclusion masking.
694          * The same masks are used for the 4xN and 2x(N+N) kernels.
695          * The masks are read either into epi32 SIMD registers or into
696          * real SIMD registers (together with a cast).
697          * In single precision this means the real and epi32 SIMD registers
698          * are of equal size.
699          * In double precision the epi32 registers can be smaller than
700          * the real registers, so depending on the architecture, we might
701          * need to use two, identical, 32-bit masks per real.
702          */
703         simd_interaction_size = NBNXN_CPU_CLUSTER_I_SIZE*simd_width;
704         snew_aligned(nbat->simd_exclusion_filter1, simd_interaction_size,   NBNXN_MEM_ALIGN);
705         snew_aligned(nbat->simd_exclusion_filter2, simd_interaction_size*2, NBNXN_MEM_ALIGN);
706         
707         for (j = 0; j < simd_interaction_size; j++)
708         {
709             /* Set the consecutive bits for filters pair exclusions masks */
710             nbat->simd_exclusion_filter1[j]       = (1U << j);
711             nbat->simd_exclusion_filter2[j*2 + 0] = (1U << j);
712             nbat->simd_exclusion_filter2[j*2 + 1] = (1U << j);
713         }
714     }
715 #endif
716
717     /* Initialize the output data structures */
718     nbat->nout    = nout;
719     snew(nbat->out, nbat->nout);
720     nbat->nalloc  = 0;
721     for (i = 0; i < nbat->nout; i++)
722     {
723         nbnxn_atomdata_output_init(&nbat->out[i],
724                                    nb_kernel_type,
725                                    nbat->nenergrp, 1<<nbat->neg_2log,
726                                    nbat->alloc);
727     }
728     nbat->buffer_flags.flag        = NULL;
729     nbat->buffer_flags.flag_nalloc = 0;
730 }
731
732 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
733                                        const int *type, int na,
734                                        real *ljparam_at)
735 {
736     int is, k, i;
737
738     /* The LJ params follow the combination rule:
739      * copy the params for the type array to the atom array.
740      */
741     for (is = 0; is < na; is += PACK_X4)
742     {
743         for (k = 0; k < PACK_X4; k++)
744         {
745             i = is + k;
746             ljparam_at[is*2        +k] = ljparam_type[type[i]*2  ];
747             ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
748         }
749     }
750 }
751
752 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
753                                        const int *type, int na,
754                                        real *ljparam_at)
755 {
756     int is, k, i;
757
758     /* The LJ params follow the combination rule:
759      * copy the params for the type array to the atom array.
760      */
761     for (is = 0; is < na; is += PACK_X8)
762     {
763         for (k = 0; k < PACK_X8; k++)
764         {
765             i = is + k;
766             ljparam_at[is*2        +k] = ljparam_type[type[i]*2  ];
767             ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
768         }
769     }
770 }
771
772 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
773 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t    *nbat,
774                                          int                  ngrid,
775                                          const nbnxn_search_t nbs,
776                                          const int           *type)
777 {
778     int                 g, i, ncz, ash;
779     const nbnxn_grid_t *grid;
780
781     for (g = 0; g < ngrid; g++)
782     {
783         grid = &nbs->grid[g];
784
785         /* Loop over all columns and copy and fill */
786         for (i = 0; i < grid->ncx*grid->ncy; i++)
787         {
788             ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
789             ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
790
791             copy_int_to_nbat_int(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
792                                  type, nbat->ntype-1, nbat->type+ash);
793
794             if (nbat->comb_rule != ljcrNONE)
795             {
796                 if (nbat->XFormat == nbatX4)
797                 {
798                     copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
799                                                nbat->type+ash, ncz*grid->na_sc,
800                                                nbat->lj_comb+ash*2);
801                 }
802                 else if (nbat->XFormat == nbatX8)
803                 {
804                     copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
805                                                nbat->type+ash, ncz*grid->na_sc,
806                                                nbat->lj_comb+ash*2);
807                 }
808             }
809         }
810     }
811 }
812
813 /* Sets the charges in nbnxn_atomdata_t *nbat */
814 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t    *nbat,
815                                        int                  ngrid,
816                                        const nbnxn_search_t nbs,
817                                        const real          *charge)
818 {
819     int                 g, cxy, ncz, ash, na, na_round, i, j;
820     real               *q;
821     const nbnxn_grid_t *grid;
822
823     for (g = 0; g < ngrid; g++)
824     {
825         grid = &nbs->grid[g];
826
827         /* Loop over all columns and copy and fill */
828         for (cxy = 0; cxy < grid->ncx*grid->ncy; cxy++)
829         {
830             ash      = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
831             na       = grid->cxy_na[cxy];
832             na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
833
834             if (nbat->XFormat == nbatXYZQ)
835             {
836                 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
837                 for (i = 0; i < na; i++)
838                 {
839                     *q = charge[nbs->a[ash+i]];
840                     q += STRIDE_XYZQ;
841                 }
842                 /* Complete the partially filled last cell with zeros */
843                 for (; i < na_round; i++)
844                 {
845                     *q = 0;
846                     q += STRIDE_XYZQ;
847                 }
848             }
849             else
850             {
851                 q = nbat->q + ash;
852                 for (i = 0; i < na; i++)
853                 {
854                     *q = charge[nbs->a[ash+i]];
855                     q++;
856                 }
857                 /* Complete the partially filled last cell with zeros */
858                 for (; i < na_round; i++)
859                 {
860                     *q = 0;
861                     q++;
862                 }
863             }
864         }
865     }
866 }
867
868 /* Copies the energy group indices to a reordered and packed array */
869 static void copy_egp_to_nbat_egps(const int *a, int na, int na_round,
870                                   int na_c, int bit_shift,
871                                   const int *in, int *innb)
872 {
873     int i, j, sa, at;
874     int comb;
875
876     j = 0;
877     for (i = 0; i < na; i += na_c)
878     {
879         /* Store na_c energy group numbers into one int */
880         comb = 0;
881         for (sa = 0; sa < na_c; sa++)
882         {
883             at = a[i+sa];
884             if (at >= 0)
885             {
886                 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
887             }
888         }
889         innb[j++] = comb;
890     }
891     /* Complete the partially filled last cell with fill */
892     for (; i < na_round; i += na_c)
893     {
894         innb[j++] = 0;
895     }
896 }
897
898 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
899 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t    *nbat,
900                                             int                  ngrid,
901                                             const nbnxn_search_t nbs,
902                                             const int           *atinfo)
903 {
904     int                 g, i, ncz, ash;
905     const nbnxn_grid_t *grid;
906
907     for (g = 0; g < ngrid; g++)
908     {
909         grid = &nbs->grid[g];
910
911         /* Loop over all columns and copy and fill */
912         for (i = 0; i < grid->ncx*grid->ncy; i++)
913         {
914             ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
915             ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
916
917             copy_egp_to_nbat_egps(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
918                                   nbat->na_c, nbat->neg_2log,
919                                   atinfo, nbat->energrp+(ash>>grid->na_c_2log));
920         }
921     }
922 }
923
924 /* Sets all required atom parameter data in nbnxn_atomdata_t */
925 void nbnxn_atomdata_set(nbnxn_atomdata_t    *nbat,
926                         int                  locality,
927                         const nbnxn_search_t nbs,
928                         const t_mdatoms     *mdatoms,
929                         const int           *atinfo)
930 {
931     int ngrid;
932
933     if (locality == eatLocal)
934     {
935         ngrid = 1;
936     }
937     else
938     {
939         ngrid = nbs->ngrid;
940     }
941
942     nbnxn_atomdata_set_atomtypes(nbat, ngrid, nbs, mdatoms->typeA);
943
944     nbnxn_atomdata_set_charges(nbat, ngrid, nbs, mdatoms->chargeA);
945
946     if (nbat->nenergrp > 1)
947     {
948         nbnxn_atomdata_set_energygroups(nbat, ngrid, nbs, atinfo);
949     }
950 }
951
952 /* Copies the shift vector array to nbnxn_atomdata_t */
953 void nbnxn_atomdata_copy_shiftvec(gmx_bool          bDynamicBox,
954                                   rvec             *shift_vec,
955                                   nbnxn_atomdata_t *nbat)
956 {
957     int i;
958
959     nbat->bDynamicBox = bDynamicBox;
960     for (i = 0; i < SHIFTS; i++)
961     {
962         copy_rvec(shift_vec[i], nbat->shift_vec[i]);
963     }
964 }
965
966 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
967 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
968                                      int                  locality,
969                                      gmx_bool             FillLocal,
970                                      rvec                *x,
971                                      nbnxn_atomdata_t    *nbat)
972 {
973     int g0 = 0, g1 = 0;
974     int nth, th;
975
976     switch (locality)
977     {
978         case eatAll:
979             g0 = 0;
980             g1 = nbs->ngrid;
981             break;
982         case eatLocal:
983             g0 = 0;
984             g1 = 1;
985             break;
986         case eatNonlocal:
987             g0 = 1;
988             g1 = nbs->ngrid;
989             break;
990     }
991
992     if (FillLocal)
993     {
994         nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
995     }
996
997     nth = gmx_omp_nthreads_get(emntPairsearch);
998
999 #pragma omp parallel for num_threads(nth) schedule(static)
1000     for (th = 0; th < nth; th++)
1001     {
1002         int g;
1003
1004         for (g = g0; g < g1; g++)
1005         {
1006             const nbnxn_grid_t *grid;
1007             int                 cxy0, cxy1, cxy;
1008
1009             grid = &nbs->grid[g];
1010
1011             cxy0 = (grid->ncx*grid->ncy* th   +nth-1)/nth;
1012             cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
1013
1014             for (cxy = cxy0; cxy < cxy1; cxy++)
1015             {
1016                 int na, ash, na_fill;
1017
1018                 na  = grid->cxy_na[cxy];
1019                 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1020
1021                 if (g == 0 && FillLocal)
1022                 {
1023                     na_fill =
1024                         (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
1025                 }
1026                 else
1027                 {
1028                     /* We fill only the real particle locations.
1029                      * We assume the filling entries at the end have been
1030                      * properly set before during ns.
1031                      */
1032                     na_fill = na;
1033                 }
1034                 copy_rvec_to_nbat_real(nbs->a+ash, na, na_fill, x,
1035                                        nbat->XFormat, nbat->x, ash,
1036                                        0, 0, 0);
1037             }
1038         }
1039     }
1040 }
1041
1042 static void
1043 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
1044                            int i0, int i1)
1045 {
1046     int i;
1047
1048     for (i = i0; i < i1; i++)
1049     {
1050         dest[i] = 0;
1051     }
1052 }
1053
1054 static void
1055 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1056                             gmx_bool bDestSet,
1057                             real ** gmx_restrict src,
1058                             int nsrc,
1059                             int i0, int i1)
1060 {
1061     int i, s;
1062
1063     if (bDestSet)
1064     {
1065         /* The destination buffer contains data, add to it */
1066         for (i = i0; i < i1; i++)
1067         {
1068             for (s = 0; s < nsrc; s++)
1069             {
1070                 dest[i] += src[s][i];
1071             }
1072         }
1073     }
1074     else
1075     {
1076         /* The destination buffer is unitialized, set it first */
1077         for (i = i0; i < i1; i++)
1078         {
1079             dest[i] = src[0][i];
1080             for (s = 1; s < nsrc; s++)
1081             {
1082                 dest[i] += src[s][i];
1083             }
1084         }
1085     }
1086 }
1087
1088 static void
1089 nbnxn_atomdata_reduce_reals_simd(real * gmx_restrict dest,
1090                                  gmx_bool bDestSet,
1091                                  real ** gmx_restrict src,
1092                                  int nsrc,
1093                                  int i0, int i1)
1094 {
1095 #ifdef GMX_NBNXN_SIMD
1096 /* The SIMD width here is actually independent of that in the kernels,
1097  * but we use the same width for simplicity (usually optimal anyhow).
1098  */
1099     int       i, s;
1100     gmx_mm_pr dest_SSE, src_SSE;
1101
1102     if (bDestSet)
1103     {
1104         for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1105         {
1106             dest_SSE = gmx_load_pr(dest+i);
1107             for (s = 0; s < nsrc; s++)
1108             {
1109                 src_SSE  = gmx_load_pr(src[s]+i);
1110                 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1111             }
1112             gmx_store_pr(dest+i, dest_SSE);
1113         }
1114     }
1115     else
1116     {
1117         for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1118         {
1119             dest_SSE = gmx_load_pr(src[0]+i);
1120             for (s = 1; s < nsrc; s++)
1121             {
1122                 src_SSE  = gmx_load_pr(src[s]+i);
1123                 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1124             }
1125             gmx_store_pr(dest+i, dest_SSE);
1126         }
1127     }
1128 #endif
1129 }
1130
1131 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1132 static void
1133 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1134                                     const nbnxn_atomdata_t *nbat,
1135                                     nbnxn_atomdata_output_t *out,
1136                                     int nfa,
1137                                     int a0, int a1,
1138                                     rvec *f)
1139 {
1140     int         a, i, fa;
1141     const int  *cell;
1142     const real *fnb;
1143
1144     cell = nbs->cell;
1145
1146     /* Loop over all columns and copy and fill */
1147     switch (nbat->FFormat)
1148     {
1149         case nbatXYZ:
1150         case nbatXYZQ:
1151             if (nfa == 1)
1152             {
1153                 fnb = out[0].f;
1154
1155                 for (a = a0; a < a1; a++)
1156                 {
1157                     i = cell[a]*nbat->fstride;
1158
1159                     f[a][XX] += fnb[i];
1160                     f[a][YY] += fnb[i+1];
1161                     f[a][ZZ] += fnb[i+2];
1162                 }
1163             }
1164             else
1165             {
1166                 for (a = a0; a < a1; a++)
1167                 {
1168                     i = cell[a]*nbat->fstride;
1169
1170                     for (fa = 0; fa < nfa; fa++)
1171                     {
1172                         f[a][XX] += out[fa].f[i];
1173                         f[a][YY] += out[fa].f[i+1];
1174                         f[a][ZZ] += out[fa].f[i+2];
1175                     }
1176                 }
1177             }
1178             break;
1179         case nbatX4:
1180             if (nfa == 1)
1181             {
1182                 fnb = out[0].f;
1183
1184                 for (a = a0; a < a1; a++)
1185                 {
1186                     i = X4_IND_A(cell[a]);
1187
1188                     f[a][XX] += fnb[i+XX*PACK_X4];
1189                     f[a][YY] += fnb[i+YY*PACK_X4];
1190                     f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1191                 }
1192             }
1193             else
1194             {
1195                 for (a = a0; a < a1; a++)
1196                 {
1197                     i = X4_IND_A(cell[a]);
1198
1199                     for (fa = 0; fa < nfa; fa++)
1200                     {
1201                         f[a][XX] += out[fa].f[i+XX*PACK_X4];
1202                         f[a][YY] += out[fa].f[i+YY*PACK_X4];
1203                         f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1204                     }
1205                 }
1206             }
1207             break;
1208         case nbatX8:
1209             if (nfa == 1)
1210             {
1211                 fnb = out[0].f;
1212
1213                 for (a = a0; a < a1; a++)
1214                 {
1215                     i = X8_IND_A(cell[a]);
1216
1217                     f[a][XX] += fnb[i+XX*PACK_X8];
1218                     f[a][YY] += fnb[i+YY*PACK_X8];
1219                     f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1220                 }
1221             }
1222             else
1223             {
1224                 for (a = a0; a < a1; a++)
1225                 {
1226                     i = X8_IND_A(cell[a]);
1227
1228                     for (fa = 0; fa < nfa; fa++)
1229                     {
1230                         f[a][XX] += out[fa].f[i+XX*PACK_X8];
1231                         f[a][YY] += out[fa].f[i+YY*PACK_X8];
1232                         f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1233                     }
1234                 }
1235             }
1236             break;
1237         default:
1238             gmx_incons("Unsupported nbnxn_atomdata_t format");
1239     }
1240 }
1241
1242 /* Add the force array(s) from nbnxn_atomdata_t to f */
1243 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t    nbs,
1244                                     int                     locality,
1245                                     const nbnxn_atomdata_t *nbat,
1246                                     rvec                   *f)
1247 {
1248     int a0 = 0, na = 0;
1249     int nth, th;
1250
1251     nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1252
1253     switch (locality)
1254     {
1255         case eatAll:
1256             a0 = 0;
1257             na = nbs->natoms_nonlocal;
1258             break;
1259         case eatLocal:
1260             a0 = 0;
1261             na = nbs->natoms_local;
1262             break;
1263         case eatNonlocal:
1264             a0 = nbs->natoms_local;
1265             na = nbs->natoms_nonlocal - nbs->natoms_local;
1266             break;
1267     }
1268
1269     nth = gmx_omp_nthreads_get(emntNonbonded);
1270
1271     if (nbat->nout > 1)
1272     {
1273         if (locality != eatAll)
1274         {
1275             gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1276         }
1277
1278         /* Reduce the force thread output buffers into buffer 0, before adding
1279          * them to the, differently ordered, "real" force buffer.
1280          */
1281 #pragma omp parallel for num_threads(nth) schedule(static)
1282         for (th = 0; th < nth; th++)
1283         {
1284             const nbnxn_buffer_flags_t *flags;
1285             int   b0, b1, b;
1286             int   i0, i1;
1287             int   nfptr;
1288             real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1289             int   out;
1290
1291             flags = &nbat->buffer_flags;
1292
1293             /* Calculate the cell-block range for our thread */
1294             b0 = (flags->nflag* th   )/nth;
1295             b1 = (flags->nflag*(th+1))/nth;
1296
1297             for (b = b0; b < b1; b++)
1298             {
1299                 i0 =  b   *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1300                 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1301
1302                 nfptr = 0;
1303                 for (out = 1; out < nbat->nout; out++)
1304                 {
1305                     if (flags->flag[b] & (1U<<out))
1306                     {
1307                         fptr[nfptr++] = nbat->out[out].f;
1308                     }
1309                 }
1310                 if (nfptr > 0)
1311                 {
1312 #ifdef GMX_NBNXN_SIMD
1313                     nbnxn_atomdata_reduce_reals_simd
1314 #else
1315                     nbnxn_atomdata_reduce_reals
1316 #endif
1317                         (nbat->out[0].f,
1318                         flags->flag[b] & (1U<<0),
1319                         fptr, nfptr,
1320                         i0, i1);
1321                 }
1322                 else if (!(flags->flag[b] & (1U<<0)))
1323                 {
1324                     nbnxn_atomdata_clear_reals(nbat->out[0].f,
1325                                                i0, i1);
1326                 }
1327             }
1328         }
1329     }
1330
1331 #pragma omp parallel for num_threads(nth) schedule(static)
1332     for (th = 0; th < nth; th++)
1333     {
1334         nbnxn_atomdata_add_nbat_f_to_f_part(nbs, nbat,
1335                                             nbat->out,
1336                                             1,
1337                                             a0+((th+0)*na)/nth,
1338                                             a0+((th+1)*na)/nth,
1339                                             f);
1340     }
1341
1342     nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1343 }
1344
1345 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1346 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1347                                               rvec                   *fshift)
1348 {
1349     const nbnxn_atomdata_output_t *out;
1350     int  th;
1351     int  s;
1352     rvec sum;
1353
1354     out = nbat->out;
1355
1356     for (s = 0; s < SHIFTS; s++)
1357     {
1358         clear_rvec(sum);
1359         for (th = 0; th < nbat->nout; th++)
1360         {
1361             sum[XX] += out[th].fshift[s*DIM+XX];
1362             sum[YY] += out[th].fshift[s*DIM+YY];
1363             sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1364         }
1365         rvec_inc(fshift[s], sum);
1366     }
1367 }