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 == nbk4xN_X86_SIMD128 ||
149 nb_kernel_type == nbk4xN_X86_SIMD256)
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;
595 switch (nb_kernel_type)
597 case nbk4xN_X86_SIMD128:
598 nbat->XFormat = nbatX4;
600 case nbk4xN_X86_SIMD256:
602 nbat->XFormat = nbatX8;
604 nbat->XFormat = nbatX4;
608 nbat->XFormat = nbatXYZ;
612 nbat->FFormat = nbat->XFormat;
616 nbat->XFormat = nbatXYZQ;
617 nbat->FFormat = nbatXYZ;
620 nbat->nenergrp = n_energygroups;
623 /* Energy groups not supported yet for super-sub lists */
624 if (n_energygroups > 1 && fp != NULL)
626 fprintf(fp,"\nNOTE: With GPUs, reporting energy group contributions is not supported\n\n");
630 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
631 if (nbat->nenergrp > 64)
633 gmx_fatal(FARGS,"With NxN kernels not more than 64 energy groups are supported\n");
636 while (nbat->nenergrp > (1<<nbat->neg_2log))
640 nbat->energrp = NULL;
641 nbat->alloc((void **)&nbat->shift_vec,SHIFTS*sizeof(*nbat->shift_vec));
642 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
643 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
646 snew(nbat->out,nbat->nout);
648 for(i=0; i<nbat->nout; i++)
650 nbnxn_atomdata_output_init(&nbat->out[i],
652 nbat->nenergrp,1<<nbat->neg_2log,
655 nbat->buffer_flags.flag = NULL;
656 nbat->buffer_flags.flag_nalloc = 0;
659 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
660 const int *type,int na,
665 /* The LJ params follow the combination rule:
666 * copy the params for the type array to the atom array.
668 for(is=0; is<na; is+=PACK_X4)
670 for(k=0; k<PACK_X4; k++)
673 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
674 ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
679 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
680 const int *type,int na,
685 /* The LJ params follow the combination rule:
686 * copy the params for the type array to the atom array.
688 for(is=0; is<na; is+=PACK_X8)
690 for(k=0; k<PACK_X8; k++)
693 ljparam_at[is*2 +k] = ljparam_type[type[i]*2 ];
694 ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
699 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
700 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t *nbat,
702 const nbnxn_search_t nbs,
706 const nbnxn_grid_t *grid;
708 for(g=0; g<ngrid; g++)
710 grid = &nbs->grid[g];
712 /* Loop over all columns and copy and fill */
713 for(i=0; i<grid->ncx*grid->ncy; i++)
715 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
716 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
718 copy_int_to_nbat_int(nbs->a+ash,grid->cxy_na[i],ncz*grid->na_sc,
719 type,nbat->ntype-1,nbat->type+ash);
721 if (nbat->comb_rule != ljcrNONE)
723 if (nbat->XFormat == nbatX4)
725 copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
726 nbat->type+ash,ncz*grid->na_sc,
727 nbat->lj_comb+ash*2);
729 else if (nbat->XFormat == nbatX8)
731 copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
732 nbat->type+ash,ncz*grid->na_sc,
733 nbat->lj_comb+ash*2);
740 /* Sets the charges in nbnxn_atomdata_t *nbat */
741 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t *nbat,
743 const nbnxn_search_t nbs,
746 int g,cxy,ncz,ash,na,na_round,i,j;
748 const nbnxn_grid_t *grid;
750 for(g=0; g<ngrid; g++)
752 grid = &nbs->grid[g];
754 /* Loop over all columns and copy and fill */
755 for(cxy=0; cxy<grid->ncx*grid->ncy; cxy++)
757 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
758 na = grid->cxy_na[cxy];
759 na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
761 if (nbat->XFormat == nbatXYZQ)
763 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
766 *q = charge[nbs->a[ash+i]];
769 /* Complete the partially filled last cell with zeros */
770 for(; i<na_round; i++)
781 *q = charge[nbs->a[ash+i]];
784 /* Complete the partially filled last cell with zeros */
785 for(; i<na_round; i++)
795 /* Copies the energy group indices to a reordered and packed array */
796 static void copy_egp_to_nbat_egps(const int *a,int na,int na_round,
797 int na_c,int bit_shift,
798 const int *in,int *innb)
804 for(i=0; i<na; i+=na_c)
806 /* Store na_c energy group numbers into one int */
808 for(sa=0; sa<na_c; sa++)
813 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
818 /* Complete the partially filled last cell with fill */
819 for(; i<na_round; i+=na_c)
825 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
826 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t *nbat,
828 const nbnxn_search_t nbs,
832 const nbnxn_grid_t *grid;
834 for(g=0; g<ngrid; g++)
836 grid = &nbs->grid[g];
838 /* Loop over all columns and copy and fill */
839 for(i=0; i<grid->ncx*grid->ncy; i++)
841 ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
842 ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
844 copy_egp_to_nbat_egps(nbs->a+ash,grid->cxy_na[i],ncz*grid->na_sc,
845 nbat->na_c,nbat->neg_2log,
846 atinfo,nbat->energrp+(ash>>grid->na_c_2log));
851 /* Sets all required atom parameter data in nbnxn_atomdata_t */
852 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
854 const nbnxn_search_t nbs,
855 const t_mdatoms *mdatoms,
860 if (locality == eatLocal)
869 nbnxn_atomdata_set_atomtypes(nbat,ngrid,nbs,mdatoms->typeA);
871 nbnxn_atomdata_set_charges(nbat,ngrid,nbs,mdatoms->chargeA);
873 if (nbat->nenergrp > 1)
875 nbnxn_atomdata_set_energygroups(nbat,ngrid,nbs,atinfo);
879 /* Copies the shift vector array to nbnxn_atomdata_t */
880 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox,
882 nbnxn_atomdata_t *nbat)
886 nbat->bDynamicBox = bDynamicBox;
887 for(i=0; i<SHIFTS; i++)
889 copy_rvec(shift_vec[i],nbat->shift_vec[i]);
893 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
894 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
898 nbnxn_atomdata_t *nbat)
921 nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
924 nth = gmx_omp_nthreads_get(emntPairsearch);
926 #pragma omp parallel for num_threads(nth) schedule(static)
927 for(th=0; th<nth; th++)
933 const nbnxn_grid_t *grid;
936 grid = &nbs->grid[g];
938 cxy0 = (grid->ncx*grid->ncy* th +nth-1)/nth;
939 cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
941 for(cxy=cxy0; cxy<cxy1; cxy++)
945 na = grid->cxy_na[cxy];
946 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
948 if (g == 0 && FillLocal)
951 (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
955 /* We fill only the real particle locations.
956 * We assume the filling entries at the end have been
957 * properly set before during ns.
961 copy_rvec_to_nbat_real(nbs->a+ash,na,na_fill,x,
962 nbat->XFormat,nbat->x,ash,
970 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
982 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
984 real ** gmx_restrict src,
992 /* The destination buffer contains data, add to it */
995 for(s=0; s<nsrc; s++)
997 dest[i] += src[s][i];
1003 /* The destination buffer is unitialized, set it first */
1004 for(i=i0; i<i1; i++)
1006 dest[i] = src[0][i];
1007 for(s=1; s<nsrc; s++)
1009 dest[i] += src[s][i];
1016 nbnxn_atomdata_reduce_reals_x86_simd(real * gmx_restrict dest,
1018 real ** gmx_restrict src,
1022 #ifdef NBNXN_SEARCH_SSE
1023 /* We can use AVX256 here, but not when AVX128 kernels are selected.
1024 * As this reduction is not faster with AVX256 anyway, we use 128-bit SIMD.
1026 #ifdef GMX_X86_AVX_256
1027 #define GMX_MM256_HERE
1029 #define GMX_MM128_HERE
1031 #include "gmx_x86_simd_macros.h"
1034 gmx_mm_pr dest_SSE,src_SSE;
1038 for(i=i0; i<i1; i+=GMX_X86_SIMD_WIDTH_HERE)
1040 dest_SSE = gmx_load_pr(dest+i);
1041 for(s=0; s<nsrc; s++)
1043 src_SSE = gmx_load_pr(src[s]+i);
1044 dest_SSE = gmx_add_pr(dest_SSE,src_SSE);
1046 gmx_store_pr(dest+i,dest_SSE);
1051 for(i=i0; i<i1; i+=GMX_X86_SIMD_WIDTH_HERE)
1053 dest_SSE = gmx_load_pr(src[0]+i);
1054 for(s=1; s<nsrc; s++)
1056 src_SSE = gmx_load_pr(src[s]+i);
1057 dest_SSE = gmx_add_pr(dest_SSE,src_SSE);
1059 gmx_store_pr(dest+i,dest_SSE);
1063 #undef GMX_MM128_HERE
1064 #undef GMX_MM256_HERE
1068 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1070 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1071 const nbnxn_atomdata_t *nbat,
1072 nbnxn_atomdata_output_t *out,
1083 /* Loop over all columns and copy and fill */
1084 switch (nbat->FFormat)
1092 for(a=a0; a<a1; a++)
1094 i = cell[a]*nbat->fstride;
1097 f[a][YY] += fnb[i+1];
1098 f[a][ZZ] += fnb[i+2];
1103 for(a=a0; a<a1; a++)
1105 i = cell[a]*nbat->fstride;
1107 for(fa=0; fa<nfa; fa++)
1109 f[a][XX] += out[fa].f[i];
1110 f[a][YY] += out[fa].f[i+1];
1111 f[a][ZZ] += out[fa].f[i+2];
1121 for(a=a0; a<a1; a++)
1123 i = X4_IND_A(cell[a]);
1125 f[a][XX] += fnb[i+XX*PACK_X4];
1126 f[a][YY] += fnb[i+YY*PACK_X4];
1127 f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1132 for(a=a0; a<a1; a++)
1134 i = X4_IND_A(cell[a]);
1136 for(fa=0; fa<nfa; fa++)
1138 f[a][XX] += out[fa].f[i+XX*PACK_X4];
1139 f[a][YY] += out[fa].f[i+YY*PACK_X4];
1140 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1150 for(a=a0; a<a1; a++)
1152 i = X8_IND_A(cell[a]);
1154 f[a][XX] += fnb[i+XX*PACK_X8];
1155 f[a][YY] += fnb[i+YY*PACK_X8];
1156 f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1161 for(a=a0; a<a1; a++)
1163 i = X8_IND_A(cell[a]);
1165 for(fa=0; fa<nfa; fa++)
1167 f[a][XX] += out[fa].f[i+XX*PACK_X8];
1168 f[a][YY] += out[fa].f[i+YY*PACK_X8];
1169 f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1177 /* Add the force array(s) from nbnxn_atomdata_t to f */
1178 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t nbs,
1180 const nbnxn_atomdata_t *nbat,
1186 nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1192 na = nbs->natoms_nonlocal;
1196 na = nbs->natoms_local;
1199 a0 = nbs->natoms_local;
1200 na = nbs->natoms_nonlocal - nbs->natoms_local;
1204 nth = gmx_omp_nthreads_get(emntNonbonded);
1208 if (locality != eatAll)
1210 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1213 /* Reduce the force thread output buffers into buffer 0, before adding
1214 * them to the, differently ordered, "real" force buffer.
1216 #pragma omp parallel for num_threads(nth) schedule(static)
1217 for(th=0; th<nth; th++)
1219 const nbnxn_buffer_flags_t *flags;
1223 real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1226 flags = &nbat->buffer_flags;
1228 /* Calculate the cell-block range for our thread */
1229 b0 = (flags->nflag* th )/nth;
1230 b1 = (flags->nflag*(th+1))/nth;
1232 for(b=b0; b<b1; b++)
1234 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1235 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1238 for(out=1; out<nbat->nout; out++)
1240 if (flags->flag[b] & (1U<<out))
1242 fptr[nfptr++] = nbat->out[out].f;
1247 #ifdef NBNXN_SEARCH_SSE
1248 nbnxn_atomdata_reduce_reals_x86_simd
1250 nbnxn_atomdata_reduce_reals
1253 flags->flag[b] & (1U<<0),
1257 else if (!(flags->flag[b] & (1U<<0)))
1259 nbnxn_atomdata_clear_reals(nbat->out[0].f,
1266 #pragma omp parallel for num_threads(nth) schedule(static)
1267 for(th=0; th<nth; th++)
1269 nbnxn_atomdata_add_nbat_f_to_f_part(nbs,nbat,
1277 nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1280 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1281 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1284 const nbnxn_atomdata_output_t *out;
1291 for(s=0; s<SHIFTS; s++)
1294 for(th=0; th<nbat->nout; th++)
1296 sum[XX] += out[th].fshift[s*DIM+XX];
1297 sum[YY] += out[th].fshift[s*DIM+YY];
1298 sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1300 rvec_inc(fshift[s],sum);