2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2012, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
48 #include "nbnxn_consts.h"
49 #include "nbnxn_internal.h"
50 #include "nbnxn_search.h"
51 #include "gmx_omp_nthreads.h"
53 /* Default nbnxn allocation routine, allocates NBNXN_MEM_ALIGN byte aligned */
54 void nbnxn_alloc_aligned(void **ptr, size_t nbytes)
56 *ptr = save_malloc_aligned("ptr", __FILE__, __LINE__, nbytes, 1, NBNXN_MEM_ALIGN);
59 /* Free function for memory allocated with nbnxn_alloc_aligned */
60 void nbnxn_free_aligned(void *ptr)
65 /* Reallocation wrapper function for nbnxn data structures */
66 void nbnxn_realloc_void(void **ptr,
67 int nbytes_copy, int nbytes_new,
73 ma(&ptr_new, nbytes_new);
75 if (nbytes_new > 0 && ptr_new == NULL)
77 gmx_fatal(FARGS, "Allocation of %d bytes failed", nbytes_new);
82 if (nbytes_new < nbytes_copy)
84 gmx_incons("In nbnxn_realloc_void: new size less than copy size");
86 memcpy(ptr_new, *ptr, nbytes_copy);
95 /* Reallocate the nbnxn_atomdata_t for a size of n atoms */
96 void nbnxn_atomdata_realloc(nbnxn_atomdata_t *nbat, int n)
100 nbnxn_realloc_void((void **)&nbat->type,
101 nbat->natoms*sizeof(*nbat->type),
102 n*sizeof(*nbat->type),
103 nbat->alloc, nbat->free);
104 nbnxn_realloc_void((void **)&nbat->lj_comb,
105 nbat->natoms*2*sizeof(*nbat->lj_comb),
106 n*2*sizeof(*nbat->lj_comb),
107 nbat->alloc, nbat->free);
108 if (nbat->XFormat != nbatXYZQ)
110 nbnxn_realloc_void((void **)&nbat->q,
111 nbat->natoms*sizeof(*nbat->q),
113 nbat->alloc, nbat->free);
115 if (nbat->nenergrp > 1)
117 nbnxn_realloc_void((void **)&nbat->energrp,
118 nbat->natoms/nbat->na_c*sizeof(*nbat->energrp),
119 n/nbat->na_c*sizeof(*nbat->energrp),
120 nbat->alloc, nbat->free);
122 nbnxn_realloc_void((void **)&nbat->x,
123 nbat->natoms*nbat->xstride*sizeof(*nbat->x),
124 n*nbat->xstride*sizeof(*nbat->x),
125 nbat->alloc, nbat->free);
126 for (t = 0; t < nbat->nout; t++)
128 /* Allocate one element extra for possible signaling with CUDA */
129 nbnxn_realloc_void((void **)&nbat->out[t].f,
130 nbat->natoms*nbat->fstride*sizeof(*nbat->out[t].f),
131 n*nbat->fstride*sizeof(*nbat->out[t].f),
132 nbat->alloc, nbat->free);
137 /* Initializes an nbnxn_atomdata_output_t data structure */
138 static void nbnxn_atomdata_output_init(nbnxn_atomdata_output_t *out,
140 int nenergrp, int stride,
146 ma((void **)&out->fshift, SHIFTS*DIM*sizeof(*out->fshift));
147 out->nV = nenergrp*nenergrp;
148 ma((void **)&out->Vvdw, out->nV*sizeof(*out->Vvdw));
149 ma((void **)&out->Vc, out->nV*sizeof(*out->Vc ));
151 if (nb_kernel_type == nbnxnk4xN_SIMD_4xN ||
152 nb_kernel_type == nbnxnk4xN_SIMD_2xNN)
154 cj_size = nbnxn_kernel_to_cj_size(nb_kernel_type);
155 out->nVS = nenergrp*nenergrp*stride*(cj_size>>1)*cj_size;
156 ma((void **)&out->VSvdw, out->nVS*sizeof(*out->VSvdw));
157 ma((void **)&out->VSc, out->nVS*sizeof(*out->VSc ));
165 static void copy_int_to_nbat_int(const int *a, int na, int na_round,
166 const int *in, int fill, int *innb)
171 for (i = 0; i < na; i++)
173 innb[j++] = in[a[i]];
175 /* Complete the partially filled last cell with fill */
176 for (; i < na_round; i++)
182 static void clear_nbat_real(int na, int nbatFormat, real *xnb, int a0)
189 for (a = 0; a < na; a++)
191 for (d = 0; d < DIM; d++)
193 xnb[(a0+a)*STRIDE_XYZ+d] = 0;
198 for (a = 0; a < na; a++)
200 for (d = 0; d < DIM; d++)
202 xnb[(a0+a)*STRIDE_XYZQ+d] = 0;
208 c = a0 & (PACK_X4-1);
209 for (a = 0; a < na; a++)
211 xnb[j+XX*PACK_X4] = 0;
212 xnb[j+YY*PACK_X4] = 0;
213 xnb[j+ZZ*PACK_X4] = 0;
218 j += (DIM-1)*PACK_X4;
225 c = a0 & (PACK_X8-1);
226 for (a = 0; a < na; a++)
228 xnb[j+XX*PACK_X8] = 0;
229 xnb[j+YY*PACK_X8] = 0;
230 xnb[j+ZZ*PACK_X8] = 0;
235 j += (DIM-1)*PACK_X8;
243 void copy_rvec_to_nbat_real(const int *a, int na, int na_round,
244 rvec *x, int nbatFormat, real *xnb, int a0,
245 int cx, int cy, int cz)
249 /* We might need to place filler particles to fill up the cell to na_round.
250 * The coefficients (LJ and q) for such particles are zero.
251 * But we might still get NaN as 0*NaN when distances are too small.
252 * We hope that -107 nm is far away enough from to zero
253 * to avoid accidental short distances to particles shifted down for pbc.
255 #define NBAT_FAR_AWAY 107
261 for (i = 0; i < na; i++)
263 xnb[j++] = x[a[i]][XX];
264 xnb[j++] = x[a[i]][YY];
265 xnb[j++] = x[a[i]][ZZ];
267 /* Complete the partially filled last cell with copies of the last element.
268 * This simplifies the bounding box calculation and avoid
269 * numerical issues with atoms that are coincidentally close.
271 for (; i < na_round; i++)
273 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
274 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
275 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
280 for (i = 0; i < na; i++)
282 xnb[j++] = x[a[i]][XX];
283 xnb[j++] = x[a[i]][YY];
284 xnb[j++] = x[a[i]][ZZ];
287 /* Complete the partially filled last cell with particles far apart */
288 for (; i < na_round; i++)
290 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
291 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
292 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
298 c = a0 & (PACK_X4-1);
299 for (i = 0; i < na; i++)
301 xnb[j+XX*PACK_X4] = x[a[i]][XX];
302 xnb[j+YY*PACK_X4] = x[a[i]][YY];
303 xnb[j+ZZ*PACK_X4] = x[a[i]][ZZ];
308 j += (DIM-1)*PACK_X4;
312 /* Complete the partially filled last cell with particles far apart */
313 for (; i < na_round; i++)
315 xnb[j+XX*PACK_X4] = -NBAT_FAR_AWAY*(1 + cx);
316 xnb[j+YY*PACK_X4] = -NBAT_FAR_AWAY*(1 + cy);
317 xnb[j+ZZ*PACK_X4] = -NBAT_FAR_AWAY*(1 + cz + i);
322 j += (DIM-1)*PACK_X4;
329 c = a0 & (PACK_X8 - 1);
330 for (i = 0; i < na; i++)
332 xnb[j+XX*PACK_X8] = x[a[i]][XX];
333 xnb[j+YY*PACK_X8] = x[a[i]][YY];
334 xnb[j+ZZ*PACK_X8] = x[a[i]][ZZ];
339 j += (DIM-1)*PACK_X8;
343 /* Complete the partially filled last cell with particles far apart */
344 for (; i < na_round; i++)
346 xnb[j+XX*PACK_X8] = -NBAT_FAR_AWAY*(1 + cx);
347 xnb[j+YY*PACK_X8] = -NBAT_FAR_AWAY*(1 + cy);
348 xnb[j+ZZ*PACK_X8] = -NBAT_FAR_AWAY*(1 + cz + i);
353 j += (DIM-1)*PACK_X8;
359 gmx_incons("Unsupported nbnxn_atomdata_t format");
363 /* Determines the combination rule (or none) to be used, stores it,
364 * and sets the LJ parameters required with the rule.
366 static void set_combination_rule_data(nbnxn_atomdata_t *nbat)
373 switch (nbat->comb_rule)
376 nbat->comb_rule = ljcrGEOM;
378 for (i = 0; i < nt; i++)
380 /* Copy the diagonal from the nbfp matrix */
381 nbat->nbfp_comb[i*2 ] = sqrt(nbat->nbfp[(i*nt+i)*2 ]);
382 nbat->nbfp_comb[i*2+1] = sqrt(nbat->nbfp[(i*nt+i)*2+1]);
386 for (i = 0; i < nt; i++)
388 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
389 c6 = nbat->nbfp[(i*nt+i)*2 ];
390 c12 = nbat->nbfp[(i*nt+i)*2+1];
391 if (c6 > 0 && c12 > 0)
393 /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
394 * so we get 6*C6 and 12*C12 after combining.
396 nbat->nbfp_comb[i*2 ] = 0.5*pow(c12/c6, 1.0/6.0);
397 nbat->nbfp_comb[i*2+1] = sqrt(c6*c6/c12);
401 nbat->nbfp_comb[i*2 ] = 0;
402 nbat->nbfp_comb[i*2+1] = 0;
407 /* nbfp_s4 stores two parameters using a stride of 4,
408 * because this would suit x86 SIMD single-precision
409 * quad-load intrinsics. There's a slight inefficiency in
410 * allocating and initializing nbfp_s4 when it might not
411 * be used, but introducing the conditional code is not
412 * really worth it. */
413 nbat->alloc((void **)&nbat->nbfp_s4, nt*nt*4*sizeof(*nbat->nbfp_s4));
414 for (i = 0; i < nt; i++)
416 for (j = 0; j < nt; j++)
418 nbat->nbfp_s4[(i*nt+j)*4+0] = nbat->nbfp[(i*nt+j)*2+0];
419 nbat->nbfp_s4[(i*nt+j)*4+1] = nbat->nbfp[(i*nt+j)*2+1];
420 nbat->nbfp_s4[(i*nt+j)*4+2] = 0;
421 nbat->nbfp_s4[(i*nt+j)*4+3] = 0;
426 gmx_incons("Unknown combination rule");
431 #ifdef GMX_NBNXN_SIMD
433 nbnxn_atomdata_init_simple_exclusion_masks(nbnxn_atomdata_t *nbat)
436 const int simd_width = GMX_SIMD_WIDTH_HERE;
438 /* Set the diagonal cluster pair exclusion mask setup data.
439 * In the kernel we check 0 < j - i to generate the masks.
440 * Here we store j - i for generating the mask for the first i,
441 * we substract 0.5 to avoid rounding issues.
442 * In the kernel we can subtract 1 to generate the subsequent mask.
444 int simd_4xn_diag_size;
445 const real simdFalse = -1, simdTrue = 1;
446 real *simd_interaction_array;
448 simd_4xn_diag_size = max(NBNXN_CPU_CLUSTER_I_SIZE, simd_width);
449 snew_aligned(nbat->simd_4xn_diagonal_j_minus_i, simd_4xn_diag_size, NBNXN_MEM_ALIGN);
450 for (j = 0; j < simd_4xn_diag_size; j++)
452 nbat->simd_4xn_diagonal_j_minus_i[j] = j - 0.5;
455 snew_aligned(nbat->simd_2xnn_diagonal_j_minus_i, simd_width, NBNXN_MEM_ALIGN);
456 for (j = 0; j < simd_width/2; j++)
458 /* The j-cluster size is half the SIMD width */
459 nbat->simd_2xnn_diagonal_j_minus_i[j] = j - 0.5;
460 /* The next half of the SIMD width is for i + 1 */
461 nbat->simd_2xnn_diagonal_j_minus_i[simd_width/2+j] = j - 1 - 0.5;
464 /* We use up to 32 bits for exclusion masking.
465 * The same masks are used for the 4xN and 2x(N+N) kernels.
466 * The masks are read either into epi32 SIMD registers or into
467 * real SIMD registers (together with a cast).
468 * In single precision this means the real and epi32 SIMD registers
470 * In double precision the epi32 registers can be smaller than
471 * the real registers, so depending on the architecture, we might
472 * need to use two, identical, 32-bit masks per real.
474 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*simd_width;
475 snew_aligned(nbat->simd_exclusion_filter1, simd_excl_size, NBNXN_MEM_ALIGN);
476 snew_aligned(nbat->simd_exclusion_filter2, simd_excl_size*2, NBNXN_MEM_ALIGN);
478 for (j = 0; j < simd_excl_size; j++)
480 /* Set the consecutive bits for masking pair exclusions */
481 nbat->simd_exclusion_filter1[j] = (1U << j);
482 nbat->simd_exclusion_filter2[j*2 + 0] = (1U << j);
483 nbat->simd_exclusion_filter2[j*2 + 1] = (1U << j);
486 #if (defined GMX_CPU_ACCELERATION_IBM_QPX)
487 /* The QPX kernels shouldn't do the bit masking that is done on
488 * x86, because the SIMD units lack bit-wise operations. Instead,
489 * we generate a vector of all 2^4 possible ways an i atom
490 * interacts with its 4 j atoms. Each array entry contains
491 * simd_width signed ints that are read in a single SIMD
492 * load. These ints must contain values that will be interpreted
493 * as true and false when loaded in the SIMD floating-point
494 * registers, ie. any positive or any negative value,
495 * respectively. Each array entry encodes how this i atom will
496 * interact with the 4 j atoms. Matching code exists in
497 * set_ci_top_excls() to generate indices into this array. Those
498 * indices are used in the kernels. */
500 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
501 const int qpx_simd_width = GMX_SIMD_WIDTH_HERE;
502 snew_aligned(simd_interaction_array, simd_excl_size * qpx_simd_width, NBNXN_MEM_ALIGN);
503 for (j = 0; j < simd_excl_size; j++)
505 int index = j * qpx_simd_width;
506 for (i = 0; i < qpx_simd_width; i++)
508 simd_interaction_array[index + i] = (j & (1 << i)) ? simdTrue : simdFalse;
511 nbat->simd_interaction_array = simd_interaction_array;
516 /* Initializes an nbnxn_atomdata_t data structure */
517 void nbnxn_atomdata_init(FILE *fp,
518 nbnxn_atomdata_t *nbat,
520 int ntype, const real *nbfp,
523 nbnxn_alloc_t *alloc,
529 gmx_bool simple, bCombGeom, bCombLB;
533 nbat->alloc = nbnxn_alloc_aligned;
541 nbat->free = nbnxn_free_aligned;
550 fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
552 nbat->ntype = ntype + 1;
553 nbat->alloc((void **)&nbat->nbfp,
554 nbat->ntype*nbat->ntype*2*sizeof(*nbat->nbfp));
555 nbat->alloc((void **)&nbat->nbfp_comb, nbat->ntype*2*sizeof(*nbat->nbfp_comb));
557 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
558 * force-field floating point parameters.
561 ptr = getenv("GMX_LJCOMB_TOL");
566 sscanf(ptr, "%lf", &dbl);
572 /* Temporarily fill nbat->nbfp_comb with sigma and epsilon
573 * to check for the LB rule.
575 for (i = 0; i < ntype; i++)
577 c6 = nbfp[(i*ntype+i)*2 ]/6.0;
578 c12 = nbfp[(i*ntype+i)*2+1]/12.0;
579 if (c6 > 0 && c12 > 0)
581 nbat->nbfp_comb[i*2 ] = pow(c12/c6, 1.0/6.0);
582 nbat->nbfp_comb[i*2+1] = 0.25*c6*c6/c12;
584 else if (c6 == 0 && c12 == 0)
586 nbat->nbfp_comb[i*2 ] = 0;
587 nbat->nbfp_comb[i*2+1] = 0;
591 /* Can not use LB rule with only dispersion or repulsion */
596 for (i = 0; i < nbat->ntype; i++)
598 for (j = 0; j < nbat->ntype; j++)
600 if (i < ntype && j < ntype)
602 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
603 * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
605 c6 = nbfp[(i*ntype+j)*2 ];
606 c12 = nbfp[(i*ntype+j)*2+1];
607 nbat->nbfp[(i*nbat->ntype+j)*2 ] = c6;
608 nbat->nbfp[(i*nbat->ntype+j)*2+1] = c12;
610 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
611 bCombGeom = bCombGeom &&
612 gmx_within_tol(c6*c6, nbfp[(i*ntype+i)*2 ]*nbfp[(j*ntype+j)*2 ], tol) &&
613 gmx_within_tol(c12*c12, nbfp[(i*ntype+i)*2+1]*nbfp[(j*ntype+j)*2+1], tol);
615 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
619 ((c6 == 0 && c12 == 0 &&
620 (nbat->nbfp_comb[i*2+1] == 0 || nbat->nbfp_comb[j*2+1] == 0)) ||
621 (c6 > 0 && c12 > 0 &&
622 gmx_within_tol(pow(c12/c6, 1.0/6.0), 0.5*(nbat->nbfp_comb[i*2]+nbat->nbfp_comb[j*2]), tol) &&
623 gmx_within_tol(0.25*c6*c6/c12, sqrt(nbat->nbfp_comb[i*2+1]*nbat->nbfp_comb[j*2+1]), tol)));
627 /* Add zero parameters for the additional dummy atom type */
628 nbat->nbfp[(i*nbat->ntype+j)*2 ] = 0;
629 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 0;
635 fprintf(debug, "Combination rules: geometric %d Lorentz-Berthelot %d\n",
639 simple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
643 /* We prefer the geometic combination rule,
644 * as that gives a slightly faster kernel than the LB rule.
648 nbat->comb_rule = ljcrGEOM;
652 nbat->comb_rule = ljcrLB;
656 nbat->comb_rule = ljcrNONE;
658 nbat->free(nbat->nbfp_comb);
663 if (nbat->comb_rule == ljcrNONE)
665 fprintf(fp, "Using full Lennard-Jones parameter combination matrix\n\n");
669 fprintf(fp, "Using %s Lennard-Jones combination rule\n\n",
670 nbat->comb_rule == ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
674 set_combination_rule_data(nbat);
678 nbat->comb_rule = ljcrNONE;
680 nbat->free(nbat->nbfp_comb);
685 nbat->lj_comb = NULL;
690 switch (nb_kernel_type)
692 case nbnxnk4xN_SIMD_4xN:
693 case nbnxnk4xN_SIMD_2xNN:
694 pack_x = max(NBNXN_CPU_CLUSTER_I_SIZE,
695 nbnxn_kernel_to_cj_size(nb_kernel_type));
699 nbat->XFormat = nbatX4;
702 nbat->XFormat = nbatX8;
705 gmx_incons("Unsupported packing width");
709 nbat->XFormat = nbatXYZ;
713 nbat->FFormat = nbat->XFormat;
717 nbat->XFormat = nbatXYZQ;
718 nbat->FFormat = nbatXYZ;
721 nbat->nenergrp = n_energygroups;
724 /* Energy groups not supported yet for super-sub lists */
725 if (n_energygroups > 1 && fp != NULL)
727 fprintf(fp, "\nNOTE: With GPUs, reporting energy group contributions is not supported\n\n");
731 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
732 if (nbat->nenergrp > 64)
734 gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
737 while (nbat->nenergrp > (1<<nbat->neg_2log))
741 nbat->energrp = NULL;
742 nbat->alloc((void **)&nbat->shift_vec, SHIFTS*sizeof(*nbat->shift_vec));
743 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
744 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
747 #ifdef GMX_NBNXN_SIMD
750 nbnxn_atomdata_init_simple_exclusion_masks(nbat);
754 /* Initialize the output data structures */
756 snew(nbat->out, nbat->nout);
758 for (i = 0; i < nbat->nout; i++)
760 nbnxn_atomdata_output_init(&nbat->out[i],
762 nbat->nenergrp, 1<<nbat->neg_2log,
765 nbat->buffer_flags.flag = NULL;
766 nbat->buffer_flags.flag_nalloc = 0;
769 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
770 const int *type, int na,
775 /* The LJ params follow the combination rule:
776 * copy the params for the type array to the atom array.
778 for (is = 0; is < na; is += PACK_X4)
780 for (k = 0; k < PACK_X4; k++)
783 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
784 ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
789 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
790 const int *type, int na,
795 /* The LJ params follow the combination rule:
796 * copy the params for the type array to the atom array.
798 for (is = 0; is < na; is += PACK_X8)
800 for (k = 0; k < PACK_X8; k++)
803 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
804 ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
809 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
810 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t *nbat,
812 const nbnxn_search_t nbs,
816 const nbnxn_grid_t *grid;
818 for (g = 0; g < ngrid; g++)
820 grid = &nbs->grid[g];
822 /* Loop over all columns and copy and fill */
823 for (i = 0; i < grid->ncx*grid->ncy; i++)
825 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
826 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
828 copy_int_to_nbat_int(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
829 type, nbat->ntype-1, nbat->type+ash);
831 if (nbat->comb_rule != ljcrNONE)
833 if (nbat->XFormat == nbatX4)
835 copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
836 nbat->type+ash, ncz*grid->na_sc,
837 nbat->lj_comb+ash*2);
839 else if (nbat->XFormat == nbatX8)
841 copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
842 nbat->type+ash, ncz*grid->na_sc,
843 nbat->lj_comb+ash*2);
850 /* Sets the charges in nbnxn_atomdata_t *nbat */
851 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t *nbat,
853 const nbnxn_search_t nbs,
856 int g, cxy, ncz, ash, na, na_round, i, j;
858 const nbnxn_grid_t *grid;
860 for (g = 0; g < ngrid; g++)
862 grid = &nbs->grid[g];
864 /* Loop over all columns and copy and fill */
865 for (cxy = 0; cxy < grid->ncx*grid->ncy; cxy++)
867 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
868 na = grid->cxy_na[cxy];
869 na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
871 if (nbat->XFormat == nbatXYZQ)
873 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
874 for (i = 0; i < na; i++)
876 *q = charge[nbs->a[ash+i]];
879 /* Complete the partially filled last cell with zeros */
880 for (; i < na_round; i++)
889 for (i = 0; i < na; i++)
891 *q = charge[nbs->a[ash+i]];
894 /* Complete the partially filled last cell with zeros */
895 for (; i < na_round; i++)
905 /* Copies the energy group indices to a reordered and packed array */
906 static void copy_egp_to_nbat_egps(const int *a, int na, int na_round,
907 int na_c, int bit_shift,
908 const int *in, int *innb)
914 for (i = 0; i < na; i += na_c)
916 /* Store na_c energy group numbers into one int */
918 for (sa = 0; sa < na_c; sa++)
923 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
928 /* Complete the partially filled last cell with fill */
929 for (; i < na_round; i += na_c)
935 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
936 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t *nbat,
938 const nbnxn_search_t nbs,
942 const nbnxn_grid_t *grid;
944 for (g = 0; g < ngrid; g++)
946 grid = &nbs->grid[g];
948 /* Loop over all columns and copy and fill */
949 for (i = 0; i < grid->ncx*grid->ncy; i++)
951 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
952 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
954 copy_egp_to_nbat_egps(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
955 nbat->na_c, nbat->neg_2log,
956 atinfo, nbat->energrp+(ash>>grid->na_c_2log));
961 /* Sets all required atom parameter data in nbnxn_atomdata_t */
962 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
964 const nbnxn_search_t nbs,
965 const t_mdatoms *mdatoms,
970 if (locality == eatLocal)
979 nbnxn_atomdata_set_atomtypes(nbat, ngrid, nbs, mdatoms->typeA);
981 nbnxn_atomdata_set_charges(nbat, ngrid, nbs, mdatoms->chargeA);
983 if (nbat->nenergrp > 1)
985 nbnxn_atomdata_set_energygroups(nbat, ngrid, nbs, atinfo);
989 /* Copies the shift vector array to nbnxn_atomdata_t */
990 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox,
992 nbnxn_atomdata_t *nbat)
996 nbat->bDynamicBox = bDynamicBox;
997 for (i = 0; i < SHIFTS; i++)
999 copy_rvec(shift_vec[i], nbat->shift_vec[i]);
1003 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1004 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
1008 nbnxn_atomdata_t *nbat)
1031 nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
1034 nth = gmx_omp_nthreads_get(emntPairsearch);
1036 #pragma omp parallel for num_threads(nth) schedule(static)
1037 for (th = 0; th < nth; th++)
1041 for (g = g0; g < g1; g++)
1043 const nbnxn_grid_t *grid;
1044 int cxy0, cxy1, cxy;
1046 grid = &nbs->grid[g];
1048 cxy0 = (grid->ncx*grid->ncy* th +nth-1)/nth;
1049 cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
1051 for (cxy = cxy0; cxy < cxy1; cxy++)
1053 int na, ash, na_fill;
1055 na = grid->cxy_na[cxy];
1056 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1058 if (g == 0 && FillLocal)
1061 (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
1065 /* We fill only the real particle locations.
1066 * We assume the filling entries at the end have been
1067 * properly set before during ns.
1071 copy_rvec_to_nbat_real(nbs->a+ash, na, na_fill, x,
1072 nbat->XFormat, nbat->x, ash,
1080 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
1085 for (i = i0; i < i1; i++)
1092 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1094 real ** gmx_restrict src,
1102 /* The destination buffer contains data, add to it */
1103 for (i = i0; i < i1; i++)
1105 for (s = 0; s < nsrc; s++)
1107 dest[i] += src[s][i];
1113 /* The destination buffer is unitialized, set it first */
1114 for (i = i0; i < i1; i++)
1116 dest[i] = src[0][i];
1117 for (s = 1; s < nsrc; s++)
1119 dest[i] += src[s][i];
1126 nbnxn_atomdata_reduce_reals_simd(real * gmx_restrict dest,
1128 real ** gmx_restrict src,
1132 #ifdef GMX_NBNXN_SIMD
1133 /* The SIMD width here is actually independent of that in the kernels,
1134 * but we use the same width for simplicity (usually optimal anyhow).
1137 gmx_mm_pr dest_SSE, src_SSE;
1141 for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1143 dest_SSE = gmx_load_pr(dest+i);
1144 for (s = 0; s < nsrc; s++)
1146 src_SSE = gmx_load_pr(src[s]+i);
1147 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1149 gmx_store_pr(dest+i, dest_SSE);
1154 for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1156 dest_SSE = gmx_load_pr(src[0]+i);
1157 for (s = 1; s < nsrc; s++)
1159 src_SSE = gmx_load_pr(src[s]+i);
1160 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1162 gmx_store_pr(dest+i, dest_SSE);
1168 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1170 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1171 const nbnxn_atomdata_t *nbat,
1172 nbnxn_atomdata_output_t *out,
1183 /* Loop over all columns and copy and fill */
1184 switch (nbat->FFormat)
1192 for (a = a0; a < a1; a++)
1194 i = cell[a]*nbat->fstride;
1197 f[a][YY] += fnb[i+1];
1198 f[a][ZZ] += fnb[i+2];
1203 for (a = a0; a < a1; a++)
1205 i = cell[a]*nbat->fstride;
1207 for (fa = 0; fa < nfa; fa++)
1209 f[a][XX] += out[fa].f[i];
1210 f[a][YY] += out[fa].f[i+1];
1211 f[a][ZZ] += out[fa].f[i+2];
1221 for (a = a0; a < a1; a++)
1223 i = X4_IND_A(cell[a]);
1225 f[a][XX] += fnb[i+XX*PACK_X4];
1226 f[a][YY] += fnb[i+YY*PACK_X4];
1227 f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1232 for (a = a0; a < a1; a++)
1234 i = X4_IND_A(cell[a]);
1236 for (fa = 0; fa < nfa; fa++)
1238 f[a][XX] += out[fa].f[i+XX*PACK_X4];
1239 f[a][YY] += out[fa].f[i+YY*PACK_X4];
1240 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1250 for (a = a0; a < a1; a++)
1252 i = X8_IND_A(cell[a]);
1254 f[a][XX] += fnb[i+XX*PACK_X8];
1255 f[a][YY] += fnb[i+YY*PACK_X8];
1256 f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1261 for (a = a0; a < a1; a++)
1263 i = X8_IND_A(cell[a]);
1265 for (fa = 0; fa < nfa; fa++)
1267 f[a][XX] += out[fa].f[i+XX*PACK_X8];
1268 f[a][YY] += out[fa].f[i+YY*PACK_X8];
1269 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1275 gmx_incons("Unsupported nbnxn_atomdata_t format");
1279 /* Add the force array(s) from nbnxn_atomdata_t to f */
1280 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t nbs,
1282 const nbnxn_atomdata_t *nbat,
1288 nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1294 na = nbs->natoms_nonlocal;
1298 na = nbs->natoms_local;
1301 a0 = nbs->natoms_local;
1302 na = nbs->natoms_nonlocal - nbs->natoms_local;
1306 nth = gmx_omp_nthreads_get(emntNonbonded);
1310 if (locality != eatAll)
1312 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1315 /* Reduce the force thread output buffers into buffer 0, before adding
1316 * them to the, differently ordered, "real" force buffer.
1318 #pragma omp parallel for num_threads(nth) schedule(static)
1319 for (th = 0; th < nth; th++)
1321 const nbnxn_buffer_flags_t *flags;
1325 real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1328 flags = &nbat->buffer_flags;
1330 /* Calculate the cell-block range for our thread */
1331 b0 = (flags->nflag* th )/nth;
1332 b1 = (flags->nflag*(th+1))/nth;
1334 for (b = b0; b < b1; b++)
1336 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1337 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1340 for (out = 1; out < nbat->nout; out++)
1342 if (flags->flag[b] & (1U<<out))
1344 fptr[nfptr++] = nbat->out[out].f;
1349 #ifdef GMX_NBNXN_SIMD
1350 nbnxn_atomdata_reduce_reals_simd
1352 nbnxn_atomdata_reduce_reals
1355 flags->flag[b] & (1U<<0),
1359 else if (!(flags->flag[b] & (1U<<0)))
1361 nbnxn_atomdata_clear_reals(nbat->out[0].f,
1368 #pragma omp parallel for num_threads(nth) schedule(static)
1369 for (th = 0; th < nth; th++)
1371 nbnxn_atomdata_add_nbat_f_to_f_part(nbs, nbat,
1379 nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1382 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1383 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1386 const nbnxn_atomdata_output_t *out;
1393 for (s = 0; s < SHIFTS; s++)
1396 for (th = 0; th < nbat->nout; th++)
1398 sum[XX] += out[th].fshift[s*DIM+XX];
1399 sum[YY] += out[th].fshift[s*DIM+YY];
1400 sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1402 rvec_inc(fshift[s], sum);