2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
46 #include "nbnxn_consts.h"
47 #include "nbnxn_internal.h"
48 #include "nbnxn_atomdata.h"
49 #include "nbnxn_search.h"
50 #include "gromacs/utility/gmxomp.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 /* Stores the LJ parameter data in a format convenient for the SIMD kernels */
364 static void set_ljparam_simd_data(nbnxn_atomdata_t *nbat)
371 /* nbfp_s4 stores two parameters using a stride of 4,
372 * because this would suit x86 SIMD single-precision
373 * quad-load intrinsics. There's a slight inefficiency in
374 * allocating and initializing nbfp_s4 when it might not
375 * be used, but introducing the conditional code is not
376 * really worth it. */
377 nbat->alloc((void **)&nbat->nbfp_s4, nt*nt*4*sizeof(*nbat->nbfp_s4));
378 for (i = 0; i < nt; i++)
380 for (j = 0; j < nt; j++)
382 nbat->nbfp_s4[(i*nt+j)*4+0] = nbat->nbfp[(i*nt+j)*2+0];
383 nbat->nbfp_s4[(i*nt+j)*4+1] = nbat->nbfp[(i*nt+j)*2+1];
384 nbat->nbfp_s4[(i*nt+j)*4+2] = 0;
385 nbat->nbfp_s4[(i*nt+j)*4+3] = 0;
389 switch (nbat->comb_rule)
392 nbat->comb_rule = ljcrGEOM;
394 for (i = 0; i < nt; i++)
396 /* Copy the diagonal from the nbfp matrix */
397 nbat->nbfp_comb[i*2 ] = sqrt(nbat->nbfp[(i*nt+i)*2 ]);
398 nbat->nbfp_comb[i*2+1] = sqrt(nbat->nbfp[(i*nt+i)*2+1]);
402 for (i = 0; i < nt; i++)
404 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
405 c6 = nbat->nbfp[(i*nt+i)*2 ];
406 c12 = nbat->nbfp[(i*nt+i)*2+1];
407 if (c6 > 0 && c12 > 0)
409 /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
410 * so we get 6*C6 and 12*C12 after combining.
412 nbat->nbfp_comb[i*2 ] = 0.5*pow(c12/c6, 1.0/6.0);
413 nbat->nbfp_comb[i*2+1] = sqrt(c6*c6/c12);
417 nbat->nbfp_comb[i*2 ] = 0;
418 nbat->nbfp_comb[i*2+1] = 0;
423 /* We always store the full matrix (see code above) */
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_REAL_WIDTH;
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_SIMD_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_REAL_WIDTH;
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 enbnxninitcombrule,
521 int ntype, const real *nbfp,
524 nbnxn_alloc_t *alloc,
530 gmx_bool simple, bCombGeom, bCombLB;
534 nbat->alloc = nbnxn_alloc_aligned;
542 nbat->free = nbnxn_free_aligned;
551 fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
553 nbat->ntype = ntype + 1;
554 nbat->alloc((void **)&nbat->nbfp,
555 nbat->ntype*nbat->ntype*2*sizeof(*nbat->nbfp));
556 nbat->alloc((void **)&nbat->nbfp_comb, nbat->ntype*2*sizeof(*nbat->nbfp_comb));
558 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
559 * force-field floating point parameters.
562 ptr = getenv("GMX_LJCOMB_TOL");
567 sscanf(ptr, "%lf", &dbl);
573 /* Temporarily fill nbat->nbfp_comb with sigma and epsilon
574 * to check for the LB rule.
576 for (i = 0; i < ntype; i++)
578 c6 = nbfp[(i*ntype+i)*2 ]/6.0;
579 c12 = nbfp[(i*ntype+i)*2+1]/12.0;
580 if (c6 > 0 && c12 > 0)
582 nbat->nbfp_comb[i*2 ] = pow(c12/c6, 1.0/6.0);
583 nbat->nbfp_comb[i*2+1] = 0.25*c6*c6/c12;
585 else if (c6 == 0 && c12 == 0)
587 nbat->nbfp_comb[i*2 ] = 0;
588 nbat->nbfp_comb[i*2+1] = 0;
592 /* Can not use LB rule with only dispersion or repulsion */
597 for (i = 0; i < nbat->ntype; i++)
599 for (j = 0; j < nbat->ntype; j++)
601 if (i < ntype && j < ntype)
603 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
604 * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
606 c6 = nbfp[(i*ntype+j)*2 ];
607 c12 = nbfp[(i*ntype+j)*2+1];
608 nbat->nbfp[(i*nbat->ntype+j)*2 ] = c6;
609 nbat->nbfp[(i*nbat->ntype+j)*2+1] = c12;
611 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
612 bCombGeom = bCombGeom &&
613 gmx_within_tol(c6*c6, nbfp[(i*ntype+i)*2 ]*nbfp[(j*ntype+j)*2 ], tol) &&
614 gmx_within_tol(c12*c12, nbfp[(i*ntype+i)*2+1]*nbfp[(j*ntype+j)*2+1], tol);
616 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
620 ((c6 == 0 && c12 == 0 &&
621 (nbat->nbfp_comb[i*2+1] == 0 || nbat->nbfp_comb[j*2+1] == 0)) ||
622 (c6 > 0 && c12 > 0 &&
623 gmx_within_tol(pow(c12/c6, 1.0/6.0), 0.5*(nbat->nbfp_comb[i*2]+nbat->nbfp_comb[j*2]), tol) &&
624 gmx_within_tol(0.25*c6*c6/c12, sqrt(nbat->nbfp_comb[i*2+1]*nbat->nbfp_comb[j*2+1]), tol)));
628 /* Add zero parameters for the additional dummy atom type */
629 nbat->nbfp[(i*nbat->ntype+j)*2 ] = 0;
630 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 0;
636 fprintf(debug, "Combination rules: geometric %d Lorentz-Berthelot %d\n",
640 simple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
642 switch (enbnxninitcombrule)
644 case enbnxninitcombruleDETECT:
645 /* We prefer the geometic combination rule,
646 * as that gives a slightly faster kernel than the LB rule.
650 nbat->comb_rule = ljcrGEOM;
654 nbat->comb_rule = ljcrLB;
658 nbat->comb_rule = ljcrNONE;
660 nbat->free(nbat->nbfp_comb);
665 if (nbat->comb_rule == ljcrNONE)
667 fprintf(fp, "Using full Lennard-Jones parameter combination matrix\n\n");
671 fprintf(fp, "Using %s Lennard-Jones combination rule\n\n",
672 nbat->comb_rule == ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
676 case enbnxninitcombruleGEOM:
677 nbat->comb_rule = ljcrGEOM;
679 case enbnxninitcombruleLB:
680 nbat->comb_rule = ljcrLB;
682 case enbnxninitcombruleNONE:
683 nbat->comb_rule = ljcrNONE;
685 nbat->free(nbat->nbfp_comb);
688 gmx_incons("Unknown enbnxninitcombrule");
693 set_ljparam_simd_data(nbat);
698 nbat->lj_comb = NULL;
703 switch (nb_kernel_type)
705 case nbnxnk4xN_SIMD_4xN:
706 case nbnxnk4xN_SIMD_2xNN:
707 pack_x = max(NBNXN_CPU_CLUSTER_I_SIZE,
708 nbnxn_kernel_to_cj_size(nb_kernel_type));
712 nbat->XFormat = nbatX4;
715 nbat->XFormat = nbatX8;
718 gmx_incons("Unsupported packing width");
722 nbat->XFormat = nbatXYZ;
726 nbat->FFormat = nbat->XFormat;
730 nbat->XFormat = nbatXYZQ;
731 nbat->FFormat = nbatXYZ;
734 nbat->nenergrp = n_energygroups;
737 /* Energy groups not supported yet for super-sub lists */
738 if (n_energygroups > 1 && fp != NULL)
740 fprintf(fp, "\nNOTE: With GPUs, reporting energy group contributions is not supported\n\n");
744 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
745 if (nbat->nenergrp > 64)
747 gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
750 while (nbat->nenergrp > (1<<nbat->neg_2log))
754 nbat->energrp = NULL;
755 nbat->alloc((void **)&nbat->shift_vec, SHIFTS*sizeof(*nbat->shift_vec));
756 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
757 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
760 #ifdef GMX_NBNXN_SIMD
763 nbnxn_atomdata_init_simple_exclusion_masks(nbat);
767 /* Initialize the output data structures */
769 snew(nbat->out, nbat->nout);
771 for (i = 0; i < nbat->nout; i++)
773 nbnxn_atomdata_output_init(&nbat->out[i],
775 nbat->nenergrp, 1<<nbat->neg_2log,
778 nbat->buffer_flags.flag = NULL;
779 nbat->buffer_flags.flag_nalloc = 0;
781 nth = gmx_omp_nthreads_get(emntNonbonded);
783 ptr = getenv("GMX_USE_TREEREDUCE");
786 nbat->bUseTreeReduce = strtol(ptr, 0, 10);
789 else if (nth > 8) /*on the CPU we currently don't benefit even at 32*/
791 nbat->bUseTreeReduce = 1;
796 nbat->bUseTreeReduce = 0;
798 if (nbat->bUseTreeReduce)
802 fprintf(fp, "Using tree force reduction\n\n");
804 snew(nbat->syncStep, nth);
808 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
809 const int *type, int na,
814 /* The LJ params follow the combination rule:
815 * copy the params for the type array to the atom array.
817 for (is = 0; is < na; is += PACK_X4)
819 for (k = 0; k < PACK_X4; k++)
822 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
823 ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
828 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
829 const int *type, int na,
834 /* The LJ params follow the combination rule:
835 * copy the params for the type array to the atom array.
837 for (is = 0; is < na; is += PACK_X8)
839 for (k = 0; k < PACK_X8; k++)
842 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
843 ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
848 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
849 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t *nbat,
851 const nbnxn_search_t nbs,
855 const nbnxn_grid_t *grid;
857 for (g = 0; g < ngrid; g++)
859 grid = &nbs->grid[g];
861 /* Loop over all columns and copy and fill */
862 for (i = 0; i < grid->ncx*grid->ncy; i++)
864 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
865 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
867 copy_int_to_nbat_int(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
868 type, nbat->ntype-1, nbat->type+ash);
870 if (nbat->comb_rule != ljcrNONE)
872 if (nbat->XFormat == nbatX4)
874 copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
875 nbat->type+ash, ncz*grid->na_sc,
876 nbat->lj_comb+ash*2);
878 else if (nbat->XFormat == nbatX8)
880 copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
881 nbat->type+ash, ncz*grid->na_sc,
882 nbat->lj_comb+ash*2);
889 /* Sets the charges in nbnxn_atomdata_t *nbat */
890 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t *nbat,
892 const nbnxn_search_t nbs,
895 int g, cxy, ncz, ash, na, na_round, i, j;
897 const nbnxn_grid_t *grid;
899 for (g = 0; g < ngrid; g++)
901 grid = &nbs->grid[g];
903 /* Loop over all columns and copy and fill */
904 for (cxy = 0; cxy < grid->ncx*grid->ncy; cxy++)
906 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
907 na = grid->cxy_na[cxy];
908 na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
910 if (nbat->XFormat == nbatXYZQ)
912 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
913 for (i = 0; i < na; i++)
915 *q = charge[nbs->a[ash+i]];
918 /* Complete the partially filled last cell with zeros */
919 for (; i < na_round; i++)
928 for (i = 0; i < na; i++)
930 *q = charge[nbs->a[ash+i]];
933 /* Complete the partially filled last cell with zeros */
934 for (; i < na_round; i++)
944 /* Copies the energy group indices to a reordered and packed array */
945 static void copy_egp_to_nbat_egps(const int *a, int na, int na_round,
946 int na_c, int bit_shift,
947 const int *in, int *innb)
953 for (i = 0; i < na; i += na_c)
955 /* Store na_c energy group numbers into one int */
957 for (sa = 0; sa < na_c; sa++)
962 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
967 /* Complete the partially filled last cell with fill */
968 for (; i < na_round; i += na_c)
974 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
975 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t *nbat,
977 const nbnxn_search_t nbs,
981 const nbnxn_grid_t *grid;
983 for (g = 0; g < ngrid; g++)
985 grid = &nbs->grid[g];
987 /* Loop over all columns and copy and fill */
988 for (i = 0; i < grid->ncx*grid->ncy; i++)
990 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
991 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
993 copy_egp_to_nbat_egps(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
994 nbat->na_c, nbat->neg_2log,
995 atinfo, nbat->energrp+(ash>>grid->na_c_2log));
1000 /* Sets all required atom parameter data in nbnxn_atomdata_t */
1001 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
1003 const nbnxn_search_t nbs,
1004 const t_mdatoms *mdatoms,
1009 if (locality == eatLocal)
1018 nbnxn_atomdata_set_atomtypes(nbat, ngrid, nbs, mdatoms->typeA);
1020 nbnxn_atomdata_set_charges(nbat, ngrid, nbs, mdatoms->chargeA);
1022 if (nbat->nenergrp > 1)
1024 nbnxn_atomdata_set_energygroups(nbat, ngrid, nbs, atinfo);
1028 /* Copies the shift vector array to nbnxn_atomdata_t */
1029 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox,
1031 nbnxn_atomdata_t *nbat)
1035 nbat->bDynamicBox = bDynamicBox;
1036 for (i = 0; i < SHIFTS; i++)
1038 copy_rvec(shift_vec[i], nbat->shift_vec[i]);
1042 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1043 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
1047 nbnxn_atomdata_t *nbat)
1070 nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
1073 nth = gmx_omp_nthreads_get(emntPairsearch);
1075 #pragma omp parallel for num_threads(nth) schedule(static)
1076 for (th = 0; th < nth; th++)
1080 for (g = g0; g < g1; g++)
1082 const nbnxn_grid_t *grid;
1083 int cxy0, cxy1, cxy;
1085 grid = &nbs->grid[g];
1087 cxy0 = (grid->ncx*grid->ncy* th +nth-1)/nth;
1088 cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
1090 for (cxy = cxy0; cxy < cxy1; cxy++)
1092 int na, ash, na_fill;
1094 na = grid->cxy_na[cxy];
1095 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1097 if (g == 0 && FillLocal)
1100 (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
1104 /* We fill only the real particle locations.
1105 * We assume the filling entries at the end have been
1106 * properly set before during ns.
1110 copy_rvec_to_nbat_real(nbs->a+ash, na, na_fill, x,
1111 nbat->XFormat, nbat->x, ash,
1119 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
1124 for (i = i0; i < i1; i++)
1131 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1133 real ** gmx_restrict src,
1141 /* The destination buffer contains data, add to it */
1142 for (i = i0; i < i1; i++)
1144 for (s = 0; s < nsrc; s++)
1146 dest[i] += src[s][i];
1152 /* The destination buffer is unitialized, set it first */
1153 for (i = i0; i < i1; i++)
1155 dest[i] = src[0][i];
1156 for (s = 1; s < nsrc; s++)
1158 dest[i] += src[s][i];
1165 nbnxn_atomdata_reduce_reals_simd(real gmx_unused * gmx_restrict dest,
1166 gmx_bool gmx_unused bDestSet,
1167 real gmx_unused ** gmx_restrict src,
1168 int gmx_unused nsrc,
1169 int gmx_unused i0, int gmx_unused i1)
1171 #ifdef GMX_NBNXN_SIMD
1172 /* The SIMD width here is actually independent of that in the kernels,
1173 * but we use the same width for simplicity (usually optimal anyhow).
1176 gmx_simd_real_t dest_SSE, src_SSE;
1180 for (i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1182 dest_SSE = gmx_simd_load_r(dest+i);
1183 for (s = 0; s < nsrc; s++)
1185 src_SSE = gmx_simd_load_r(src[s]+i);
1186 dest_SSE = gmx_simd_add_r(dest_SSE, src_SSE);
1188 gmx_simd_store_r(dest+i, dest_SSE);
1193 for (i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1195 dest_SSE = gmx_simd_load_r(src[0]+i);
1196 for (s = 1; s < nsrc; s++)
1198 src_SSE = gmx_simd_load_r(src[s]+i);
1199 dest_SSE = gmx_simd_add_r(dest_SSE, src_SSE);
1201 gmx_simd_store_r(dest+i, dest_SSE);
1207 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1209 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1210 const nbnxn_atomdata_t *nbat,
1211 nbnxn_atomdata_output_t *out,
1222 /* Loop over all columns and copy and fill */
1223 switch (nbat->FFormat)
1231 for (a = a0; a < a1; a++)
1233 i = cell[a]*nbat->fstride;
1236 f[a][YY] += fnb[i+1];
1237 f[a][ZZ] += fnb[i+2];
1242 for (a = a0; a < a1; a++)
1244 i = cell[a]*nbat->fstride;
1246 for (fa = 0; fa < nfa; fa++)
1248 f[a][XX] += out[fa].f[i];
1249 f[a][YY] += out[fa].f[i+1];
1250 f[a][ZZ] += out[fa].f[i+2];
1260 for (a = a0; a < a1; a++)
1262 i = X4_IND_A(cell[a]);
1264 f[a][XX] += fnb[i+XX*PACK_X4];
1265 f[a][YY] += fnb[i+YY*PACK_X4];
1266 f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1271 for (a = a0; a < a1; a++)
1273 i = X4_IND_A(cell[a]);
1275 for (fa = 0; fa < nfa; fa++)
1277 f[a][XX] += out[fa].f[i+XX*PACK_X4];
1278 f[a][YY] += out[fa].f[i+YY*PACK_X4];
1279 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1289 for (a = a0; a < a1; a++)
1291 i = X8_IND_A(cell[a]);
1293 f[a][XX] += fnb[i+XX*PACK_X8];
1294 f[a][YY] += fnb[i+YY*PACK_X8];
1295 f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1300 for (a = a0; a < a1; a++)
1302 i = X8_IND_A(cell[a]);
1304 for (fa = 0; fa < nfa; fa++)
1306 f[a][XX] += out[fa].f[i+XX*PACK_X8];
1307 f[a][YY] += out[fa].f[i+YY*PACK_X8];
1308 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1314 gmx_incons("Unsupported nbnxn_atomdata_t format");
1318 static gmx_inline unsigned char reverse_bits(unsigned char b)
1320 /* http://graphics.stanford.edu/~seander/bithacks.html#ReverseByteWith64BitsDiv */
1321 return (b * 0x0202020202ULL & 0x010884422010ULL) % 1023;
1324 static void nbnxn_atomdata_add_nbat_f_to_f_treereduce(const nbnxn_atomdata_t *nbat,
1327 const nbnxn_buffer_flags_t *flags = &nbat->buffer_flags;
1329 int next_pow2 = 1<<(gmx_log2i(nth-1)+1);
1331 assert(nbat->nout == nth); /* tree-reduce currently only works for nout==nth */
1333 memset(nbat->syncStep, 0, sizeof(*(nbat->syncStep))*nth);
1335 #pragma omp parallel num_threads(nth)
1341 th = gmx_omp_get_thread_num();
1343 for (group_size = 2; group_size < 2*next_pow2; group_size *= 2)
1345 int index[2], group_pos, partner_pos, wu;
1346 int partner_th = th ^ (group_size/2);
1351 /* wait on partner thread - replaces full barrier */
1352 int sync_th, sync_group_size;
1354 tMPI_Atomic_memory_barrier(); /* gurantee data is saved before marking work as done */
1355 tMPI_Atomic_set(&(nbat->syncStep[th]), group_size/2); /* mark previous step as completed */
1357 /* find thread to sync with. Equal to partner_th unless nth is not a power of two. */
1358 for (sync_th = partner_th, sync_group_size = group_size; sync_th >= nth && sync_group_size > 2; sync_group_size /= 2)
1360 sync_th &= ~(sync_group_size/4);
1362 if (sync_th < nth) /* otherwise nothing to sync index[1] will be >=nout */
1364 /* wait on the thread which computed input data in previous step */
1365 while (tMPI_Atomic_get((volatile tMPI_Atomic_t*)&(nbat->syncStep[sync_th])) < group_size/2)
1369 /* guarantee that no later load happens before wait loop is finisehd */
1370 tMPI_Atomic_memory_barrier();
1372 #else /* TMPI_ATOMICS */
1377 /* Calculate buffers to sum (result goes into first buffer) */
1378 group_pos = th % group_size;
1379 index[0] = th - group_pos;
1380 index[1] = index[0] + group_size/2;
1382 /* If no second buffer, nothing to do */
1383 if (index[1] >= nbat->nout && group_size > 2)
1388 #if NBNXN_BUFFERFLAG_MAX_THREADS > 256
1389 #error reverse_bits assumes max 256 threads
1391 /* Position is permuted so that one of the 2 vectors being added was computed on the same thread in the previous step.
1392 This improves locality and enables to sync with just a single thread between steps (=the levels in the btree).
1393 The permutation which allows this corresponds to reversing the bits of the group position.
1395 group_pos = reverse_bits(group_pos)/(256/group_size);
1397 partner_pos = group_pos ^ 1;
1399 /* loop over two work-units (own and partner) */
1400 for (wu = 0; wu < 2; wu++)
1404 if (partner_th < nth)
1406 break; /* partner exists we don't have to do his work */
1410 group_pos = partner_pos;
1414 /* Calculate the cell-block range for our thread */
1415 b0 = (flags->nflag* group_pos )/group_size;
1416 b1 = (flags->nflag*(group_pos+1))/group_size;
1418 for (b = b0; b < b1; b++)
1420 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1421 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1423 if ((flags->flag[b] & (1ULL<<index[1])) || group_size > 2)
1425 #ifdef GMX_NBNXN_SIMD
1426 nbnxn_atomdata_reduce_reals_simd
1428 nbnxn_atomdata_reduce_reals
1430 (nbat->out[index[0]].f,
1431 (flags->flag[b] & (1ULL<<index[0])) || group_size > 2,
1432 &(nbat->out[index[1]].f), 1, i0, i1);
1435 else if (!(flags->flag[b] & (1ULL<<index[0])))
1437 nbnxn_atomdata_clear_reals(nbat->out[index[0]].f,
1447 static void nbnxn_atomdata_add_nbat_f_to_f_stdreduce(const nbnxn_atomdata_t *nbat,
1451 #pragma omp parallel for num_threads(nth) schedule(static)
1452 for (th = 0; th < nth; th++)
1454 const nbnxn_buffer_flags_t *flags;
1458 real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1461 flags = &nbat->buffer_flags;
1463 /* Calculate the cell-block range for our thread */
1464 b0 = (flags->nflag* th )/nth;
1465 b1 = (flags->nflag*(th+1))/nth;
1467 for (b = b0; b < b1; b++)
1469 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1470 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1473 for (out = 1; out < nbat->nout; out++)
1475 if (flags->flag[b] & (1U<<out))
1477 fptr[nfptr++] = nbat->out[out].f;
1482 #ifdef GMX_NBNXN_SIMD
1483 nbnxn_atomdata_reduce_reals_simd
1485 nbnxn_atomdata_reduce_reals
1488 flags->flag[b] & (1U<<0),
1492 else if (!(flags->flag[b] & (1U<<0)))
1494 nbnxn_atomdata_clear_reals(nbat->out[0].f,
1501 /* Add the force array(s) from nbnxn_atomdata_t to f */
1502 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t nbs,
1504 const nbnxn_atomdata_t *nbat,
1510 nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1516 na = nbs->natoms_nonlocal;
1520 na = nbs->natoms_local;
1523 a0 = nbs->natoms_local;
1524 na = nbs->natoms_nonlocal - nbs->natoms_local;
1528 nth = gmx_omp_nthreads_get(emntNonbonded);
1532 if (locality != eatAll)
1534 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1537 /* Reduce the force thread output buffers into buffer 0, before adding
1538 * them to the, differently ordered, "real" force buffer.
1540 if (nbat->bUseTreeReduce)
1542 nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbat, nth);
1546 nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbat, nth);
1549 #pragma omp parallel for num_threads(nth) schedule(static)
1550 for (th = 0; th < nth; th++)
1552 nbnxn_atomdata_add_nbat_f_to_f_part(nbs, nbat,
1560 nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1563 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1564 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1567 const nbnxn_atomdata_output_t *out;
1574 for (s = 0; s < SHIFTS; s++)
1577 for (th = 0; th < nbat->nout; th++)
1579 sum[XX] += out[th].fshift[s*DIM+XX];
1580 sum[YY] += out[th].fshift[s*DIM+YY];
1581 sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1583 rvec_inc(fshift[s], sum);