Merge branch 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 nbnxn_atomdata_t format");
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             /* nbfp_s4 stores two parameters using a stride of 4,
403              * because this would suit x86 SIMD single-precision
404              * quad-load intrinsics. There's a slight inefficiency in
405              * allocating and initializing nbfp_s4 when it might not
406              * be used, but introducing the conditional code is not
407              * really worth it. */
408             nbat->alloc((void **)&nbat->nbfp_s4, nt*nt*4*sizeof(*nbat->nbfp_s4));
409             for (i = 0; i < nt; i++)
410             {
411                 for (j = 0; j < nt; j++)
412                 {
413                     nbat->nbfp_s4[(i*nt+j)*4+0] = nbat->nbfp[(i*nt+j)*2+0];
414                     nbat->nbfp_s4[(i*nt+j)*4+1] = nbat->nbfp[(i*nt+j)*2+1];
415                     nbat->nbfp_s4[(i*nt+j)*4+2] = 0;
416                     nbat->nbfp_s4[(i*nt+j)*4+3] = 0;
417                 }
418             }
419             break;
420         default:
421             gmx_incons("Unknown combination rule");
422             break;
423     }
424 }
425
426 /* Initializes an nbnxn_atomdata_t data structure */
427 void nbnxn_atomdata_init(FILE *fp,
428                          nbnxn_atomdata_t *nbat,
429                          int nb_kernel_type,
430                          int ntype, const real *nbfp,
431                          int n_energygroups,
432                          int nout,
433                          nbnxn_alloc_t *alloc,
434                          nbnxn_free_t  *free)
435 {
436     int      i, j;
437     real     c6, c12, tol;
438     char    *ptr;
439     gmx_bool simple, bCombGeom, bCombLB;
440
441     if (alloc == NULL)
442     {
443         nbat->alloc = nbnxn_alloc_aligned;
444     }
445     else
446     {
447         nbat->alloc = alloc;
448     }
449     if (free == NULL)
450     {
451         nbat->free = nbnxn_free_aligned;
452     }
453     else
454     {
455         nbat->free = free;
456     }
457
458     if (debug)
459     {
460         fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
461     }
462     nbat->ntype = ntype + 1;
463     nbat->alloc((void **)&nbat->nbfp,
464                 nbat->ntype*nbat->ntype*2*sizeof(*nbat->nbfp));
465     nbat->alloc((void **)&nbat->nbfp_comb, nbat->ntype*2*sizeof(*nbat->nbfp_comb));
466
467     /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
468      * force-field floating point parameters.
469      */
470     tol = 1e-5;
471     ptr = getenv("GMX_LJCOMB_TOL");
472     if (ptr != NULL)
473     {
474         double dbl;
475
476         sscanf(ptr, "%lf", &dbl);
477         tol = dbl;
478     }
479     bCombGeom = TRUE;
480     bCombLB   = TRUE;
481
482     /* Temporarily fill nbat->nbfp_comb with sigma and epsilon
483      * to check for the LB rule.
484      */
485     for (i = 0; i < ntype; i++)
486     {
487         c6  = nbfp[(i*ntype+i)*2  ]/6.0;
488         c12 = nbfp[(i*ntype+i)*2+1]/12.0;
489         if (c6 > 0 && c12 > 0)
490         {
491             nbat->nbfp_comb[i*2  ] = pow(c12/c6, 1.0/6.0);
492             nbat->nbfp_comb[i*2+1] = 0.25*c6*c6/c12;
493         }
494         else if (c6 == 0 && c12 == 0)
495         {
496             nbat->nbfp_comb[i*2  ] = 0;
497             nbat->nbfp_comb[i*2+1] = 0;
498         }
499         else
500         {
501             /* Can not use LB rule with only dispersion or repulsion */
502             bCombLB = FALSE;
503         }
504     }
505
506     for (i = 0; i < nbat->ntype; i++)
507     {
508         for (j = 0; j < nbat->ntype; j++)
509         {
510             if (i < ntype && j < ntype)
511             {
512                 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
513                  * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
514                  */
515                 c6  = nbfp[(i*ntype+j)*2  ];
516                 c12 = nbfp[(i*ntype+j)*2+1];
517                 nbat->nbfp[(i*nbat->ntype+j)*2  ] = c6;
518                 nbat->nbfp[(i*nbat->ntype+j)*2+1] = c12;
519
520                 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
521                 bCombGeom = bCombGeom &&
522                     gmx_within_tol(c6*c6, nbfp[(i*ntype+i)*2  ]*nbfp[(j*ntype+j)*2  ], tol) &&
523                     gmx_within_tol(c12*c12, nbfp[(i*ntype+i)*2+1]*nbfp[(j*ntype+j)*2+1], tol);
524
525                 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
526                 c6     /= 6.0;
527                 c12    /= 12.0;
528                 bCombLB = bCombLB &&
529                     ((c6 == 0 && c12 == 0 &&
530                       (nbat->nbfp_comb[i*2+1] == 0 || nbat->nbfp_comb[j*2+1] == 0)) ||
531                      (c6 > 0 && c12 > 0 &&
532                       gmx_within_tol(pow(c12/c6, 1.0/6.0), 0.5*(nbat->nbfp_comb[i*2]+nbat->nbfp_comb[j*2]), tol) &&
533                       gmx_within_tol(0.25*c6*c6/c12, sqrt(nbat->nbfp_comb[i*2+1]*nbat->nbfp_comb[j*2+1]), tol)));
534             }
535             else
536             {
537                 /* Add zero parameters for the additional dummy atom type */
538                 nbat->nbfp[(i*nbat->ntype+j)*2  ] = 0;
539                 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 0;
540             }
541         }
542     }
543     if (debug)
544     {
545         fprintf(debug, "Combination rules: geometric %d Lorentz-Berthelot %d\n",
546                 bCombGeom, bCombLB);
547     }
548
549     simple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
550
551     if (simple)
552     {
553         /* We prefer the geometic combination rule,
554          * as that gives a slightly faster kernel than the LB rule.
555          */
556         if (bCombGeom)
557         {
558             nbat->comb_rule = ljcrGEOM;
559         }
560         else if (bCombLB)
561         {
562             nbat->comb_rule = ljcrLB;
563         }
564         else
565         {
566             nbat->comb_rule = ljcrNONE;
567
568             nbat->free(nbat->nbfp_comb);
569         }
570
571         if (fp)
572         {
573             if (nbat->comb_rule == ljcrNONE)
574             {
575                 fprintf(fp, "Using full Lennard-Jones parameter combination matrix\n\n");
576             }
577             else
578             {
579                 fprintf(fp, "Using %s Lennard-Jones combination rule\n\n",
580                         nbat->comb_rule == ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
581             }
582         }
583
584         set_combination_rule_data(nbat);
585     }
586     else
587     {
588         nbat->comb_rule = ljcrNONE;
589
590         nbat->free(nbat->nbfp_comb);
591     }
592
593     nbat->natoms  = 0;
594     nbat->type    = NULL;
595     nbat->lj_comb = NULL;
596     if (simple)
597     {
598         int pack_x;
599
600         switch (nb_kernel_type)
601         {
602             case nbnxnk4xN_SIMD_4xN:
603             case nbnxnk4xN_SIMD_2xNN:
604                 pack_x = max(NBNXN_CPU_CLUSTER_I_SIZE,
605                              nbnxn_kernel_to_cj_size(nb_kernel_type));
606                 switch (pack_x)
607                 {
608                     case 4:
609                         nbat->XFormat = nbatX4;
610                         break;
611                     case 8:
612                         nbat->XFormat = nbatX8;
613                         break;
614                     default:
615                         gmx_incons("Unsupported packing width");
616                 }
617                 break;
618             default:
619                 nbat->XFormat = nbatXYZ;
620                 break;
621         }
622
623         nbat->FFormat = nbat->XFormat;
624     }
625     else
626     {
627         nbat->XFormat = nbatXYZQ;
628         nbat->FFormat = nbatXYZ;
629     }
630     nbat->q        = NULL;
631     nbat->nenergrp = n_energygroups;
632     if (!simple)
633     {
634         /* Energy groups not supported yet for super-sub lists */
635         if (n_energygroups > 1 && fp != NULL)
636         {
637             fprintf(fp, "\nNOTE: With GPUs, reporting energy group contributions is not supported\n\n");
638         }
639         nbat->nenergrp = 1;
640     }
641     /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
642     if (nbat->nenergrp > 64)
643     {
644         gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
645     }
646     nbat->neg_2log = 1;
647     while (nbat->nenergrp > (1<<nbat->neg_2log))
648     {
649         nbat->neg_2log++;
650     }
651     nbat->energrp = NULL;
652     nbat->alloc((void **)&nbat->shift_vec, SHIFTS*sizeof(*nbat->shift_vec));
653     nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
654     nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
655     nbat->x       = NULL;
656
657 #ifdef GMX_NBNXN_SIMD
658     if (simple)
659     {
660         /* Set the diagonal cluster pair exclusion mask setup data.
661          * In the kernel we check 0 < j - i to generate the masks.
662          * Here we store j - i for generating the mask for the first i,
663          * we substract 0.5 to avoid rounding issues.
664          * In the kernel we can subtract 1 to generate the subsequent mask.
665          */
666         const int simd_width = GMX_NBNXN_SIMD_BITWIDTH/(sizeof(real)*8);
667         int       simd_4xn_diag_size, real_excl, simd_excl_size, j, s;
668
669         simd_4xn_diag_size = max(NBNXN_CPU_CLUSTER_I_SIZE, simd_width);
670         snew_aligned(nbat->simd_4xn_diag, simd_4xn_diag_size, NBNXN_MEM_ALIGN);
671         for (j = 0; j < simd_4xn_diag_size; j++)
672         {
673             nbat->simd_4xn_diag[j] = j - 0.5;
674         }
675
676         snew_aligned(nbat->simd_2xnn_diag, simd_width, NBNXN_MEM_ALIGN);
677         for (j = 0; j < simd_width/2; j++)
678         {
679             /* The j-cluster size is half the SIMD width */
680             nbat->simd_2xnn_diag[j]              = j - 0.5;
681             /* The next half of the SIMD width is for i + 1 */
682             nbat->simd_2xnn_diag[simd_width/2+j] = j - 1 - 0.5;
683         }
684
685         /* We always use 32-bit integer exclusion masks. When we use
686          * double precision, we fit two integers in a double SIMD register.
687          */
688         real_excl = sizeof(real)/sizeof(*nbat->simd_excl_mask);
689         /* Set bits for use with both 4xN and 2x(N+N) kernels */
690         simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*simd_width*real_excl;
691         snew_aligned(nbat->simd_excl_mask, simd_excl_size*real_excl, NBNXN_MEM_ALIGN);
692         for (j = 0; j < simd_excl_size; j++)
693         {
694             /* Set the consecutive bits for masking pair exclusions.
695              * For double a single-bit mask would be enough.
696              * But using two bits avoids endianness issues.
697              */
698             for (s = 0; s < real_excl; s++)
699             {
700                 /* Set the consecutive bits for masking pair exclusions */
701                 nbat->simd_excl_mask[j*real_excl + s] = (1U << j);
702             }
703         }
704     }
705 #endif
706
707     /* Initialize the output data structures */
708     nbat->nout    = nout;
709     snew(nbat->out, nbat->nout);
710     nbat->nalloc  = 0;
711     for (i = 0; i < nbat->nout; i++)
712     {
713         nbnxn_atomdata_output_init(&nbat->out[i],
714                                    nb_kernel_type,
715                                    nbat->nenergrp, 1<<nbat->neg_2log,
716                                    nbat->alloc);
717     }
718     nbat->buffer_flags.flag        = NULL;
719     nbat->buffer_flags.flag_nalloc = 0;
720 }
721
722 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
723                                        const int *type, int na,
724                                        real *ljparam_at)
725 {
726     int is, k, i;
727
728     /* The LJ params follow the combination rule:
729      * copy the params for the type array to the atom array.
730      */
731     for (is = 0; is < na; is += PACK_X4)
732     {
733         for (k = 0; k < PACK_X4; k++)
734         {
735             i = is + k;
736             ljparam_at[is*2        +k] = ljparam_type[type[i]*2  ];
737             ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
738         }
739     }
740 }
741
742 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
743                                        const int *type, int na,
744                                        real *ljparam_at)
745 {
746     int is, k, i;
747
748     /* The LJ params follow the combination rule:
749      * copy the params for the type array to the atom array.
750      */
751     for (is = 0; is < na; is += PACK_X8)
752     {
753         for (k = 0; k < PACK_X8; k++)
754         {
755             i = is + k;
756             ljparam_at[is*2        +k] = ljparam_type[type[i]*2  ];
757             ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
758         }
759     }
760 }
761
762 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
763 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t    *nbat,
764                                          int                  ngrid,
765                                          const nbnxn_search_t nbs,
766                                          const int           *type)
767 {
768     int                 g, i, ncz, ash;
769     const nbnxn_grid_t *grid;
770
771     for (g = 0; g < ngrid; g++)
772     {
773         grid = &nbs->grid[g];
774
775         /* Loop over all columns and copy and fill */
776         for (i = 0; i < grid->ncx*grid->ncy; i++)
777         {
778             ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
779             ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
780
781             copy_int_to_nbat_int(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
782                                  type, nbat->ntype-1, nbat->type+ash);
783
784             if (nbat->comb_rule != ljcrNONE)
785             {
786                 if (nbat->XFormat == nbatX4)
787                 {
788                     copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
789                                                nbat->type+ash, ncz*grid->na_sc,
790                                                nbat->lj_comb+ash*2);
791                 }
792                 else if (nbat->XFormat == nbatX8)
793                 {
794                     copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
795                                                nbat->type+ash, ncz*grid->na_sc,
796                                                nbat->lj_comb+ash*2);
797                 }
798             }
799         }
800     }
801 }
802
803 /* Sets the charges in nbnxn_atomdata_t *nbat */
804 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t    *nbat,
805                                        int                  ngrid,
806                                        const nbnxn_search_t nbs,
807                                        const real          *charge)
808 {
809     int                 g, cxy, ncz, ash, na, na_round, i, j;
810     real               *q;
811     const nbnxn_grid_t *grid;
812
813     for (g = 0; g < ngrid; g++)
814     {
815         grid = &nbs->grid[g];
816
817         /* Loop over all columns and copy and fill */
818         for (cxy = 0; cxy < grid->ncx*grid->ncy; cxy++)
819         {
820             ash      = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
821             na       = grid->cxy_na[cxy];
822             na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
823
824             if (nbat->XFormat == nbatXYZQ)
825             {
826                 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
827                 for (i = 0; i < na; i++)
828                 {
829                     *q = charge[nbs->a[ash+i]];
830                     q += STRIDE_XYZQ;
831                 }
832                 /* Complete the partially filled last cell with zeros */
833                 for (; i < na_round; i++)
834                 {
835                     *q = 0;
836                     q += STRIDE_XYZQ;
837                 }
838             }
839             else
840             {
841                 q = nbat->q + ash;
842                 for (i = 0; i < na; i++)
843                 {
844                     *q = charge[nbs->a[ash+i]];
845                     q++;
846                 }
847                 /* Complete the partially filled last cell with zeros */
848                 for (; i < na_round; i++)
849                 {
850                     *q = 0;
851                     q++;
852                 }
853             }
854         }
855     }
856 }
857
858 /* Copies the energy group indices to a reordered and packed array */
859 static void copy_egp_to_nbat_egps(const int *a, int na, int na_round,
860                                   int na_c, int bit_shift,
861                                   const int *in, int *innb)
862 {
863     int i, j, sa, at;
864     int comb;
865
866     j = 0;
867     for (i = 0; i < na; i += na_c)
868     {
869         /* Store na_c energy group numbers into one int */
870         comb = 0;
871         for (sa = 0; sa < na_c; sa++)
872         {
873             at = a[i+sa];
874             if (at >= 0)
875             {
876                 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
877             }
878         }
879         innb[j++] = comb;
880     }
881     /* Complete the partially filled last cell with fill */
882     for (; i < na_round; i += na_c)
883     {
884         innb[j++] = 0;
885     }
886 }
887
888 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
889 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t    *nbat,
890                                             int                  ngrid,
891                                             const nbnxn_search_t nbs,
892                                             const int           *atinfo)
893 {
894     int                 g, i, ncz, ash;
895     const nbnxn_grid_t *grid;
896
897     for (g = 0; g < ngrid; g++)
898     {
899         grid = &nbs->grid[g];
900
901         /* Loop over all columns and copy and fill */
902         for (i = 0; i < grid->ncx*grid->ncy; i++)
903         {
904             ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
905             ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
906
907             copy_egp_to_nbat_egps(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
908                                   nbat->na_c, nbat->neg_2log,
909                                   atinfo, nbat->energrp+(ash>>grid->na_c_2log));
910         }
911     }
912 }
913
914 /* Sets all required atom parameter data in nbnxn_atomdata_t */
915 void nbnxn_atomdata_set(nbnxn_atomdata_t    *nbat,
916                         int                  locality,
917                         const nbnxn_search_t nbs,
918                         const t_mdatoms     *mdatoms,
919                         const int           *atinfo)
920 {
921     int ngrid;
922
923     if (locality == eatLocal)
924     {
925         ngrid = 1;
926     }
927     else
928     {
929         ngrid = nbs->ngrid;
930     }
931
932     nbnxn_atomdata_set_atomtypes(nbat, ngrid, nbs, mdatoms->typeA);
933
934     nbnxn_atomdata_set_charges(nbat, ngrid, nbs, mdatoms->chargeA);
935
936     if (nbat->nenergrp > 1)
937     {
938         nbnxn_atomdata_set_energygroups(nbat, ngrid, nbs, atinfo);
939     }
940 }
941
942 /* Copies the shift vector array to nbnxn_atomdata_t */
943 void nbnxn_atomdata_copy_shiftvec(gmx_bool          bDynamicBox,
944                                   rvec             *shift_vec,
945                                   nbnxn_atomdata_t *nbat)
946 {
947     int i;
948
949     nbat->bDynamicBox = bDynamicBox;
950     for (i = 0; i < SHIFTS; i++)
951     {
952         copy_rvec(shift_vec[i], nbat->shift_vec[i]);
953     }
954 }
955
956 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
957 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
958                                      int                  locality,
959                                      gmx_bool             FillLocal,
960                                      rvec                *x,
961                                      nbnxn_atomdata_t    *nbat)
962 {
963     int g0 = 0, g1 = 0;
964     int nth, th;
965
966     switch (locality)
967     {
968         case eatAll:
969             g0 = 0;
970             g1 = nbs->ngrid;
971             break;
972         case eatLocal:
973             g0 = 0;
974             g1 = 1;
975             break;
976         case eatNonlocal:
977             g0 = 1;
978             g1 = nbs->ngrid;
979             break;
980     }
981
982     if (FillLocal)
983     {
984         nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
985     }
986
987     nth = gmx_omp_nthreads_get(emntPairsearch);
988
989 #pragma omp parallel for num_threads(nth) schedule(static)
990     for (th = 0; th < nth; th++)
991     {
992         int g;
993
994         for (g = g0; g < g1; g++)
995         {
996             const nbnxn_grid_t *grid;
997             int                 cxy0, cxy1, cxy;
998
999             grid = &nbs->grid[g];
1000
1001             cxy0 = (grid->ncx*grid->ncy* th   +nth-1)/nth;
1002             cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
1003
1004             for (cxy = cxy0; cxy < cxy1; cxy++)
1005             {
1006                 int na, ash, na_fill;
1007
1008                 na  = grid->cxy_na[cxy];
1009                 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1010
1011                 if (g == 0 && FillLocal)
1012                 {
1013                     na_fill =
1014                         (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
1015                 }
1016                 else
1017                 {
1018                     /* We fill only the real particle locations.
1019                      * We assume the filling entries at the end have been
1020                      * properly set before during ns.
1021                      */
1022                     na_fill = na;
1023                 }
1024                 copy_rvec_to_nbat_real(nbs->a+ash, na, na_fill, x,
1025                                        nbat->XFormat, nbat->x, ash,
1026                                        0, 0, 0);
1027             }
1028         }
1029     }
1030 }
1031
1032 static void
1033 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
1034                            int i0, int i1)
1035 {
1036     int i;
1037
1038     for (i = i0; i < i1; i++)
1039     {
1040         dest[i] = 0;
1041     }
1042 }
1043
1044 static void
1045 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1046                             gmx_bool bDestSet,
1047                             real ** gmx_restrict src,
1048                             int nsrc,
1049                             int i0, int i1)
1050 {
1051     int i, s;
1052
1053     if (bDestSet)
1054     {
1055         /* The destination buffer contains data, add to it */
1056         for (i = i0; i < i1; i++)
1057         {
1058             for (s = 0; s < nsrc; s++)
1059             {
1060                 dest[i] += src[s][i];
1061             }
1062         }
1063     }
1064     else
1065     {
1066         /* The destination buffer is unitialized, set it first */
1067         for (i = i0; i < i1; i++)
1068         {
1069             dest[i] = src[0][i];
1070             for (s = 1; s < nsrc; s++)
1071             {
1072                 dest[i] += src[s][i];
1073             }
1074         }
1075     }
1076 }
1077
1078 static void
1079 nbnxn_atomdata_reduce_reals_simd(real * gmx_restrict dest,
1080                                  gmx_bool bDestSet,
1081                                  real ** gmx_restrict src,
1082                                  int nsrc,
1083                                  int i0, int i1)
1084 {
1085 #ifdef GMX_NBNXN_SIMD
1086 /* The SIMD width here is actually independent of that in the kernels,
1087  * but we use the same width for simplicity (usually optimal anyhow).
1088  */
1089 #ifdef GMX_NBNXN_HALF_WIDTH_SIMD
1090 #define GMX_USE_HALF_WIDTH_SIMD_HERE
1091 #endif
1092 #include "gmx_simd_macros.h"
1093
1094     int       i, s;
1095     gmx_mm_pr dest_SSE, src_SSE;
1096
1097     if (bDestSet)
1098     {
1099         for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1100         {
1101             dest_SSE = gmx_load_pr(dest+i);
1102             for (s = 0; s < nsrc; s++)
1103             {
1104                 src_SSE  = gmx_load_pr(src[s]+i);
1105                 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1106             }
1107             gmx_store_pr(dest+i, dest_SSE);
1108         }
1109     }
1110     else
1111     {
1112         for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1113         {
1114             dest_SSE = gmx_load_pr(src[0]+i);
1115             for (s = 1; s < nsrc; s++)
1116             {
1117                 src_SSE  = gmx_load_pr(src[s]+i);
1118                 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1119             }
1120             gmx_store_pr(dest+i, dest_SSE);
1121         }
1122     }
1123 #endif
1124 }
1125
1126 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1127 static void
1128 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1129                                     const nbnxn_atomdata_t *nbat,
1130                                     nbnxn_atomdata_output_t *out,
1131                                     int nfa,
1132                                     int a0, int a1,
1133                                     rvec *f)
1134 {
1135     int         a, i, fa;
1136     const int  *cell;
1137     const real *fnb;
1138
1139     cell = nbs->cell;
1140
1141     /* Loop over all columns and copy and fill */
1142     switch (nbat->FFormat)
1143     {
1144         case nbatXYZ:
1145         case nbatXYZQ:
1146             if (nfa == 1)
1147             {
1148                 fnb = out[0].f;
1149
1150                 for (a = a0; a < a1; a++)
1151                 {
1152                     i = cell[a]*nbat->fstride;
1153
1154                     f[a][XX] += fnb[i];
1155                     f[a][YY] += fnb[i+1];
1156                     f[a][ZZ] += fnb[i+2];
1157                 }
1158             }
1159             else
1160             {
1161                 for (a = a0; a < a1; a++)
1162                 {
1163                     i = cell[a]*nbat->fstride;
1164
1165                     for (fa = 0; fa < nfa; fa++)
1166                     {
1167                         f[a][XX] += out[fa].f[i];
1168                         f[a][YY] += out[fa].f[i+1];
1169                         f[a][ZZ] += out[fa].f[i+2];
1170                     }
1171                 }
1172             }
1173             break;
1174         case nbatX4:
1175             if (nfa == 1)
1176             {
1177                 fnb = out[0].f;
1178
1179                 for (a = a0; a < a1; a++)
1180                 {
1181                     i = X4_IND_A(cell[a]);
1182
1183                     f[a][XX] += fnb[i+XX*PACK_X4];
1184                     f[a][YY] += fnb[i+YY*PACK_X4];
1185                     f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1186                 }
1187             }
1188             else
1189             {
1190                 for (a = a0; a < a1; a++)
1191                 {
1192                     i = X4_IND_A(cell[a]);
1193
1194                     for (fa = 0; fa < nfa; fa++)
1195                     {
1196                         f[a][XX] += out[fa].f[i+XX*PACK_X4];
1197                         f[a][YY] += out[fa].f[i+YY*PACK_X4];
1198                         f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1199                     }
1200                 }
1201             }
1202             break;
1203         case nbatX8:
1204             if (nfa == 1)
1205             {
1206                 fnb = out[0].f;
1207
1208                 for (a = a0; a < a1; a++)
1209                 {
1210                     i = X8_IND_A(cell[a]);
1211
1212                     f[a][XX] += fnb[i+XX*PACK_X8];
1213                     f[a][YY] += fnb[i+YY*PACK_X8];
1214                     f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1215                 }
1216             }
1217             else
1218             {
1219                 for (a = a0; a < a1; a++)
1220                 {
1221                     i = X8_IND_A(cell[a]);
1222
1223                     for (fa = 0; fa < nfa; fa++)
1224                     {
1225                         f[a][XX] += out[fa].f[i+XX*PACK_X8];
1226                         f[a][YY] += out[fa].f[i+YY*PACK_X8];
1227                         f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1228                     }
1229                 }
1230             }
1231             break;
1232         default:
1233             gmx_incons("Unsupported nbnxn_atomdata_t format");
1234     }
1235 }
1236
1237 /* Add the force array(s) from nbnxn_atomdata_t to f */
1238 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t    nbs,
1239                                     int                     locality,
1240                                     const nbnxn_atomdata_t *nbat,
1241                                     rvec                   *f)
1242 {
1243     int a0 = 0, na = 0;
1244     int nth, th;
1245
1246     nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1247
1248     switch (locality)
1249     {
1250         case eatAll:
1251             a0 = 0;
1252             na = nbs->natoms_nonlocal;
1253             break;
1254         case eatLocal:
1255             a0 = 0;
1256             na = nbs->natoms_local;
1257             break;
1258         case eatNonlocal:
1259             a0 = nbs->natoms_local;
1260             na = nbs->natoms_nonlocal - nbs->natoms_local;
1261             break;
1262     }
1263
1264     nth = gmx_omp_nthreads_get(emntNonbonded);
1265
1266     if (nbat->nout > 1)
1267     {
1268         if (locality != eatAll)
1269         {
1270             gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1271         }
1272
1273         /* Reduce the force thread output buffers into buffer 0, before adding
1274          * them to the, differently ordered, "real" force buffer.
1275          */
1276 #pragma omp parallel for num_threads(nth) schedule(static)
1277         for (th = 0; th < nth; th++)
1278         {
1279             const nbnxn_buffer_flags_t *flags;
1280             int   b0, b1, b;
1281             int   i0, i1;
1282             int   nfptr;
1283             real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1284             int   out;
1285
1286             flags = &nbat->buffer_flags;
1287
1288             /* Calculate the cell-block range for our thread */
1289             b0 = (flags->nflag* th   )/nth;
1290             b1 = (flags->nflag*(th+1))/nth;
1291
1292             for (b = b0; b < b1; b++)
1293             {
1294                 i0 =  b   *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1295                 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1296
1297                 nfptr = 0;
1298                 for (out = 1; out < nbat->nout; out++)
1299                 {
1300                     if (flags->flag[b] & (1U<<out))
1301                     {
1302                         fptr[nfptr++] = nbat->out[out].f;
1303                     }
1304                 }
1305                 if (nfptr > 0)
1306                 {
1307 #ifdef GMX_NBNXN_SIMD
1308                     nbnxn_atomdata_reduce_reals_simd
1309 #else
1310                     nbnxn_atomdata_reduce_reals
1311 #endif
1312                         (nbat->out[0].f,
1313                         flags->flag[b] & (1U<<0),
1314                         fptr, nfptr,
1315                         i0, i1);
1316                 }
1317                 else if (!(flags->flag[b] & (1U<<0)))
1318                 {
1319                     nbnxn_atomdata_clear_reals(nbat->out[0].f,
1320                                                i0, i1);
1321                 }
1322             }
1323         }
1324     }
1325
1326 #pragma omp parallel for num_threads(nth) schedule(static)
1327     for (th = 0; th < nth; th++)
1328     {
1329         nbnxn_atomdata_add_nbat_f_to_f_part(nbs, nbat,
1330                                             nbat->out,
1331                                             1,
1332                                             a0+((th+0)*na)/nth,
1333                                             a0+((th+1)*na)/nth,
1334                                             f);
1335     }
1336
1337     nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1338 }
1339
1340 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1341 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1342                                               rvec                   *fshift)
1343 {
1344     const nbnxn_atomdata_output_t *out;
1345     int  th;
1346     int  s;
1347     rvec sum;
1348
1349     out = nbat->out;
1350
1351     for (s = 0; s < SHIFTS; s++)
1352     {
1353         clear_rvec(sum);
1354         for (th = 0; th < nbat->nout; th++)
1355         {
1356             sum[XX] += out[th].fshift[s*DIM+XX];
1357             sum[YY] += out[th].fshift[s*DIM+YY];
1358             sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1359         }
1360         rvec_inc(fshift[s], sum);
1361     }
1362 }