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_search.h"
49 #include "gromacs/utility/gmxomp.h"
50 #include "gmx_omp_nthreads.h"
52 /* Default nbnxn allocation routine, allocates NBNXN_MEM_ALIGN byte aligned */
53 void nbnxn_alloc_aligned(void **ptr, size_t nbytes)
55 *ptr = save_malloc_aligned("ptr", __FILE__, __LINE__, nbytes, 1, NBNXN_MEM_ALIGN);
58 /* Free function for memory allocated with nbnxn_alloc_aligned */
59 void nbnxn_free_aligned(void *ptr)
64 /* Reallocation wrapper function for nbnxn data structures */
65 void nbnxn_realloc_void(void **ptr,
66 int nbytes_copy, int nbytes_new,
72 ma(&ptr_new, nbytes_new);
74 if (nbytes_new > 0 && ptr_new == NULL)
76 gmx_fatal(FARGS, "Allocation of %d bytes failed", nbytes_new);
81 if (nbytes_new < nbytes_copy)
83 gmx_incons("In nbnxn_realloc_void: new size less than copy size");
85 memcpy(ptr_new, *ptr, nbytes_copy);
94 /* Reallocate the nbnxn_atomdata_t for a size of n atoms */
95 void nbnxn_atomdata_realloc(nbnxn_atomdata_t *nbat, int n)
99 nbnxn_realloc_void((void **)&nbat->type,
100 nbat->natoms*sizeof(*nbat->type),
101 n*sizeof(*nbat->type),
102 nbat->alloc, nbat->free);
103 nbnxn_realloc_void((void **)&nbat->lj_comb,
104 nbat->natoms*2*sizeof(*nbat->lj_comb),
105 n*2*sizeof(*nbat->lj_comb),
106 nbat->alloc, nbat->free);
107 if (nbat->XFormat != nbatXYZQ)
109 nbnxn_realloc_void((void **)&nbat->q,
110 nbat->natoms*sizeof(*nbat->q),
112 nbat->alloc, nbat->free);
114 if (nbat->nenergrp > 1)
116 nbnxn_realloc_void((void **)&nbat->energrp,
117 nbat->natoms/nbat->na_c*sizeof(*nbat->energrp),
118 n/nbat->na_c*sizeof(*nbat->energrp),
119 nbat->alloc, nbat->free);
121 nbnxn_realloc_void((void **)&nbat->x,
122 nbat->natoms*nbat->xstride*sizeof(*nbat->x),
123 n*nbat->xstride*sizeof(*nbat->x),
124 nbat->alloc, nbat->free);
125 for (t = 0; t < nbat->nout; t++)
127 /* Allocate one element extra for possible signaling with CUDA */
128 nbnxn_realloc_void((void **)&nbat->out[t].f,
129 nbat->natoms*nbat->fstride*sizeof(*nbat->out[t].f),
130 n*nbat->fstride*sizeof(*nbat->out[t].f),
131 nbat->alloc, nbat->free);
136 /* Initializes an nbnxn_atomdata_output_t data structure */
137 static void nbnxn_atomdata_output_init(nbnxn_atomdata_output_t *out,
139 int nenergrp, int stride,
145 ma((void **)&out->fshift, SHIFTS*DIM*sizeof(*out->fshift));
146 out->nV = nenergrp*nenergrp;
147 ma((void **)&out->Vvdw, out->nV*sizeof(*out->Vvdw));
148 ma((void **)&out->Vc, out->nV*sizeof(*out->Vc ));
150 if (nb_kernel_type == nbnxnk4xN_SIMD_4xN ||
151 nb_kernel_type == nbnxnk4xN_SIMD_2xNN)
153 cj_size = nbnxn_kernel_to_cj_size(nb_kernel_type);
154 out->nVS = nenergrp*nenergrp*stride*(cj_size>>1)*cj_size;
155 ma((void **)&out->VSvdw, out->nVS*sizeof(*out->VSvdw));
156 ma((void **)&out->VSc, out->nVS*sizeof(*out->VSc ));
164 static void copy_int_to_nbat_int(const int *a, int na, int na_round,
165 const int *in, int fill, int *innb)
170 for (i = 0; i < na; i++)
172 innb[j++] = in[a[i]];
174 /* Complete the partially filled last cell with fill */
175 for (; i < na_round; i++)
181 static void clear_nbat_real(int na, int nbatFormat, real *xnb, int a0)
188 for (a = 0; a < na; a++)
190 for (d = 0; d < DIM; d++)
192 xnb[(a0+a)*STRIDE_XYZ+d] = 0;
197 for (a = 0; a < na; a++)
199 for (d = 0; d < DIM; d++)
201 xnb[(a0+a)*STRIDE_XYZQ+d] = 0;
207 c = a0 & (PACK_X4-1);
208 for (a = 0; a < na; a++)
210 xnb[j+XX*PACK_X4] = 0;
211 xnb[j+YY*PACK_X4] = 0;
212 xnb[j+ZZ*PACK_X4] = 0;
217 j += (DIM-1)*PACK_X4;
224 c = a0 & (PACK_X8-1);
225 for (a = 0; a < na; a++)
227 xnb[j+XX*PACK_X8] = 0;
228 xnb[j+YY*PACK_X8] = 0;
229 xnb[j+ZZ*PACK_X8] = 0;
234 j += (DIM-1)*PACK_X8;
242 void copy_rvec_to_nbat_real(const int *a, int na, int na_round,
243 rvec *x, int nbatFormat, real *xnb, int a0,
244 int cx, int cy, int cz)
248 /* We might need to place filler particles to fill up the cell to na_round.
249 * The coefficients (LJ and q) for such particles are zero.
250 * But we might still get NaN as 0*NaN when distances are too small.
251 * We hope that -107 nm is far away enough from to zero
252 * to avoid accidental short distances to particles shifted down for pbc.
254 #define NBAT_FAR_AWAY 107
260 for (i = 0; i < na; i++)
262 xnb[j++] = x[a[i]][XX];
263 xnb[j++] = x[a[i]][YY];
264 xnb[j++] = x[a[i]][ZZ];
266 /* Complete the partially filled last cell with copies of the last element.
267 * This simplifies the bounding box calculation and avoid
268 * numerical issues with atoms that are coincidentally close.
270 for (; i < na_round; i++)
272 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
273 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
274 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
279 for (i = 0; i < na; i++)
281 xnb[j++] = x[a[i]][XX];
282 xnb[j++] = x[a[i]][YY];
283 xnb[j++] = x[a[i]][ZZ];
286 /* Complete the partially filled last cell with particles far apart */
287 for (; i < na_round; i++)
289 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
290 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
291 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
297 c = a0 & (PACK_X4-1);
298 for (i = 0; i < na; i++)
300 xnb[j+XX*PACK_X4] = x[a[i]][XX];
301 xnb[j+YY*PACK_X4] = x[a[i]][YY];
302 xnb[j+ZZ*PACK_X4] = x[a[i]][ZZ];
307 j += (DIM-1)*PACK_X4;
311 /* Complete the partially filled last cell with particles far apart */
312 for (; i < na_round; i++)
314 xnb[j+XX*PACK_X4] = -NBAT_FAR_AWAY*(1 + cx);
315 xnb[j+YY*PACK_X4] = -NBAT_FAR_AWAY*(1 + cy);
316 xnb[j+ZZ*PACK_X4] = -NBAT_FAR_AWAY*(1 + cz + i);
321 j += (DIM-1)*PACK_X4;
328 c = a0 & (PACK_X8 - 1);
329 for (i = 0; i < na; i++)
331 xnb[j+XX*PACK_X8] = x[a[i]][XX];
332 xnb[j+YY*PACK_X8] = x[a[i]][YY];
333 xnb[j+ZZ*PACK_X8] = x[a[i]][ZZ];
338 j += (DIM-1)*PACK_X8;
342 /* Complete the partially filled last cell with particles far apart */
343 for (; i < na_round; i++)
345 xnb[j+XX*PACK_X8] = -NBAT_FAR_AWAY*(1 + cx);
346 xnb[j+YY*PACK_X8] = -NBAT_FAR_AWAY*(1 + cy);
347 xnb[j+ZZ*PACK_X8] = -NBAT_FAR_AWAY*(1 + cz + i);
352 j += (DIM-1)*PACK_X8;
358 gmx_incons("Unsupported nbnxn_atomdata_t format");
362 /* Determines the combination rule (or none) to be used, stores it,
363 * and sets the LJ parameters required with the rule.
365 static void set_combination_rule_data(nbnxn_atomdata_t *nbat)
372 switch (nbat->comb_rule)
375 nbat->comb_rule = ljcrGEOM;
377 for (i = 0; i < nt; i++)
379 /* Copy the diagonal from the nbfp matrix */
380 nbat->nbfp_comb[i*2 ] = sqrt(nbat->nbfp[(i*nt+i)*2 ]);
381 nbat->nbfp_comb[i*2+1] = sqrt(nbat->nbfp[(i*nt+i)*2+1]);
385 for (i = 0; i < nt; i++)
387 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
388 c6 = nbat->nbfp[(i*nt+i)*2 ];
389 c12 = nbat->nbfp[(i*nt+i)*2+1];
390 if (c6 > 0 && c12 > 0)
392 /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
393 * so we get 6*C6 and 12*C12 after combining.
395 nbat->nbfp_comb[i*2 ] = 0.5*pow(c12/c6, 1.0/6.0);
396 nbat->nbfp_comb[i*2+1] = sqrt(c6*c6/c12);
400 nbat->nbfp_comb[i*2 ] = 0;
401 nbat->nbfp_comb[i*2+1] = 0;
406 /* nbfp_s4 stores two parameters using a stride of 4,
407 * because this would suit x86 SIMD single-precision
408 * quad-load intrinsics. There's a slight inefficiency in
409 * allocating and initializing nbfp_s4 when it might not
410 * be used, but introducing the conditional code is not
411 * really worth it. */
412 nbat->alloc((void **)&nbat->nbfp_s4, nt*nt*4*sizeof(*nbat->nbfp_s4));
413 for (i = 0; i < nt; i++)
415 for (j = 0; j < nt; j++)
417 nbat->nbfp_s4[(i*nt+j)*4+0] = nbat->nbfp[(i*nt+j)*2+0];
418 nbat->nbfp_s4[(i*nt+j)*4+1] = nbat->nbfp[(i*nt+j)*2+1];
419 nbat->nbfp_s4[(i*nt+j)*4+2] = 0;
420 nbat->nbfp_s4[(i*nt+j)*4+3] = 0;
425 gmx_incons("Unknown combination rule");
430 #ifdef GMX_NBNXN_SIMD
432 nbnxn_atomdata_init_simple_exclusion_masks(nbnxn_atomdata_t *nbat)
435 const int simd_width = GMX_SIMD_WIDTH_HERE;
437 /* Set the diagonal cluster pair exclusion mask setup data.
438 * In the kernel we check 0 < j - i to generate the masks.
439 * Here we store j - i for generating the mask for the first i,
440 * we substract 0.5 to avoid rounding issues.
441 * In the kernel we can subtract 1 to generate the subsequent mask.
443 int simd_4xn_diag_size;
444 const real simdFalse = -1, simdTrue = 1;
445 real *simd_interaction_array;
447 simd_4xn_diag_size = max(NBNXN_CPU_CLUSTER_I_SIZE, simd_width);
448 snew_aligned(nbat->simd_4xn_diagonal_j_minus_i, simd_4xn_diag_size, NBNXN_MEM_ALIGN);
449 for (j = 0; j < simd_4xn_diag_size; j++)
451 nbat->simd_4xn_diagonal_j_minus_i[j] = j - 0.5;
454 snew_aligned(nbat->simd_2xnn_diagonal_j_minus_i, simd_width, NBNXN_MEM_ALIGN);
455 for (j = 0; j < simd_width/2; j++)
457 /* The j-cluster size is half the SIMD width */
458 nbat->simd_2xnn_diagonal_j_minus_i[j] = j - 0.5;
459 /* The next half of the SIMD width is for i + 1 */
460 nbat->simd_2xnn_diagonal_j_minus_i[simd_width/2+j] = j - 1 - 0.5;
463 /* We use up to 32 bits for exclusion masking.
464 * The same masks are used for the 4xN and 2x(N+N) kernels.
465 * The masks are read either into epi32 SIMD registers or into
466 * real SIMD registers (together with a cast).
467 * In single precision this means the real and epi32 SIMD registers
469 * In double precision the epi32 registers can be smaller than
470 * the real registers, so depending on the architecture, we might
471 * need to use two, identical, 32-bit masks per real.
473 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*simd_width;
474 snew_aligned(nbat->simd_exclusion_filter1, simd_excl_size, NBNXN_MEM_ALIGN);
475 snew_aligned(nbat->simd_exclusion_filter2, simd_excl_size*2, NBNXN_MEM_ALIGN);
477 for (j = 0; j < simd_excl_size; j++)
479 /* Set the consecutive bits for masking pair exclusions */
480 nbat->simd_exclusion_filter1[j] = (1U << j);
481 nbat->simd_exclusion_filter2[j*2 + 0] = (1U << j);
482 nbat->simd_exclusion_filter2[j*2 + 1] = (1U << j);
485 #if (defined GMX_CPU_ACCELERATION_IBM_QPX)
486 /* The QPX kernels shouldn't do the bit masking that is done on
487 * x86, because the SIMD units lack bit-wise operations. Instead,
488 * we generate a vector of all 2^4 possible ways an i atom
489 * interacts with its 4 j atoms. Each array entry contains
490 * simd_width signed ints that are read in a single SIMD
491 * load. These ints must contain values that will be interpreted
492 * as true and false when loaded in the SIMD floating-point
493 * registers, ie. any positive or any negative value,
494 * respectively. Each array entry encodes how this i atom will
495 * interact with the 4 j atoms. Matching code exists in
496 * set_ci_top_excls() to generate indices into this array. Those
497 * indices are used in the kernels. */
499 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
500 const int qpx_simd_width = GMX_SIMD_WIDTH_HERE;
501 snew_aligned(simd_interaction_array, simd_excl_size * qpx_simd_width, NBNXN_MEM_ALIGN);
502 for (j = 0; j < simd_excl_size; j++)
504 int index = j * qpx_simd_width;
505 for (i = 0; i < qpx_simd_width; i++)
507 simd_interaction_array[index + i] = (j & (1 << i)) ? simdTrue : simdFalse;
510 nbat->simd_interaction_array = simd_interaction_array;
515 /* Initializes an nbnxn_atomdata_t data structure */
516 void nbnxn_atomdata_init(FILE *fp,
517 nbnxn_atomdata_t *nbat,
519 int ntype, const real *nbfp,
522 nbnxn_alloc_t *alloc,
528 gmx_bool simple, bCombGeom, bCombLB;
532 nbat->alloc = nbnxn_alloc_aligned;
540 nbat->free = nbnxn_free_aligned;
549 fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
551 nbat->ntype = ntype + 1;
552 nbat->alloc((void **)&nbat->nbfp,
553 nbat->ntype*nbat->ntype*2*sizeof(*nbat->nbfp));
554 nbat->alloc((void **)&nbat->nbfp_comb, nbat->ntype*2*sizeof(*nbat->nbfp_comb));
556 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
557 * force-field floating point parameters.
560 ptr = getenv("GMX_LJCOMB_TOL");
565 sscanf(ptr, "%lf", &dbl);
571 /* Temporarily fill nbat->nbfp_comb with sigma and epsilon
572 * to check for the LB rule.
574 for (i = 0; i < ntype; i++)
576 c6 = nbfp[(i*ntype+i)*2 ]/6.0;
577 c12 = nbfp[(i*ntype+i)*2+1]/12.0;
578 if (c6 > 0 && c12 > 0)
580 nbat->nbfp_comb[i*2 ] = pow(c12/c6, 1.0/6.0);
581 nbat->nbfp_comb[i*2+1] = 0.25*c6*c6/c12;
583 else if (c6 == 0 && c12 == 0)
585 nbat->nbfp_comb[i*2 ] = 0;
586 nbat->nbfp_comb[i*2+1] = 0;
590 /* Can not use LB rule with only dispersion or repulsion */
595 for (i = 0; i < nbat->ntype; i++)
597 for (j = 0; j < nbat->ntype; j++)
599 if (i < ntype && j < ntype)
601 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
602 * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
604 c6 = nbfp[(i*ntype+j)*2 ];
605 c12 = nbfp[(i*ntype+j)*2+1];
606 nbat->nbfp[(i*nbat->ntype+j)*2 ] = c6;
607 nbat->nbfp[(i*nbat->ntype+j)*2+1] = c12;
609 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
610 bCombGeom = bCombGeom &&
611 gmx_within_tol(c6*c6, nbfp[(i*ntype+i)*2 ]*nbfp[(j*ntype+j)*2 ], tol) &&
612 gmx_within_tol(c12*c12, nbfp[(i*ntype+i)*2+1]*nbfp[(j*ntype+j)*2+1], tol);
614 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
618 ((c6 == 0 && c12 == 0 &&
619 (nbat->nbfp_comb[i*2+1] == 0 || nbat->nbfp_comb[j*2+1] == 0)) ||
620 (c6 > 0 && c12 > 0 &&
621 gmx_within_tol(pow(c12/c6, 1.0/6.0), 0.5*(nbat->nbfp_comb[i*2]+nbat->nbfp_comb[j*2]), tol) &&
622 gmx_within_tol(0.25*c6*c6/c12, sqrt(nbat->nbfp_comb[i*2+1]*nbat->nbfp_comb[j*2+1]), tol)));
626 /* Add zero parameters for the additional dummy atom type */
627 nbat->nbfp[(i*nbat->ntype+j)*2 ] = 0;
628 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 0;
634 fprintf(debug, "Combination rules: geometric %d Lorentz-Berthelot %d\n",
638 simple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
642 /* We prefer the geometic combination rule,
643 * as that gives a slightly faster kernel than the LB rule.
647 nbat->comb_rule = ljcrGEOM;
651 nbat->comb_rule = ljcrLB;
655 nbat->comb_rule = ljcrNONE;
657 nbat->free(nbat->nbfp_comb);
662 if (nbat->comb_rule == ljcrNONE)
664 fprintf(fp, "Using full Lennard-Jones parameter combination matrix\n\n");
668 fprintf(fp, "Using %s Lennard-Jones combination rule\n\n",
669 nbat->comb_rule == ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
673 set_combination_rule_data(nbat);
677 nbat->comb_rule = ljcrNONE;
679 nbat->free(nbat->nbfp_comb);
684 nbat->lj_comb = NULL;
689 switch (nb_kernel_type)
691 case nbnxnk4xN_SIMD_4xN:
692 case nbnxnk4xN_SIMD_2xNN:
693 pack_x = max(NBNXN_CPU_CLUSTER_I_SIZE,
694 nbnxn_kernel_to_cj_size(nb_kernel_type));
698 nbat->XFormat = nbatX4;
701 nbat->XFormat = nbatX8;
704 gmx_incons("Unsupported packing width");
708 nbat->XFormat = nbatXYZ;
712 nbat->FFormat = nbat->XFormat;
716 nbat->XFormat = nbatXYZQ;
717 nbat->FFormat = nbatXYZ;
720 nbat->nenergrp = n_energygroups;
723 /* Energy groups not supported yet for super-sub lists */
724 if (n_energygroups > 1 && fp != NULL)
726 fprintf(fp, "\nNOTE: With GPUs, reporting energy group contributions is not supported\n\n");
730 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
731 if (nbat->nenergrp > 64)
733 gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
736 while (nbat->nenergrp > (1<<nbat->neg_2log))
740 nbat->energrp = NULL;
741 nbat->alloc((void **)&nbat->shift_vec, SHIFTS*sizeof(*nbat->shift_vec));
742 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
743 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
746 #ifdef GMX_NBNXN_SIMD
749 nbnxn_atomdata_init_simple_exclusion_masks(nbat);
753 /* Initialize the output data structures */
755 snew(nbat->out, nbat->nout);
757 for (i = 0; i < nbat->nout; i++)
759 nbnxn_atomdata_output_init(&nbat->out[i],
761 nbat->nenergrp, 1<<nbat->neg_2log,
764 nbat->buffer_flags.flag = NULL;
765 nbat->buffer_flags.flag_nalloc = 0;
767 nth = gmx_omp_nthreads_get(emntNonbonded);
769 ptr = getenv("GMX_USE_TREEREDUCE");
772 nbat->bUseTreeReduce = strtol(ptr, 0, 10);
775 else if (nth > 8) /*on the CPU we currently don't benefit even at 32*/
777 nbat->bUseTreeReduce = 1;
782 nbat->bUseTreeReduce = 0;
784 if (nbat->bUseTreeReduce)
788 fprintf(fp, "Using tree force reduction\n\n");
790 snew(nbat->syncStep, nth);
794 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
795 const int *type, int na,
800 /* The LJ params follow the combination rule:
801 * copy the params for the type array to the atom array.
803 for (is = 0; is < na; is += PACK_X4)
805 for (k = 0; k < PACK_X4; k++)
808 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
809 ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
814 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
815 const int *type, int na,
820 /* The LJ params follow the combination rule:
821 * copy the params for the type array to the atom array.
823 for (is = 0; is < na; is += PACK_X8)
825 for (k = 0; k < PACK_X8; k++)
828 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
829 ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
834 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
835 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t *nbat,
837 const nbnxn_search_t nbs,
841 const nbnxn_grid_t *grid;
843 for (g = 0; g < ngrid; g++)
845 grid = &nbs->grid[g];
847 /* Loop over all columns and copy and fill */
848 for (i = 0; i < grid->ncx*grid->ncy; i++)
850 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
851 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
853 copy_int_to_nbat_int(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
854 type, nbat->ntype-1, nbat->type+ash);
856 if (nbat->comb_rule != ljcrNONE)
858 if (nbat->XFormat == nbatX4)
860 copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
861 nbat->type+ash, ncz*grid->na_sc,
862 nbat->lj_comb+ash*2);
864 else if (nbat->XFormat == nbatX8)
866 copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
867 nbat->type+ash, ncz*grid->na_sc,
868 nbat->lj_comb+ash*2);
875 /* Sets the charges in nbnxn_atomdata_t *nbat */
876 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t *nbat,
878 const nbnxn_search_t nbs,
881 int g, cxy, ncz, ash, na, na_round, i, j;
883 const nbnxn_grid_t *grid;
885 for (g = 0; g < ngrid; g++)
887 grid = &nbs->grid[g];
889 /* Loop over all columns and copy and fill */
890 for (cxy = 0; cxy < grid->ncx*grid->ncy; cxy++)
892 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
893 na = grid->cxy_na[cxy];
894 na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
896 if (nbat->XFormat == nbatXYZQ)
898 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
899 for (i = 0; i < na; i++)
901 *q = charge[nbs->a[ash+i]];
904 /* Complete the partially filled last cell with zeros */
905 for (; i < na_round; i++)
914 for (i = 0; i < na; i++)
916 *q = charge[nbs->a[ash+i]];
919 /* Complete the partially filled last cell with zeros */
920 for (; i < na_round; i++)
930 /* Copies the energy group indices to a reordered and packed array */
931 static void copy_egp_to_nbat_egps(const int *a, int na, int na_round,
932 int na_c, int bit_shift,
933 const int *in, int *innb)
939 for (i = 0; i < na; i += na_c)
941 /* Store na_c energy group numbers into one int */
943 for (sa = 0; sa < na_c; sa++)
948 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
953 /* Complete the partially filled last cell with fill */
954 for (; i < na_round; i += na_c)
960 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
961 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t *nbat,
963 const nbnxn_search_t nbs,
967 const nbnxn_grid_t *grid;
969 for (g = 0; g < ngrid; g++)
971 grid = &nbs->grid[g];
973 /* Loop over all columns and copy and fill */
974 for (i = 0; i < grid->ncx*grid->ncy; i++)
976 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
977 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
979 copy_egp_to_nbat_egps(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
980 nbat->na_c, nbat->neg_2log,
981 atinfo, nbat->energrp+(ash>>grid->na_c_2log));
986 /* Sets all required atom parameter data in nbnxn_atomdata_t */
987 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
989 const nbnxn_search_t nbs,
990 const t_mdatoms *mdatoms,
995 if (locality == eatLocal)
1004 nbnxn_atomdata_set_atomtypes(nbat, ngrid, nbs, mdatoms->typeA);
1006 nbnxn_atomdata_set_charges(nbat, ngrid, nbs, mdatoms->chargeA);
1008 if (nbat->nenergrp > 1)
1010 nbnxn_atomdata_set_energygroups(nbat, ngrid, nbs, atinfo);
1014 /* Copies the shift vector array to nbnxn_atomdata_t */
1015 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox,
1017 nbnxn_atomdata_t *nbat)
1021 nbat->bDynamicBox = bDynamicBox;
1022 for (i = 0; i < SHIFTS; i++)
1024 copy_rvec(shift_vec[i], nbat->shift_vec[i]);
1028 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1029 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
1033 nbnxn_atomdata_t *nbat)
1056 nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
1059 nth = gmx_omp_nthreads_get(emntPairsearch);
1061 #pragma omp parallel for num_threads(nth) schedule(static)
1062 for (th = 0; th < nth; th++)
1066 for (g = g0; g < g1; g++)
1068 const nbnxn_grid_t *grid;
1069 int cxy0, cxy1, cxy;
1071 grid = &nbs->grid[g];
1073 cxy0 = (grid->ncx*grid->ncy* th +nth-1)/nth;
1074 cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
1076 for (cxy = cxy0; cxy < cxy1; cxy++)
1078 int na, ash, na_fill;
1080 na = grid->cxy_na[cxy];
1081 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1083 if (g == 0 && FillLocal)
1086 (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
1090 /* We fill only the real particle locations.
1091 * We assume the filling entries at the end have been
1092 * properly set before during ns.
1096 copy_rvec_to_nbat_real(nbs->a+ash, na, na_fill, x,
1097 nbat->XFormat, nbat->x, ash,
1105 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
1110 for (i = i0; i < i1; i++)
1117 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1119 real ** gmx_restrict src,
1127 /* The destination buffer contains data, add to it */
1128 for (i = i0; i < i1; i++)
1130 for (s = 0; s < nsrc; s++)
1132 dest[i] += src[s][i];
1138 /* The destination buffer is unitialized, set it first */
1139 for (i = i0; i < i1; i++)
1141 dest[i] = src[0][i];
1142 for (s = 1; s < nsrc; s++)
1144 dest[i] += src[s][i];
1151 nbnxn_atomdata_reduce_reals_simd(real gmx_unused * gmx_restrict dest,
1152 gmx_bool gmx_unused bDestSet,
1153 real gmx_unused ** gmx_restrict src,
1154 int gmx_unused nsrc,
1155 int gmx_unused i0, int gmx_unused i1)
1157 #ifdef GMX_NBNXN_SIMD
1158 /* The SIMD width here is actually independent of that in the kernels,
1159 * but we use the same width for simplicity (usually optimal anyhow).
1162 gmx_mm_pr dest_SSE, src_SSE;
1166 for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1168 dest_SSE = gmx_load_pr(dest+i);
1169 for (s = 0; s < nsrc; s++)
1171 src_SSE = gmx_load_pr(src[s]+i);
1172 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1174 gmx_store_pr(dest+i, dest_SSE);
1179 for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1181 dest_SSE = gmx_load_pr(src[0]+i);
1182 for (s = 1; s < nsrc; s++)
1184 src_SSE = gmx_load_pr(src[s]+i);
1185 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1187 gmx_store_pr(dest+i, dest_SSE);
1193 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1195 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1196 const nbnxn_atomdata_t *nbat,
1197 nbnxn_atomdata_output_t *out,
1208 /* Loop over all columns and copy and fill */
1209 switch (nbat->FFormat)
1217 for (a = a0; a < a1; a++)
1219 i = cell[a]*nbat->fstride;
1222 f[a][YY] += fnb[i+1];
1223 f[a][ZZ] += fnb[i+2];
1228 for (a = a0; a < a1; a++)
1230 i = cell[a]*nbat->fstride;
1232 for (fa = 0; fa < nfa; fa++)
1234 f[a][XX] += out[fa].f[i];
1235 f[a][YY] += out[fa].f[i+1];
1236 f[a][ZZ] += out[fa].f[i+2];
1246 for (a = a0; a < a1; a++)
1248 i = X4_IND_A(cell[a]);
1250 f[a][XX] += fnb[i+XX*PACK_X4];
1251 f[a][YY] += fnb[i+YY*PACK_X4];
1252 f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1257 for (a = a0; a < a1; a++)
1259 i = X4_IND_A(cell[a]);
1261 for (fa = 0; fa < nfa; fa++)
1263 f[a][XX] += out[fa].f[i+XX*PACK_X4];
1264 f[a][YY] += out[fa].f[i+YY*PACK_X4];
1265 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1275 for (a = a0; a < a1; a++)
1277 i = X8_IND_A(cell[a]);
1279 f[a][XX] += fnb[i+XX*PACK_X8];
1280 f[a][YY] += fnb[i+YY*PACK_X8];
1281 f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1286 for (a = a0; a < a1; a++)
1288 i = X8_IND_A(cell[a]);
1290 for (fa = 0; fa < nfa; fa++)
1292 f[a][XX] += out[fa].f[i+XX*PACK_X8];
1293 f[a][YY] += out[fa].f[i+YY*PACK_X8];
1294 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1300 gmx_incons("Unsupported nbnxn_atomdata_t format");
1304 static gmx_inline unsigned char reverse_bits(unsigned char b)
1306 /* http://graphics.stanford.edu/~seander/bithacks.html#ReverseByteWith64BitsDiv */
1307 return (b * 0x0202020202ULL & 0x010884422010ULL) % 1023;
1310 static void nbnxn_atomdata_add_nbat_f_to_f_treereduce(const nbnxn_atomdata_t *nbat,
1313 const nbnxn_buffer_flags_t *flags = &nbat->buffer_flags;
1315 int next_pow2 = 1<<(gmx_log2i(nth-1)+1);
1317 assert(nbat->nout == nth); /* tree-reduce currently only works for nout==nth */
1319 memset(nbat->syncStep, 0, sizeof(*(nbat->syncStep))*nth);
1321 #pragma omp parallel num_threads(nth)
1327 th = gmx_omp_get_thread_num();
1329 for (group_size = 2; group_size < 2*next_pow2; group_size *= 2)
1331 int index[2], group_pos, partner_pos, wu;
1332 int partner_th = th ^ (group_size/2);
1337 /* wait on partner thread - replaces full barrier */
1338 int sync_th, sync_group_size;
1340 tMPI_Atomic_memory_barrier(); /* gurantee data is saved before marking work as done */
1341 tMPI_Atomic_set(&(nbat->syncStep[th]), group_size/2); /* mark previous step as completed */
1343 /* find thread to sync with. Equal to partner_th unless nth is not a power of two. */
1344 for (sync_th = partner_th, sync_group_size = group_size; sync_th >= nth && sync_group_size > 2; sync_group_size /= 2)
1346 sync_th &= ~(sync_group_size/4);
1348 if (sync_th < nth) /* otherwise nothing to sync index[1] will be >=nout */
1350 /* wait on the thread which computed input data in previous step */
1351 while (tMPI_Atomic_get((volatile tMPI_Atomic_t*)&(nbat->syncStep[sync_th])) < group_size/2)
1355 /* guarantee that no later load happens before wait loop is finisehd */
1356 tMPI_Atomic_memory_barrier();
1358 #else /* TMPI_ATOMICS */
1363 /* Calculate buffers to sum (result goes into first buffer) */
1364 group_pos = th % group_size;
1365 index[0] = th - group_pos;
1366 index[1] = index[0] + group_size/2;
1368 /* If no second buffer, nothing to do */
1369 if (index[1] >= nbat->nout && group_size > 2)
1374 #if NBNXN_BUFFERFLAG_MAX_THREADS > 256
1375 #error reverse_bits assumes max 256 threads
1377 /* Position is permuted so that one of the 2 vectors being added was computed on the same thread in the previous step.
1378 This improves locality and enables to sync with just a single thread between steps (=the levels in the btree).
1379 The permutation which allows this corresponds to reversing the bits of the group position.
1381 group_pos = reverse_bits(group_pos)/(256/group_size);
1383 partner_pos = group_pos ^ 1;
1385 /* loop over two work-units (own and partner) */
1386 for (wu = 0; wu < 2; wu++)
1390 if (partner_th < nth)
1392 break; /* partner exists we don't have to do his work */
1396 group_pos = partner_pos;
1400 /* Calculate the cell-block range for our thread */
1401 b0 = (flags->nflag* group_pos )/group_size;
1402 b1 = (flags->nflag*(group_pos+1))/group_size;
1404 for (b = b0; b < b1; b++)
1406 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1407 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1409 if ((flags->flag[b] & (1ULL<<index[1])) || group_size > 2)
1411 #ifdef GMX_NBNXN_SIMD
1412 nbnxn_atomdata_reduce_reals_simd
1414 nbnxn_atomdata_reduce_reals
1416 (nbat->out[index[0]].f,
1417 (flags->flag[b] & (1ULL<<index[0])) || group_size > 2,
1418 &(nbat->out[index[1]].f), 1, i0, i1);
1421 else if (!(flags->flag[b] & (1ULL<<index[0])))
1423 nbnxn_atomdata_clear_reals(nbat->out[index[0]].f,
1433 static void nbnxn_atomdata_add_nbat_f_to_f_stdreduce(const nbnxn_atomdata_t *nbat,
1437 #pragma omp parallel for num_threads(nth) schedule(static)
1438 for (th = 0; th < nth; th++)
1440 const nbnxn_buffer_flags_t *flags;
1444 real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1447 flags = &nbat->buffer_flags;
1449 /* Calculate the cell-block range for our thread */
1450 b0 = (flags->nflag* th )/nth;
1451 b1 = (flags->nflag*(th+1))/nth;
1453 for (b = b0; b < b1; b++)
1455 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1456 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1459 for (out = 1; out < nbat->nout; out++)
1461 if (flags->flag[b] & (1U<<out))
1463 fptr[nfptr++] = nbat->out[out].f;
1468 #ifdef GMX_NBNXN_SIMD
1469 nbnxn_atomdata_reduce_reals_simd
1471 nbnxn_atomdata_reduce_reals
1474 flags->flag[b] & (1U<<0),
1478 else if (!(flags->flag[b] & (1U<<0)))
1480 nbnxn_atomdata_clear_reals(nbat->out[0].f,
1487 /* Add the force array(s) from nbnxn_atomdata_t to f */
1488 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t nbs,
1490 const nbnxn_atomdata_t *nbat,
1496 nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1502 na = nbs->natoms_nonlocal;
1506 na = nbs->natoms_local;
1509 a0 = nbs->natoms_local;
1510 na = nbs->natoms_nonlocal - nbs->natoms_local;
1514 nth = gmx_omp_nthreads_get(emntNonbonded);
1518 if (locality != eatAll)
1520 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1523 /* Reduce the force thread output buffers into buffer 0, before adding
1524 * them to the, differently ordered, "real" force buffer.
1526 if (nbat->bUseTreeReduce)
1528 nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbat, nth);
1532 nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbat, nth);
1535 #pragma omp parallel for num_threads(nth) schedule(static)
1536 for (th = 0; th < nth; th++)
1538 nbnxn_atomdata_add_nbat_f_to_f_part(nbs, nbat,
1546 nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1549 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1550 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1553 const nbnxn_atomdata_output_t *out;
1560 for (s = 0; s < SHIFTS; s++)
1563 for (th = 0; th < nbat->nout; th++)
1565 sum[XX] += out[th].fshift[s*DIM+XX];
1566 sum[YY] += out[th].fshift[s*DIM+YY];
1567 sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1569 rvec_inc(fshift[s], sum);