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