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.
41 #include "types/simple.h"
42 #include "gromacs/math/utilities.h"
46 #include "nbnxn_kernel_gpu_ref.h"
47 #include "../nbnxn_consts.h"
48 #include "nbnxn_kernel_common.h"
50 #define NCL_PER_SUPERCL (NBNXN_GPU_NCLUSTER_PER_SUPERCLUSTER)
51 #define CL_SIZE (NBNXN_GPU_CLUSTER_SIZE)
54 nbnxn_kernel_gpu_ref(const nbnxn_pairlist_t *nbl,
55 const nbnxn_atomdata_t *nbat,
56 const interaction_const_t *iconst,
65 const nbnxn_sci_t *nbln;
69 const real *Ftab = NULL;
70 real rcut2, rvdw2, rlist2;
76 int cj4_ind0, cj4_ind1, cj4_ind;
78 int ic, jc, ia, ja, is, ifs, js, jfs, im, jm;
82 real fscal, tx, ty, tz;
85 real qq, vcoul = 0, krsq, vctot;
91 real Vvdw_rep, Vvdw_disp;
92 real ix, iy, iz, fix, fiy, fiz;
94 real dx, dy, dz, rsq, rinv;
97 real c6, c12, cexp1, cexp2, br;
98 const real * shiftvec;
102 const nbnxn_excl_t *excl[2];
104 int npair_tot, npair;
105 int nhwu, nhwu_pruned;
107 if (nbl->na_ci != CL_SIZE)
109 gmx_fatal(FARGS, "The neighborlist cluster size in the GPU reference kernel is %d, expected it to be %d", nbl->na_ci, CL_SIZE);
112 if (clearF == enbvClearFYes)
117 bEner = (force_flags & GMX_FORCE_ENERGY);
119 bEwald = EEL_FULL(iconst->eeltype);
122 Ftab = iconst->tabq_coul_F;
125 rcut2 = iconst->rcoulomb*iconst->rcoulomb;
126 rvdw2 = iconst->rvdw*iconst->rvdw;
128 rlist2 = nbl->rlist*nbl->rlist;
131 facel = iconst->epsfac;
132 shiftvec = shift_vec[0];
133 vdwparam = nbat->nbfp;
142 for (n = 0; n < nbl->nsci; n++)
146 ish3 = 3*nbln->shift;
147 shX = shiftvec[ish3];
148 shY = shiftvec[ish3+1];
149 shZ = shiftvec[ish3+2];
150 cj4_ind0 = nbln->cj4_ind_start;
151 cj4_ind1 = nbln->cj4_ind_end;
156 if (nbln->shift == CENTRAL &&
157 nbl->cj4[cj4_ind0].cj[0] == sci*NCL_PER_SUPERCL)
159 /* we have the diagonal:
160 * add the charge self interaction energy term
162 for (im = 0; im < NCL_PER_SUPERCL; im++)
164 ci = sci*NCL_PER_SUPERCL + im;
165 for (ic = 0; ic < CL_SIZE; ic++)
167 ia = ci*CL_SIZE + ic;
168 iq = x[ia*nbat->xstride+3];
174 vctot *= -facel*0.5*iconst->c_rf;
178 /* last factor 1/sqrt(pi) */
179 vctot *= -facel*iconst->ewaldcoeff_q*M_1_SQRTPI;
183 for (cj4_ind = cj4_ind0; (cj4_ind < cj4_ind1); cj4_ind++)
185 excl[0] = &nbl->excl[nbl->cj4[cj4_ind].imei[0].excl_ind];
186 excl[1] = &nbl->excl[nbl->cj4[cj4_ind].imei[1].excl_ind];
188 for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
190 cj = nbl->cj4[cj4_ind].cj[jm];
192 for (im = 0; im < NCL_PER_SUPERCL; im++)
194 /* We're only using the first imask,
195 * but here imei[1].imask is identical.
197 if ((nbl->cj4[cj4_ind].imei[0].imask >> (jm*NCL_PER_SUPERCL+im)) & 1)
199 gmx_bool within_rlist;
201 ci = sci*NCL_PER_SUPERCL + im;
203 within_rlist = FALSE;
205 for (ic = 0; ic < CL_SIZE; ic++)
207 ia = ci*CL_SIZE + ic;
209 is = ia*nbat->xstride;
210 ifs = ia*nbat->fstride;
215 nti = ntype*2*type[ia];
221 for (jc = 0; jc < CL_SIZE; jc++)
223 ja = cj*CL_SIZE + jc;
225 if (nbln->shift == CENTRAL &&
226 ci == cj && ja <= ia)
231 int_bit = ((excl[jc>>2]->pair[(jc & 3)*CL_SIZE+ic] >> (jm*NCL_PER_SUPERCL+im)) & 1);
233 js = ja*nbat->xstride;
234 jfs = ja*nbat->fstride;
241 rsq = dx*dx + dy*dy + dz*dz;
251 if (type[ia] != ntype-1 && type[ja] != ntype-1)
256 /* avoid NaN for excluded pairs at r=0 */
257 rsq += (1.0 - int_bit)*NBNXN_AVOID_SING_R2_INC;
259 rinv = gmx_invsqrt(rsq);
267 krsq = iconst->k_rf*rsq;
268 fscal = qq*(int_bit*rinv - 2*krsq)*rinvsq;
271 vcoul = qq*(int_bit*rinv + krsq - iconst->c_rf);
277 rt = r*iconst->tabq_scale;
281 fexcl = (1 - eps)*Ftab[n0] + eps*Ftab[n0+1];
283 fscal = qq*(int_bit*rinvsq - fexcl)*rinv;
287 vcoul = qq*((int_bit - gmx_erf(iconst->ewaldcoeff_q*r))*rinv - int_bit*iconst->sh_ewald);
293 tj = nti + 2*type[ja];
295 /* Vanilla Lennard-Jones cutoff */
297 c12 = vdwparam[tj+1];
299 rinvsix = int_bit*rinvsq*rinvsq*rinvsq;
300 Vvdw_disp = c6*rinvsix;
301 Vvdw_rep = c12*rinvsix*rinvsix;
302 fscal += (Vvdw_rep - Vvdw_disp)*rinvsq;
309 (Vvdw_rep - int_bit*c12*iconst->sh_invrc6*iconst->sh_invrc6)/12 -
310 (Vvdw_disp - int_bit*c6*iconst->sh_invrc6)/6;
328 fshift[ish3] = fshift[ish3] + fix;
329 fshift[ish3+1] = fshift[ish3+1] + fiy;
330 fshift[ish3+2] = fshift[ish3+2] + fiz;
332 /* Count in half work-units.
333 * In CUDA one work-unit is 2 warps.
335 if ((ic+1) % (CL_SIZE/2) == 0)
345 within_rlist = FALSE;
357 Vc[ggid] = Vc[ggid] + vctot;
358 Vvdw[ggid] = Vvdw[ggid] + Vvdwtot;
364 fprintf(debug, "number of half %dx%d atom pairs: %d after pruning: %d fraction %4.2f\n",
365 nbl->na_ci, nbl->na_ci,
366 nhwu, nhwu_pruned, nhwu_pruned/(double)nhwu);
367 fprintf(debug, "generic kernel pair interactions: %d\n",
368 nhwu*nbl->na_ci/2*nbl->na_ci);
369 fprintf(debug, "generic kernel post-prune pair interactions: %d\n",
370 nhwu_pruned*nbl->na_ci/2*nbl->na_ci);
371 fprintf(debug, "generic kernel non-zero pair interactions: %d\n",
373 fprintf(debug, "ratio non-zero/post-prune pair interactions: %4.2f\n",
374 npair_tot/(double)(nhwu_pruned*nbl->na_ci/2*nbl->na_ci));