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 "gromacs/math/vec.h"
47 #include "nbnxn_consts.h"
48 #include "nbnxn_internal.h"
49 #include "nbnxn_atomdata.h"
50 #include "nbnxn_search.h"
51 #include "gmx_omp_nthreads.h"
52 #include "thread_mpi/atomic.h"
54 #include "gromacs/pbcutil/ishift.h"
55 #include "gromacs/utility/gmxomp.h"
56 #include "gromacs/utility/smalloc.h"
58 /* Default nbnxn allocation routine, allocates NBNXN_MEM_ALIGN byte aligned */
59 void nbnxn_alloc_aligned(void **ptr, size_t nbytes)
61 *ptr = save_malloc_aligned("ptr", __FILE__, __LINE__, nbytes, 1, NBNXN_MEM_ALIGN);
64 /* Free function for memory allocated with nbnxn_alloc_aligned */
65 void nbnxn_free_aligned(void *ptr)
70 /* Reallocation wrapper function for nbnxn data structures */
71 void nbnxn_realloc_void(void **ptr,
72 int nbytes_copy, int nbytes_new,
78 ma(&ptr_new, nbytes_new);
80 if (nbytes_new > 0 && ptr_new == NULL)
82 gmx_fatal(FARGS, "Allocation of %d bytes failed", nbytes_new);
87 if (nbytes_new < nbytes_copy)
89 gmx_incons("In nbnxn_realloc_void: new size less than copy size");
91 memcpy(ptr_new, *ptr, nbytes_copy);
100 /* Reallocate the nbnxn_atomdata_t for a size of n atoms */
101 void nbnxn_atomdata_realloc(nbnxn_atomdata_t *nbat, int n)
105 nbnxn_realloc_void((void **)&nbat->type,
106 nbat->natoms*sizeof(*nbat->type),
107 n*sizeof(*nbat->type),
108 nbat->alloc, nbat->free);
109 nbnxn_realloc_void((void **)&nbat->lj_comb,
110 nbat->natoms*2*sizeof(*nbat->lj_comb),
111 n*2*sizeof(*nbat->lj_comb),
112 nbat->alloc, nbat->free);
113 if (nbat->XFormat != nbatXYZQ)
115 nbnxn_realloc_void((void **)&nbat->q,
116 nbat->natoms*sizeof(*nbat->q),
118 nbat->alloc, nbat->free);
120 if (nbat->nenergrp > 1)
122 nbnxn_realloc_void((void **)&nbat->energrp,
123 nbat->natoms/nbat->na_c*sizeof(*nbat->energrp),
124 n/nbat->na_c*sizeof(*nbat->energrp),
125 nbat->alloc, nbat->free);
127 nbnxn_realloc_void((void **)&nbat->x,
128 nbat->natoms*nbat->xstride*sizeof(*nbat->x),
129 n*nbat->xstride*sizeof(*nbat->x),
130 nbat->alloc, nbat->free);
131 for (t = 0; t < nbat->nout; t++)
133 /* Allocate one element extra for possible signaling with CUDA */
134 nbnxn_realloc_void((void **)&nbat->out[t].f,
135 nbat->natoms*nbat->fstride*sizeof(*nbat->out[t].f),
136 n*nbat->fstride*sizeof(*nbat->out[t].f),
137 nbat->alloc, nbat->free);
142 /* Initializes an nbnxn_atomdata_output_t data structure */
143 static void nbnxn_atomdata_output_init(nbnxn_atomdata_output_t *out,
145 int nenergrp, int stride,
151 ma((void **)&out->fshift, SHIFTS*DIM*sizeof(*out->fshift));
152 out->nV = nenergrp*nenergrp;
153 ma((void **)&out->Vvdw, out->nV*sizeof(*out->Vvdw));
154 ma((void **)&out->Vc, out->nV*sizeof(*out->Vc ));
156 if (nb_kernel_type == nbnxnk4xN_SIMD_4xN ||
157 nb_kernel_type == nbnxnk4xN_SIMD_2xNN)
159 cj_size = nbnxn_kernel_to_cj_size(nb_kernel_type);
160 out->nVS = nenergrp*nenergrp*stride*(cj_size>>1)*cj_size;
161 ma((void **)&out->VSvdw, out->nVS*sizeof(*out->VSvdw));
162 ma((void **)&out->VSc, out->nVS*sizeof(*out->VSc ));
170 static void copy_int_to_nbat_int(const int *a, int na, int na_round,
171 const int *in, int fill, int *innb)
176 for (i = 0; i < na; i++)
178 innb[j++] = in[a[i]];
180 /* Complete the partially filled last cell with fill */
181 for (; i < na_round; i++)
187 static void clear_nbat_real(int na, int nbatFormat, real *xnb, int a0)
194 for (a = 0; a < na; a++)
196 for (d = 0; d < DIM; d++)
198 xnb[(a0+a)*STRIDE_XYZ+d] = 0;
203 for (a = 0; a < na; a++)
205 for (d = 0; d < DIM; d++)
207 xnb[(a0+a)*STRIDE_XYZQ+d] = 0;
213 c = a0 & (PACK_X4-1);
214 for (a = 0; a < na; a++)
216 xnb[j+XX*PACK_X4] = 0;
217 xnb[j+YY*PACK_X4] = 0;
218 xnb[j+ZZ*PACK_X4] = 0;
223 j += (DIM-1)*PACK_X4;
230 c = a0 & (PACK_X8-1);
231 for (a = 0; a < na; a++)
233 xnb[j+XX*PACK_X8] = 0;
234 xnb[j+YY*PACK_X8] = 0;
235 xnb[j+ZZ*PACK_X8] = 0;
240 j += (DIM-1)*PACK_X8;
248 void copy_rvec_to_nbat_real(const int *a, int na, int na_round,
249 rvec *x, int nbatFormat, real *xnb, int a0,
250 int cx, int cy, int cz)
254 /* We might need to place filler particles to fill up the cell to na_round.
255 * The coefficients (LJ and q) for such particles are zero.
256 * But we might still get NaN as 0*NaN when distances are too small.
257 * We hope that -107 nm is far away enough from to zero
258 * to avoid accidental short distances to particles shifted down for pbc.
260 #define NBAT_FAR_AWAY 107
266 for (i = 0; i < na; i++)
268 xnb[j++] = x[a[i]][XX];
269 xnb[j++] = x[a[i]][YY];
270 xnb[j++] = x[a[i]][ZZ];
272 /* Complete the partially filled last cell with copies of the last element.
273 * This simplifies the bounding box calculation and avoid
274 * numerical issues with atoms that are coincidentally close.
276 for (; i < na_round; i++)
278 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
279 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
280 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
285 for (i = 0; i < na; i++)
287 xnb[j++] = x[a[i]][XX];
288 xnb[j++] = x[a[i]][YY];
289 xnb[j++] = x[a[i]][ZZ];
292 /* Complete the partially filled last cell with particles far apart */
293 for (; i < na_round; i++)
295 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
296 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
297 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
303 c = a0 & (PACK_X4-1);
304 for (i = 0; i < na; i++)
306 xnb[j+XX*PACK_X4] = x[a[i]][XX];
307 xnb[j+YY*PACK_X4] = x[a[i]][YY];
308 xnb[j+ZZ*PACK_X4] = x[a[i]][ZZ];
313 j += (DIM-1)*PACK_X4;
317 /* Complete the partially filled last cell with particles far apart */
318 for (; i < na_round; i++)
320 xnb[j+XX*PACK_X4] = -NBAT_FAR_AWAY*(1 + cx);
321 xnb[j+YY*PACK_X4] = -NBAT_FAR_AWAY*(1 + cy);
322 xnb[j+ZZ*PACK_X4] = -NBAT_FAR_AWAY*(1 + cz + i);
327 j += (DIM-1)*PACK_X4;
334 c = a0 & (PACK_X8 - 1);
335 for (i = 0; i < na; i++)
337 xnb[j+XX*PACK_X8] = x[a[i]][XX];
338 xnb[j+YY*PACK_X8] = x[a[i]][YY];
339 xnb[j+ZZ*PACK_X8] = x[a[i]][ZZ];
344 j += (DIM-1)*PACK_X8;
348 /* Complete the partially filled last cell with particles far apart */
349 for (; i < na_round; i++)
351 xnb[j+XX*PACK_X8] = -NBAT_FAR_AWAY*(1 + cx);
352 xnb[j+YY*PACK_X8] = -NBAT_FAR_AWAY*(1 + cy);
353 xnb[j+ZZ*PACK_X8] = -NBAT_FAR_AWAY*(1 + cz + i);
358 j += (DIM-1)*PACK_X8;
364 gmx_incons("Unsupported nbnxn_atomdata_t format");
368 /* Stores the LJ parameter data in a format convenient for different kernels */
369 static void set_lj_parameter_data(nbnxn_atomdata_t *nbat, gmx_bool bSIMD)
378 /* nbfp_s4 stores two parameters using a stride of 4,
379 * because this would suit x86 SIMD single-precision
380 * quad-load intrinsics. There's a slight inefficiency in
381 * allocating and initializing nbfp_s4 when it might not
382 * be used, but introducing the conditional code is not
383 * really worth it. */
384 nbat->alloc((void **)&nbat->nbfp_s4, nt*nt*4*sizeof(*nbat->nbfp_s4));
385 for (i = 0; i < nt; i++)
387 for (j = 0; j < nt; j++)
389 nbat->nbfp_s4[(i*nt+j)*4+0] = nbat->nbfp[(i*nt+j)*2+0];
390 nbat->nbfp_s4[(i*nt+j)*4+1] = nbat->nbfp[(i*nt+j)*2+1];
391 nbat->nbfp_s4[(i*nt+j)*4+2] = 0;
392 nbat->nbfp_s4[(i*nt+j)*4+3] = 0;
397 /* We use combination rule data for SIMD combination rule kernels
398 * and with LJ-PME kernels. We then only need parameters per atom type,
399 * not per pair of atom types.
401 switch (nbat->comb_rule)
404 nbat->comb_rule = ljcrGEOM;
406 for (i = 0; i < nt; i++)
408 /* Store the sqrt of the diagonal from the nbfp matrix */
409 nbat->nbfp_comb[i*2 ] = sqrt(nbat->nbfp[(i*nt+i)*2 ]);
410 nbat->nbfp_comb[i*2+1] = sqrt(nbat->nbfp[(i*nt+i)*2+1]);
414 for (i = 0; i < nt; i++)
416 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
417 c6 = nbat->nbfp[(i*nt+i)*2 ];
418 c12 = nbat->nbfp[(i*nt+i)*2+1];
419 if (c6 > 0 && c12 > 0)
421 /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
422 * so we get 6*C6 and 12*C12 after combining.
424 nbat->nbfp_comb[i*2 ] = 0.5*pow(c12/c6, 1.0/6.0);
425 nbat->nbfp_comb[i*2+1] = sqrt(c6*c6/c12);
429 nbat->nbfp_comb[i*2 ] = 0;
430 nbat->nbfp_comb[i*2+1] = 0;
435 /* We always store the full matrix (see code above) */
438 gmx_incons("Unknown combination rule");
443 #ifdef GMX_NBNXN_SIMD
445 nbnxn_atomdata_init_simple_exclusion_masks(nbnxn_atomdata_t *nbat)
448 const int simd_width = GMX_SIMD_REAL_WIDTH;
450 /* Set the diagonal cluster pair exclusion mask setup data.
451 * In the kernel we check 0 < j - i to generate the masks.
452 * Here we store j - i for generating the mask for the first i,
453 * we substract 0.5 to avoid rounding issues.
454 * In the kernel we can subtract 1 to generate the subsequent mask.
456 int simd_4xn_diag_size;
457 const real simdFalse = -1, simdTrue = 1;
458 real *simd_interaction_array;
460 simd_4xn_diag_size = max(NBNXN_CPU_CLUSTER_I_SIZE, simd_width);
461 snew_aligned(nbat->simd_4xn_diagonal_j_minus_i, simd_4xn_diag_size, NBNXN_MEM_ALIGN);
462 for (j = 0; j < simd_4xn_diag_size; j++)
464 nbat->simd_4xn_diagonal_j_minus_i[j] = j - 0.5;
467 snew_aligned(nbat->simd_2xnn_diagonal_j_minus_i, simd_width, NBNXN_MEM_ALIGN);
468 for (j = 0; j < simd_width/2; j++)
470 /* The j-cluster size is half the SIMD width */
471 nbat->simd_2xnn_diagonal_j_minus_i[j] = j - 0.5;
472 /* The next half of the SIMD width is for i + 1 */
473 nbat->simd_2xnn_diagonal_j_minus_i[simd_width/2+j] = j - 1 - 0.5;
476 /* We use up to 32 bits for exclusion masking.
477 * The same masks are used for the 4xN and 2x(N+N) kernels.
478 * The masks are read either into epi32 SIMD registers or into
479 * real SIMD registers (together with a cast).
480 * In single precision this means the real and epi32 SIMD registers
482 * In double precision the epi32 registers can be smaller than
483 * the real registers, so depending on the architecture, we might
484 * need to use two, identical, 32-bit masks per real.
486 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*simd_width;
487 snew_aligned(nbat->simd_exclusion_filter1, simd_excl_size, NBNXN_MEM_ALIGN);
488 snew_aligned(nbat->simd_exclusion_filter2, simd_excl_size*2, NBNXN_MEM_ALIGN);
490 for (j = 0; j < simd_excl_size; j++)
492 /* Set the consecutive bits for masking pair exclusions */
493 nbat->simd_exclusion_filter1[j] = (1U << j);
494 nbat->simd_exclusion_filter2[j*2 + 0] = (1U << j);
495 nbat->simd_exclusion_filter2[j*2 + 1] = (1U << j);
498 #if (defined GMX_SIMD_IBM_QPX)
499 /* The QPX kernels shouldn't do the bit masking that is done on
500 * x86, because the SIMD units lack bit-wise operations. Instead,
501 * we generate a vector of all 2^4 possible ways an i atom
502 * interacts with its 4 j atoms. Each array entry contains
503 * simd_width signed ints that are read in a single SIMD
504 * load. These ints must contain values that will be interpreted
505 * as true and false when loaded in the SIMD floating-point
506 * registers, ie. any positive or any negative value,
507 * respectively. Each array entry encodes how this i atom will
508 * interact with the 4 j atoms. Matching code exists in
509 * set_ci_top_excls() to generate indices into this array. Those
510 * indices are used in the kernels. */
512 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
513 const int qpx_simd_width = GMX_SIMD_REAL_WIDTH;
514 snew_aligned(simd_interaction_array, simd_excl_size * qpx_simd_width, NBNXN_MEM_ALIGN);
515 for (j = 0; j < simd_excl_size; j++)
517 int index = j * qpx_simd_width;
518 for (i = 0; i < qpx_simd_width; i++)
520 simd_interaction_array[index + i] = (j & (1 << i)) ? simdTrue : simdFalse;
523 nbat->simd_interaction_array = simd_interaction_array;
528 /* Initializes an nbnxn_atomdata_t data structure */
529 void nbnxn_atomdata_init(FILE *fp,
530 nbnxn_atomdata_t *nbat,
532 int enbnxninitcombrule,
533 int ntype, const real *nbfp,
536 nbnxn_alloc_t *alloc,
542 gmx_bool simple, bCombGeom, bCombLB, bSIMD;
546 nbat->alloc = nbnxn_alloc_aligned;
554 nbat->free = nbnxn_free_aligned;
563 fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
565 nbat->ntype = ntype + 1;
566 nbat->alloc((void **)&nbat->nbfp,
567 nbat->ntype*nbat->ntype*2*sizeof(*nbat->nbfp));
568 nbat->alloc((void **)&nbat->nbfp_comb, nbat->ntype*2*sizeof(*nbat->nbfp_comb));
570 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
571 * force-field floating point parameters.
574 ptr = getenv("GMX_LJCOMB_TOL");
579 sscanf(ptr, "%lf", &dbl);
585 /* Temporarily fill nbat->nbfp_comb with sigma and epsilon
586 * to check for the LB rule.
588 for (i = 0; i < ntype; i++)
590 c6 = nbfp[(i*ntype+i)*2 ]/6.0;
591 c12 = nbfp[(i*ntype+i)*2+1]/12.0;
592 if (c6 > 0 && c12 > 0)
594 nbat->nbfp_comb[i*2 ] = pow(c12/c6, 1.0/6.0);
595 nbat->nbfp_comb[i*2+1] = 0.25*c6*c6/c12;
597 else if (c6 == 0 && c12 == 0)
599 nbat->nbfp_comb[i*2 ] = 0;
600 nbat->nbfp_comb[i*2+1] = 0;
604 /* Can not use LB rule with only dispersion or repulsion */
609 for (i = 0; i < nbat->ntype; i++)
611 for (j = 0; j < nbat->ntype; j++)
613 if (i < ntype && j < ntype)
615 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
616 * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
618 c6 = nbfp[(i*ntype+j)*2 ];
619 c12 = nbfp[(i*ntype+j)*2+1];
620 nbat->nbfp[(i*nbat->ntype+j)*2 ] = c6;
621 nbat->nbfp[(i*nbat->ntype+j)*2+1] = c12;
623 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
624 bCombGeom = bCombGeom &&
625 gmx_within_tol(c6*c6, nbfp[(i*ntype+i)*2 ]*nbfp[(j*ntype+j)*2 ], tol) &&
626 gmx_within_tol(c12*c12, nbfp[(i*ntype+i)*2+1]*nbfp[(j*ntype+j)*2+1], tol);
628 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
632 ((c6 == 0 && c12 == 0 &&
633 (nbat->nbfp_comb[i*2+1] == 0 || nbat->nbfp_comb[j*2+1] == 0)) ||
634 (c6 > 0 && c12 > 0 &&
635 gmx_within_tol(pow(c12/c6, 1.0/6.0), 0.5*(nbat->nbfp_comb[i*2]+nbat->nbfp_comb[j*2]), tol) &&
636 gmx_within_tol(0.25*c6*c6/c12, sqrt(nbat->nbfp_comb[i*2+1]*nbat->nbfp_comb[j*2+1]), tol)));
640 /* Add zero parameters for the additional dummy atom type */
641 nbat->nbfp[(i*nbat->ntype+j)*2 ] = 0;
642 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 0;
648 fprintf(debug, "Combination rules: geometric %d Lorentz-Berthelot %d\n",
652 simple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
654 switch (enbnxninitcombrule)
656 case enbnxninitcombruleDETECT:
657 /* We prefer the geometic combination rule,
658 * as that gives a slightly faster kernel than the LB rule.
662 nbat->comb_rule = ljcrGEOM;
666 nbat->comb_rule = ljcrLB;
670 nbat->comb_rule = ljcrNONE;
672 nbat->free(nbat->nbfp_comb);
677 if (nbat->comb_rule == ljcrNONE)
679 fprintf(fp, "Using full Lennard-Jones parameter combination matrix\n\n");
683 fprintf(fp, "Using %s Lennard-Jones combination rule\n\n",
684 nbat->comb_rule == ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
688 case enbnxninitcombruleGEOM:
689 nbat->comb_rule = ljcrGEOM;
691 case enbnxninitcombruleLB:
692 nbat->comb_rule = ljcrLB;
694 case enbnxninitcombruleNONE:
695 nbat->comb_rule = ljcrNONE;
697 nbat->free(nbat->nbfp_comb);
700 gmx_incons("Unknown enbnxninitcombrule");
703 bSIMD = (nb_kernel_type == nbnxnk4xN_SIMD_4xN ||
704 nb_kernel_type == nbnxnk4xN_SIMD_2xNN);
706 set_lj_parameter_data(nbat, bSIMD);
710 nbat->lj_comb = NULL;
717 pack_x = max(NBNXN_CPU_CLUSTER_I_SIZE,
718 nbnxn_kernel_to_cj_size(nb_kernel_type));
722 nbat->XFormat = nbatX4;
725 nbat->XFormat = nbatX8;
728 gmx_incons("Unsupported packing width");
733 nbat->XFormat = nbatXYZ;
736 nbat->FFormat = nbat->XFormat;
740 nbat->XFormat = nbatXYZQ;
741 nbat->FFormat = nbatXYZ;
744 nbat->nenergrp = n_energygroups;
747 /* Energy groups not supported yet for super-sub lists */
748 if (n_energygroups > 1 && fp != NULL)
750 fprintf(fp, "\nNOTE: With GPUs, reporting energy group contributions is not supported\n\n");
754 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
755 if (nbat->nenergrp > 64)
757 gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
760 while (nbat->nenergrp > (1<<nbat->neg_2log))
764 nbat->energrp = NULL;
765 nbat->alloc((void **)&nbat->shift_vec, SHIFTS*sizeof(*nbat->shift_vec));
766 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
767 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
770 #ifdef GMX_NBNXN_SIMD
773 nbnxn_atomdata_init_simple_exclusion_masks(nbat);
777 /* Initialize the output data structures */
779 snew(nbat->out, nbat->nout);
781 for (i = 0; i < nbat->nout; i++)
783 nbnxn_atomdata_output_init(&nbat->out[i],
785 nbat->nenergrp, 1<<nbat->neg_2log,
788 nbat->buffer_flags.flag = NULL;
789 nbat->buffer_flags.flag_nalloc = 0;
791 nth = gmx_omp_nthreads_get(emntNonbonded);
793 ptr = getenv("GMX_USE_TREEREDUCE");
796 nbat->bUseTreeReduce = strtol(ptr, 0, 10);
799 else if (nth > 8) /*on the CPU we currently don't benefit even at 32*/
801 nbat->bUseTreeReduce = 1;
806 nbat->bUseTreeReduce = 0;
808 if (nbat->bUseTreeReduce)
812 fprintf(fp, "Using tree force reduction\n\n");
814 snew(nbat->syncStep, nth);
818 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
819 const int *type, int na,
824 /* The LJ params follow the combination rule:
825 * copy the params for the type array to the atom array.
827 for (is = 0; is < na; is += PACK_X4)
829 for (k = 0; k < PACK_X4; k++)
832 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
833 ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
838 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
839 const int *type, int na,
844 /* The LJ params follow the combination rule:
845 * copy the params for the type array to the atom array.
847 for (is = 0; is < na; is += PACK_X8)
849 for (k = 0; k < PACK_X8; k++)
852 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
853 ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
858 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
859 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t *nbat,
861 const nbnxn_search_t nbs,
865 const nbnxn_grid_t *grid;
867 for (g = 0; g < ngrid; g++)
869 grid = &nbs->grid[g];
871 /* Loop over all columns and copy and fill */
872 for (i = 0; i < grid->ncx*grid->ncy; i++)
874 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
875 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
877 copy_int_to_nbat_int(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
878 type, nbat->ntype-1, nbat->type+ash);
880 if (nbat->comb_rule != ljcrNONE)
882 if (nbat->XFormat == nbatX4)
884 copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
885 nbat->type+ash, ncz*grid->na_sc,
886 nbat->lj_comb+ash*2);
888 else if (nbat->XFormat == nbatX8)
890 copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
891 nbat->type+ash, ncz*grid->na_sc,
892 nbat->lj_comb+ash*2);
899 /* Sets the charges in nbnxn_atomdata_t *nbat */
900 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t *nbat,
902 const nbnxn_search_t nbs,
905 int g, cxy, ncz, ash, na, na_round, i, j;
907 const nbnxn_grid_t *grid;
909 for (g = 0; g < ngrid; g++)
911 grid = &nbs->grid[g];
913 /* Loop over all columns and copy and fill */
914 for (cxy = 0; cxy < grid->ncx*grid->ncy; cxy++)
916 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
917 na = grid->cxy_na[cxy];
918 na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
920 if (nbat->XFormat == nbatXYZQ)
922 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
923 for (i = 0; i < na; i++)
925 *q = charge[nbs->a[ash+i]];
928 /* Complete the partially filled last cell with zeros */
929 for (; i < na_round; i++)
938 for (i = 0; i < na; i++)
940 *q = charge[nbs->a[ash+i]];
943 /* Complete the partially filled last cell with zeros */
944 for (; i < na_round; i++)
954 /* Set the charges of perturbed atoms in nbnxn_atomdata_t to 0.
955 * This is to automatically remove the RF/PME self term in the nbnxn kernels.
956 * Part of the zero interactions are still calculated in the normal kernels.
957 * All perturbed interactions are calculated in the free energy kernel,
958 * using the original charge and LJ data, not nbnxn_atomdata_t.
960 static void nbnxn_atomdata_mask_fep(nbnxn_atomdata_t *nbat,
962 const nbnxn_search_t nbs)
965 int stride_q, g, nsubc, c_offset, c, subc, i, ind;
966 const nbnxn_grid_t *grid;
968 if (nbat->XFormat == nbatXYZQ)
970 q = nbat->x + ZZ + 1;
971 stride_q = STRIDE_XYZQ;
979 for (g = 0; g < ngrid; g++)
981 grid = &nbs->grid[g];
988 nsubc = GPU_NSUBCELL;
991 c_offset = grid->cell0*grid->na_sc;
993 /* Loop over all columns and copy and fill */
994 for (c = 0; c < grid->nc*nsubc; c++)
996 /* Does this cluster contain perturbed particles? */
997 if (grid->fep[c] != 0)
999 for (i = 0; i < grid->na_c; i++)
1001 /* Is this a perturbed particle? */
1002 if (grid->fep[c] & (1 << i))
1004 ind = c_offset + c*grid->na_c + i;
1005 /* Set atom type and charge to non-interacting */
1006 nbat->type[ind] = nbat->ntype - 1;
1007 q[ind*stride_q] = 0;
1015 /* Copies the energy group indices to a reordered and packed array */
1016 static void copy_egp_to_nbat_egps(const int *a, int na, int na_round,
1017 int na_c, int bit_shift,
1018 const int *in, int *innb)
1024 for (i = 0; i < na; i += na_c)
1026 /* Store na_c energy group numbers into one int */
1028 for (sa = 0; sa < na_c; sa++)
1033 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
1038 /* Complete the partially filled last cell with fill */
1039 for (; i < na_round; i += na_c)
1045 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
1046 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t *nbat,
1048 const nbnxn_search_t nbs,
1052 const nbnxn_grid_t *grid;
1054 if (nbat->nenergrp == 1)
1059 for (g = 0; g < ngrid; g++)
1061 grid = &nbs->grid[g];
1063 /* Loop over all columns and copy and fill */
1064 for (i = 0; i < grid->ncx*grid->ncy; i++)
1066 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
1067 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
1069 copy_egp_to_nbat_egps(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
1070 nbat->na_c, nbat->neg_2log,
1071 atinfo, nbat->energrp+(ash>>grid->na_c_2log));
1076 /* Sets all required atom parameter data in nbnxn_atomdata_t */
1077 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
1079 const nbnxn_search_t nbs,
1080 const t_mdatoms *mdatoms,
1085 if (locality == eatLocal)
1094 nbnxn_atomdata_set_atomtypes(nbat, ngrid, nbs, mdatoms->typeA);
1096 nbnxn_atomdata_set_charges(nbat, ngrid, nbs, mdatoms->chargeA);
1100 nbnxn_atomdata_mask_fep(nbat, ngrid, nbs);
1103 nbnxn_atomdata_set_energygroups(nbat, ngrid, nbs, atinfo);
1106 /* Copies the shift vector array to nbnxn_atomdata_t */
1107 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox,
1109 nbnxn_atomdata_t *nbat)
1113 nbat->bDynamicBox = bDynamicBox;
1114 for (i = 0; i < SHIFTS; i++)
1116 copy_rvec(shift_vec[i], nbat->shift_vec[i]);
1120 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1121 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
1125 nbnxn_atomdata_t *nbat)
1148 nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
1151 nth = gmx_omp_nthreads_get(emntPairsearch);
1153 #pragma omp parallel for num_threads(nth) schedule(static)
1154 for (th = 0; th < nth; th++)
1158 for (g = g0; g < g1; g++)
1160 const nbnxn_grid_t *grid;
1161 int cxy0, cxy1, cxy;
1163 grid = &nbs->grid[g];
1165 cxy0 = (grid->ncx*grid->ncy* th +nth-1)/nth;
1166 cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
1168 for (cxy = cxy0; cxy < cxy1; cxy++)
1170 int na, ash, na_fill;
1172 na = grid->cxy_na[cxy];
1173 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1175 if (g == 0 && FillLocal)
1178 (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
1182 /* We fill only the real particle locations.
1183 * We assume the filling entries at the end have been
1184 * properly set before during ns.
1188 copy_rvec_to_nbat_real(nbs->a+ash, na, na_fill, x,
1189 nbat->XFormat, nbat->x, ash,
1197 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
1202 for (i = i0; i < i1; i++)
1209 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1211 real ** gmx_restrict src,
1219 /* The destination buffer contains data, add to it */
1220 for (i = i0; i < i1; i++)
1222 for (s = 0; s < nsrc; s++)
1224 dest[i] += src[s][i];
1230 /* The destination buffer is unitialized, set it first */
1231 for (i = i0; i < i1; i++)
1233 dest[i] = src[0][i];
1234 for (s = 1; s < nsrc; s++)
1236 dest[i] += src[s][i];
1243 nbnxn_atomdata_reduce_reals_simd(real gmx_unused * gmx_restrict dest,
1244 gmx_bool gmx_unused bDestSet,
1245 real gmx_unused ** gmx_restrict src,
1246 int gmx_unused nsrc,
1247 int gmx_unused i0, int gmx_unused i1)
1249 #ifdef GMX_NBNXN_SIMD
1250 /* The SIMD width here is actually independent of that in the kernels,
1251 * but we use the same width for simplicity (usually optimal anyhow).
1254 gmx_simd_real_t dest_SSE, src_SSE;
1258 for (i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1260 dest_SSE = gmx_simd_load_r(dest+i);
1261 for (s = 0; s < nsrc; s++)
1263 src_SSE = gmx_simd_load_r(src[s]+i);
1264 dest_SSE = gmx_simd_add_r(dest_SSE, src_SSE);
1266 gmx_simd_store_r(dest+i, dest_SSE);
1271 for (i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1273 dest_SSE = gmx_simd_load_r(src[0]+i);
1274 for (s = 1; s < nsrc; s++)
1276 src_SSE = gmx_simd_load_r(src[s]+i);
1277 dest_SSE = gmx_simd_add_r(dest_SSE, src_SSE);
1279 gmx_simd_store_r(dest+i, dest_SSE);
1285 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1287 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1288 const nbnxn_atomdata_t *nbat,
1289 nbnxn_atomdata_output_t *out,
1300 /* Loop over all columns and copy and fill */
1301 switch (nbat->FFormat)
1309 for (a = a0; a < a1; a++)
1311 i = cell[a]*nbat->fstride;
1314 f[a][YY] += fnb[i+1];
1315 f[a][ZZ] += fnb[i+2];
1320 for (a = a0; a < a1; a++)
1322 i = cell[a]*nbat->fstride;
1324 for (fa = 0; fa < nfa; fa++)
1326 f[a][XX] += out[fa].f[i];
1327 f[a][YY] += out[fa].f[i+1];
1328 f[a][ZZ] += out[fa].f[i+2];
1338 for (a = a0; a < a1; a++)
1340 i = X4_IND_A(cell[a]);
1342 f[a][XX] += fnb[i+XX*PACK_X4];
1343 f[a][YY] += fnb[i+YY*PACK_X4];
1344 f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1349 for (a = a0; a < a1; a++)
1351 i = X4_IND_A(cell[a]);
1353 for (fa = 0; fa < nfa; fa++)
1355 f[a][XX] += out[fa].f[i+XX*PACK_X4];
1356 f[a][YY] += out[fa].f[i+YY*PACK_X4];
1357 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1367 for (a = a0; a < a1; a++)
1369 i = X8_IND_A(cell[a]);
1371 f[a][XX] += fnb[i+XX*PACK_X8];
1372 f[a][YY] += fnb[i+YY*PACK_X8];
1373 f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1378 for (a = a0; a < a1; a++)
1380 i = X8_IND_A(cell[a]);
1382 for (fa = 0; fa < nfa; fa++)
1384 f[a][XX] += out[fa].f[i+XX*PACK_X8];
1385 f[a][YY] += out[fa].f[i+YY*PACK_X8];
1386 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1392 gmx_incons("Unsupported nbnxn_atomdata_t format");
1396 static gmx_inline unsigned char reverse_bits(unsigned char b)
1398 /* http://graphics.stanford.edu/~seander/bithacks.html#ReverseByteWith64BitsDiv */
1399 return (b * 0x0202020202ULL & 0x010884422010ULL) % 1023;
1402 static void nbnxn_atomdata_add_nbat_f_to_f_treereduce(const nbnxn_atomdata_t *nbat,
1405 const nbnxn_buffer_flags_t *flags = &nbat->buffer_flags;
1407 int next_pow2 = 1<<(gmx_log2i(nth-1)+1);
1409 assert(nbat->nout == nth); /* tree-reduce currently only works for nout==nth */
1411 memset(nbat->syncStep, 0, sizeof(*(nbat->syncStep))*nth);
1413 #pragma omp parallel num_threads(nth)
1419 th = gmx_omp_get_thread_num();
1421 for (group_size = 2; group_size < 2*next_pow2; group_size *= 2)
1423 int index[2], group_pos, partner_pos, wu;
1424 int partner_th = th ^ (group_size/2);
1429 /* wait on partner thread - replaces full barrier */
1430 int sync_th, sync_group_size;
1432 tMPI_Atomic_memory_barrier(); /* gurantee data is saved before marking work as done */
1433 tMPI_Atomic_set(&(nbat->syncStep[th]), group_size/2); /* mark previous step as completed */
1435 /* find thread to sync with. Equal to partner_th unless nth is not a power of two. */
1436 for (sync_th = partner_th, sync_group_size = group_size; sync_th >= nth && sync_group_size > 2; sync_group_size /= 2)
1438 sync_th &= ~(sync_group_size/4);
1440 if (sync_th < nth) /* otherwise nothing to sync index[1] will be >=nout */
1442 /* wait on the thread which computed input data in previous step */
1443 while (tMPI_Atomic_get((volatile tMPI_Atomic_t*)&(nbat->syncStep[sync_th])) < group_size/2)
1447 /* guarantee that no later load happens before wait loop is finisehd */
1448 tMPI_Atomic_memory_barrier();
1450 #else /* TMPI_ATOMICS */
1455 /* Calculate buffers to sum (result goes into first buffer) */
1456 group_pos = th % group_size;
1457 index[0] = th - group_pos;
1458 index[1] = index[0] + group_size/2;
1460 /* If no second buffer, nothing to do */
1461 if (index[1] >= nbat->nout && group_size > 2)
1466 #if NBNXN_BUFFERFLAG_MAX_THREADS > 256
1467 #error reverse_bits assumes max 256 threads
1469 /* Position is permuted so that one of the 2 vectors being added was computed on the same thread in the previous step.
1470 This improves locality and enables to sync with just a single thread between steps (=the levels in the btree).
1471 The permutation which allows this corresponds to reversing the bits of the group position.
1473 group_pos = reverse_bits(group_pos)/(256/group_size);
1475 partner_pos = group_pos ^ 1;
1477 /* loop over two work-units (own and partner) */
1478 for (wu = 0; wu < 2; wu++)
1482 if (partner_th < nth)
1484 break; /* partner exists we don't have to do his work */
1488 group_pos = partner_pos;
1492 /* Calculate the cell-block range for our thread */
1493 b0 = (flags->nflag* group_pos )/group_size;
1494 b1 = (flags->nflag*(group_pos+1))/group_size;
1496 for (b = b0; b < b1; b++)
1498 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1499 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1501 if ((flags->flag[b] & (1ULL<<index[1])) || group_size > 2)
1503 #ifdef GMX_NBNXN_SIMD
1504 nbnxn_atomdata_reduce_reals_simd
1506 nbnxn_atomdata_reduce_reals
1508 (nbat->out[index[0]].f,
1509 (flags->flag[b] & (1ULL<<index[0])) || group_size > 2,
1510 &(nbat->out[index[1]].f), 1, i0, i1);
1513 else if (!(flags->flag[b] & (1ULL<<index[0])))
1515 nbnxn_atomdata_clear_reals(nbat->out[index[0]].f,
1525 static void nbnxn_atomdata_add_nbat_f_to_f_stdreduce(const nbnxn_atomdata_t *nbat,
1529 #pragma omp parallel for num_threads(nth) schedule(static)
1530 for (th = 0; th < nth; th++)
1532 const nbnxn_buffer_flags_t *flags;
1536 real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1539 flags = &nbat->buffer_flags;
1541 /* Calculate the cell-block range for our thread */
1542 b0 = (flags->nflag* th )/nth;
1543 b1 = (flags->nflag*(th+1))/nth;
1545 for (b = b0; b < b1; b++)
1547 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1548 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1551 for (out = 1; out < nbat->nout; out++)
1553 if (flags->flag[b] & (1U<<out))
1555 fptr[nfptr++] = nbat->out[out].f;
1560 #ifdef GMX_NBNXN_SIMD
1561 nbnxn_atomdata_reduce_reals_simd
1563 nbnxn_atomdata_reduce_reals
1566 flags->flag[b] & (1U<<0),
1570 else if (!(flags->flag[b] & (1U<<0)))
1572 nbnxn_atomdata_clear_reals(nbat->out[0].f,
1579 /* Add the force array(s) from nbnxn_atomdata_t to f */
1580 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t nbs,
1582 const nbnxn_atomdata_t *nbat,
1588 nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1594 na = nbs->natoms_nonlocal;
1598 na = nbs->natoms_local;
1601 a0 = nbs->natoms_local;
1602 na = nbs->natoms_nonlocal - nbs->natoms_local;
1606 nth = gmx_omp_nthreads_get(emntNonbonded);
1610 if (locality != eatAll)
1612 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1615 /* Reduce the force thread output buffers into buffer 0, before adding
1616 * them to the, differently ordered, "real" force buffer.
1618 if (nbat->bUseTreeReduce)
1620 nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbat, nth);
1624 nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbat, nth);
1627 #pragma omp parallel for num_threads(nth) schedule(static)
1628 for (th = 0; th < nth; th++)
1630 nbnxn_atomdata_add_nbat_f_to_f_part(nbs, nbat,
1638 nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1641 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1642 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1645 const nbnxn_atomdata_output_t *out;
1652 for (s = 0; s < SHIFTS; s++)
1655 for (th = 0; th < nbat->nout; th++)
1657 sum[XX] += out[th].fshift[s*DIM+XX];
1658 sum[YY] += out[th].fshift[s*DIM+YY];
1659 sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1661 rvec_inc(fshift[s], sum);