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/utility/gmxomp.h"
55 #include "gromacs/utility/smalloc.h"
57 /* Default nbnxn allocation routine, allocates NBNXN_MEM_ALIGN byte aligned */
58 void nbnxn_alloc_aligned(void **ptr, size_t nbytes)
60 *ptr = save_malloc_aligned("ptr", __FILE__, __LINE__, nbytes, 1, NBNXN_MEM_ALIGN);
63 /* Free function for memory allocated with nbnxn_alloc_aligned */
64 void nbnxn_free_aligned(void *ptr)
69 /* Reallocation wrapper function for nbnxn data structures */
70 void nbnxn_realloc_void(void **ptr,
71 int nbytes_copy, int nbytes_new,
77 ma(&ptr_new, nbytes_new);
79 if (nbytes_new > 0 && ptr_new == NULL)
81 gmx_fatal(FARGS, "Allocation of %d bytes failed", nbytes_new);
86 if (nbytes_new < nbytes_copy)
88 gmx_incons("In nbnxn_realloc_void: new size less than copy size");
90 memcpy(ptr_new, *ptr, nbytes_copy);
99 /* Reallocate the nbnxn_atomdata_t for a size of n atoms */
100 void nbnxn_atomdata_realloc(nbnxn_atomdata_t *nbat, int n)
104 nbnxn_realloc_void((void **)&nbat->type,
105 nbat->natoms*sizeof(*nbat->type),
106 n*sizeof(*nbat->type),
107 nbat->alloc, nbat->free);
108 nbnxn_realloc_void((void **)&nbat->lj_comb,
109 nbat->natoms*2*sizeof(*nbat->lj_comb),
110 n*2*sizeof(*nbat->lj_comb),
111 nbat->alloc, nbat->free);
112 if (nbat->XFormat != nbatXYZQ)
114 nbnxn_realloc_void((void **)&nbat->q,
115 nbat->natoms*sizeof(*nbat->q),
117 nbat->alloc, nbat->free);
119 if (nbat->nenergrp > 1)
121 nbnxn_realloc_void((void **)&nbat->energrp,
122 nbat->natoms/nbat->na_c*sizeof(*nbat->energrp),
123 n/nbat->na_c*sizeof(*nbat->energrp),
124 nbat->alloc, nbat->free);
126 nbnxn_realloc_void((void **)&nbat->x,
127 nbat->natoms*nbat->xstride*sizeof(*nbat->x),
128 n*nbat->xstride*sizeof(*nbat->x),
129 nbat->alloc, nbat->free);
130 for (t = 0; t < nbat->nout; t++)
132 /* Allocate one element extra for possible signaling with CUDA */
133 nbnxn_realloc_void((void **)&nbat->out[t].f,
134 nbat->natoms*nbat->fstride*sizeof(*nbat->out[t].f),
135 n*nbat->fstride*sizeof(*nbat->out[t].f),
136 nbat->alloc, nbat->free);
141 /* Initializes an nbnxn_atomdata_output_t data structure */
142 static void nbnxn_atomdata_output_init(nbnxn_atomdata_output_t *out,
144 int nenergrp, int stride,
150 ma((void **)&out->fshift, SHIFTS*DIM*sizeof(*out->fshift));
151 out->nV = nenergrp*nenergrp;
152 ma((void **)&out->Vvdw, out->nV*sizeof(*out->Vvdw));
153 ma((void **)&out->Vc, out->nV*sizeof(*out->Vc ));
155 if (nb_kernel_type == nbnxnk4xN_SIMD_4xN ||
156 nb_kernel_type == nbnxnk4xN_SIMD_2xNN)
158 cj_size = nbnxn_kernel_to_cj_size(nb_kernel_type);
159 out->nVS = nenergrp*nenergrp*stride*(cj_size>>1)*cj_size;
160 ma((void **)&out->VSvdw, out->nVS*sizeof(*out->VSvdw));
161 ma((void **)&out->VSc, out->nVS*sizeof(*out->VSc ));
169 static void copy_int_to_nbat_int(const int *a, int na, int na_round,
170 const int *in, int fill, int *innb)
175 for (i = 0; i < na; i++)
177 innb[j++] = in[a[i]];
179 /* Complete the partially filled last cell with fill */
180 for (; i < na_round; i++)
186 static void clear_nbat_real(int na, int nbatFormat, real *xnb, int a0)
193 for (a = 0; a < na; a++)
195 for (d = 0; d < DIM; d++)
197 xnb[(a0+a)*STRIDE_XYZ+d] = 0;
202 for (a = 0; a < na; a++)
204 for (d = 0; d < DIM; d++)
206 xnb[(a0+a)*STRIDE_XYZQ+d] = 0;
212 c = a0 & (PACK_X4-1);
213 for (a = 0; a < na; a++)
215 xnb[j+XX*PACK_X4] = 0;
216 xnb[j+YY*PACK_X4] = 0;
217 xnb[j+ZZ*PACK_X4] = 0;
222 j += (DIM-1)*PACK_X4;
229 c = a0 & (PACK_X8-1);
230 for (a = 0; a < na; a++)
232 xnb[j+XX*PACK_X8] = 0;
233 xnb[j+YY*PACK_X8] = 0;
234 xnb[j+ZZ*PACK_X8] = 0;
239 j += (DIM-1)*PACK_X8;
247 void copy_rvec_to_nbat_real(const int *a, int na, int na_round,
248 rvec *x, int nbatFormat, real *xnb, int a0,
249 int cx, int cy, int cz)
253 /* We might need to place filler particles to fill up the cell to na_round.
254 * The coefficients (LJ and q) for such particles are zero.
255 * But we might still get NaN as 0*NaN when distances are too small.
256 * We hope that -107 nm is far away enough from to zero
257 * to avoid accidental short distances to particles shifted down for pbc.
259 #define NBAT_FAR_AWAY 107
265 for (i = 0; i < na; i++)
267 xnb[j++] = x[a[i]][XX];
268 xnb[j++] = x[a[i]][YY];
269 xnb[j++] = x[a[i]][ZZ];
271 /* Complete the partially filled last cell with copies of the last element.
272 * This simplifies the bounding box calculation and avoid
273 * numerical issues with atoms that are coincidentally close.
275 for (; i < na_round; i++)
277 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
278 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
279 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
284 for (i = 0; i < na; i++)
286 xnb[j++] = x[a[i]][XX];
287 xnb[j++] = x[a[i]][YY];
288 xnb[j++] = x[a[i]][ZZ];
291 /* Complete the partially filled last cell with particles far apart */
292 for (; i < na_round; i++)
294 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
295 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
296 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
302 c = a0 & (PACK_X4-1);
303 for (i = 0; i < na; i++)
305 xnb[j+XX*PACK_X4] = x[a[i]][XX];
306 xnb[j+YY*PACK_X4] = x[a[i]][YY];
307 xnb[j+ZZ*PACK_X4] = x[a[i]][ZZ];
312 j += (DIM-1)*PACK_X4;
316 /* Complete the partially filled last cell with particles far apart */
317 for (; i < na_round; i++)
319 xnb[j+XX*PACK_X4] = -NBAT_FAR_AWAY*(1 + cx);
320 xnb[j+YY*PACK_X4] = -NBAT_FAR_AWAY*(1 + cy);
321 xnb[j+ZZ*PACK_X4] = -NBAT_FAR_AWAY*(1 + cz + i);
326 j += (DIM-1)*PACK_X4;
333 c = a0 & (PACK_X8 - 1);
334 for (i = 0; i < na; i++)
336 xnb[j+XX*PACK_X8] = x[a[i]][XX];
337 xnb[j+YY*PACK_X8] = x[a[i]][YY];
338 xnb[j+ZZ*PACK_X8] = x[a[i]][ZZ];
343 j += (DIM-1)*PACK_X8;
347 /* Complete the partially filled last cell with particles far apart */
348 for (; i < na_round; i++)
350 xnb[j+XX*PACK_X8] = -NBAT_FAR_AWAY*(1 + cx);
351 xnb[j+YY*PACK_X8] = -NBAT_FAR_AWAY*(1 + cy);
352 xnb[j+ZZ*PACK_X8] = -NBAT_FAR_AWAY*(1 + cz + i);
357 j += (DIM-1)*PACK_X8;
363 gmx_incons("Unsupported nbnxn_atomdata_t format");
367 /* Stores the LJ parameter data in a format convenient for different kernels */
368 static void set_lj_parameter_data(nbnxn_atomdata_t *nbat, gmx_bool bSIMD)
377 /* nbfp_s4 stores two parameters using a stride of 4,
378 * because this would suit x86 SIMD single-precision
379 * quad-load intrinsics. There's a slight inefficiency in
380 * allocating and initializing nbfp_s4 when it might not
381 * be used, but introducing the conditional code is not
382 * really worth it. */
383 nbat->alloc((void **)&nbat->nbfp_s4, nt*nt*4*sizeof(*nbat->nbfp_s4));
384 for (i = 0; i < nt; i++)
386 for (j = 0; j < nt; j++)
388 nbat->nbfp_s4[(i*nt+j)*4+0] = nbat->nbfp[(i*nt+j)*2+0];
389 nbat->nbfp_s4[(i*nt+j)*4+1] = nbat->nbfp[(i*nt+j)*2+1];
390 nbat->nbfp_s4[(i*nt+j)*4+2] = 0;
391 nbat->nbfp_s4[(i*nt+j)*4+3] = 0;
396 /* We use combination rule data for SIMD combination rule kernels
397 * and with LJ-PME kernels. We then only need parameters per atom type,
398 * not per pair of atom types.
400 switch (nbat->comb_rule)
403 nbat->comb_rule = ljcrGEOM;
405 for (i = 0; i < nt; i++)
407 /* Store the sqrt of the diagonal from the nbfp matrix */
408 nbat->nbfp_comb[i*2 ] = sqrt(nbat->nbfp[(i*nt+i)*2 ]);
409 nbat->nbfp_comb[i*2+1] = sqrt(nbat->nbfp[(i*nt+i)*2+1]);
413 for (i = 0; i < nt; i++)
415 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
416 c6 = nbat->nbfp[(i*nt+i)*2 ];
417 c12 = nbat->nbfp[(i*nt+i)*2+1];
418 if (c6 > 0 && c12 > 0)
420 /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
421 * so we get 6*C6 and 12*C12 after combining.
423 nbat->nbfp_comb[i*2 ] = 0.5*pow(c12/c6, 1.0/6.0);
424 nbat->nbfp_comb[i*2+1] = sqrt(c6*c6/c12);
428 nbat->nbfp_comb[i*2 ] = 0;
429 nbat->nbfp_comb[i*2+1] = 0;
434 /* We always store the full matrix (see code above) */
437 gmx_incons("Unknown combination rule");
442 #ifdef GMX_NBNXN_SIMD
444 nbnxn_atomdata_init_simple_exclusion_masks(nbnxn_atomdata_t *nbat)
447 const int simd_width = GMX_SIMD_REAL_WIDTH;
449 /* Set the diagonal cluster pair exclusion mask setup data.
450 * In the kernel we check 0 < j - i to generate the masks.
451 * Here we store j - i for generating the mask for the first i,
452 * we substract 0.5 to avoid rounding issues.
453 * In the kernel we can subtract 1 to generate the subsequent mask.
455 int simd_4xn_diag_size;
456 const real simdFalse = -1, simdTrue = 1;
457 real *simd_interaction_array;
459 simd_4xn_diag_size = max(NBNXN_CPU_CLUSTER_I_SIZE, simd_width);
460 snew_aligned(nbat->simd_4xn_diagonal_j_minus_i, simd_4xn_diag_size, NBNXN_MEM_ALIGN);
461 for (j = 0; j < simd_4xn_diag_size; j++)
463 nbat->simd_4xn_diagonal_j_minus_i[j] = j - 0.5;
466 snew_aligned(nbat->simd_2xnn_diagonal_j_minus_i, simd_width, NBNXN_MEM_ALIGN);
467 for (j = 0; j < simd_width/2; j++)
469 /* The j-cluster size is half the SIMD width */
470 nbat->simd_2xnn_diagonal_j_minus_i[j] = j - 0.5;
471 /* The next half of the SIMD width is for i + 1 */
472 nbat->simd_2xnn_diagonal_j_minus_i[simd_width/2+j] = j - 1 - 0.5;
475 /* We use up to 32 bits for exclusion masking.
476 * The same masks are used for the 4xN and 2x(N+N) kernels.
477 * The masks are read either into epi32 SIMD registers or into
478 * real SIMD registers (together with a cast).
479 * In single precision this means the real and epi32 SIMD registers
481 * In double precision the epi32 registers can be smaller than
482 * the real registers, so depending on the architecture, we might
483 * need to use two, identical, 32-bit masks per real.
485 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*simd_width;
486 snew_aligned(nbat->simd_exclusion_filter1, simd_excl_size, NBNXN_MEM_ALIGN);
487 snew_aligned(nbat->simd_exclusion_filter2, simd_excl_size*2, NBNXN_MEM_ALIGN);
489 for (j = 0; j < simd_excl_size; j++)
491 /* Set the consecutive bits for masking pair exclusions */
492 nbat->simd_exclusion_filter1[j] = (1U << j);
493 nbat->simd_exclusion_filter2[j*2 + 0] = (1U << j);
494 nbat->simd_exclusion_filter2[j*2 + 1] = (1U << j);
497 #if (defined GMX_SIMD_IBM_QPX)
498 /* The QPX kernels shouldn't do the bit masking that is done on
499 * x86, because the SIMD units lack bit-wise operations. Instead,
500 * we generate a vector of all 2^4 possible ways an i atom
501 * interacts with its 4 j atoms. Each array entry contains
502 * simd_width signed ints that are read in a single SIMD
503 * load. These ints must contain values that will be interpreted
504 * as true and false when loaded in the SIMD floating-point
505 * registers, ie. any positive or any negative value,
506 * respectively. Each array entry encodes how this i atom will
507 * interact with the 4 j atoms. Matching code exists in
508 * set_ci_top_excls() to generate indices into this array. Those
509 * indices are used in the kernels. */
511 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
512 const int qpx_simd_width = GMX_SIMD_REAL_WIDTH;
513 snew_aligned(simd_interaction_array, simd_excl_size * qpx_simd_width, NBNXN_MEM_ALIGN);
514 for (j = 0; j < simd_excl_size; j++)
516 int index = j * qpx_simd_width;
517 for (i = 0; i < qpx_simd_width; i++)
519 simd_interaction_array[index + i] = (j & (1 << i)) ? simdTrue : simdFalse;
522 nbat->simd_interaction_array = simd_interaction_array;
527 /* Initializes an nbnxn_atomdata_t data structure */
528 void nbnxn_atomdata_init(FILE *fp,
529 nbnxn_atomdata_t *nbat,
531 int enbnxninitcombrule,
532 int ntype, const real *nbfp,
535 nbnxn_alloc_t *alloc,
541 gmx_bool simple, bCombGeom, bCombLB, bSIMD;
545 nbat->alloc = nbnxn_alloc_aligned;
553 nbat->free = nbnxn_free_aligned;
562 fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
564 nbat->ntype = ntype + 1;
565 nbat->alloc((void **)&nbat->nbfp,
566 nbat->ntype*nbat->ntype*2*sizeof(*nbat->nbfp));
567 nbat->alloc((void **)&nbat->nbfp_comb, nbat->ntype*2*sizeof(*nbat->nbfp_comb));
569 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
570 * force-field floating point parameters.
573 ptr = getenv("GMX_LJCOMB_TOL");
578 sscanf(ptr, "%lf", &dbl);
584 /* Temporarily fill nbat->nbfp_comb with sigma and epsilon
585 * to check for the LB rule.
587 for (i = 0; i < ntype; i++)
589 c6 = nbfp[(i*ntype+i)*2 ]/6.0;
590 c12 = nbfp[(i*ntype+i)*2+1]/12.0;
591 if (c6 > 0 && c12 > 0)
593 nbat->nbfp_comb[i*2 ] = pow(c12/c6, 1.0/6.0);
594 nbat->nbfp_comb[i*2+1] = 0.25*c6*c6/c12;
596 else if (c6 == 0 && c12 == 0)
598 nbat->nbfp_comb[i*2 ] = 0;
599 nbat->nbfp_comb[i*2+1] = 0;
603 /* Can not use LB rule with only dispersion or repulsion */
608 for (i = 0; i < nbat->ntype; i++)
610 for (j = 0; j < nbat->ntype; j++)
612 if (i < ntype && j < ntype)
614 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
615 * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
617 c6 = nbfp[(i*ntype+j)*2 ];
618 c12 = nbfp[(i*ntype+j)*2+1];
619 nbat->nbfp[(i*nbat->ntype+j)*2 ] = c6;
620 nbat->nbfp[(i*nbat->ntype+j)*2+1] = c12;
622 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
623 bCombGeom = bCombGeom &&
624 gmx_within_tol(c6*c6, nbfp[(i*ntype+i)*2 ]*nbfp[(j*ntype+j)*2 ], tol) &&
625 gmx_within_tol(c12*c12, nbfp[(i*ntype+i)*2+1]*nbfp[(j*ntype+j)*2+1], tol);
627 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
631 ((c6 == 0 && c12 == 0 &&
632 (nbat->nbfp_comb[i*2+1] == 0 || nbat->nbfp_comb[j*2+1] == 0)) ||
633 (c6 > 0 && c12 > 0 &&
634 gmx_within_tol(pow(c12/c6, 1.0/6.0), 0.5*(nbat->nbfp_comb[i*2]+nbat->nbfp_comb[j*2]), tol) &&
635 gmx_within_tol(0.25*c6*c6/c12, sqrt(nbat->nbfp_comb[i*2+1]*nbat->nbfp_comb[j*2+1]), tol)));
639 /* Add zero parameters for the additional dummy atom type */
640 nbat->nbfp[(i*nbat->ntype+j)*2 ] = 0;
641 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 0;
647 fprintf(debug, "Combination rules: geometric %d Lorentz-Berthelot %d\n",
651 simple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
653 switch (enbnxninitcombrule)
655 case enbnxninitcombruleDETECT:
656 /* We prefer the geometic combination rule,
657 * as that gives a slightly faster kernel than the LB rule.
661 nbat->comb_rule = ljcrGEOM;
665 nbat->comb_rule = ljcrLB;
669 nbat->comb_rule = ljcrNONE;
671 nbat->free(nbat->nbfp_comb);
676 if (nbat->comb_rule == ljcrNONE)
678 fprintf(fp, "Using full Lennard-Jones parameter combination matrix\n\n");
682 fprintf(fp, "Using %s Lennard-Jones combination rule\n\n",
683 nbat->comb_rule == ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
687 case enbnxninitcombruleGEOM:
688 nbat->comb_rule = ljcrGEOM;
690 case enbnxninitcombruleLB:
691 nbat->comb_rule = ljcrLB;
693 case enbnxninitcombruleNONE:
694 nbat->comb_rule = ljcrNONE;
696 nbat->free(nbat->nbfp_comb);
699 gmx_incons("Unknown enbnxninitcombrule");
702 bSIMD = (nb_kernel_type == nbnxnk4xN_SIMD_4xN ||
703 nb_kernel_type == nbnxnk4xN_SIMD_2xNN);
705 set_lj_parameter_data(nbat, bSIMD);
709 nbat->lj_comb = NULL;
716 pack_x = max(NBNXN_CPU_CLUSTER_I_SIZE,
717 nbnxn_kernel_to_cj_size(nb_kernel_type));
721 nbat->XFormat = nbatX4;
724 nbat->XFormat = nbatX8;
727 gmx_incons("Unsupported packing width");
732 nbat->XFormat = nbatXYZ;
735 nbat->FFormat = nbat->XFormat;
739 nbat->XFormat = nbatXYZQ;
740 nbat->FFormat = nbatXYZ;
743 nbat->nenergrp = n_energygroups;
746 /* Energy groups not supported yet for super-sub lists */
747 if (n_energygroups > 1 && fp != NULL)
749 fprintf(fp, "\nNOTE: With GPUs, reporting energy group contributions is not supported\n\n");
753 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
754 if (nbat->nenergrp > 64)
756 gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
759 while (nbat->nenergrp > (1<<nbat->neg_2log))
763 nbat->energrp = NULL;
764 nbat->alloc((void **)&nbat->shift_vec, SHIFTS*sizeof(*nbat->shift_vec));
765 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
766 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
769 #ifdef GMX_NBNXN_SIMD
772 nbnxn_atomdata_init_simple_exclusion_masks(nbat);
776 /* Initialize the output data structures */
778 snew(nbat->out, nbat->nout);
780 for (i = 0; i < nbat->nout; i++)
782 nbnxn_atomdata_output_init(&nbat->out[i],
784 nbat->nenergrp, 1<<nbat->neg_2log,
787 nbat->buffer_flags.flag = NULL;
788 nbat->buffer_flags.flag_nalloc = 0;
790 nth = gmx_omp_nthreads_get(emntNonbonded);
792 ptr = getenv("GMX_USE_TREEREDUCE");
795 nbat->bUseTreeReduce = strtol(ptr, 0, 10);
798 else if (nth > 8) /*on the CPU we currently don't benefit even at 32*/
800 nbat->bUseTreeReduce = 1;
805 nbat->bUseTreeReduce = 0;
807 if (nbat->bUseTreeReduce)
811 fprintf(fp, "Using tree force reduction\n\n");
813 snew(nbat->syncStep, nth);
817 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
818 const int *type, int na,
823 /* The LJ params follow the combination rule:
824 * copy the params for the type array to the atom array.
826 for (is = 0; is < na; is += PACK_X4)
828 for (k = 0; k < PACK_X4; k++)
831 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
832 ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
837 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
838 const int *type, int na,
843 /* The LJ params follow the combination rule:
844 * copy the params for the type array to the atom array.
846 for (is = 0; is < na; is += PACK_X8)
848 for (k = 0; k < PACK_X8; k++)
851 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
852 ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
857 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
858 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t *nbat,
860 const nbnxn_search_t nbs,
864 const nbnxn_grid_t *grid;
866 for (g = 0; g < ngrid; g++)
868 grid = &nbs->grid[g];
870 /* Loop over all columns and copy and fill */
871 for (i = 0; i < grid->ncx*grid->ncy; i++)
873 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
874 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
876 copy_int_to_nbat_int(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
877 type, nbat->ntype-1, nbat->type+ash);
879 if (nbat->comb_rule != ljcrNONE)
881 if (nbat->XFormat == nbatX4)
883 copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
884 nbat->type+ash, ncz*grid->na_sc,
885 nbat->lj_comb+ash*2);
887 else if (nbat->XFormat == nbatX8)
889 copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
890 nbat->type+ash, ncz*grid->na_sc,
891 nbat->lj_comb+ash*2);
898 /* Sets the charges in nbnxn_atomdata_t *nbat */
899 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t *nbat,
901 const nbnxn_search_t nbs,
904 int g, cxy, ncz, ash, na, na_round, i, j;
906 const nbnxn_grid_t *grid;
908 for (g = 0; g < ngrid; g++)
910 grid = &nbs->grid[g];
912 /* Loop over all columns and copy and fill */
913 for (cxy = 0; cxy < grid->ncx*grid->ncy; cxy++)
915 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
916 na = grid->cxy_na[cxy];
917 na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
919 if (nbat->XFormat == nbatXYZQ)
921 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
922 for (i = 0; i < na; i++)
924 *q = charge[nbs->a[ash+i]];
927 /* Complete the partially filled last cell with zeros */
928 for (; i < na_round; i++)
937 for (i = 0; i < na; i++)
939 *q = charge[nbs->a[ash+i]];
942 /* Complete the partially filled last cell with zeros */
943 for (; i < na_round; i++)
953 /* Set the charges of perturbed atoms in nbnxn_atomdata_t to 0.
954 * This is to automatically remove the RF/PME self term in the nbnxn kernels.
955 * Part of the zero interactions are still calculated in the normal kernels.
956 * All perturbed interactions are calculated in the free energy kernel,
957 * using the original charge and LJ data, not nbnxn_atomdata_t.
959 static void nbnxn_atomdata_mask_fep(nbnxn_atomdata_t *nbat,
961 const nbnxn_search_t nbs)
964 int stride_q, g, nsubc, c_offset, c, subc, i, ind;
965 const nbnxn_grid_t *grid;
967 if (nbat->XFormat == nbatXYZQ)
969 q = nbat->x + ZZ + 1;
970 stride_q = STRIDE_XYZQ;
978 for (g = 0; g < ngrid; g++)
980 grid = &nbs->grid[g];
987 nsubc = GPU_NSUBCELL;
990 c_offset = grid->cell0*grid->na_sc;
992 /* Loop over all columns and copy and fill */
993 for (c = 0; c < grid->nc*nsubc; c++)
995 /* Does this cluster contain perturbed particles? */
996 if (grid->fep[c] != 0)
998 for (i = 0; i < grid->na_c; i++)
1000 /* Is this a perturbed particle? */
1001 if (grid->fep[c] & (1 << i))
1003 ind = c_offset + c*grid->na_c + i;
1004 /* Set atom type and charge to non-interacting */
1005 nbat->type[ind] = nbat->ntype - 1;
1006 q[ind*stride_q] = 0;
1014 /* Copies the energy group indices to a reordered and packed array */
1015 static void copy_egp_to_nbat_egps(const int *a, int na, int na_round,
1016 int na_c, int bit_shift,
1017 const int *in, int *innb)
1023 for (i = 0; i < na; i += na_c)
1025 /* Store na_c energy group numbers into one int */
1027 for (sa = 0; sa < na_c; sa++)
1032 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
1037 /* Complete the partially filled last cell with fill */
1038 for (; i < na_round; i += na_c)
1044 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
1045 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t *nbat,
1047 const nbnxn_search_t nbs,
1051 const nbnxn_grid_t *grid;
1053 if (nbat->nenergrp == 1)
1058 for (g = 0; g < ngrid; g++)
1060 grid = &nbs->grid[g];
1062 /* Loop over all columns and copy and fill */
1063 for (i = 0; i < grid->ncx*grid->ncy; i++)
1065 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
1066 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
1068 copy_egp_to_nbat_egps(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
1069 nbat->na_c, nbat->neg_2log,
1070 atinfo, nbat->energrp+(ash>>grid->na_c_2log));
1075 /* Sets all required atom parameter data in nbnxn_atomdata_t */
1076 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
1078 const nbnxn_search_t nbs,
1079 const t_mdatoms *mdatoms,
1084 if (locality == eatLocal)
1093 nbnxn_atomdata_set_atomtypes(nbat, ngrid, nbs, mdatoms->typeA);
1095 nbnxn_atomdata_set_charges(nbat, ngrid, nbs, mdatoms->chargeA);
1099 nbnxn_atomdata_mask_fep(nbat, ngrid, nbs);
1102 nbnxn_atomdata_set_energygroups(nbat, ngrid, nbs, atinfo);
1105 /* Copies the shift vector array to nbnxn_atomdata_t */
1106 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox,
1108 nbnxn_atomdata_t *nbat)
1112 nbat->bDynamicBox = bDynamicBox;
1113 for (i = 0; i < SHIFTS; i++)
1115 copy_rvec(shift_vec[i], nbat->shift_vec[i]);
1119 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1120 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
1124 nbnxn_atomdata_t *nbat)
1147 nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
1150 nth = gmx_omp_nthreads_get(emntPairsearch);
1152 #pragma omp parallel for num_threads(nth) schedule(static)
1153 for (th = 0; th < nth; th++)
1157 for (g = g0; g < g1; g++)
1159 const nbnxn_grid_t *grid;
1160 int cxy0, cxy1, cxy;
1162 grid = &nbs->grid[g];
1164 cxy0 = (grid->ncx*grid->ncy* th +nth-1)/nth;
1165 cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
1167 for (cxy = cxy0; cxy < cxy1; cxy++)
1169 int na, ash, na_fill;
1171 na = grid->cxy_na[cxy];
1172 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1174 if (g == 0 && FillLocal)
1177 (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
1181 /* We fill only the real particle locations.
1182 * We assume the filling entries at the end have been
1183 * properly set before during ns.
1187 copy_rvec_to_nbat_real(nbs->a+ash, na, na_fill, x,
1188 nbat->XFormat, nbat->x, ash,
1196 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
1201 for (i = i0; i < i1; i++)
1208 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1210 real ** gmx_restrict src,
1218 /* The destination buffer contains data, add to it */
1219 for (i = i0; i < i1; i++)
1221 for (s = 0; s < nsrc; s++)
1223 dest[i] += src[s][i];
1229 /* The destination buffer is unitialized, set it first */
1230 for (i = i0; i < i1; i++)
1232 dest[i] = src[0][i];
1233 for (s = 1; s < nsrc; s++)
1235 dest[i] += src[s][i];
1242 nbnxn_atomdata_reduce_reals_simd(real gmx_unused * gmx_restrict dest,
1243 gmx_bool gmx_unused bDestSet,
1244 real gmx_unused ** gmx_restrict src,
1245 int gmx_unused nsrc,
1246 int gmx_unused i0, int gmx_unused i1)
1248 #ifdef GMX_NBNXN_SIMD
1249 /* The SIMD width here is actually independent of that in the kernels,
1250 * but we use the same width for simplicity (usually optimal anyhow).
1253 gmx_simd_real_t dest_SSE, src_SSE;
1257 for (i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1259 dest_SSE = gmx_simd_load_r(dest+i);
1260 for (s = 0; s < nsrc; s++)
1262 src_SSE = gmx_simd_load_r(src[s]+i);
1263 dest_SSE = gmx_simd_add_r(dest_SSE, src_SSE);
1265 gmx_simd_store_r(dest+i, dest_SSE);
1270 for (i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1272 dest_SSE = gmx_simd_load_r(src[0]+i);
1273 for (s = 1; s < nsrc; s++)
1275 src_SSE = gmx_simd_load_r(src[s]+i);
1276 dest_SSE = gmx_simd_add_r(dest_SSE, src_SSE);
1278 gmx_simd_store_r(dest+i, dest_SSE);
1284 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1286 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1287 const nbnxn_atomdata_t *nbat,
1288 nbnxn_atomdata_output_t *out,
1299 /* Loop over all columns and copy and fill */
1300 switch (nbat->FFormat)
1308 for (a = a0; a < a1; a++)
1310 i = cell[a]*nbat->fstride;
1313 f[a][YY] += fnb[i+1];
1314 f[a][ZZ] += fnb[i+2];
1319 for (a = a0; a < a1; a++)
1321 i = cell[a]*nbat->fstride;
1323 for (fa = 0; fa < nfa; fa++)
1325 f[a][XX] += out[fa].f[i];
1326 f[a][YY] += out[fa].f[i+1];
1327 f[a][ZZ] += out[fa].f[i+2];
1337 for (a = a0; a < a1; a++)
1339 i = X4_IND_A(cell[a]);
1341 f[a][XX] += fnb[i+XX*PACK_X4];
1342 f[a][YY] += fnb[i+YY*PACK_X4];
1343 f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1348 for (a = a0; a < a1; a++)
1350 i = X4_IND_A(cell[a]);
1352 for (fa = 0; fa < nfa; fa++)
1354 f[a][XX] += out[fa].f[i+XX*PACK_X4];
1355 f[a][YY] += out[fa].f[i+YY*PACK_X4];
1356 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1366 for (a = a0; a < a1; a++)
1368 i = X8_IND_A(cell[a]);
1370 f[a][XX] += fnb[i+XX*PACK_X8];
1371 f[a][YY] += fnb[i+YY*PACK_X8];
1372 f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1377 for (a = a0; a < a1; a++)
1379 i = X8_IND_A(cell[a]);
1381 for (fa = 0; fa < nfa; fa++)
1383 f[a][XX] += out[fa].f[i+XX*PACK_X8];
1384 f[a][YY] += out[fa].f[i+YY*PACK_X8];
1385 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1391 gmx_incons("Unsupported nbnxn_atomdata_t format");
1395 static gmx_inline unsigned char reverse_bits(unsigned char b)
1397 /* http://graphics.stanford.edu/~seander/bithacks.html#ReverseByteWith64BitsDiv */
1398 return (b * 0x0202020202ULL & 0x010884422010ULL) % 1023;
1401 static void nbnxn_atomdata_add_nbat_f_to_f_treereduce(const nbnxn_atomdata_t *nbat,
1404 const nbnxn_buffer_flags_t *flags = &nbat->buffer_flags;
1406 int next_pow2 = 1<<(gmx_log2i(nth-1)+1);
1408 assert(nbat->nout == nth); /* tree-reduce currently only works for nout==nth */
1410 memset(nbat->syncStep, 0, sizeof(*(nbat->syncStep))*nth);
1412 #pragma omp parallel num_threads(nth)
1418 th = gmx_omp_get_thread_num();
1420 for (group_size = 2; group_size < 2*next_pow2; group_size *= 2)
1422 int index[2], group_pos, partner_pos, wu;
1423 int partner_th = th ^ (group_size/2);
1428 /* wait on partner thread - replaces full barrier */
1429 int sync_th, sync_group_size;
1431 tMPI_Atomic_memory_barrier(); /* gurantee data is saved before marking work as done */
1432 tMPI_Atomic_set(&(nbat->syncStep[th]), group_size/2); /* mark previous step as completed */
1434 /* find thread to sync with. Equal to partner_th unless nth is not a power of two. */
1435 for (sync_th = partner_th, sync_group_size = group_size; sync_th >= nth && sync_group_size > 2; sync_group_size /= 2)
1437 sync_th &= ~(sync_group_size/4);
1439 if (sync_th < nth) /* otherwise nothing to sync index[1] will be >=nout */
1441 /* wait on the thread which computed input data in previous step */
1442 while (tMPI_Atomic_get((volatile tMPI_Atomic_t*)&(nbat->syncStep[sync_th])) < group_size/2)
1446 /* guarantee that no later load happens before wait loop is finisehd */
1447 tMPI_Atomic_memory_barrier();
1449 #else /* TMPI_ATOMICS */
1454 /* Calculate buffers to sum (result goes into first buffer) */
1455 group_pos = th % group_size;
1456 index[0] = th - group_pos;
1457 index[1] = index[0] + group_size/2;
1459 /* If no second buffer, nothing to do */
1460 if (index[1] >= nbat->nout && group_size > 2)
1465 #if NBNXN_BUFFERFLAG_MAX_THREADS > 256
1466 #error reverse_bits assumes max 256 threads
1468 /* Position is permuted so that one of the 2 vectors being added was computed on the same thread in the previous step.
1469 This improves locality and enables to sync with just a single thread between steps (=the levels in the btree).
1470 The permutation which allows this corresponds to reversing the bits of the group position.
1472 group_pos = reverse_bits(group_pos)/(256/group_size);
1474 partner_pos = group_pos ^ 1;
1476 /* loop over two work-units (own and partner) */
1477 for (wu = 0; wu < 2; wu++)
1481 if (partner_th < nth)
1483 break; /* partner exists we don't have to do his work */
1487 group_pos = partner_pos;
1491 /* Calculate the cell-block range for our thread */
1492 b0 = (flags->nflag* group_pos )/group_size;
1493 b1 = (flags->nflag*(group_pos+1))/group_size;
1495 for (b = b0; b < b1; b++)
1497 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1498 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1500 if ((flags->flag[b] & (1ULL<<index[1])) || group_size > 2)
1502 #ifdef GMX_NBNXN_SIMD
1503 nbnxn_atomdata_reduce_reals_simd
1505 nbnxn_atomdata_reduce_reals
1507 (nbat->out[index[0]].f,
1508 (flags->flag[b] & (1ULL<<index[0])) || group_size > 2,
1509 &(nbat->out[index[1]].f), 1, i0, i1);
1512 else if (!(flags->flag[b] & (1ULL<<index[0])))
1514 nbnxn_atomdata_clear_reals(nbat->out[index[0]].f,
1524 static void nbnxn_atomdata_add_nbat_f_to_f_stdreduce(const nbnxn_atomdata_t *nbat,
1528 #pragma omp parallel for num_threads(nth) schedule(static)
1529 for (th = 0; th < nth; th++)
1531 const nbnxn_buffer_flags_t *flags;
1535 real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1538 flags = &nbat->buffer_flags;
1540 /* Calculate the cell-block range for our thread */
1541 b0 = (flags->nflag* th )/nth;
1542 b1 = (flags->nflag*(th+1))/nth;
1544 for (b = b0; b < b1; b++)
1546 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1547 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1550 for (out = 1; out < nbat->nout; out++)
1552 if (flags->flag[b] & (1U<<out))
1554 fptr[nfptr++] = nbat->out[out].f;
1559 #ifdef GMX_NBNXN_SIMD
1560 nbnxn_atomdata_reduce_reals_simd
1562 nbnxn_atomdata_reduce_reals
1565 flags->flag[b] & (1U<<0),
1569 else if (!(flags->flag[b] & (1U<<0)))
1571 nbnxn_atomdata_clear_reals(nbat->out[0].f,
1578 /* Add the force array(s) from nbnxn_atomdata_t to f */
1579 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t nbs,
1581 const nbnxn_atomdata_t *nbat,
1587 nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1593 na = nbs->natoms_nonlocal;
1597 na = nbs->natoms_local;
1600 a0 = nbs->natoms_local;
1601 na = nbs->natoms_nonlocal - nbs->natoms_local;
1605 nth = gmx_omp_nthreads_get(emntNonbonded);
1609 if (locality != eatAll)
1611 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1614 /* Reduce the force thread output buffers into buffer 0, before adding
1615 * them to the, differently ordered, "real" force buffer.
1617 if (nbat->bUseTreeReduce)
1619 nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbat, nth);
1623 nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbat, nth);
1626 #pragma omp parallel for num_threads(nth) schedule(static)
1627 for (th = 0; th < nth; th++)
1629 nbnxn_atomdata_add_nbat_f_to_f_part(nbs, nbat,
1637 nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1640 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1641 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1644 const nbnxn_atomdata_output_t *out;
1651 for (s = 0; s < SHIFTS; s++)
1654 for (th = 0; th < nbat->nout; th++)
1656 sum[XX] += out[th].fshift[s*DIM+XX];
1657 sum[YY] += out[th].fshift[s*DIM+YY];
1658 sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1660 rvec_inc(fshift[s], sum);