Merge remote-tracking branch 'origin/release-4-6' into HEAD
[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  ];
485         c12 = nbfp[(i*ntype+i)*2+1];
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                 /* We store the prefactor in the derivative of the potential
510                  * in the parameter to avoid multiplications in the inner loop.
511                  */
512                 c6  = nbfp[(i*ntype+j)*2  ];
513                 c12 = nbfp[(i*ntype+j)*2+1];
514                 nbat->nbfp[(i*nbat->ntype+j)*2  ] =  6.0*c6;
515                 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 12.0*c12;
516
517                 bCombGeom = bCombGeom &&
518                     gmx_within_tol(c6*c6  ,nbfp[(i*ntype+i)*2  ]*nbfp[(j*ntype+j)*2  ],tol) &&
519                     gmx_within_tol(c12*c12,nbfp[(i*ntype+i)*2+1]*nbfp[(j*ntype+j)*2+1],tol);
520
521                 bCombLB = bCombLB &&
522                     ((c6 == 0 && c12 == 0 &&
523                       (nbat->nbfp_comb[i*2+1] == 0 || nbat->nbfp_comb[j*2+1] == 0)) ||
524                      (c6 > 0 && c12 > 0 &&
525                       gmx_within_tol(pow(c12/c6,1.0/6.0),0.5*(nbat->nbfp_comb[i*2]+nbat->nbfp_comb[j*2]),tol) &&
526                       gmx_within_tol(0.25*c6*c6/c12,sqrt(nbat->nbfp_comb[i*2+1]*nbat->nbfp_comb[j*2+1]),tol)));
527             }
528             else
529             {
530                 /* Add zero parameters for the additional dummy atom type */
531                 nbat->nbfp[(i*nbat->ntype+j)*2  ] = 0;
532                 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 0;
533             }
534         }
535     }
536     if (debug)
537     {
538         fprintf(debug,"Combination rules: geometric %d Lorentz-Berthelot %d\n",
539                 bCombGeom,bCombLB);
540     }
541
542     simple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
543
544     if (simple)
545     {
546         /* We prefer the geometic combination rule,
547          * as that gives a slightly faster kernel than the LB rule.
548          */
549         if (bCombGeom)
550         {
551             nbat->comb_rule = ljcrGEOM;
552         }
553         else if (bCombLB)
554         {
555             nbat->comb_rule = ljcrLB;
556         }
557         else
558         {
559             nbat->comb_rule = ljcrNONE;
560
561             nbat->free(nbat->nbfp_comb);
562         }
563
564         if (fp)
565         {
566             if (nbat->comb_rule == ljcrNONE)
567             {
568                 fprintf(fp,"Using full Lennard-Jones parameter combination matrix\n\n");
569             }
570             else
571             {
572                 fprintf(fp,"Using %s Lennard-Jones combination rule\n\n",
573                         nbat->comb_rule==ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
574             }
575         }
576
577         set_combination_rule_data(nbat);
578     }
579     else
580     {
581         nbat->comb_rule = ljcrNONE;
582
583         nbat->free(nbat->nbfp_comb);
584     }
585
586     nbat->natoms  = 0;
587     nbat->type    = NULL;
588     nbat->lj_comb = NULL;
589     if (simple)
590     {
591         switch (nb_kernel_type)
592         {
593         case nbk4xN_X86_SIMD128:
594             nbat->XFormat = nbatX4;
595             break;
596         case nbk4xN_X86_SIMD256:
597 #ifndef GMX_DOUBLE
598             nbat->XFormat = nbatX8;
599 #else
600             nbat->XFormat = nbatX4;
601 #endif
602             break;
603         default:
604             nbat->XFormat = nbatXYZ;
605             break;
606         }
607
608         nbat->FFormat = nbat->XFormat;
609     }
610     else
611     {
612         nbat->XFormat = nbatXYZQ;
613         nbat->FFormat = nbatXYZ;
614     }
615     nbat->q       = NULL;
616     nbat->nenergrp = n_energygroups;
617     if (!simple)
618     {
619         /* Energy groups not supported yet for super-sub lists */
620         nbat->nenergrp = 1;
621     }
622     /* Temporary storage goes is #grp^3*8 real, so limit to 64 */
623     if (nbat->nenergrp > 64)
624     {
625         gmx_fatal(FARGS,"With NxN kernels not more than 64 energy groups are supported\n");
626     }
627     nbat->neg_2log = 1;
628     while (nbat->nenergrp > (1<<nbat->neg_2log))
629     {
630         nbat->neg_2log++;
631     }
632     nbat->energrp = NULL;
633     nbat->alloc((void **)&nbat->shift_vec,SHIFTS*sizeof(*nbat->shift_vec));
634     nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
635     nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
636     nbat->x       = NULL;
637     nbat->nout    = nout;
638     snew(nbat->out,nbat->nout);
639     nbat->nalloc  = 0;
640     for(i=0; i<nbat->nout; i++)
641     {
642         nbnxn_atomdata_output_init(&nbat->out[i],
643                                    nb_kernel_type,
644                                    nbat->nenergrp,1<<nbat->neg_2log,
645                                    nbat->alloc);
646     }
647 }
648
649 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
650                                        const int *type,int na,
651                                        real *ljparam_at)
652 {
653     int is,k,i;
654
655     /* The LJ params follow the combination rule:
656      * copy the params for the type array to the atom array.
657      */
658     for(is=0; is<na; is+=PACK_X4)
659     {
660         for(k=0; k<PACK_X4; k++)
661         {
662             i = is + k;
663             ljparam_at[is*2        +k] = ljparam_type[type[i]*2  ];
664             ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
665         }
666     }
667 }
668
669 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
670                                        const int *type,int na,
671                                        real *ljparam_at)
672 {
673     int is,k,i;
674
675     /* The LJ params follow the combination rule:
676      * copy the params for the type array to the atom array.
677      */
678     for(is=0; is<na; is+=PACK_X8)
679     {
680         for(k=0; k<PACK_X8; k++)
681         {
682             i = is + k;
683             ljparam_at[is*2        +k] = ljparam_type[type[i]*2  ];
684             ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
685         }
686     }
687 }
688
689 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
690 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t *nbat,
691                                          int ngrid,
692                                          const nbnxn_search_t nbs,
693                                          const int *type)
694 {
695     int g,i,ncz,ash;
696     const nbnxn_grid_t *grid;
697
698     for(g=0; g<ngrid; g++)
699     {
700         grid = &nbs->grid[g];
701
702         /* Loop over all columns and copy and fill */
703         for(i=0; i<grid->ncx*grid->ncy; i++)
704         {
705             ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
706             ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
707
708             copy_int_to_nbat_int(nbs->a+ash,grid->cxy_na[i],ncz*grid->na_sc,
709                                  type,nbat->ntype-1,nbat->type+ash);
710
711             if (nbat->comb_rule != ljcrNONE)
712             {
713                 if (nbat->XFormat == nbatX4)
714                 {
715                     copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
716                                                nbat->type+ash,ncz*grid->na_sc,
717                                                nbat->lj_comb+ash*2);
718                 }
719                 else if (nbat->XFormat == nbatX8)
720                 {
721                     copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
722                                                nbat->type+ash,ncz*grid->na_sc,
723                                                nbat->lj_comb+ash*2);
724                 }
725             }
726         }
727     }
728 }
729
730 /* Sets the charges in nbnxn_atomdata_t *nbat */
731 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t *nbat,
732                                        int ngrid,
733                                        const nbnxn_search_t nbs,
734                                        const real *charge)
735 {
736     int  g,cxy,ncz,ash,na,na_round,i,j;
737     real *q;
738     const nbnxn_grid_t *grid;
739
740     for(g=0; g<ngrid; g++)
741     {
742         grid = &nbs->grid[g];
743
744         /* Loop over all columns and copy and fill */
745         for(cxy=0; cxy<grid->ncx*grid->ncy; cxy++)
746         {
747             ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
748             na  = grid->cxy_na[cxy];
749             na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
750
751             if (nbat->XFormat == nbatXYZQ)
752             {
753                 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
754                 for(i=0; i<na; i++)
755                 {
756                     *q = charge[nbs->a[ash+i]];
757                     q += STRIDE_XYZQ;
758                 }
759                 /* Complete the partially filled last cell with zeros */
760                 for(; i<na_round; i++)
761                 {
762                     *q = 0;
763                     q += STRIDE_XYZQ;
764                 }
765             }
766             else
767             {
768                 q = nbat->q + ash;
769                 for(i=0; i<na; i++)
770                 {
771                     *q = charge[nbs->a[ash+i]];
772                     q++;
773                 }
774                 /* Complete the partially filled last cell with zeros */
775                 for(; i<na_round; i++)
776                 {
777                     *q = 0;
778                     q++;
779                 }
780             }
781         }
782     }
783 }
784
785 /* Copies the energy group indices to a reordered and packed array */
786 static void copy_egp_to_nbat_egps(const int *a,int na,int na_round,
787                                   int na_c,int bit_shift,
788                                   const int *in,int *innb)
789 {
790     int i,j,sa,at;
791     int comb;
792
793     j = 0;
794     for(i=0; i<na; i+=na_c)
795     {
796         /* Store na_c energy groups number into one int */
797         comb = 0;
798         for(sa=0; sa<na_c; sa++)
799         {
800             at = a[i+sa];
801             if (at >= 0)
802             {
803                 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
804             }
805         }
806         innb[j++] = comb;
807     }
808     /* Complete the partially filled last cell with fill */
809     for(; i<na_round; i+=na_c)
810     {
811         innb[j++] = 0;
812     }
813 }
814
815 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
816 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t *nbat,
817                                             int ngrid,
818                                             const nbnxn_search_t nbs,
819                                             const int *atinfo)
820 {
821     int g,i,ncz,ash;
822     const nbnxn_grid_t *grid;
823
824     for(g=0; g<ngrid; g++)
825     {
826         grid = &nbs->grid[g];
827
828         /* Loop over all columns and copy and fill */
829         for(i=0; i<grid->ncx*grid->ncy; i++)
830         {
831             ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
832             ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
833
834             copy_egp_to_nbat_egps(nbs->a+ash,grid->cxy_na[i],ncz*grid->na_sc,
835                                   nbat->na_c,nbat->neg_2log,
836                                   atinfo,nbat->energrp+(ash>>grid->na_c_2log));
837         }
838     }
839 }
840
841 /* Sets all required atom parameter data in nbnxn_atomdata_t */
842 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
843                         int locality,
844                         const nbnxn_search_t nbs,
845                         const t_mdatoms *mdatoms,
846                         const int *atinfo)
847 {
848     int ngrid;
849
850     if (locality == eatLocal)
851     {
852         ngrid = 1;
853     }
854     else
855     {
856         ngrid = nbs->ngrid;
857     }
858
859     nbnxn_atomdata_set_atomtypes(nbat,ngrid,nbs,mdatoms->typeA);
860
861     nbnxn_atomdata_set_charges(nbat,ngrid,nbs,mdatoms->chargeA);
862
863     if (nbat->nenergrp > 1)
864     {
865         nbnxn_atomdata_set_energygroups(nbat,ngrid,nbs,atinfo);
866     }
867 }
868
869 /* Copies the shift vector array to nbnxn_atomdata_t */
870 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox,
871                                    rvec *shift_vec,
872                                    nbnxn_atomdata_t *nbat)
873 {
874     int i;
875
876     nbat->bDynamicBox = bDynamicBox;
877     for(i=0; i<SHIFTS; i++)
878     {
879         copy_rvec(shift_vec[i],nbat->shift_vec[i]);
880     }
881 }
882
883 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
884 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
885                                       int locality,
886                                       gmx_bool FillLocal,
887                                       rvec *x,
888                                       nbnxn_atomdata_t *nbat)
889 {
890     int g0=0,g1=0;
891     int nth,th;
892
893     switch (locality)
894     {
895     case eatAll:
896         g0 = 0;
897         g1 = nbs->ngrid;
898         break;
899     case eatLocal:
900         g0 = 0;
901         g1 = 1;
902         break;
903     case eatNonlocal:
904         g0 = 1;
905         g1 = nbs->ngrid;
906         break;
907     }
908
909     if (FillLocal)
910     {
911         nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
912     }
913
914     nth = gmx_omp_nthreads_get(emntPairsearch);
915
916 #pragma omp parallel for num_threads(nth) schedule(static)
917     for(th=0; th<nth; th++)
918     {
919         int g;
920
921         for(g=g0; g<g1; g++)
922         {
923             const nbnxn_grid_t *grid;
924             int cxy0,cxy1,cxy;
925
926             grid = &nbs->grid[g];
927
928             cxy0 = (grid->ncx*grid->ncy* th   +nth-1)/nth;
929             cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
930
931             for(cxy=cxy0; cxy<cxy1; cxy++)
932             {
933                 int na,ash,na_fill;
934
935                 na  = grid->cxy_na[cxy];
936                 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
937
938                 if (g == 0 && FillLocal)
939                 {
940                     na_fill =
941                         (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
942                 }
943                 else
944                 {
945                     /* We fill only the real particle locations.
946                      * We assume the filling entries at the end have been
947                      * properly set before during ns.
948                      */
949                     na_fill = na;
950                 }
951                 copy_rvec_to_nbat_real(nbs->a+ash,na,na_fill,x,
952                                        nbat->XFormat,nbat->x,ash,
953                                        0,0,0);
954             }
955         }
956     }
957 }
958
959 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
960 static void
961 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
962                                     const nbnxn_atomdata_t *nbat,
963                                     nbnxn_atomdata_output_t *out,
964                                     int nfa,
965                                     int a0,int a1,
966                                     rvec *f)
967 {
968     int  a,i,fa;
969     const int  *cell;
970     const real *fnb;
971
972     cell = nbs->cell;
973
974     /* Loop over all columns and copy and fill */
975     switch (nbat->FFormat)
976     {
977     case nbatXYZ:
978     case nbatXYZQ:
979         if (nfa == 1)
980         {
981             fnb = out[0].f;
982
983             for(a=a0; a<a1; a++)
984             {
985                 i = cell[a]*nbat->fstride;
986
987                 f[a][XX] += fnb[i];
988                 f[a][YY] += fnb[i+1];
989                 f[a][ZZ] += fnb[i+2];
990             }
991         }
992         else
993         {
994             for(a=a0; a<a1; a++)
995             {
996                 i = cell[a]*nbat->fstride;
997
998                 for(fa=0; fa<nfa; fa++)
999                 {
1000                     f[a][XX] += out[fa].f[i];
1001                     f[a][YY] += out[fa].f[i+1];
1002                     f[a][ZZ] += out[fa].f[i+2];
1003                 }
1004             }
1005         }
1006         break;
1007     case nbatX4:
1008         if (nfa == 1)
1009         {
1010             fnb = out[0].f;
1011
1012             for(a=a0; a<a1; a++)
1013             {
1014                 i = X4_IND_A(cell[a]);
1015
1016                 f[a][XX] += fnb[i+XX*PACK_X4];
1017                 f[a][YY] += fnb[i+YY*PACK_X4];
1018                 f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1019             }
1020         }
1021         else
1022         {
1023             for(a=a0; a<a1; a++)
1024             {
1025                 i = X4_IND_A(cell[a]);
1026                 
1027                 for(fa=0; fa<nfa; fa++)
1028                 {
1029                     f[a][XX] += out[fa].f[i+XX*PACK_X4];
1030                     f[a][YY] += out[fa].f[i+YY*PACK_X4];
1031                     f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1032                 }
1033             }
1034         }
1035         break;
1036     case nbatX8:
1037         if (nfa == 1)
1038         {
1039             fnb = out[0].f;
1040
1041             for(a=a0; a<a1; a++)
1042             {
1043                 i = X8_IND_A(cell[a]);
1044
1045                 f[a][XX] += fnb[i+XX*PACK_X8];
1046                 f[a][YY] += fnb[i+YY*PACK_X8];
1047                 f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1048             }
1049         }
1050         else
1051         {
1052             for(a=a0; a<a1; a++)
1053             {
1054                 i = X8_IND_A(cell[a]);
1055                 
1056                 for(fa=0; fa<nfa; fa++)
1057                 {
1058                     f[a][XX] += out[fa].f[i+XX*PACK_X8];
1059                     f[a][YY] += out[fa].f[i+YY*PACK_X8];
1060                     f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1061                 }
1062             }
1063         }
1064         break;
1065     }
1066 }
1067
1068 /* Add the force array(s) from nbnxn_atomdata_t to f */
1069 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t nbs,
1070                                     int locality,
1071                                     const nbnxn_atomdata_t *nbat,
1072                                     rvec *f)
1073 {
1074     int a0=0,na=0;
1075     int nth,th;
1076
1077     nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1078
1079     switch (locality)
1080     {
1081     case eatAll:
1082         a0 = 0;
1083         na = nbs->natoms_nonlocal;
1084         break;
1085     case eatLocal:
1086         a0 = 0;
1087         na = nbs->natoms_local;
1088         break;
1089     case eatNonlocal:
1090         a0 = nbs->natoms_local;
1091         na = nbs->natoms_nonlocal - nbs->natoms_local;
1092         break;
1093     }
1094
1095     nth = gmx_omp_nthreads_get(emntNonbonded);
1096 #pragma omp parallel for num_threads(nth) schedule(static)
1097     for(th=0; th<nth; th++)
1098     {
1099         nbnxn_atomdata_add_nbat_f_to_f_part(nbs,nbat,
1100                                              nbat->out,
1101                                              nbat->nout,
1102                                              a0+((th+0)*na)/nth,
1103                                              a0+((th+1)*na)/nth,
1104                                              f);
1105     }
1106
1107     nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1108 }
1109
1110 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1111 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1112                                               rvec *fshift)
1113 {
1114     const nbnxn_atomdata_output_t *out;
1115     int  th;
1116     int  s;
1117     rvec sum;
1118
1119     out = nbat->out;
1120     
1121     for(s=0; s<SHIFTS; s++)
1122     {
1123         clear_rvec(sum);
1124         for(th=0; th<nbat->nout; th++)
1125         {
1126             sum[XX] += out[th].fshift[s*DIM+XX];
1127             sum[YY] += out[th].fshift[s*DIM+YY];
1128             sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1129         }
1130         rvec_inc(fshift[s],sum);
1131     }
1132 }