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.
43 #include "gromacs/utility/smalloc.h"
46 #include "nbnxn_consts.h"
47 #include "nbnxn_internal.h"
48 #include "nbnxn_atomdata.h"
49 #include "nbnxn_search.h"
50 #include "gromacs/utility/gmxomp.h"
51 #include "gmx_omp_nthreads.h"
52 #include "thread_mpi/atomic.h"
54 /* Default nbnxn allocation routine, allocates NBNXN_MEM_ALIGN byte aligned */
55 void nbnxn_alloc_aligned(void **ptr, size_t nbytes)
57 *ptr = save_malloc_aligned("ptr", __FILE__, __LINE__, nbytes, 1, NBNXN_MEM_ALIGN);
60 /* Free function for memory allocated with nbnxn_alloc_aligned */
61 void nbnxn_free_aligned(void *ptr)
66 /* Reallocation wrapper function for nbnxn data structures */
67 void nbnxn_realloc_void(void **ptr,
68 int nbytes_copy, int nbytes_new,
74 ma(&ptr_new, nbytes_new);
76 if (nbytes_new > 0 && ptr_new == NULL)
78 gmx_fatal(FARGS, "Allocation of %d bytes failed", nbytes_new);
83 if (nbytes_new < nbytes_copy)
85 gmx_incons("In nbnxn_realloc_void: new size less than copy size");
87 memcpy(ptr_new, *ptr, nbytes_copy);
96 /* Reallocate the nbnxn_atomdata_t for a size of n atoms */
97 void nbnxn_atomdata_realloc(nbnxn_atomdata_t *nbat, int n)
101 nbnxn_realloc_void((void **)&nbat->type,
102 nbat->natoms*sizeof(*nbat->type),
103 n*sizeof(*nbat->type),
104 nbat->alloc, nbat->free);
105 nbnxn_realloc_void((void **)&nbat->lj_comb,
106 nbat->natoms*2*sizeof(*nbat->lj_comb),
107 n*2*sizeof(*nbat->lj_comb),
108 nbat->alloc, nbat->free);
109 if (nbat->XFormat != nbatXYZQ)
111 nbnxn_realloc_void((void **)&nbat->q,
112 nbat->natoms*sizeof(*nbat->q),
114 nbat->alloc, nbat->free);
116 if (nbat->nenergrp > 1)
118 nbnxn_realloc_void((void **)&nbat->energrp,
119 nbat->natoms/nbat->na_c*sizeof(*nbat->energrp),
120 n/nbat->na_c*sizeof(*nbat->energrp),
121 nbat->alloc, nbat->free);
123 nbnxn_realloc_void((void **)&nbat->x,
124 nbat->natoms*nbat->xstride*sizeof(*nbat->x),
125 n*nbat->xstride*sizeof(*nbat->x),
126 nbat->alloc, nbat->free);
127 for (t = 0; t < nbat->nout; t++)
129 /* Allocate one element extra for possible signaling with CUDA */
130 nbnxn_realloc_void((void **)&nbat->out[t].f,
131 nbat->natoms*nbat->fstride*sizeof(*nbat->out[t].f),
132 n*nbat->fstride*sizeof(*nbat->out[t].f),
133 nbat->alloc, nbat->free);
138 /* Initializes an nbnxn_atomdata_output_t data structure */
139 static void nbnxn_atomdata_output_init(nbnxn_atomdata_output_t *out,
141 int nenergrp, int stride,
147 ma((void **)&out->fshift, SHIFTS*DIM*sizeof(*out->fshift));
148 out->nV = nenergrp*nenergrp;
149 ma((void **)&out->Vvdw, out->nV*sizeof(*out->Vvdw));
150 ma((void **)&out->Vc, out->nV*sizeof(*out->Vc ));
152 if (nb_kernel_type == nbnxnk4xN_SIMD_4xN ||
153 nb_kernel_type == nbnxnk4xN_SIMD_2xNN)
155 cj_size = nbnxn_kernel_to_cj_size(nb_kernel_type);
156 out->nVS = nenergrp*nenergrp*stride*(cj_size>>1)*cj_size;
157 ma((void **)&out->VSvdw, out->nVS*sizeof(*out->VSvdw));
158 ma((void **)&out->VSc, out->nVS*sizeof(*out->VSc ));
166 static void copy_int_to_nbat_int(const int *a, int na, int na_round,
167 const int *in, int fill, int *innb)
172 for (i = 0; i < na; i++)
174 innb[j++] = in[a[i]];
176 /* Complete the partially filled last cell with fill */
177 for (; i < na_round; i++)
183 static void clear_nbat_real(int na, int nbatFormat, real *xnb, int a0)
190 for (a = 0; a < na; a++)
192 for (d = 0; d < DIM; d++)
194 xnb[(a0+a)*STRIDE_XYZ+d] = 0;
199 for (a = 0; a < na; a++)
201 for (d = 0; d < DIM; d++)
203 xnb[(a0+a)*STRIDE_XYZQ+d] = 0;
209 c = a0 & (PACK_X4-1);
210 for (a = 0; a < na; a++)
212 xnb[j+XX*PACK_X4] = 0;
213 xnb[j+YY*PACK_X4] = 0;
214 xnb[j+ZZ*PACK_X4] = 0;
219 j += (DIM-1)*PACK_X4;
226 c = a0 & (PACK_X8-1);
227 for (a = 0; a < na; a++)
229 xnb[j+XX*PACK_X8] = 0;
230 xnb[j+YY*PACK_X8] = 0;
231 xnb[j+ZZ*PACK_X8] = 0;
236 j += (DIM-1)*PACK_X8;
244 void copy_rvec_to_nbat_real(const int *a, int na, int na_round,
245 rvec *x, int nbatFormat, real *xnb, int a0,
246 int cx, int cy, int cz)
250 /* We might need to place filler particles to fill up the cell to na_round.
251 * The coefficients (LJ and q) for such particles are zero.
252 * But we might still get NaN as 0*NaN when distances are too small.
253 * We hope that -107 nm is far away enough from to zero
254 * to avoid accidental short distances to particles shifted down for pbc.
256 #define NBAT_FAR_AWAY 107
262 for (i = 0; i < na; i++)
264 xnb[j++] = x[a[i]][XX];
265 xnb[j++] = x[a[i]][YY];
266 xnb[j++] = x[a[i]][ZZ];
268 /* Complete the partially filled last cell with copies of the last element.
269 * This simplifies the bounding box calculation and avoid
270 * numerical issues with atoms that are coincidentally close.
272 for (; i < na_round; i++)
274 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
275 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
276 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
281 for (i = 0; i < na; i++)
283 xnb[j++] = x[a[i]][XX];
284 xnb[j++] = x[a[i]][YY];
285 xnb[j++] = x[a[i]][ZZ];
288 /* Complete the partially filled last cell with particles far apart */
289 for (; i < na_round; i++)
291 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
292 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
293 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
299 c = a0 & (PACK_X4-1);
300 for (i = 0; i < na; i++)
302 xnb[j+XX*PACK_X4] = x[a[i]][XX];
303 xnb[j+YY*PACK_X4] = x[a[i]][YY];
304 xnb[j+ZZ*PACK_X4] = x[a[i]][ZZ];
309 j += (DIM-1)*PACK_X4;
313 /* Complete the partially filled last cell with particles far apart */
314 for (; i < na_round; i++)
316 xnb[j+XX*PACK_X4] = -NBAT_FAR_AWAY*(1 + cx);
317 xnb[j+YY*PACK_X4] = -NBAT_FAR_AWAY*(1 + cy);
318 xnb[j+ZZ*PACK_X4] = -NBAT_FAR_AWAY*(1 + cz + i);
323 j += (DIM-1)*PACK_X4;
330 c = a0 & (PACK_X8 - 1);
331 for (i = 0; i < na; i++)
333 xnb[j+XX*PACK_X8] = x[a[i]][XX];
334 xnb[j+YY*PACK_X8] = x[a[i]][YY];
335 xnb[j+ZZ*PACK_X8] = x[a[i]][ZZ];
340 j += (DIM-1)*PACK_X8;
344 /* Complete the partially filled last cell with particles far apart */
345 for (; i < na_round; i++)
347 xnb[j+XX*PACK_X8] = -NBAT_FAR_AWAY*(1 + cx);
348 xnb[j+YY*PACK_X8] = -NBAT_FAR_AWAY*(1 + cy);
349 xnb[j+ZZ*PACK_X8] = -NBAT_FAR_AWAY*(1 + cz + i);
354 j += (DIM-1)*PACK_X8;
360 gmx_incons("Unsupported nbnxn_atomdata_t format");
364 /* Stores the LJ parameter data in a format convenient for different kernels */
365 static void set_lj_parameter_data(nbnxn_atomdata_t *nbat, gmx_bool bSIMD)
374 /* nbfp_s4 stores two parameters using a stride of 4,
375 * because this would suit x86 SIMD single-precision
376 * quad-load intrinsics. There's a slight inefficiency in
377 * allocating and initializing nbfp_s4 when it might not
378 * be used, but introducing the conditional code is not
379 * really worth it. */
380 nbat->alloc((void **)&nbat->nbfp_s4, nt*nt*4*sizeof(*nbat->nbfp_s4));
381 for (i = 0; i < nt; i++)
383 for (j = 0; j < nt; j++)
385 nbat->nbfp_s4[(i*nt+j)*4+0] = nbat->nbfp[(i*nt+j)*2+0];
386 nbat->nbfp_s4[(i*nt+j)*4+1] = nbat->nbfp[(i*nt+j)*2+1];
387 nbat->nbfp_s4[(i*nt+j)*4+2] = 0;
388 nbat->nbfp_s4[(i*nt+j)*4+3] = 0;
393 /* We use combination rule data for SIMD combination rule kernels
394 * and with LJ-PME kernels. We then only need parameters per atom type,
395 * not per pair of atom types.
397 switch (nbat->comb_rule)
400 nbat->comb_rule = ljcrGEOM;
402 for (i = 0; i < nt; i++)
404 /* Store the sqrt of the diagonal from the nbfp matrix */
405 nbat->nbfp_comb[i*2 ] = sqrt(nbat->nbfp[(i*nt+i)*2 ]);
406 nbat->nbfp_comb[i*2+1] = sqrt(nbat->nbfp[(i*nt+i)*2+1]);
410 for (i = 0; i < nt; i++)
412 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
413 c6 = nbat->nbfp[(i*nt+i)*2 ];
414 c12 = nbat->nbfp[(i*nt+i)*2+1];
415 if (c6 > 0 && c12 > 0)
417 /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
418 * so we get 6*C6 and 12*C12 after combining.
420 nbat->nbfp_comb[i*2 ] = 0.5*pow(c12/c6, 1.0/6.0);
421 nbat->nbfp_comb[i*2+1] = sqrt(c6*c6/c12);
425 nbat->nbfp_comb[i*2 ] = 0;
426 nbat->nbfp_comb[i*2+1] = 0;
431 /* We always store the full matrix (see code above) */
434 gmx_incons("Unknown combination rule");
439 #ifdef GMX_NBNXN_SIMD
441 nbnxn_atomdata_init_simple_exclusion_masks(nbnxn_atomdata_t *nbat)
444 const int simd_width = GMX_SIMD_REAL_WIDTH;
446 /* Set the diagonal cluster pair exclusion mask setup data.
447 * In the kernel we check 0 < j - i to generate the masks.
448 * Here we store j - i for generating the mask for the first i,
449 * we substract 0.5 to avoid rounding issues.
450 * In the kernel we can subtract 1 to generate the subsequent mask.
452 int simd_4xn_diag_size;
453 const real simdFalse = -1, simdTrue = 1;
454 real *simd_interaction_array;
456 simd_4xn_diag_size = max(NBNXN_CPU_CLUSTER_I_SIZE, simd_width);
457 snew_aligned(nbat->simd_4xn_diagonal_j_minus_i, simd_4xn_diag_size, NBNXN_MEM_ALIGN);
458 for (j = 0; j < simd_4xn_diag_size; j++)
460 nbat->simd_4xn_diagonal_j_minus_i[j] = j - 0.5;
463 snew_aligned(nbat->simd_2xnn_diagonal_j_minus_i, simd_width, NBNXN_MEM_ALIGN);
464 for (j = 0; j < simd_width/2; j++)
466 /* The j-cluster size is half the SIMD width */
467 nbat->simd_2xnn_diagonal_j_minus_i[j] = j - 0.5;
468 /* The next half of the SIMD width is for i + 1 */
469 nbat->simd_2xnn_diagonal_j_minus_i[simd_width/2+j] = j - 1 - 0.5;
472 /* We use up to 32 bits for exclusion masking.
473 * The same masks are used for the 4xN and 2x(N+N) kernels.
474 * The masks are read either into epi32 SIMD registers or into
475 * real SIMD registers (together with a cast).
476 * In single precision this means the real and epi32 SIMD registers
478 * In double precision the epi32 registers can be smaller than
479 * the real registers, so depending on the architecture, we might
480 * need to use two, identical, 32-bit masks per real.
482 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*simd_width;
483 snew_aligned(nbat->simd_exclusion_filter1, simd_excl_size, NBNXN_MEM_ALIGN);
484 snew_aligned(nbat->simd_exclusion_filter2, simd_excl_size*2, NBNXN_MEM_ALIGN);
486 for (j = 0; j < simd_excl_size; j++)
488 /* Set the consecutive bits for masking pair exclusions */
489 nbat->simd_exclusion_filter1[j] = (1U << j);
490 nbat->simd_exclusion_filter2[j*2 + 0] = (1U << j);
491 nbat->simd_exclusion_filter2[j*2 + 1] = (1U << j);
494 #if (defined GMX_SIMD_IBM_QPX)
495 /* The QPX kernels shouldn't do the bit masking that is done on
496 * x86, because the SIMD units lack bit-wise operations. Instead,
497 * we generate a vector of all 2^4 possible ways an i atom
498 * interacts with its 4 j atoms. Each array entry contains
499 * simd_width signed ints that are read in a single SIMD
500 * load. These ints must contain values that will be interpreted
501 * as true and false when loaded in the SIMD floating-point
502 * registers, ie. any positive or any negative value,
503 * respectively. Each array entry encodes how this i atom will
504 * interact with the 4 j atoms. Matching code exists in
505 * set_ci_top_excls() to generate indices into this array. Those
506 * indices are used in the kernels. */
508 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
509 const int qpx_simd_width = GMX_SIMD_REAL_WIDTH;
510 snew_aligned(simd_interaction_array, simd_excl_size * qpx_simd_width, NBNXN_MEM_ALIGN);
511 for (j = 0; j < simd_excl_size; j++)
513 int index = j * qpx_simd_width;
514 for (i = 0; i < qpx_simd_width; i++)
516 simd_interaction_array[index + i] = (j & (1 << i)) ? simdTrue : simdFalse;
519 nbat->simd_interaction_array = simd_interaction_array;
524 /* Initializes an nbnxn_atomdata_t data structure */
525 void nbnxn_atomdata_init(FILE *fp,
526 nbnxn_atomdata_t *nbat,
528 int enbnxninitcombrule,
529 int ntype, const real *nbfp,
532 nbnxn_alloc_t *alloc,
538 gmx_bool simple, bCombGeom, bCombLB, bSIMD;
542 nbat->alloc = nbnxn_alloc_aligned;
550 nbat->free = nbnxn_free_aligned;
559 fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
561 nbat->ntype = ntype + 1;
562 nbat->alloc((void **)&nbat->nbfp,
563 nbat->ntype*nbat->ntype*2*sizeof(*nbat->nbfp));
564 nbat->alloc((void **)&nbat->nbfp_comb, nbat->ntype*2*sizeof(*nbat->nbfp_comb));
566 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
567 * force-field floating point parameters.
570 ptr = getenv("GMX_LJCOMB_TOL");
575 sscanf(ptr, "%lf", &dbl);
581 /* Temporarily fill nbat->nbfp_comb with sigma and epsilon
582 * to check for the LB rule.
584 for (i = 0; i < ntype; i++)
586 c6 = nbfp[(i*ntype+i)*2 ]/6.0;
587 c12 = nbfp[(i*ntype+i)*2+1]/12.0;
588 if (c6 > 0 && c12 > 0)
590 nbat->nbfp_comb[i*2 ] = pow(c12/c6, 1.0/6.0);
591 nbat->nbfp_comb[i*2+1] = 0.25*c6*c6/c12;
593 else if (c6 == 0 && c12 == 0)
595 nbat->nbfp_comb[i*2 ] = 0;
596 nbat->nbfp_comb[i*2+1] = 0;
600 /* Can not use LB rule with only dispersion or repulsion */
605 for (i = 0; i < nbat->ntype; i++)
607 for (j = 0; j < nbat->ntype; j++)
609 if (i < ntype && j < ntype)
611 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
612 * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
614 c6 = nbfp[(i*ntype+j)*2 ];
615 c12 = nbfp[(i*ntype+j)*2+1];
616 nbat->nbfp[(i*nbat->ntype+j)*2 ] = c6;
617 nbat->nbfp[(i*nbat->ntype+j)*2+1] = c12;
619 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
620 bCombGeom = bCombGeom &&
621 gmx_within_tol(c6*c6, nbfp[(i*ntype+i)*2 ]*nbfp[(j*ntype+j)*2 ], tol) &&
622 gmx_within_tol(c12*c12, nbfp[(i*ntype+i)*2+1]*nbfp[(j*ntype+j)*2+1], tol);
624 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
628 ((c6 == 0 && c12 == 0 &&
629 (nbat->nbfp_comb[i*2+1] == 0 || nbat->nbfp_comb[j*2+1] == 0)) ||
630 (c6 > 0 && c12 > 0 &&
631 gmx_within_tol(pow(c12/c6, 1.0/6.0), 0.5*(nbat->nbfp_comb[i*2]+nbat->nbfp_comb[j*2]), tol) &&
632 gmx_within_tol(0.25*c6*c6/c12, sqrt(nbat->nbfp_comb[i*2+1]*nbat->nbfp_comb[j*2+1]), tol)));
636 /* Add zero parameters for the additional dummy atom type */
637 nbat->nbfp[(i*nbat->ntype+j)*2 ] = 0;
638 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 0;
644 fprintf(debug, "Combination rules: geometric %d Lorentz-Berthelot %d\n",
648 simple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
650 switch (enbnxninitcombrule)
652 case enbnxninitcombruleDETECT:
653 /* We prefer the geometic combination rule,
654 * as that gives a slightly faster kernel than the LB rule.
658 nbat->comb_rule = ljcrGEOM;
662 nbat->comb_rule = ljcrLB;
666 nbat->comb_rule = ljcrNONE;
668 nbat->free(nbat->nbfp_comb);
673 if (nbat->comb_rule == ljcrNONE)
675 fprintf(fp, "Using full Lennard-Jones parameter combination matrix\n\n");
679 fprintf(fp, "Using %s Lennard-Jones combination rule\n\n",
680 nbat->comb_rule == ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
684 case enbnxninitcombruleGEOM:
685 nbat->comb_rule = ljcrGEOM;
687 case enbnxninitcombruleLB:
688 nbat->comb_rule = ljcrLB;
690 case enbnxninitcombruleNONE:
691 nbat->comb_rule = ljcrNONE;
693 nbat->free(nbat->nbfp_comb);
696 gmx_incons("Unknown enbnxninitcombrule");
699 bSIMD = (nb_kernel_type == nbnxnk4xN_SIMD_4xN ||
700 nb_kernel_type == nbnxnk4xN_SIMD_2xNN);
702 set_lj_parameter_data(nbat, bSIMD);
706 nbat->lj_comb = NULL;
713 pack_x = max(NBNXN_CPU_CLUSTER_I_SIZE,
714 nbnxn_kernel_to_cj_size(nb_kernel_type));
718 nbat->XFormat = nbatX4;
721 nbat->XFormat = nbatX8;
724 gmx_incons("Unsupported packing width");
729 nbat->XFormat = nbatXYZ;
732 nbat->FFormat = nbat->XFormat;
736 nbat->XFormat = nbatXYZQ;
737 nbat->FFormat = nbatXYZ;
740 nbat->nenergrp = n_energygroups;
743 /* Energy groups not supported yet for super-sub lists */
744 if (n_energygroups > 1 && fp != NULL)
746 fprintf(fp, "\nNOTE: With GPUs, reporting energy group contributions is not supported\n\n");
750 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
751 if (nbat->nenergrp > 64)
753 gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
756 while (nbat->nenergrp > (1<<nbat->neg_2log))
760 nbat->energrp = NULL;
761 nbat->alloc((void **)&nbat->shift_vec, SHIFTS*sizeof(*nbat->shift_vec));
762 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
763 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
766 #ifdef GMX_NBNXN_SIMD
769 nbnxn_atomdata_init_simple_exclusion_masks(nbat);
773 /* Initialize the output data structures */
775 snew(nbat->out, nbat->nout);
777 for (i = 0; i < nbat->nout; i++)
779 nbnxn_atomdata_output_init(&nbat->out[i],
781 nbat->nenergrp, 1<<nbat->neg_2log,
784 nbat->buffer_flags.flag = NULL;
785 nbat->buffer_flags.flag_nalloc = 0;
787 nth = gmx_omp_nthreads_get(emntNonbonded);
789 ptr = getenv("GMX_USE_TREEREDUCE");
792 nbat->bUseTreeReduce = strtol(ptr, 0, 10);
795 else if (nth > 8) /*on the CPU we currently don't benefit even at 32*/
797 nbat->bUseTreeReduce = 1;
802 nbat->bUseTreeReduce = 0;
804 if (nbat->bUseTreeReduce)
808 fprintf(fp, "Using tree force reduction\n\n");
810 snew(nbat->syncStep, nth);
814 static void copy_lj_to_nbat_lj_comb_x4(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_X4)
825 for (k = 0; k < PACK_X4; k++)
828 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
829 ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
834 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
835 const int *type, int na,
840 /* The LJ params follow the combination rule:
841 * copy the params for the type array to the atom array.
843 for (is = 0; is < na; is += PACK_X8)
845 for (k = 0; k < PACK_X8; k++)
848 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
849 ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
854 /* Sets the atom type in nbnxn_atomdata_t */
855 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t *nbat,
857 const nbnxn_search_t nbs,
861 const nbnxn_grid_t *grid;
863 for (g = 0; g < ngrid; g++)
865 grid = &nbs->grid[g];
867 /* Loop over all columns and copy and fill */
868 for (i = 0; i < grid->ncx*grid->ncy; i++)
870 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
871 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
873 copy_int_to_nbat_int(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
874 type, nbat->ntype-1, nbat->type+ash);
879 /* Sets the LJ combination rule parameters in nbnxn_atomdata_t */
880 static void nbnxn_atomdata_set_ljcombparams(nbnxn_atomdata_t *nbat,
882 const nbnxn_search_t nbs)
885 const nbnxn_grid_t *grid;
887 if (nbat->comb_rule != ljcrNONE)
889 for (g = 0; g < ngrid; g++)
891 grid = &nbs->grid[g];
893 /* Loop over all columns and copy and fill */
894 for (i = 0; i < grid->ncx*grid->ncy; i++)
896 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
897 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
899 if (nbat->XFormat == nbatX4)
901 copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
902 nbat->type+ash, ncz*grid->na_sc,
903 nbat->lj_comb+ash*2);
905 else if (nbat->XFormat == nbatX8)
907 copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
908 nbat->type+ash, ncz*grid->na_sc,
909 nbat->lj_comb+ash*2);
916 /* Sets the charges in nbnxn_atomdata_t *nbat */
917 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t *nbat,
919 const nbnxn_search_t nbs,
922 int g, cxy, ncz, ash, na, na_round, i, j;
924 const nbnxn_grid_t *grid;
926 for (g = 0; g < ngrid; g++)
928 grid = &nbs->grid[g];
930 /* Loop over all columns and copy and fill */
931 for (cxy = 0; cxy < grid->ncx*grid->ncy; cxy++)
933 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
934 na = grid->cxy_na[cxy];
935 na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
937 if (nbat->XFormat == nbatXYZQ)
939 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
940 for (i = 0; i < na; i++)
942 *q = charge[nbs->a[ash+i]];
945 /* Complete the partially filled last cell with zeros */
946 for (; i < na_round; i++)
955 for (i = 0; i < na; i++)
957 *q = charge[nbs->a[ash+i]];
960 /* Complete the partially filled last cell with zeros */
961 for (; i < na_round; i++)
971 /* Set the charges of perturbed atoms in nbnxn_atomdata_t to 0.
972 * This is to automatically remove the RF/PME self term in the nbnxn kernels.
973 * Part of the zero interactions are still calculated in the normal kernels.
974 * All perturbed interactions are calculated in the free energy kernel,
975 * using the original charge and LJ data, not nbnxn_atomdata_t.
977 static void nbnxn_atomdata_mask_fep(nbnxn_atomdata_t *nbat,
979 const nbnxn_search_t nbs)
982 int stride_q, g, nsubc, c_offset, c, subc, i, ind;
983 const nbnxn_grid_t *grid;
985 if (nbat->XFormat == nbatXYZQ)
987 q = nbat->x + ZZ + 1;
988 stride_q = STRIDE_XYZQ;
996 for (g = 0; g < ngrid; g++)
998 grid = &nbs->grid[g];
1005 nsubc = GPU_NSUBCELL;
1008 c_offset = grid->cell0*grid->na_sc;
1010 /* Loop over all columns and copy and fill */
1011 for (c = 0; c < grid->nc*nsubc; c++)
1013 /* Does this cluster contain perturbed particles? */
1014 if (grid->fep[c] != 0)
1016 for (i = 0; i < grid->na_c; i++)
1018 /* Is this a perturbed particle? */
1019 if (grid->fep[c] & (1 << i))
1021 ind = c_offset + c*grid->na_c + i;
1022 /* Set atom type and charge to non-interacting */
1023 nbat->type[ind] = nbat->ntype - 1;
1024 q[ind*stride_q] = 0;
1032 /* Copies the energy group indices to a reordered and packed array */
1033 static void copy_egp_to_nbat_egps(const int *a, int na, int na_round,
1034 int na_c, int bit_shift,
1035 const int *in, int *innb)
1041 for (i = 0; i < na; i += na_c)
1043 /* Store na_c energy group numbers into one int */
1045 for (sa = 0; sa < na_c; sa++)
1050 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
1055 /* Complete the partially filled last cell with fill */
1056 for (; i < na_round; i += na_c)
1062 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
1063 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t *nbat,
1065 const nbnxn_search_t nbs,
1069 const nbnxn_grid_t *grid;
1071 if (nbat->nenergrp == 1)
1076 for (g = 0; g < ngrid; g++)
1078 grid = &nbs->grid[g];
1080 /* Loop over all columns and copy and fill */
1081 for (i = 0; i < grid->ncx*grid->ncy; i++)
1083 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
1084 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
1086 copy_egp_to_nbat_egps(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
1087 nbat->na_c, nbat->neg_2log,
1088 atinfo, nbat->energrp+(ash>>grid->na_c_2log));
1093 /* Sets all required atom parameter data in nbnxn_atomdata_t */
1094 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
1096 const nbnxn_search_t nbs,
1097 const t_mdatoms *mdatoms,
1102 if (locality == eatLocal)
1111 nbnxn_atomdata_set_atomtypes(nbat, ngrid, nbs, mdatoms->typeA);
1113 nbnxn_atomdata_set_charges(nbat, ngrid, nbs, mdatoms->chargeA);
1117 nbnxn_atomdata_mask_fep(nbat, ngrid, nbs);
1120 /* This must be done after masking types for FEP */
1121 nbnxn_atomdata_set_ljcombparams(nbat, ngrid, nbs);
1123 nbnxn_atomdata_set_energygroups(nbat, ngrid, nbs, atinfo);
1126 /* Copies the shift vector array to nbnxn_atomdata_t */
1127 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox,
1129 nbnxn_atomdata_t *nbat)
1133 nbat->bDynamicBox = bDynamicBox;
1134 for (i = 0; i < SHIFTS; i++)
1136 copy_rvec(shift_vec[i], nbat->shift_vec[i]);
1140 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1141 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
1145 nbnxn_atomdata_t *nbat)
1168 nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
1171 nth = gmx_omp_nthreads_get(emntPairsearch);
1173 #pragma omp parallel for num_threads(nth) schedule(static)
1174 for (th = 0; th < nth; th++)
1178 for (g = g0; g < g1; g++)
1180 const nbnxn_grid_t *grid;
1181 int cxy0, cxy1, cxy;
1183 grid = &nbs->grid[g];
1185 cxy0 = (grid->ncx*grid->ncy* th +nth-1)/nth;
1186 cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
1188 for (cxy = cxy0; cxy < cxy1; cxy++)
1190 int na, ash, na_fill;
1192 na = grid->cxy_na[cxy];
1193 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1195 if (g == 0 && FillLocal)
1198 (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
1202 /* We fill only the real particle locations.
1203 * We assume the filling entries at the end have been
1204 * properly set before during ns.
1208 copy_rvec_to_nbat_real(nbs->a+ash, na, na_fill, x,
1209 nbat->XFormat, nbat->x, ash,
1217 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
1222 for (i = i0; i < i1; i++)
1229 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1231 real ** gmx_restrict src,
1239 /* The destination buffer contains data, add to it */
1240 for (i = i0; i < i1; i++)
1242 for (s = 0; s < nsrc; s++)
1244 dest[i] += src[s][i];
1250 /* The destination buffer is unitialized, set it first */
1251 for (i = i0; i < i1; i++)
1253 dest[i] = src[0][i];
1254 for (s = 1; s < nsrc; s++)
1256 dest[i] += src[s][i];
1263 nbnxn_atomdata_reduce_reals_simd(real gmx_unused * gmx_restrict dest,
1264 gmx_bool gmx_unused bDestSet,
1265 real gmx_unused ** gmx_restrict src,
1266 int gmx_unused nsrc,
1267 int gmx_unused i0, int gmx_unused i1)
1269 #ifdef GMX_NBNXN_SIMD
1270 /* The SIMD width here is actually independent of that in the kernels,
1271 * but we use the same width for simplicity (usually optimal anyhow).
1274 gmx_simd_real_t dest_SSE, src_SSE;
1278 for (i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1280 dest_SSE = gmx_simd_load_r(dest+i);
1281 for (s = 0; s < nsrc; s++)
1283 src_SSE = gmx_simd_load_r(src[s]+i);
1284 dest_SSE = gmx_simd_add_r(dest_SSE, src_SSE);
1286 gmx_simd_store_r(dest+i, dest_SSE);
1291 for (i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1293 dest_SSE = gmx_simd_load_r(src[0]+i);
1294 for (s = 1; s < nsrc; s++)
1296 src_SSE = gmx_simd_load_r(src[s]+i);
1297 dest_SSE = gmx_simd_add_r(dest_SSE, src_SSE);
1299 gmx_simd_store_r(dest+i, dest_SSE);
1305 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1307 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1308 const nbnxn_atomdata_t *nbat,
1309 nbnxn_atomdata_output_t *out,
1320 /* Loop over all columns and copy and fill */
1321 switch (nbat->FFormat)
1329 for (a = a0; a < a1; a++)
1331 i = cell[a]*nbat->fstride;
1334 f[a][YY] += fnb[i+1];
1335 f[a][ZZ] += fnb[i+2];
1340 for (a = a0; a < a1; a++)
1342 i = cell[a]*nbat->fstride;
1344 for (fa = 0; fa < nfa; fa++)
1346 f[a][XX] += out[fa].f[i];
1347 f[a][YY] += out[fa].f[i+1];
1348 f[a][ZZ] += out[fa].f[i+2];
1358 for (a = a0; a < a1; a++)
1360 i = X4_IND_A(cell[a]);
1362 f[a][XX] += fnb[i+XX*PACK_X4];
1363 f[a][YY] += fnb[i+YY*PACK_X4];
1364 f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1369 for (a = a0; a < a1; a++)
1371 i = X4_IND_A(cell[a]);
1373 for (fa = 0; fa < nfa; fa++)
1375 f[a][XX] += out[fa].f[i+XX*PACK_X4];
1376 f[a][YY] += out[fa].f[i+YY*PACK_X4];
1377 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1387 for (a = a0; a < a1; a++)
1389 i = X8_IND_A(cell[a]);
1391 f[a][XX] += fnb[i+XX*PACK_X8];
1392 f[a][YY] += fnb[i+YY*PACK_X8];
1393 f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1398 for (a = a0; a < a1; a++)
1400 i = X8_IND_A(cell[a]);
1402 for (fa = 0; fa < nfa; fa++)
1404 f[a][XX] += out[fa].f[i+XX*PACK_X8];
1405 f[a][YY] += out[fa].f[i+YY*PACK_X8];
1406 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1412 gmx_incons("Unsupported nbnxn_atomdata_t format");
1416 static gmx_inline unsigned char reverse_bits(unsigned char b)
1418 /* http://graphics.stanford.edu/~seander/bithacks.html#ReverseByteWith64BitsDiv */
1419 return (b * 0x0202020202ULL & 0x010884422010ULL) % 1023;
1422 static void nbnxn_atomdata_add_nbat_f_to_f_treereduce(const nbnxn_atomdata_t *nbat,
1425 const nbnxn_buffer_flags_t *flags = &nbat->buffer_flags;
1427 int next_pow2 = 1<<(gmx_log2i(nth-1)+1);
1429 assert(nbat->nout == nth); /* tree-reduce currently only works for nout==nth */
1431 memset(nbat->syncStep, 0, sizeof(*(nbat->syncStep))*nth);
1433 #pragma omp parallel num_threads(nth)
1439 th = gmx_omp_get_thread_num();
1441 for (group_size = 2; group_size < 2*next_pow2; group_size *= 2)
1443 int index[2], group_pos, partner_pos, wu;
1444 int partner_th = th ^ (group_size/2);
1449 /* wait on partner thread - replaces full barrier */
1450 int sync_th, sync_group_size;
1452 tMPI_Atomic_memory_barrier(); /* gurantee data is saved before marking work as done */
1453 tMPI_Atomic_set(&(nbat->syncStep[th]), group_size/2); /* mark previous step as completed */
1455 /* find thread to sync with. Equal to partner_th unless nth is not a power of two. */
1456 for (sync_th = partner_th, sync_group_size = group_size; sync_th >= nth && sync_group_size > 2; sync_group_size /= 2)
1458 sync_th &= ~(sync_group_size/4);
1460 if (sync_th < nth) /* otherwise nothing to sync index[1] will be >=nout */
1462 /* wait on the thread which computed input data in previous step */
1463 while (tMPI_Atomic_get((volatile tMPI_Atomic_t*)&(nbat->syncStep[sync_th])) < group_size/2)
1467 /* guarantee that no later load happens before wait loop is finisehd */
1468 tMPI_Atomic_memory_barrier();
1470 #else /* TMPI_ATOMICS */
1475 /* Calculate buffers to sum (result goes into first buffer) */
1476 group_pos = th % group_size;
1477 index[0] = th - group_pos;
1478 index[1] = index[0] + group_size/2;
1480 /* If no second buffer, nothing to do */
1481 if (index[1] >= nbat->nout && group_size > 2)
1486 #if NBNXN_BUFFERFLAG_MAX_THREADS > 256
1487 #error reverse_bits assumes max 256 threads
1489 /* Position is permuted so that one of the 2 vectors being added was computed on the same thread in the previous step.
1490 This improves locality and enables to sync with just a single thread between steps (=the levels in the btree).
1491 The permutation which allows this corresponds to reversing the bits of the group position.
1493 group_pos = reverse_bits(group_pos)/(256/group_size);
1495 partner_pos = group_pos ^ 1;
1497 /* loop over two work-units (own and partner) */
1498 for (wu = 0; wu < 2; wu++)
1502 if (partner_th < nth)
1504 break; /* partner exists we don't have to do his work */
1508 group_pos = partner_pos;
1512 /* Calculate the cell-block range for our thread */
1513 b0 = (flags->nflag* group_pos )/group_size;
1514 b1 = (flags->nflag*(group_pos+1))/group_size;
1516 for (b = b0; b < b1; b++)
1518 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1519 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1521 if ((flags->flag[b] & (1ULL<<index[1])) || group_size > 2)
1523 #ifdef GMX_NBNXN_SIMD
1524 nbnxn_atomdata_reduce_reals_simd
1526 nbnxn_atomdata_reduce_reals
1528 (nbat->out[index[0]].f,
1529 (flags->flag[b] & (1ULL<<index[0])) || group_size > 2,
1530 &(nbat->out[index[1]].f), 1, i0, i1);
1533 else if (!(flags->flag[b] & (1ULL<<index[0])))
1535 nbnxn_atomdata_clear_reals(nbat->out[index[0]].f,
1545 static void nbnxn_atomdata_add_nbat_f_to_f_stdreduce(const nbnxn_atomdata_t *nbat,
1549 #pragma omp parallel for num_threads(nth) schedule(static)
1550 for (th = 0; th < nth; th++)
1552 const nbnxn_buffer_flags_t *flags;
1556 real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1559 flags = &nbat->buffer_flags;
1561 /* Calculate the cell-block range for our thread */
1562 b0 = (flags->nflag* th )/nth;
1563 b1 = (flags->nflag*(th+1))/nth;
1565 for (b = b0; b < b1; b++)
1567 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1568 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1571 for (out = 1; out < nbat->nout; out++)
1573 if (flags->flag[b] & (1U<<out))
1575 fptr[nfptr++] = nbat->out[out].f;
1580 #ifdef GMX_NBNXN_SIMD
1581 nbnxn_atomdata_reduce_reals_simd
1583 nbnxn_atomdata_reduce_reals
1586 flags->flag[b] & (1U<<0),
1590 else if (!(flags->flag[b] & (1U<<0)))
1592 nbnxn_atomdata_clear_reals(nbat->out[0].f,
1599 /* Add the force array(s) from nbnxn_atomdata_t to f */
1600 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t nbs,
1602 const nbnxn_atomdata_t *nbat,
1608 nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1614 na = nbs->natoms_nonlocal;
1618 na = nbs->natoms_local;
1621 a0 = nbs->natoms_local;
1622 na = nbs->natoms_nonlocal - nbs->natoms_local;
1626 nth = gmx_omp_nthreads_get(emntNonbonded);
1630 if (locality != eatAll)
1632 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1635 /* Reduce the force thread output buffers into buffer 0, before adding
1636 * them to the, differently ordered, "real" force buffer.
1638 if (nbat->bUseTreeReduce)
1640 nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbat, nth);
1644 nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbat, nth);
1647 #pragma omp parallel for num_threads(nth) schedule(static)
1648 for (th = 0; th < nth; th++)
1650 nbnxn_atomdata_add_nbat_f_to_f_part(nbs, nbat,
1658 nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1661 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1662 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1665 const nbnxn_atomdata_output_t *out;
1672 for (s = 0; s < SHIFTS; s++)
1675 for (th = 0; th < nbat->nout; th++)
1677 sum[XX] += out[th].fshift[s*DIM+XX];
1678 sum[YY] += out[th].fshift[s*DIM+YY];
1679 sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1681 rvec_inc(fshift[s], sum);