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
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROningen Mixture of Alchemy and Childrens' Stories
42 #include "types/simple.h"
45 #include "nb_generic_cg.h"
48 gmx_nb_generic_cg_kernel(t_nblist * nlist,
61 int nri,ntype,table_nelements,icoul,ivdw;
62 real facel,gbtabscale;
63 int n,is3,i3,k,nj0,nj1,j3,ggid,nnn,n0;
64 int ai0,ai1,ai,aj0,aj1,aj;
69 real qq,vcoul,krsq,vctot;
72 real rt,r,eps,eps2,Y,F,Geps,Heps2,VV,FF,Fp,fijD,fijR;
75 real Vvdw_rep,Vvdw_disp;
76 real ix,iy,iz,fix,fiy,fiz;
78 real dx,dy,dz,rsq,rinv;
79 real c6,c12,cexp1,cexp2,br;
90 /* avoid compiler warnings for cases that cannot happen */
96 /* 3 VdW parameters for buckingham, otherwise 2 */
97 nvdwparam = (nlist->ivdw==2) ? 3 : 2;
98 table_nelements = (icoul==3) ? 4 : 0;
99 table_nelements += (ivdw==3) ? 8 : 0;
101 charge = mdatoms->chargeA;
102 type = mdatoms->typeA;
104 shiftvec = fr->shift_vec[0];
108 for(n=0; (n<nlist->nri); n++)
110 is3 = 3*nlist->shift[n];
112 shY = shiftvec[is3+1];
113 shZ = shiftvec[is3+2];
114 nj0 = nlist->jindex[n];
115 nj1 = nlist->jindex[n+1];
116 ai0 = nlist->iinr[n];
117 ai1 = nlist->iinr_end[n];
124 for(k=nj0; (k<nj1); k++)
126 aj0 = nlist->jjnr[k];
127 aj1 = nlist->jjnr_end[k];
128 excl = &nlist->excl[k*MAX_CGCGSIZE];
130 for(ai=ai0; (ai<ai1); ai++)
136 iq = facel*charge[ai];
137 nti = nvdwparam*ntype*type[ai];
139 /* Note that this code currently calculates
140 * all LJ and Coulomb interactions,
141 * even if the LJ parameters or charges are zero.
142 * If required, this can be optimized.
145 for(aj=aj0; (aj<aj1); aj++)
147 /* Check if this interaction is excluded */
148 if (excl[aj-aj0] & (1<<(ai-ai0)))
160 rsq = dx*dx+dy*dy+dz*dz;
161 rinv = gmx_invsqrt(rsq);
165 if (icoul==3 || ivdw==3)
172 nnn = table_nelements*n0;
175 /* Coulomb interaction. icoul==0 means no interaction */
183 /* Vanilla cutoff coulomb */
185 fscal = vcoul*rinvsq;
191 vcoul = qq*(rinv+krsq-fr->c_rf);
192 fscal = qq*(rinv-2.0*krsq)*rinvsq;
196 /* Tabulated coulomb */
199 Geps = eps*VFtab[nnn+2];
200 Heps2 = eps2*VFtab[nnn+3];
204 FF = Fp+Geps+2.0*Heps2;
206 fscal = -qq*FF*tabscale*rinv;
211 gmx_fatal(FARGS,"Death & horror! GB generic interaction not implemented.\n");
215 gmx_fatal(FARGS,"Death & horror! No generic coulomb interaction for icoul=%d.\n",icoul);
219 } /* End of coulomb interactions */
222 /* VdW interaction. ivdw==0 means no interaction */
225 tj = nti+nvdwparam*type[aj];
230 /* Vanilla Lennard-Jones cutoff */
232 c12 = vdwparam[tj+1];
234 rinvsix = rinvsq*rinvsq*rinvsq;
235 Vvdw_disp = c6*rinvsix;
236 Vvdw_rep = c12*rinvsix*rinvsix;
237 fscal += (12.0*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
238 Vvdwtot = Vvdwtot+Vvdw_rep-Vvdw_disp;
244 cexp1 = vdwparam[tj+1];
245 cexp2 = vdwparam[tj+2];
247 rinvsix = rinvsq*rinvsq*rinvsq;
248 Vvdw_disp = c6*rinvsix;
250 Vvdw_rep = cexp1*exp(-br);
251 fscal += (br*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
252 Vvdwtot = Vvdwtot+Vvdw_rep-Vvdw_disp;
258 c12 = vdwparam[tj+1];
262 Geps = eps*VFtab[nnn+2];
263 Heps2 = eps2*VFtab[nnn+3];
266 FF = Fp+Geps+2.0*Heps2;
272 Geps = eps*VFtab[nnn+2];
273 Heps2 = eps2*VFtab[nnn+3];
276 FF = Fp+Geps+2.0*Heps2;
279 fscal += -(fijD+fijR)*tabscale*rinv;
280 Vvdwtot = Vvdwtot + Vvdw_disp + Vvdw_rep;
284 gmx_fatal(FARGS,"Death & horror! No generic VdW interaction for ivdw=%d.\n",ivdw);
287 } /* end VdW interactions */
307 fshift[is3+1] += fiy;
308 fshift[is3+2] += fiz;
309 ggid = nlist->gid[n];
311 Vvdw[ggid] += Vvdwtot;
314 *outeriter = nlist->nri;
315 *inneriter = nlist->jindex[n];