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 == nbk4xN_X86_SIMD128 ||
149         nb_kernel_type == nbk4xN_X86_SIMD256)
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         switch (nb_kernel_type)
596         {
597         case nbk4xN_X86_SIMD128:
598             nbat->XFormat = nbatX4;
599             break;
600         case nbk4xN_X86_SIMD256:
601 #ifndef GMX_DOUBLE
602             nbat->XFormat = nbatX8;
603 #else
604             nbat->XFormat = nbatX4;
605 #endif
606             break;
607         default:
608             nbat->XFormat = nbatXYZ;
609             break;
610         }
611
612         nbat->FFormat = nbat->XFormat;
613     }
614     else
615     {
616         nbat->XFormat = nbatXYZQ;
617         nbat->FFormat = nbatXYZ;
618     }
619     nbat->q       = NULL;
620     nbat->nenergrp = n_energygroups;
621     if (!simple)
622     {
623         /* Energy groups not supported yet for super-sub lists */
624         if (n_energygroups > 1 && fp != NULL)
625         {
626             fprintf(fp,"\nNOTE: With GPUs, reporting energy group contributions is not supported\n\n");
627         }
628         nbat->nenergrp = 1;
629     }
630     /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
631     if (nbat->nenergrp > 64)
632     {
633         gmx_fatal(FARGS,"With NxN kernels not more than 64 energy groups are supported\n");
634     }
635     nbat->neg_2log = 1;
636     while (nbat->nenergrp > (1<<nbat->neg_2log))
637     {
638         nbat->neg_2log++;
639     }
640     nbat->energrp = NULL;
641     nbat->alloc((void **)&nbat->shift_vec,SHIFTS*sizeof(*nbat->shift_vec));
642     nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
643     nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
644     nbat->x       = NULL;
645     nbat->nout    = nout;
646     snew(nbat->out,nbat->nout);
647     nbat->nalloc  = 0;
648     for(i=0; i<nbat->nout; i++)
649     {
650         nbnxn_atomdata_output_init(&nbat->out[i],
651                                    nb_kernel_type,
652                                    nbat->nenergrp,1<<nbat->neg_2log,
653                                    nbat->alloc);
654     }
655 }
656
657 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
658                                        const int *type,int na,
659                                        real *ljparam_at)
660 {
661     int is,k,i;
662
663     /* The LJ params follow the combination rule:
664      * copy the params for the type array to the atom array.
665      */
666     for(is=0; is<na; is+=PACK_X4)
667     {
668         for(k=0; k<PACK_X4; k++)
669         {
670             i = is + k;
671             ljparam_at[is*2        +k] = ljparam_type[type[i]*2  ];
672             ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
673         }
674     }
675 }
676
677 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
678                                        const int *type,int na,
679                                        real *ljparam_at)
680 {
681     int is,k,i;
682
683     /* The LJ params follow the combination rule:
684      * copy the params for the type array to the atom array.
685      */
686     for(is=0; is<na; is+=PACK_X8)
687     {
688         for(k=0; k<PACK_X8; k++)
689         {
690             i = is + k;
691             ljparam_at[is*2        +k] = ljparam_type[type[i]*2  ];
692             ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
693         }
694     }
695 }
696
697 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
698 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t *nbat,
699                                          int ngrid,
700                                          const nbnxn_search_t nbs,
701                                          const int *type)
702 {
703     int g,i,ncz,ash;
704     const nbnxn_grid_t *grid;
705
706     for(g=0; g<ngrid; g++)
707     {
708         grid = &nbs->grid[g];
709
710         /* Loop over all columns and copy and fill */
711         for(i=0; i<grid->ncx*grid->ncy; i++)
712         {
713             ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
714             ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
715
716             copy_int_to_nbat_int(nbs->a+ash,grid->cxy_na[i],ncz*grid->na_sc,
717                                  type,nbat->ntype-1,nbat->type+ash);
718
719             if (nbat->comb_rule != ljcrNONE)
720             {
721                 if (nbat->XFormat == nbatX4)
722                 {
723                     copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
724                                                nbat->type+ash,ncz*grid->na_sc,
725                                                nbat->lj_comb+ash*2);
726                 }
727                 else if (nbat->XFormat == nbatX8)
728                 {
729                     copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
730                                                nbat->type+ash,ncz*grid->na_sc,
731                                                nbat->lj_comb+ash*2);
732                 }
733             }
734         }
735     }
736 }
737
738 /* Sets the charges in nbnxn_atomdata_t *nbat */
739 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t *nbat,
740                                        int ngrid,
741                                        const nbnxn_search_t nbs,
742                                        const real *charge)
743 {
744     int  g,cxy,ncz,ash,na,na_round,i,j;
745     real *q;
746     const nbnxn_grid_t *grid;
747
748     for(g=0; g<ngrid; g++)
749     {
750         grid = &nbs->grid[g];
751
752         /* Loop over all columns and copy and fill */
753         for(cxy=0; cxy<grid->ncx*grid->ncy; cxy++)
754         {
755             ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
756             na  = grid->cxy_na[cxy];
757             na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
758
759             if (nbat->XFormat == nbatXYZQ)
760             {
761                 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
762                 for(i=0; i<na; i++)
763                 {
764                     *q = charge[nbs->a[ash+i]];
765                     q += STRIDE_XYZQ;
766                 }
767                 /* Complete the partially filled last cell with zeros */
768                 for(; i<na_round; i++)
769                 {
770                     *q = 0;
771                     q += STRIDE_XYZQ;
772                 }
773             }
774             else
775             {
776                 q = nbat->q + ash;
777                 for(i=0; i<na; i++)
778                 {
779                     *q = charge[nbs->a[ash+i]];
780                     q++;
781                 }
782                 /* Complete the partially filled last cell with zeros */
783                 for(; i<na_round; i++)
784                 {
785                     *q = 0;
786                     q++;
787                 }
788             }
789         }
790     }
791 }
792
793 /* Copies the energy group indices to a reordered and packed array */
794 static void copy_egp_to_nbat_egps(const int *a,int na,int na_round,
795                                   int na_c,int bit_shift,
796                                   const int *in,int *innb)
797 {
798     int i,j,sa,at;
799     int comb;
800
801     j = 0;
802     for(i=0; i<na; i+=na_c)
803     {
804         /* Store na_c energy group numbers into one int */
805         comb = 0;
806         for(sa=0; sa<na_c; sa++)
807         {
808             at = a[i+sa];
809             if (at >= 0)
810             {
811                 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
812             }
813         }
814         innb[j++] = comb;
815     }
816     /* Complete the partially filled last cell with fill */
817     for(; i<na_round; i+=na_c)
818     {
819         innb[j++] = 0;
820     }
821 }
822
823 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
824 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t *nbat,
825                                             int ngrid,
826                                             const nbnxn_search_t nbs,
827                                             const int *atinfo)
828 {
829     int g,i,ncz,ash;
830     const nbnxn_grid_t *grid;
831
832     for(g=0; g<ngrid; g++)
833     {
834         grid = &nbs->grid[g];
835
836         /* Loop over all columns and copy and fill */
837         for(i=0; i<grid->ncx*grid->ncy; i++)
838         {
839             ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
840             ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
841
842             copy_egp_to_nbat_egps(nbs->a+ash,grid->cxy_na[i],ncz*grid->na_sc,
843                                   nbat->na_c,nbat->neg_2log,
844                                   atinfo,nbat->energrp+(ash>>grid->na_c_2log));
845         }
846     }
847 }
848
849 /* Sets all required atom parameter data in nbnxn_atomdata_t */
850 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
851                         int locality,
852                         const nbnxn_search_t nbs,
853                         const t_mdatoms *mdatoms,
854                         const int *atinfo)
855 {
856     int ngrid;
857
858     if (locality == eatLocal)
859     {
860         ngrid = 1;
861     }
862     else
863     {
864         ngrid = nbs->ngrid;
865     }
866
867     nbnxn_atomdata_set_atomtypes(nbat,ngrid,nbs,mdatoms->typeA);
868
869     nbnxn_atomdata_set_charges(nbat,ngrid,nbs,mdatoms->chargeA);
870
871     if (nbat->nenergrp > 1)
872     {
873         nbnxn_atomdata_set_energygroups(nbat,ngrid,nbs,atinfo);
874     }
875 }
876
877 /* Copies the shift vector array to nbnxn_atomdata_t */
878 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox,
879                                    rvec *shift_vec,
880                                    nbnxn_atomdata_t *nbat)
881 {
882     int i;
883
884     nbat->bDynamicBox = bDynamicBox;
885     for(i=0; i<SHIFTS; i++)
886     {
887         copy_rvec(shift_vec[i],nbat->shift_vec[i]);
888     }
889 }
890
891 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
892 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
893                                       int locality,
894                                       gmx_bool FillLocal,
895                                       rvec *x,
896                                       nbnxn_atomdata_t *nbat)
897 {
898     int g0=0,g1=0;
899     int nth,th;
900
901     switch (locality)
902     {
903     case eatAll:
904         g0 = 0;
905         g1 = nbs->ngrid;
906         break;
907     case eatLocal:
908         g0 = 0;
909         g1 = 1;
910         break;
911     case eatNonlocal:
912         g0 = 1;
913         g1 = nbs->ngrid;
914         break;
915     }
916
917     if (FillLocal)
918     {
919         nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
920     }
921
922     nth = gmx_omp_nthreads_get(emntPairsearch);
923
924 #pragma omp parallel for num_threads(nth) schedule(static)
925     for(th=0; th<nth; th++)
926     {
927         int g;
928
929         for(g=g0; g<g1; g++)
930         {
931             const nbnxn_grid_t *grid;
932             int cxy0,cxy1,cxy;
933
934             grid = &nbs->grid[g];
935
936             cxy0 = (grid->ncx*grid->ncy* th   +nth-1)/nth;
937             cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
938
939             for(cxy=cxy0; cxy<cxy1; cxy++)
940             {
941                 int na,ash,na_fill;
942
943                 na  = grid->cxy_na[cxy];
944                 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
945
946                 if (g == 0 && FillLocal)
947                 {
948                     na_fill =
949                         (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
950                 }
951                 else
952                 {
953                     /* We fill only the real particle locations.
954                      * We assume the filling entries at the end have been
955                      * properly set before during ns.
956                      */
957                     na_fill = na;
958                 }
959                 copy_rvec_to_nbat_real(nbs->a+ash,na,na_fill,x,
960                                        nbat->XFormat,nbat->x,ash,
961                                        0,0,0);
962             }
963         }
964     }
965 }
966
967 static void
968 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
969                             real ** gmx_restrict src,
970                             int nsrc,
971                             int i0, int i1)
972 {
973     int i,s;
974
975     for(i=i0; i<i1; i++)
976     {
977         for(s=0; s<nsrc; s++)
978         {
979             dest[i] += src[s][i];
980         }
981     }
982 }
983
984 static void
985 nbnxn_atomdata_reduce_reals_x86_simd(real * gmx_restrict dest,
986                                      real ** gmx_restrict src,
987                                      int nsrc,
988                                      int i0, int i1)
989 {
990 #ifdef NBNXN_SEARCH_SSE
991 /* We can use AVX256 here, but not when AVX128 kernels are selected.
992  * As this reduction is not faster with AVX256 anyway, we use 128-bit SIMD.
993  */
994 #define GMX_MM128_HERE
995 #include "gmx_x86_simd_macros.h"
996
997     int       i,s;
998     gmx_mm_pr dest_SSE,src_SSE;
999
1000     if ((i0 & (GMX_X86_SIMD_WIDTH_HERE-1)) ||
1001         (i1 & (GMX_X86_SIMD_WIDTH_HERE-1)))
1002     {
1003         gmx_incons("bounds not a multiple of GMX_X86_SIMD_WIDTH_HERE in nbnxn_atomdata_reduce_reals_x86_simd");
1004     }
1005
1006     for(i=i0; i<i1; i+=GMX_X86_SIMD_WIDTH_HERE)
1007     {
1008         dest_SSE = gmx_load_pr(dest+i);
1009         for(s=0; s<nsrc; s++)
1010         {
1011             src_SSE  = gmx_load_pr(src[s]+i);
1012             dest_SSE = gmx_add_pr(dest_SSE,src_SSE);
1013         }
1014         gmx_store_pr(dest+i,dest_SSE);
1015     }
1016
1017 #undef GMX_MM128_HERE
1018 #undef GMX_MM256_HERE
1019 #endif
1020 }
1021
1022 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1023 static void
1024 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1025                                     const nbnxn_atomdata_t *nbat,
1026                                     nbnxn_atomdata_output_t *out,
1027                                     int nfa,
1028                                     int a0,int a1,
1029                                     rvec *f)
1030 {
1031     int  a,i,fa;
1032     const int  *cell;
1033     const real *fnb;
1034
1035     cell = nbs->cell;
1036
1037     /* Loop over all columns and copy and fill */
1038     switch (nbat->FFormat)
1039     {
1040     case nbatXYZ:
1041     case nbatXYZQ:
1042         if (nfa == 1)
1043         {
1044             fnb = out[0].f;
1045
1046             for(a=a0; a<a1; a++)
1047             {
1048                 i = cell[a]*nbat->fstride;
1049
1050                 f[a][XX] += fnb[i];
1051                 f[a][YY] += fnb[i+1];
1052                 f[a][ZZ] += fnb[i+2];
1053             }
1054         }
1055         else
1056         {
1057             for(a=a0; a<a1; a++)
1058             {
1059                 i = cell[a]*nbat->fstride;
1060
1061                 for(fa=0; fa<nfa; fa++)
1062                 {
1063                     f[a][XX] += out[fa].f[i];
1064                     f[a][YY] += out[fa].f[i+1];
1065                     f[a][ZZ] += out[fa].f[i+2];
1066                 }
1067             }
1068         }
1069         break;
1070     case nbatX4:
1071         if (nfa == 1)
1072         {
1073             fnb = out[0].f;
1074
1075             for(a=a0; a<a1; a++)
1076             {
1077                 i = X4_IND_A(cell[a]);
1078
1079                 f[a][XX] += fnb[i+XX*PACK_X4];
1080                 f[a][YY] += fnb[i+YY*PACK_X4];
1081                 f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1082             }
1083         }
1084         else
1085         {
1086             for(a=a0; a<a1; a++)
1087             {
1088                 i = X4_IND_A(cell[a]);
1089                 
1090                 for(fa=0; fa<nfa; fa++)
1091                 {
1092                     f[a][XX] += out[fa].f[i+XX*PACK_X4];
1093                     f[a][YY] += out[fa].f[i+YY*PACK_X4];
1094                     f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1095                 }
1096             }
1097         }
1098         break;
1099     case nbatX8:
1100         if (nfa == 1)
1101         {
1102             fnb = out[0].f;
1103
1104             for(a=a0; a<a1; a++)
1105             {
1106                 i = X8_IND_A(cell[a]);
1107
1108                 f[a][XX] += fnb[i+XX*PACK_X8];
1109                 f[a][YY] += fnb[i+YY*PACK_X8];
1110                 f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1111             }
1112         }
1113         else
1114         {
1115             for(a=a0; a<a1; a++)
1116             {
1117                 i = X8_IND_A(cell[a]);
1118                 
1119                 for(fa=0; fa<nfa; fa++)
1120                 {
1121                     f[a][XX] += out[fa].f[i+XX*PACK_X8];
1122                     f[a][YY] += out[fa].f[i+YY*PACK_X8];
1123                     f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1124                 }
1125             }
1126         }
1127         break;
1128     }
1129 }
1130
1131 /* Add the force array(s) from nbnxn_atomdata_t to f */
1132 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t nbs,
1133                                     int locality,
1134                                     const nbnxn_atomdata_t *nbat,
1135                                     rvec *f)
1136 {
1137     int a0=0,na=0;
1138     int nth,th;
1139     gmx_bool bStreamingReduce;
1140
1141     nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1142
1143     switch (locality)
1144     {
1145     case eatAll:
1146         a0 = 0;
1147         na = nbs->natoms_nonlocal;
1148         break;
1149     case eatLocal:
1150         a0 = 0;
1151         na = nbs->natoms_local;
1152         break;
1153     case eatNonlocal:
1154         a0 = nbs->natoms_local;
1155         na = nbs->natoms_nonlocal - nbs->natoms_local;
1156         break;
1157     }
1158
1159     nth = gmx_omp_nthreads_get(emntNonbonded);
1160
1161     /* Using the two-step streaming reduction is probably always faster */
1162     bStreamingReduce = (nbat->nout > 1);
1163
1164     if (bStreamingReduce)
1165     {
1166         /* Reduce the force thread output buffers into buffer 0, before adding
1167          * them to the, differently ordered, "real" force buffer.
1168          */
1169 #pragma omp parallel for num_threads(nth) schedule(static)
1170         for(th=0; th<nth; th++)
1171         {
1172             int g0,g1,g;
1173
1174             /* For which grids should we reduce the force output? */
1175             g0 = ((locality==eatLocal || locality==eatAll) ? 0 : 1);
1176             g1 = (locality==eatLocal ? 1 : nbs->ngrid);
1177
1178             for(g=g0; g<g1; g++)
1179             {
1180                 nbnxn_grid_t *grid;
1181                 int b0,b1,b;
1182                 int c0,c1,i0,i1;
1183                 int nfptr;
1184                 real *fptr[NBNXN_CELLBLOCK_MAX_THREADS];
1185                 int out;
1186
1187                 grid = &nbs->grid[g];
1188
1189                 /* Calculate the cell-block range for our thread */
1190                 b0 = (grid->cellblock_flags.ncb* th   )/nth;
1191                 b1 = (grid->cellblock_flags.ncb*(th+1))/nth;
1192
1193                 if (grid->cellblock_flags.bUse)
1194                 {
1195                     for(b=b0; b<b1; b++)
1196                     {
1197                         c0 = b*NBNXN_CELLBLOCK_SIZE;
1198                         c1 = min(c0 + NBNXN_CELLBLOCK_SIZE,grid->nc);
1199                         i0 = (grid->cell0 + c0)*grid->na_c*nbat->fstride;
1200                         i1 = (grid->cell0 + c1)*grid->na_c*nbat->fstride;
1201
1202                         nfptr = 0;
1203                         for(out=1; out<nbat->nout; out++)
1204                         {
1205                             if (grid->cellblock_flags.flag[b] & (1U<<out))
1206                             {
1207                                 fptr[nfptr++] = nbat->out[out].f;
1208                             }
1209                         }
1210                         if (nfptr > 0)
1211                         {
1212 #ifdef NBNXN_SEARCH_SSE
1213                             nbnxn_atomdata_reduce_reals_x86_simd
1214 #else
1215                             nbnxn_atomdata_reduce_reals
1216 #endif
1217                                                        (nbat->out[0].f,
1218                                                         fptr,nfptr,
1219                                                         i0,i1);
1220                         }
1221                     }
1222                 }
1223                 else
1224                 {
1225                     c0 = b0*NBNXN_CELLBLOCK_SIZE;
1226                     c1 = min(b1*NBNXN_CELLBLOCK_SIZE,grid->nc);
1227                     i0 = (grid->cell0 + c0)*grid->na_c*nbat->fstride;
1228                     i1 = (grid->cell0 + c1)*grid->na_c*nbat->fstride;
1229
1230                     nfptr = 0;
1231                     for(out=1; out<nbat->nout; out++)
1232                     {
1233                         fptr[nfptr++] = nbat->out[out].f;
1234                     }
1235
1236 #ifdef NBNXN_SEARCH_SSE
1237                     nbnxn_atomdata_reduce_reals_x86_simd
1238 #else
1239                     nbnxn_atomdata_reduce_reals
1240 #endif
1241                                                (nbat->out[0].f,
1242                                                 fptr,nfptr,
1243                                                 i0,i1);
1244                 }
1245             }
1246         }
1247     }
1248
1249 #pragma omp parallel for num_threads(nth) schedule(static)
1250     for(th=0; th<nth; th++)
1251     {
1252         nbnxn_atomdata_add_nbat_f_to_f_part(nbs,nbat,
1253                                             nbat->out,
1254                                             bStreamingReduce ? 1 : nbat->nout,
1255                                             a0+((th+0)*na)/nth,
1256                                             a0+((th+1)*na)/nth,
1257                                             f);
1258     }
1259
1260     nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1261 }
1262
1263 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1264 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1265                                               rvec *fshift)
1266 {
1267     const nbnxn_atomdata_output_t *out;
1268     int  th;
1269     int  s;
1270     rvec sum;
1271
1272     out = nbat->out;
1273     
1274     for(s=0; s<SHIFTS; s++)
1275     {
1276         clear_rvec(sum);
1277         for(th=0; th<nbat->nout; th++)
1278         {
1279             sum[XX] += out[th].fshift[s*DIM+XX];
1280             sum[YY] += out[th].fshift[s*DIM+YY];
1281             sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1282         }
1283         rvec_inc(fshift[s],sum);
1284     }
1285 }