1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2012, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
42 #include "nbnxn_consts.h"
43 #include "nbnxn_internal.h"
44 #include "nbnxn_search.h"
45 #include "nbnxn_atomdata.h"
46 #include "gmx_omp_nthreads.h"
48 /* Default nbnxn allocation routine, allocates 32 byte aligned,
49 * which works for plain C and aligned SSE and AVX loads/stores.
51 void nbnxn_alloc_aligned(void **ptr,size_t nbytes)
53 *ptr = save_malloc_aligned("ptr",__FILE__,__LINE__,nbytes,1,32);
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)
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)
190 xnb[(a0+a)*STRIDE_XYZ+d] = 0;
199 xnb[(a0+a)*STRIDE_XYZQ+d] = 0;
205 c = a0 & (PACK_X4-1);
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);
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
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);
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);
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);
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 stride");
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;
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]);
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 /* In nbfp_s4 we use a stride of 4 for storing two parameters */
405 nbat->alloc((void **)&nbat->nbfp_s4,nt*nt*4*sizeof(*nbat->nbfp_s4));
410 nbat->nbfp_s4[(i*nt+j)*4+0] = nbat->nbfp[(i*nt+j)*2+0];
411 nbat->nbfp_s4[(i*nt+j)*4+1] = nbat->nbfp[(i*nt+j)*2+1];
412 nbat->nbfp_s4[(i*nt+j)*4+2] = 0;
413 nbat->nbfp_s4[(i*nt+j)*4+3] = 0;
418 gmx_incons("Unknown combination rule");
423 /* Initializes an nbnxn_atomdata_t data structure */
424 void nbnxn_atomdata_init(FILE *fp,
425 nbnxn_atomdata_t *nbat,
427 int ntype,const real *nbfp,
430 nbnxn_alloc_t *alloc,
436 gmx_bool simple,bCombGeom,bCombLB;
440 nbat->alloc = nbnxn_alloc_aligned;
448 nbat->free = nbnxn_free_aligned;
457 fprintf(debug,"There are %d atom types in the system, adding one for nbnxn_atomdata_t\n",ntype);
459 nbat->ntype = ntype + 1;
460 nbat->alloc((void **)&nbat->nbfp,
461 nbat->ntype*nbat->ntype*2*sizeof(*nbat->nbfp));
462 nbat->alloc((void **)&nbat->nbfp_comb,nbat->ntype*2*sizeof(*nbat->nbfp_comb));
464 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
465 * force-field floating point parameters.
468 ptr = getenv("GMX_LJCOMB_TOL");
473 sscanf(ptr,"%lf",&dbl);
479 /* Temporarily fill nbat->nbfp_comb with sigma and epsilon
480 * to check for the LB rule.
482 for(i=0; i<ntype; i++)
484 c6 = nbfp[(i*ntype+i)*2 ]/6.0;
485 c12 = nbfp[(i*ntype+i)*2+1]/12.0;
486 if (c6 > 0 && c12 > 0)
488 nbat->nbfp_comb[i*2 ] = pow(c12/c6,1.0/6.0);
489 nbat->nbfp_comb[i*2+1] = 0.25*c6*c6/c12;
491 else if (c6 == 0 && c12 == 0)
493 nbat->nbfp_comb[i*2 ] = 0;
494 nbat->nbfp_comb[i*2+1] = 0;
498 /* Can not use LB rule with only dispersion or repulsion */
503 for(i=0; i<nbat->ntype; i++)
505 for(j=0; j<nbat->ntype; j++)
507 if (i < ntype && j < ntype)
509 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
510 * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
512 c6 = nbfp[(i*ntype+j)*2 ];
513 c12 = nbfp[(i*ntype+j)*2+1];
514 nbat->nbfp[(i*nbat->ntype+j)*2 ] = c6;
515 nbat->nbfp[(i*nbat->ntype+j)*2+1] = c12;
517 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
518 bCombGeom = bCombGeom &&
519 gmx_within_tol(c6*c6 ,nbfp[(i*ntype+i)*2 ]*nbfp[(j*ntype+j)*2 ],tol) &&
520 gmx_within_tol(c12*c12,nbfp[(i*ntype+i)*2+1]*nbfp[(j*ntype+j)*2+1],tol);
522 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
526 ((c6 == 0 && c12 == 0 &&
527 (nbat->nbfp_comb[i*2+1] == 0 || nbat->nbfp_comb[j*2+1] == 0)) ||
528 (c6 > 0 && c12 > 0 &&
529 gmx_within_tol(pow(c12/c6,1.0/6.0),0.5*(nbat->nbfp_comb[i*2]+nbat->nbfp_comb[j*2]),tol) &&
530 gmx_within_tol(0.25*c6*c6/c12,sqrt(nbat->nbfp_comb[i*2+1]*nbat->nbfp_comb[j*2+1]),tol)));
534 /* Add zero parameters for the additional dummy atom type */
535 nbat->nbfp[(i*nbat->ntype+j)*2 ] = 0;
536 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 0;
542 fprintf(debug,"Combination rules: geometric %d Lorentz-Berthelot %d\n",
546 simple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
550 /* We prefer the geometic combination rule,
551 * as that gives a slightly faster kernel than the LB rule.
555 nbat->comb_rule = ljcrGEOM;
559 nbat->comb_rule = ljcrLB;
563 nbat->comb_rule = ljcrNONE;
565 nbat->free(nbat->nbfp_comb);
570 if (nbat->comb_rule == ljcrNONE)
572 fprintf(fp,"Using full Lennard-Jones parameter combination matrix\n\n");
576 fprintf(fp,"Using %s Lennard-Jones combination rule\n\n",
577 nbat->comb_rule==ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
581 set_combination_rule_data(nbat);
585 nbat->comb_rule = ljcrNONE;
587 nbat->free(nbat->nbfp_comb);
592 nbat->lj_comb = NULL;
597 switch (nb_kernel_type)
599 case nbnxnk4xN_SIMD_4xN:
600 case nbnxnk4xN_SIMD_2xNN:
601 pack_x = max(NBNXN_CPU_CLUSTER_I_SIZE,
602 nbnxn_kernel_to_cj_size(nb_kernel_type));
606 nbat->XFormat = nbatX4;
609 nbat->XFormat = nbatX8;
612 gmx_incons("Unsupported packing width");
616 nbat->XFormat = nbatXYZ;
620 nbat->FFormat = nbat->XFormat;
624 nbat->XFormat = nbatXYZQ;
625 nbat->FFormat = nbatXYZ;
628 nbat->nenergrp = n_energygroups;
631 /* Energy groups not supported yet for super-sub lists */
632 if (n_energygroups > 1 && fp != NULL)
634 fprintf(fp,"\nNOTE: With GPUs, reporting energy group contributions is not supported\n\n");
638 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
639 if (nbat->nenergrp > 64)
641 gmx_fatal(FARGS,"With NxN kernels not more than 64 energy groups are supported\n");
644 while (nbat->nenergrp > (1<<nbat->neg_2log))
648 nbat->energrp = NULL;
649 nbat->alloc((void **)&nbat->shift_vec,SHIFTS*sizeof(*nbat->shift_vec));
650 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
651 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
654 snew(nbat->out,nbat->nout);
656 for(i=0; i<nbat->nout; i++)
658 nbnxn_atomdata_output_init(&nbat->out[i],
660 nbat->nenergrp,1<<nbat->neg_2log,
663 nbat->buffer_flags.flag = NULL;
664 nbat->buffer_flags.flag_nalloc = 0;
667 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
668 const int *type,int na,
673 /* The LJ params follow the combination rule:
674 * copy the params for the type array to the atom array.
676 for(is=0; is<na; is+=PACK_X4)
678 for(k=0; k<PACK_X4; k++)
681 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
682 ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
687 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
688 const int *type,int na,
693 /* The LJ params follow the combination rule:
694 * copy the params for the type array to the atom array.
696 for(is=0; is<na; is+=PACK_X8)
698 for(k=0; k<PACK_X8; k++)
701 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
702 ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
707 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
708 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t *nbat,
710 const nbnxn_search_t nbs,
714 const nbnxn_grid_t *grid;
716 for(g=0; g<ngrid; g++)
718 grid = &nbs->grid[g];
720 /* Loop over all columns and copy and fill */
721 for(i=0; i<grid->ncx*grid->ncy; i++)
723 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
724 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
726 copy_int_to_nbat_int(nbs->a+ash,grid->cxy_na[i],ncz*grid->na_sc,
727 type,nbat->ntype-1,nbat->type+ash);
729 if (nbat->comb_rule != ljcrNONE)
731 if (nbat->XFormat == nbatX4)
733 copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
734 nbat->type+ash,ncz*grid->na_sc,
735 nbat->lj_comb+ash*2);
737 else if (nbat->XFormat == nbatX8)
739 copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
740 nbat->type+ash,ncz*grid->na_sc,
741 nbat->lj_comb+ash*2);
748 /* Sets the charges in nbnxn_atomdata_t *nbat */
749 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t *nbat,
751 const nbnxn_search_t nbs,
754 int g,cxy,ncz,ash,na,na_round,i,j;
756 const nbnxn_grid_t *grid;
758 for(g=0; g<ngrid; g++)
760 grid = &nbs->grid[g];
762 /* Loop over all columns and copy and fill */
763 for(cxy=0; cxy<grid->ncx*grid->ncy; cxy++)
765 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
766 na = grid->cxy_na[cxy];
767 na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
769 if (nbat->XFormat == nbatXYZQ)
771 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
774 *q = charge[nbs->a[ash+i]];
777 /* Complete the partially filled last cell with zeros */
778 for(; i<na_round; i++)
789 *q = charge[nbs->a[ash+i]];
792 /* Complete the partially filled last cell with zeros */
793 for(; i<na_round; i++)
803 /* Copies the energy group indices to a reordered and packed array */
804 static void copy_egp_to_nbat_egps(const int *a,int na,int na_round,
805 int na_c,int bit_shift,
806 const int *in,int *innb)
812 for(i=0; i<na; i+=na_c)
814 /* Store na_c energy group numbers into one int */
816 for(sa=0; sa<na_c; sa++)
821 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
826 /* Complete the partially filled last cell with fill */
827 for(; i<na_round; i+=na_c)
833 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
834 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t *nbat,
836 const nbnxn_search_t nbs,
840 const nbnxn_grid_t *grid;
842 for(g=0; g<ngrid; g++)
844 grid = &nbs->grid[g];
846 /* Loop over all columns and copy and fill */
847 for(i=0; i<grid->ncx*grid->ncy; i++)
849 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
850 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
852 copy_egp_to_nbat_egps(nbs->a+ash,grid->cxy_na[i],ncz*grid->na_sc,
853 nbat->na_c,nbat->neg_2log,
854 atinfo,nbat->energrp+(ash>>grid->na_c_2log));
859 /* Sets all required atom parameter data in nbnxn_atomdata_t */
860 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
862 const nbnxn_search_t nbs,
863 const t_mdatoms *mdatoms,
868 if (locality == eatLocal)
877 nbnxn_atomdata_set_atomtypes(nbat,ngrid,nbs,mdatoms->typeA);
879 nbnxn_atomdata_set_charges(nbat,ngrid,nbs,mdatoms->chargeA);
881 if (nbat->nenergrp > 1)
883 nbnxn_atomdata_set_energygroups(nbat,ngrid,nbs,atinfo);
887 /* Copies the shift vector array to nbnxn_atomdata_t */
888 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox,
890 nbnxn_atomdata_t *nbat)
894 nbat->bDynamicBox = bDynamicBox;
895 for(i=0; i<SHIFTS; i++)
897 copy_rvec(shift_vec[i],nbat->shift_vec[i]);
901 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
902 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
906 nbnxn_atomdata_t *nbat)
929 nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
932 nth = gmx_omp_nthreads_get(emntPairsearch);
934 #pragma omp parallel for num_threads(nth) schedule(static)
935 for(th=0; th<nth; th++)
941 const nbnxn_grid_t *grid;
944 grid = &nbs->grid[g];
946 cxy0 = (grid->ncx*grid->ncy* th +nth-1)/nth;
947 cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
949 for(cxy=cxy0; cxy<cxy1; cxy++)
953 na = grid->cxy_na[cxy];
954 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
956 if (g == 0 && FillLocal)
959 (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
963 /* We fill only the real particle locations.
964 * We assume the filling entries at the end have been
965 * properly set before during ns.
969 copy_rvec_to_nbat_real(nbs->a+ash,na,na_fill,x,
970 nbat->XFormat,nbat->x,ash,
978 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
990 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
992 real ** gmx_restrict src,
1000 /* The destination buffer contains data, add to it */
1001 for(i=i0; i<i1; i++)
1003 for(s=0; s<nsrc; s++)
1005 dest[i] += src[s][i];
1011 /* The destination buffer is unitialized, set it first */
1012 for(i=i0; i<i1; i++)
1014 dest[i] = src[0][i];
1015 for(s=1; s<nsrc; s++)
1017 dest[i] += src[s][i];
1024 nbnxn_atomdata_reduce_reals_x86_simd(real * gmx_restrict dest,
1026 real ** gmx_restrict src,
1030 #ifdef NBNXN_SEARCH_SSE
1031 /* We can use AVX256 here, but not when AVX128 kernels are selected.
1032 * As this reduction is not faster with AVX256 anyway, we use 128-bit SIMD.
1034 #ifdef GMX_X86_AVX_256
1035 #define GMX_MM256_HERE
1037 #define GMX_MM128_HERE
1039 #include "gmx_simd_macros.h"
1042 gmx_mm_pr dest_SSE,src_SSE;
1046 for(i=i0; i<i1; i+=GMX_SIMD_WIDTH_HERE)
1048 dest_SSE = gmx_load_pr(dest+i);
1049 for(s=0; s<nsrc; s++)
1051 src_SSE = gmx_load_pr(src[s]+i);
1052 dest_SSE = gmx_add_pr(dest_SSE,src_SSE);
1054 gmx_store_pr(dest+i,dest_SSE);
1059 for(i=i0; i<i1; i+=GMX_SIMD_WIDTH_HERE)
1061 dest_SSE = gmx_load_pr(src[0]+i);
1062 for(s=1; s<nsrc; s++)
1064 src_SSE = gmx_load_pr(src[s]+i);
1065 dest_SSE = gmx_add_pr(dest_SSE,src_SSE);
1067 gmx_store_pr(dest+i,dest_SSE);
1071 #undef GMX_MM128_HERE
1072 #undef GMX_MM256_HERE
1076 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1078 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1079 const nbnxn_atomdata_t *nbat,
1080 nbnxn_atomdata_output_t *out,
1091 /* Loop over all columns and copy and fill */
1092 switch (nbat->FFormat)
1100 for(a=a0; a<a1; a++)
1102 i = cell[a]*nbat->fstride;
1105 f[a][YY] += fnb[i+1];
1106 f[a][ZZ] += fnb[i+2];
1111 for(a=a0; a<a1; a++)
1113 i = cell[a]*nbat->fstride;
1115 for(fa=0; fa<nfa; fa++)
1117 f[a][XX] += out[fa].f[i];
1118 f[a][YY] += out[fa].f[i+1];
1119 f[a][ZZ] += out[fa].f[i+2];
1129 for(a=a0; a<a1; a++)
1131 i = X4_IND_A(cell[a]);
1133 f[a][XX] += fnb[i+XX*PACK_X4];
1134 f[a][YY] += fnb[i+YY*PACK_X4];
1135 f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1140 for(a=a0; a<a1; a++)
1142 i = X4_IND_A(cell[a]);
1144 for(fa=0; fa<nfa; fa++)
1146 f[a][XX] += out[fa].f[i+XX*PACK_X4];
1147 f[a][YY] += out[fa].f[i+YY*PACK_X4];
1148 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1158 for(a=a0; a<a1; a++)
1160 i = X8_IND_A(cell[a]);
1162 f[a][XX] += fnb[i+XX*PACK_X8];
1163 f[a][YY] += fnb[i+YY*PACK_X8];
1164 f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1169 for(a=a0; a<a1; a++)
1171 i = X8_IND_A(cell[a]);
1173 for(fa=0; fa<nfa; fa++)
1175 f[a][XX] += out[fa].f[i+XX*PACK_X8];
1176 f[a][YY] += out[fa].f[i+YY*PACK_X8];
1177 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1185 /* Add the force array(s) from nbnxn_atomdata_t to f */
1186 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t nbs,
1188 const nbnxn_atomdata_t *nbat,
1194 nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1200 na = nbs->natoms_nonlocal;
1204 na = nbs->natoms_local;
1207 a0 = nbs->natoms_local;
1208 na = nbs->natoms_nonlocal - nbs->natoms_local;
1212 nth = gmx_omp_nthreads_get(emntNonbonded);
1216 if (locality != eatAll)
1218 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1221 /* Reduce the force thread output buffers into buffer 0, before adding
1222 * them to the, differently ordered, "real" force buffer.
1224 #pragma omp parallel for num_threads(nth) schedule(static)
1225 for(th=0; th<nth; th++)
1227 const nbnxn_buffer_flags_t *flags;
1231 real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1234 flags = &nbat->buffer_flags;
1236 /* Calculate the cell-block range for our thread */
1237 b0 = (flags->nflag* th )/nth;
1238 b1 = (flags->nflag*(th+1))/nth;
1240 for(b=b0; b<b1; b++)
1242 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1243 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1246 for(out=1; out<nbat->nout; out++)
1248 if (flags->flag[b] & (1U<<out))
1250 fptr[nfptr++] = nbat->out[out].f;
1255 #ifdef NBNXN_SEARCH_SSE
1256 nbnxn_atomdata_reduce_reals_x86_simd
1258 nbnxn_atomdata_reduce_reals
1261 flags->flag[b] & (1U<<0),
1265 else if (!(flags->flag[b] & (1U<<0)))
1267 nbnxn_atomdata_clear_reals(nbat->out[0].f,
1274 #pragma omp parallel for num_threads(nth) schedule(static)
1275 for(th=0; th<nth; th++)
1277 nbnxn_atomdata_add_nbat_f_to_f_part(nbs,nbat,
1285 nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1288 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1289 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1292 const nbnxn_atomdata_output_t *out;
1299 for(s=0; s<SHIFTS; s++)
1302 for(th=0; th<nbat->nout; th++)
1304 sum[XX] += out[th].fshift[s*DIM+XX];
1305 sum[YY] += out[th].fshift[s*DIM+YY];
1306 sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1308 rvec_inc(fshift[s],sum);