2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013, 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.
45 #include "nbnxn_consts.h"
46 #include "nbnxn_internal.h"
47 #include "nbnxn_search.h"
48 #include "gmx_omp_nthreads.h"
50 /* Default nbnxn allocation routine, allocates NBNXN_MEM_ALIGN byte aligned */
51 void nbnxn_alloc_aligned(void **ptr, size_t nbytes)
53 *ptr = save_malloc_aligned("ptr", __FILE__, __LINE__, nbytes, 1, NBNXN_MEM_ALIGN);
56 /* Free function for memory allocated with nbnxn_alloc_aligned */
57 void nbnxn_free_aligned(void *ptr)
62 /* Reallocation wrapper function for nbnxn data structures */
63 void nbnxn_realloc_void(void **ptr,
64 int nbytes_copy, int nbytes_new,
70 ma(&ptr_new, nbytes_new);
72 if (nbytes_new > 0 && ptr_new == NULL)
74 gmx_fatal(FARGS, "Allocation of %d bytes failed", nbytes_new);
79 if (nbytes_new < nbytes_copy)
81 gmx_incons("In nbnxn_realloc_void: new size less than copy size");
83 memcpy(ptr_new, *ptr, nbytes_copy);
92 /* Reallocate the nbnxn_atomdata_t for a size of n atoms */
93 void nbnxn_atomdata_realloc(nbnxn_atomdata_t *nbat, int n)
97 nbnxn_realloc_void((void **)&nbat->type,
98 nbat->natoms*sizeof(*nbat->type),
99 n*sizeof(*nbat->type),
100 nbat->alloc, nbat->free);
101 nbnxn_realloc_void((void **)&nbat->lj_comb,
102 nbat->natoms*2*sizeof(*nbat->lj_comb),
103 n*2*sizeof(*nbat->lj_comb),
104 nbat->alloc, nbat->free);
105 if (nbat->XFormat != nbatXYZQ)
107 nbnxn_realloc_void((void **)&nbat->q,
108 nbat->natoms*sizeof(*nbat->q),
110 nbat->alloc, nbat->free);
112 if (nbat->nenergrp > 1)
114 nbnxn_realloc_void((void **)&nbat->energrp,
115 nbat->natoms/nbat->na_c*sizeof(*nbat->energrp),
116 n/nbat->na_c*sizeof(*nbat->energrp),
117 nbat->alloc, nbat->free);
119 nbnxn_realloc_void((void **)&nbat->x,
120 nbat->natoms*nbat->xstride*sizeof(*nbat->x),
121 n*nbat->xstride*sizeof(*nbat->x),
122 nbat->alloc, nbat->free);
123 for (t = 0; t < nbat->nout; t++)
125 /* Allocate one element extra for possible signaling with CUDA */
126 nbnxn_realloc_void((void **)&nbat->out[t].f,
127 nbat->natoms*nbat->fstride*sizeof(*nbat->out[t].f),
128 n*nbat->fstride*sizeof(*nbat->out[t].f),
129 nbat->alloc, nbat->free);
134 /* Initializes an nbnxn_atomdata_output_t data structure */
135 static void nbnxn_atomdata_output_init(nbnxn_atomdata_output_t *out,
137 int nenergrp, int stride,
143 ma((void **)&out->fshift, SHIFTS*DIM*sizeof(*out->fshift));
144 out->nV = nenergrp*nenergrp;
145 ma((void **)&out->Vvdw, out->nV*sizeof(*out->Vvdw));
146 ma((void **)&out->Vc, out->nV*sizeof(*out->Vc ));
148 if (nb_kernel_type == nbnxnk4xN_SIMD_4xN ||
149 nb_kernel_type == nbnxnk4xN_SIMD_2xNN)
151 cj_size = nbnxn_kernel_to_cj_size(nb_kernel_type);
152 out->nVS = nenergrp*nenergrp*stride*(cj_size>>1)*cj_size;
153 ma((void **)&out->VSvdw, out->nVS*sizeof(*out->VSvdw));
154 ma((void **)&out->VSc, out->nVS*sizeof(*out->VSc ));
162 static void copy_int_to_nbat_int(const int *a, int na, int na_round,
163 const int *in, int fill, int *innb)
168 for (i = 0; i < na; i++)
170 innb[j++] = in[a[i]];
172 /* Complete the partially filled last cell with fill */
173 for (; i < na_round; i++)
179 static void clear_nbat_real(int na, int nbatFormat, real *xnb, int a0)
186 for (a = 0; a < na; a++)
188 for (d = 0; d < DIM; d++)
190 xnb[(a0+a)*STRIDE_XYZ+d] = 0;
195 for (a = 0; a < na; a++)
197 for (d = 0; d < DIM; d++)
199 xnb[(a0+a)*STRIDE_XYZQ+d] = 0;
205 c = a0 & (PACK_X4-1);
206 for (a = 0; a < na; a++)
208 xnb[j+XX*PACK_X4] = 0;
209 xnb[j+YY*PACK_X4] = 0;
210 xnb[j+ZZ*PACK_X4] = 0;
215 j += (DIM-1)*PACK_X4;
222 c = a0 & (PACK_X8-1);
223 for (a = 0; a < na; a++)
225 xnb[j+XX*PACK_X8] = 0;
226 xnb[j+YY*PACK_X8] = 0;
227 xnb[j+ZZ*PACK_X8] = 0;
232 j += (DIM-1)*PACK_X8;
240 void copy_rvec_to_nbat_real(const int *a, int na, int na_round,
241 rvec *x, int nbatFormat, real *xnb, int a0,
242 int cx, int cy, int cz)
246 /* We might need to place filler particles to fill up the cell to na_round.
247 * The coefficients (LJ and q) for such particles are zero.
248 * But we might still get NaN as 0*NaN when distances are too small.
249 * We hope that -107 nm is far away enough from to zero
250 * to avoid accidental short distances to particles shifted down for pbc.
252 #define NBAT_FAR_AWAY 107
258 for (i = 0; i < na; i++)
260 xnb[j++] = x[a[i]][XX];
261 xnb[j++] = x[a[i]][YY];
262 xnb[j++] = x[a[i]][ZZ];
264 /* Complete the partially filled last cell with copies of the last element.
265 * This simplifies the bounding box calculation and avoid
266 * numerical issues with atoms that are coincidentally close.
268 for (; i < na_round; i++)
270 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
271 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
272 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
277 for (i = 0; i < na; i++)
279 xnb[j++] = x[a[i]][XX];
280 xnb[j++] = x[a[i]][YY];
281 xnb[j++] = x[a[i]][ZZ];
284 /* Complete the partially filled last cell with particles far apart */
285 for (; i < na_round; i++)
287 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
288 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
289 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
295 c = a0 & (PACK_X4-1);
296 for (i = 0; i < na; i++)
298 xnb[j+XX*PACK_X4] = x[a[i]][XX];
299 xnb[j+YY*PACK_X4] = x[a[i]][YY];
300 xnb[j+ZZ*PACK_X4] = x[a[i]][ZZ];
305 j += (DIM-1)*PACK_X4;
309 /* Complete the partially filled last cell with particles far apart */
310 for (; i < na_round; i++)
312 xnb[j+XX*PACK_X4] = -NBAT_FAR_AWAY*(1 + cx);
313 xnb[j+YY*PACK_X4] = -NBAT_FAR_AWAY*(1 + cy);
314 xnb[j+ZZ*PACK_X4] = -NBAT_FAR_AWAY*(1 + cz + i);
319 j += (DIM-1)*PACK_X4;
326 c = a0 & (PACK_X8 - 1);
327 for (i = 0; i < na; i++)
329 xnb[j+XX*PACK_X8] = x[a[i]][XX];
330 xnb[j+YY*PACK_X8] = x[a[i]][YY];
331 xnb[j+ZZ*PACK_X8] = x[a[i]][ZZ];
336 j += (DIM-1)*PACK_X8;
340 /* Complete the partially filled last cell with particles far apart */
341 for (; i < na_round; i++)
343 xnb[j+XX*PACK_X8] = -NBAT_FAR_AWAY*(1 + cx);
344 xnb[j+YY*PACK_X8] = -NBAT_FAR_AWAY*(1 + cy);
345 xnb[j+ZZ*PACK_X8] = -NBAT_FAR_AWAY*(1 + cz + i);
350 j += (DIM-1)*PACK_X8;
356 gmx_incons("Unsupported nbnxn_atomdata_t format");
360 /* Determines the combination rule (or none) to be used, stores it,
361 * and sets the LJ parameters required with the rule.
363 static void set_combination_rule_data(nbnxn_atomdata_t *nbat)
370 switch (nbat->comb_rule)
373 nbat->comb_rule = ljcrGEOM;
375 for (i = 0; i < nt; i++)
377 /* Copy the diagonal from the nbfp matrix */
378 nbat->nbfp_comb[i*2 ] = sqrt(nbat->nbfp[(i*nt+i)*2 ]);
379 nbat->nbfp_comb[i*2+1] = sqrt(nbat->nbfp[(i*nt+i)*2+1]);
383 for (i = 0; i < nt; i++)
385 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
386 c6 = nbat->nbfp[(i*nt+i)*2 ];
387 c12 = nbat->nbfp[(i*nt+i)*2+1];
388 if (c6 > 0 && c12 > 0)
390 /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
391 * so we get 6*C6 and 12*C12 after combining.
393 nbat->nbfp_comb[i*2 ] = 0.5*pow(c12/c6, 1.0/6.0);
394 nbat->nbfp_comb[i*2+1] = sqrt(c6*c6/c12);
398 nbat->nbfp_comb[i*2 ] = 0;
399 nbat->nbfp_comb[i*2+1] = 0;
404 /* nbfp_s4 stores two parameters using a stride of 4,
405 * because this would suit x86 SIMD single-precision
406 * quad-load intrinsics. There's a slight inefficiency in
407 * allocating and initializing nbfp_s4 when it might not
408 * be used, but introducing the conditional code is not
409 * really worth it. */
410 nbat->alloc((void **)&nbat->nbfp_s4, nt*nt*4*sizeof(*nbat->nbfp_s4));
411 for (i = 0; i < nt; i++)
413 for (j = 0; j < nt; j++)
415 nbat->nbfp_s4[(i*nt+j)*4+0] = nbat->nbfp[(i*nt+j)*2+0];
416 nbat->nbfp_s4[(i*nt+j)*4+1] = nbat->nbfp[(i*nt+j)*2+1];
417 nbat->nbfp_s4[(i*nt+j)*4+2] = 0;
418 nbat->nbfp_s4[(i*nt+j)*4+3] = 0;
423 gmx_incons("Unknown combination rule");
428 #ifdef GMX_NBNXN_SIMD
430 nbnxn_atomdata_init_simple_exclusion_masks(nbnxn_atomdata_t *nbat)
433 const int simd_width = GMX_SIMD_WIDTH_HERE;
435 /* Set the diagonal cluster pair exclusion mask setup data.
436 * In the kernel we check 0 < j - i to generate the masks.
437 * Here we store j - i for generating the mask for the first i,
438 * we substract 0.5 to avoid rounding issues.
439 * In the kernel we can subtract 1 to generate the subsequent mask.
441 int simd_4xn_diag_size;
442 const real simdFalse = -1, simdTrue = 1;
443 real *simd_interaction_array;
445 simd_4xn_diag_size = max(NBNXN_CPU_CLUSTER_I_SIZE, simd_width);
446 snew_aligned(nbat->simd_4xn_diagonal_j_minus_i, simd_4xn_diag_size, NBNXN_MEM_ALIGN);
447 for (j = 0; j < simd_4xn_diag_size; j++)
449 nbat->simd_4xn_diagonal_j_minus_i[j] = j - 0.5;
452 snew_aligned(nbat->simd_2xnn_diagonal_j_minus_i, simd_width, NBNXN_MEM_ALIGN);
453 for (j = 0; j < simd_width/2; j++)
455 /* The j-cluster size is half the SIMD width */
456 nbat->simd_2xnn_diagonal_j_minus_i[j] = j - 0.5;
457 /* The next half of the SIMD width is for i + 1 */
458 nbat->simd_2xnn_diagonal_j_minus_i[simd_width/2+j] = j - 1 - 0.5;
461 /* We use up to 32 bits for exclusion masking.
462 * The same masks are used for the 4xN and 2x(N+N) kernels.
463 * The masks are read either into epi32 SIMD registers or into
464 * real SIMD registers (together with a cast).
465 * In single precision this means the real and epi32 SIMD registers
467 * In double precision the epi32 registers can be smaller than
468 * the real registers, so depending on the architecture, we might
469 * need to use two, identical, 32-bit masks per real.
471 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*simd_width;
472 snew_aligned(nbat->simd_exclusion_filter1, simd_excl_size, NBNXN_MEM_ALIGN);
473 snew_aligned(nbat->simd_exclusion_filter2, simd_excl_size*2, NBNXN_MEM_ALIGN);
475 for (j = 0; j < simd_excl_size; j++)
477 /* Set the consecutive bits for masking pair exclusions */
478 nbat->simd_exclusion_filter1[j] = (1U << j);
479 nbat->simd_exclusion_filter2[j*2 + 0] = (1U << j);
480 nbat->simd_exclusion_filter2[j*2 + 1] = (1U << j);
483 #if (defined GMX_CPU_ACCELERATION_IBM_QPX)
484 /* The QPX kernels shouldn't do the bit masking that is done on
485 * x86, because the SIMD units lack bit-wise operations. Instead,
486 * we generate a vector of all 2^4 possible ways an i atom
487 * interacts with its 4 j atoms. Each array entry contains
488 * simd_width signed ints that are read in a single SIMD
489 * load. These ints must contain values that will be interpreted
490 * as true and false when loaded in the SIMD floating-point
491 * registers, ie. any positive or any negative value,
492 * respectively. Each array entry encodes how this i atom will
493 * interact with the 4 j atoms. Matching code exists in
494 * set_ci_top_excls() to generate indices into this array. Those
495 * indices are used in the kernels. */
497 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
498 const int qpx_simd_width = GMX_SIMD_WIDTH_HERE;
499 snew_aligned(simd_interaction_array, simd_excl_size * qpx_simd_width, NBNXN_MEM_ALIGN);
500 for (j = 0; j < simd_excl_size; j++)
502 int index = j * qpx_simd_width;
503 for (i = 0; i < qpx_simd_width; i++)
505 simd_interaction_array[index + i] = (j & (1 << i)) ? simdTrue : simdFalse;
508 nbat->simd_interaction_array = simd_interaction_array;
513 /* Initializes an nbnxn_atomdata_t data structure */
514 void nbnxn_atomdata_init(FILE *fp,
515 nbnxn_atomdata_t *nbat,
517 int ntype, const real *nbfp,
520 nbnxn_alloc_t *alloc,
526 gmx_bool simple, bCombGeom, bCombLB;
530 nbat->alloc = nbnxn_alloc_aligned;
538 nbat->free = nbnxn_free_aligned;
547 fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
549 nbat->ntype = ntype + 1;
550 nbat->alloc((void **)&nbat->nbfp,
551 nbat->ntype*nbat->ntype*2*sizeof(*nbat->nbfp));
552 nbat->alloc((void **)&nbat->nbfp_comb, nbat->ntype*2*sizeof(*nbat->nbfp_comb));
554 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
555 * force-field floating point parameters.
558 ptr = getenv("GMX_LJCOMB_TOL");
563 sscanf(ptr, "%lf", &dbl);
569 /* Temporarily fill nbat->nbfp_comb with sigma and epsilon
570 * to check for the LB rule.
572 for (i = 0; i < ntype; i++)
574 c6 = nbfp[(i*ntype+i)*2 ]/6.0;
575 c12 = nbfp[(i*ntype+i)*2+1]/12.0;
576 if (c6 > 0 && c12 > 0)
578 nbat->nbfp_comb[i*2 ] = pow(c12/c6, 1.0/6.0);
579 nbat->nbfp_comb[i*2+1] = 0.25*c6*c6/c12;
581 else if (c6 == 0 && c12 == 0)
583 nbat->nbfp_comb[i*2 ] = 0;
584 nbat->nbfp_comb[i*2+1] = 0;
588 /* Can not use LB rule with only dispersion or repulsion */
593 for (i = 0; i < nbat->ntype; i++)
595 for (j = 0; j < nbat->ntype; j++)
597 if (i < ntype && j < ntype)
599 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
600 * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
602 c6 = nbfp[(i*ntype+j)*2 ];
603 c12 = nbfp[(i*ntype+j)*2+1];
604 nbat->nbfp[(i*nbat->ntype+j)*2 ] = c6;
605 nbat->nbfp[(i*nbat->ntype+j)*2+1] = c12;
607 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
608 bCombGeom = bCombGeom &&
609 gmx_within_tol(c6*c6, nbfp[(i*ntype+i)*2 ]*nbfp[(j*ntype+j)*2 ], tol) &&
610 gmx_within_tol(c12*c12, nbfp[(i*ntype+i)*2+1]*nbfp[(j*ntype+j)*2+1], tol);
612 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
616 ((c6 == 0 && c12 == 0 &&
617 (nbat->nbfp_comb[i*2+1] == 0 || nbat->nbfp_comb[j*2+1] == 0)) ||
618 (c6 > 0 && c12 > 0 &&
619 gmx_within_tol(pow(c12/c6, 1.0/6.0), 0.5*(nbat->nbfp_comb[i*2]+nbat->nbfp_comb[j*2]), tol) &&
620 gmx_within_tol(0.25*c6*c6/c12, sqrt(nbat->nbfp_comb[i*2+1]*nbat->nbfp_comb[j*2+1]), tol)));
624 /* Add zero parameters for the additional dummy atom type */
625 nbat->nbfp[(i*nbat->ntype+j)*2 ] = 0;
626 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 0;
632 fprintf(debug, "Combination rules: geometric %d Lorentz-Berthelot %d\n",
636 simple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
640 /* We prefer the geometic combination rule,
641 * as that gives a slightly faster kernel than the LB rule.
645 nbat->comb_rule = ljcrGEOM;
649 nbat->comb_rule = ljcrLB;
653 nbat->comb_rule = ljcrNONE;
655 nbat->free(nbat->nbfp_comb);
660 if (nbat->comb_rule == ljcrNONE)
662 fprintf(fp, "Using full Lennard-Jones parameter combination matrix\n\n");
666 fprintf(fp, "Using %s Lennard-Jones combination rule\n\n",
667 nbat->comb_rule == ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
671 set_combination_rule_data(nbat);
675 nbat->comb_rule = ljcrNONE;
677 nbat->free(nbat->nbfp_comb);
682 nbat->lj_comb = NULL;
687 switch (nb_kernel_type)
689 case nbnxnk4xN_SIMD_4xN:
690 case nbnxnk4xN_SIMD_2xNN:
691 pack_x = max(NBNXN_CPU_CLUSTER_I_SIZE,
692 nbnxn_kernel_to_cj_size(nb_kernel_type));
696 nbat->XFormat = nbatX4;
699 nbat->XFormat = nbatX8;
702 gmx_incons("Unsupported packing width");
706 nbat->XFormat = nbatXYZ;
710 nbat->FFormat = nbat->XFormat;
714 nbat->XFormat = nbatXYZQ;
715 nbat->FFormat = nbatXYZ;
718 nbat->nenergrp = n_energygroups;
721 /* Energy groups not supported yet for super-sub lists */
722 if (n_energygroups > 1 && fp != NULL)
724 fprintf(fp, "\nNOTE: With GPUs, reporting energy group contributions is not supported\n\n");
728 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
729 if (nbat->nenergrp > 64)
731 gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
734 while (nbat->nenergrp > (1<<nbat->neg_2log))
738 nbat->energrp = NULL;
739 nbat->alloc((void **)&nbat->shift_vec, SHIFTS*sizeof(*nbat->shift_vec));
740 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
741 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
744 #ifdef GMX_NBNXN_SIMD
747 nbnxn_atomdata_init_simple_exclusion_masks(nbat);
751 /* Initialize the output data structures */
753 snew(nbat->out, nbat->nout);
755 for (i = 0; i < nbat->nout; i++)
757 nbnxn_atomdata_output_init(&nbat->out[i],
759 nbat->nenergrp, 1<<nbat->neg_2log,
762 nbat->buffer_flags.flag = NULL;
763 nbat->buffer_flags.flag_nalloc = 0;
766 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
767 const int *type, int na,
772 /* The LJ params follow the combination rule:
773 * copy the params for the type array to the atom array.
775 for (is = 0; is < na; is += PACK_X4)
777 for (k = 0; k < PACK_X4; k++)
780 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
781 ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
786 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
787 const int *type, int na,
792 /* The LJ params follow the combination rule:
793 * copy the params for the type array to the atom array.
795 for (is = 0; is < na; is += PACK_X8)
797 for (k = 0; k < PACK_X8; k++)
800 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
801 ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
806 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
807 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t *nbat,
809 const nbnxn_search_t nbs,
813 const nbnxn_grid_t *grid;
815 for (g = 0; g < ngrid; g++)
817 grid = &nbs->grid[g];
819 /* Loop over all columns and copy and fill */
820 for (i = 0; i < grid->ncx*grid->ncy; i++)
822 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
823 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
825 copy_int_to_nbat_int(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
826 type, nbat->ntype-1, nbat->type+ash);
828 if (nbat->comb_rule != ljcrNONE)
830 if (nbat->XFormat == nbatX4)
832 copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
833 nbat->type+ash, ncz*grid->na_sc,
834 nbat->lj_comb+ash*2);
836 else if (nbat->XFormat == nbatX8)
838 copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
839 nbat->type+ash, ncz*grid->na_sc,
840 nbat->lj_comb+ash*2);
847 /* Sets the charges in nbnxn_atomdata_t *nbat */
848 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t *nbat,
850 const nbnxn_search_t nbs,
853 int g, cxy, ncz, ash, na, na_round, i, j;
855 const nbnxn_grid_t *grid;
857 for (g = 0; g < ngrid; g++)
859 grid = &nbs->grid[g];
861 /* Loop over all columns and copy and fill */
862 for (cxy = 0; cxy < grid->ncx*grid->ncy; cxy++)
864 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
865 na = grid->cxy_na[cxy];
866 na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
868 if (nbat->XFormat == nbatXYZQ)
870 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
871 for (i = 0; i < na; i++)
873 *q = charge[nbs->a[ash+i]];
876 /* Complete the partially filled last cell with zeros */
877 for (; i < na_round; i++)
886 for (i = 0; i < na; i++)
888 *q = charge[nbs->a[ash+i]];
891 /* Complete the partially filled last cell with zeros */
892 for (; i < na_round; i++)
902 /* Copies the energy group indices to a reordered and packed array */
903 static void copy_egp_to_nbat_egps(const int *a, int na, int na_round,
904 int na_c, int bit_shift,
905 const int *in, int *innb)
911 for (i = 0; i < na; i += na_c)
913 /* Store na_c energy group numbers into one int */
915 for (sa = 0; sa < na_c; sa++)
920 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
925 /* Complete the partially filled last cell with fill */
926 for (; i < na_round; i += na_c)
932 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
933 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t *nbat,
935 const nbnxn_search_t nbs,
939 const nbnxn_grid_t *grid;
941 for (g = 0; g < ngrid; g++)
943 grid = &nbs->grid[g];
945 /* Loop over all columns and copy and fill */
946 for (i = 0; i < grid->ncx*grid->ncy; i++)
948 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
949 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
951 copy_egp_to_nbat_egps(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
952 nbat->na_c, nbat->neg_2log,
953 atinfo, nbat->energrp+(ash>>grid->na_c_2log));
958 /* Sets all required atom parameter data in nbnxn_atomdata_t */
959 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
961 const nbnxn_search_t nbs,
962 const t_mdatoms *mdatoms,
967 if (locality == eatLocal)
976 nbnxn_atomdata_set_atomtypes(nbat, ngrid, nbs, mdatoms->typeA);
978 nbnxn_atomdata_set_charges(nbat, ngrid, nbs, mdatoms->chargeA);
980 if (nbat->nenergrp > 1)
982 nbnxn_atomdata_set_energygroups(nbat, ngrid, nbs, atinfo);
986 /* Copies the shift vector array to nbnxn_atomdata_t */
987 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox,
989 nbnxn_atomdata_t *nbat)
993 nbat->bDynamicBox = bDynamicBox;
994 for (i = 0; i < SHIFTS; i++)
996 copy_rvec(shift_vec[i], nbat->shift_vec[i]);
1000 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1001 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
1005 nbnxn_atomdata_t *nbat)
1028 nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
1031 nth = gmx_omp_nthreads_get(emntPairsearch);
1033 #pragma omp parallel for num_threads(nth) schedule(static)
1034 for (th = 0; th < nth; th++)
1038 for (g = g0; g < g1; g++)
1040 const nbnxn_grid_t *grid;
1041 int cxy0, cxy1, cxy;
1043 grid = &nbs->grid[g];
1045 cxy0 = (grid->ncx*grid->ncy* th +nth-1)/nth;
1046 cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
1048 for (cxy = cxy0; cxy < cxy1; cxy++)
1050 int na, ash, na_fill;
1052 na = grid->cxy_na[cxy];
1053 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1055 if (g == 0 && FillLocal)
1058 (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
1062 /* We fill only the real particle locations.
1063 * We assume the filling entries at the end have been
1064 * properly set before during ns.
1068 copy_rvec_to_nbat_real(nbs->a+ash, na, na_fill, x,
1069 nbat->XFormat, nbat->x, ash,
1077 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
1082 for (i = i0; i < i1; i++)
1089 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1091 real ** gmx_restrict src,
1099 /* The destination buffer contains data, add to it */
1100 for (i = i0; i < i1; i++)
1102 for (s = 0; s < nsrc; s++)
1104 dest[i] += src[s][i];
1110 /* The destination buffer is unitialized, set it first */
1111 for (i = i0; i < i1; i++)
1113 dest[i] = src[0][i];
1114 for (s = 1; s < nsrc; s++)
1116 dest[i] += src[s][i];
1123 nbnxn_atomdata_reduce_reals_simd(real gmx_unused * gmx_restrict dest,
1124 gmx_bool gmx_unused bDestSet,
1125 real gmx_unused ** gmx_restrict src,
1126 int gmx_unused nsrc,
1127 int gmx_unused i0, int gmx_unused i1)
1129 #ifdef GMX_NBNXN_SIMD
1130 /* The SIMD width here is actually independent of that in the kernels,
1131 * but we use the same width for simplicity (usually optimal anyhow).
1134 gmx_mm_pr dest_SSE, src_SSE;
1138 for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1140 dest_SSE = gmx_load_pr(dest+i);
1141 for (s = 0; s < nsrc; s++)
1143 src_SSE = gmx_load_pr(src[s]+i);
1144 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1146 gmx_store_pr(dest+i, dest_SSE);
1151 for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1153 dest_SSE = gmx_load_pr(src[0]+i);
1154 for (s = 1; s < nsrc; s++)
1156 src_SSE = gmx_load_pr(src[s]+i);
1157 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1159 gmx_store_pr(dest+i, dest_SSE);
1165 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1167 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1168 const nbnxn_atomdata_t *nbat,
1169 nbnxn_atomdata_output_t *out,
1180 /* Loop over all columns and copy and fill */
1181 switch (nbat->FFormat)
1189 for (a = a0; a < a1; a++)
1191 i = cell[a]*nbat->fstride;
1194 f[a][YY] += fnb[i+1];
1195 f[a][ZZ] += fnb[i+2];
1200 for (a = a0; a < a1; a++)
1202 i = cell[a]*nbat->fstride;
1204 for (fa = 0; fa < nfa; fa++)
1206 f[a][XX] += out[fa].f[i];
1207 f[a][YY] += out[fa].f[i+1];
1208 f[a][ZZ] += out[fa].f[i+2];
1218 for (a = a0; a < a1; a++)
1220 i = X4_IND_A(cell[a]);
1222 f[a][XX] += fnb[i+XX*PACK_X4];
1223 f[a][YY] += fnb[i+YY*PACK_X4];
1224 f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1229 for (a = a0; a < a1; a++)
1231 i = X4_IND_A(cell[a]);
1233 for (fa = 0; fa < nfa; fa++)
1235 f[a][XX] += out[fa].f[i+XX*PACK_X4];
1236 f[a][YY] += out[fa].f[i+YY*PACK_X4];
1237 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1247 for (a = a0; a < a1; a++)
1249 i = X8_IND_A(cell[a]);
1251 f[a][XX] += fnb[i+XX*PACK_X8];
1252 f[a][YY] += fnb[i+YY*PACK_X8];
1253 f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1258 for (a = a0; a < a1; a++)
1260 i = X8_IND_A(cell[a]);
1262 for (fa = 0; fa < nfa; fa++)
1264 f[a][XX] += out[fa].f[i+XX*PACK_X8];
1265 f[a][YY] += out[fa].f[i+YY*PACK_X8];
1266 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1272 gmx_incons("Unsupported nbnxn_atomdata_t format");
1276 /* Add the force array(s) from nbnxn_atomdata_t to f */
1277 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t nbs,
1279 const nbnxn_atomdata_t *nbat,
1285 nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1291 na = nbs->natoms_nonlocal;
1295 na = nbs->natoms_local;
1298 a0 = nbs->natoms_local;
1299 na = nbs->natoms_nonlocal - nbs->natoms_local;
1303 nth = gmx_omp_nthreads_get(emntNonbonded);
1307 if (locality != eatAll)
1309 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1312 /* Reduce the force thread output buffers into buffer 0, before adding
1313 * them to the, differently ordered, "real" force buffer.
1315 #pragma omp parallel for num_threads(nth) schedule(static)
1316 for (th = 0; th < nth; th++)
1318 const nbnxn_buffer_flags_t *flags;
1322 real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1325 flags = &nbat->buffer_flags;
1327 /* Calculate the cell-block range for our thread */
1328 b0 = (flags->nflag* th )/nth;
1329 b1 = (flags->nflag*(th+1))/nth;
1331 for (b = b0; b < b1; b++)
1333 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1334 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1337 for (out = 1; out < nbat->nout; out++)
1339 if (flags->flag[b] & (1U<<out))
1341 fptr[nfptr++] = nbat->out[out].f;
1346 #ifdef GMX_NBNXN_SIMD
1347 nbnxn_atomdata_reduce_reals_simd
1349 nbnxn_atomdata_reduce_reals
1352 flags->flag[b] & (1U<<0),
1356 else if (!(flags->flag[b] & (1U<<0)))
1358 nbnxn_atomdata_clear_reals(nbat->out[0].f,
1365 #pragma omp parallel for num_threads(nth) schedule(static)
1366 for (th = 0; th < nth; th++)
1368 nbnxn_atomdata_add_nbat_f_to_f_part(nbs, nbat,
1376 nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1379 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1380 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1383 const nbnxn_atomdata_output_t *out;
1390 for (s = 0; s < SHIFTS; s++)
1393 for (th = 0; th < nbat->nout; th++)
1395 sum[XX] += out[th].fshift[s*DIM+XX];
1396 sum[YY] += out[th].fshift[s*DIM+YY];
1397 sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1399 rvec_inc(fshift[s], sum);