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