3 * This source code is part of
7 * 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-2004, 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
33 * GROningen Mixture of Alchemy and Childrens' Stories
41 #include "types/simple.h"
44 #include "nb_generic.h"
47 gmx_nb_generic_kernel(t_nblist * nlist,
60 int nri,ntype,table_nelements,icoul,ivdw;
61 real facel,gbtabscale;
62 int n,ii,is3,ii3,k,nj0,nj1,jnr,j3,ggid,nnn,n0;
67 real qq,vcoul,krsq,vctot;
70 real rt,r,eps,eps2,Y,F,Geps,Heps2,VV,FF,Fp,fijD,fijR;
73 real Vvdw_rep,Vvdw_disp;
74 real ix,iy,iz,fix,fiy,fiz;
76 real dx,dy,dz,rsq,rinv;
77 real c6,c12,cexp1,cexp2,br;
87 /* avoid compiler warnings for cases that cannot happen */
93 /* 3 VdW parameters for buckingham, otherwise 2 */
94 nvdwparam = (nlist->ivdw==2) ? 3 : 2;
95 table_nelements = (icoul==3) ? 4 : 0;
96 table_nelements += (ivdw==3) ? 8 : 0;
98 charge = mdatoms->chargeA;
99 type = mdatoms->typeA;
101 shiftvec = fr->shift_vec[0];
105 for(n=0; (n<nlist->nri); n++)
107 is3 = 3*nlist->shift[n];
109 shY = shiftvec[is3+1];
110 shZ = shiftvec[is3+2];
111 nj0 = nlist->jindex[n];
112 nj1 = nlist->jindex[n+1];
118 iq = facel*charge[ii];
119 nti = nvdwparam*ntype*type[ii];
126 for(k=nj0; (k<nj1); k++)
128 jnr = nlist->jjnr[k];
136 rsq = dx*dx+dy*dy+dz*dz;
137 rinv = gmx_invsqrt(rsq);
141 if(icoul==3 || ivdw==3)
148 nnn = table_nelements*n0;
151 /* Coulomb interaction. icoul==0 means no interaction */
159 /* Vanilla cutoff coulomb */
161 fscal = vcoul*rinvsq;
167 vcoul = qq*(rinv+krsq-fr->c_rf);
168 fscal = qq*(rinv-2.0*krsq)*rinvsq;
172 /* Tabulated coulomb */
175 Geps = eps*VFtab[nnn+2];
176 Heps2 = eps2*VFtab[nnn+3];
180 FF = Fp+Geps+2.0*Heps2;
182 fscal = -qq*FF*tabscale*rinv;
187 gmx_fatal(FARGS,"Death & horror! GB generic interaction not implemented.\n");
191 gmx_fatal(FARGS,"Death & horror! No generic coulomb interaction for icoul=%d.\n",icoul);
195 } /* End of coulomb interactions */
198 /* VdW interaction. ivdw==0 means no interaction */
201 tj = nti+nvdwparam*type[jnr];
206 /* Vanilla Lennard-Jones cutoff */
208 c12 = vdwparam[tj+1];
210 rinvsix = rinvsq*rinvsq*rinvsq;
211 Vvdw_disp = c6*rinvsix;
212 Vvdw_rep = c12*rinvsix*rinvsix;
213 fscal += (12.0*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
214 Vvdwtot = Vvdwtot+Vvdw_rep-Vvdw_disp;
220 cexp1 = vdwparam[tj+1];
221 cexp2 = vdwparam[tj+2];
223 rinvsix = rinvsq*rinvsq*rinvsq;
224 Vvdw_disp = c6*rinvsix;
226 Vvdw_rep = cexp1*exp(-br);
227 fscal += (br*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
228 Vvdwtot = Vvdwtot+Vvdw_rep-Vvdw_disp;
234 c12 = vdwparam[tj+1];
238 Geps = eps*VFtab[nnn+2];
239 Heps2 = eps2*VFtab[nnn+3];
242 FF = Fp+Geps+2.0*Heps2;
248 Geps = eps*VFtab[nnn+2];
249 Heps2 = eps2*VFtab[nnn+3];
252 FF = Fp+Geps+2.0*Heps2;
255 fscal += -(fijD+fijR)*tabscale*rinv;
256 Vvdwtot = Vvdwtot + Vvdw_disp + Vvdw_rep;
260 gmx_fatal(FARGS,"Death & horror! No generic VdW interaction for ivdw=%d.\n",ivdw);
263 } /* end VdW interactions */
272 f[j3+0] = f[j3+0] - tx;
273 f[j3+1] = f[j3+1] - ty;
274 f[j3+2] = f[j3+2] - tz;
277 f[ii3+0] = f[ii3+0] + fix;
278 f[ii3+1] = f[ii3+1] + fiy;
279 f[ii3+2] = f[ii3+2] + fiz;
280 fshift[is3] = fshift[is3]+fix;
281 fshift[is3+1] = fshift[is3+1]+fiy;
282 fshift[is3+2] = fshift[is3+2]+fiz;
283 ggid = nlist->gid[n];
284 Vc[ggid] = Vc[ggid] + vctot;
285 Vvdw[ggid] = Vvdw[ggid] + Vvdwtot;
288 *outeriter = nlist->nri;
289 *inneriter = nlist->jindex[n];