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"
46 #include "nonbonded.h"
47 #include "nb_kernel.h"
51 gmx_nb_generic_cg_kernel(t_nblist * nlist,
56 nb_kernel_data_t * kernel_data,
59 int nri, ntype, table_nelements, ielec, ivdw;
60 real facel, gbtabscale;
61 int n, is3, i3, k, nj0, nj1, j3, ggid, nnn, n0;
62 int ai0, ai1, ai, aj0, aj1, aj;
64 real fscal, tx, ty, tz;
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;
97 fshift = fr->fshift[0];
98 Vc = kernel_data->energygrp_elec;
99 Vvdw = kernel_data->energygrp_vdw;
100 tabscale = kernel_data->table_elec_vdw->scale;
101 VFtab = kernel_data->table_elec_vdw->data;
103 /* avoid compiler warnings for cases that cannot happen */
109 /* 3 VdW parameters for buckingham, otherwise 2 */
110 nvdwparam = (nlist->ivdw == 2) ? 3 : 2;
111 table_nelements = (ielec == 3) ? 4 : 0;
112 table_nelements += (ivdw == 3) ? 8 : 0;
114 charge = mdatoms->chargeA;
115 type = mdatoms->typeA;
117 shiftvec = fr->shift_vec[0];
121 for (n = 0; (n < nlist->nri); n++)
123 is3 = 3*nlist->shift[n];
125 shY = shiftvec[is3+1];
126 shZ = shiftvec[is3+2];
127 nj0 = nlist->jindex[n];
128 nj1 = nlist->jindex[n+1];
129 ai0 = nlist->iinr[n];
130 ai1 = nlist->iinr_end[n];
137 for (k = nj0; (k < nj1); k++)
139 aj0 = nlist->jjnr[k];
140 aj1 = nlist->jjnr_end[k];
141 excl = &nlist->excl[k*MAX_CGCGSIZE];
143 for (ai = ai0; (ai < ai1); ai++)
149 iq = facel*charge[ai];
150 nti = nvdwparam*ntype*type[ai];
152 /* Note that this code currently calculates
153 * all LJ and Coulomb interactions,
154 * even if the LJ parameters or charges are zero.
155 * If required, this can be optimized.
158 for (aj = aj0; (aj < aj1); aj++)
160 /* Check if this interaction is excluded */
161 if (excl[aj-aj0] & (1<<(ai-ai0)))
173 rsq = dx*dx+dy*dy+dz*dz;
174 rinv = gmx_invsqrt(rsq);
178 if (ielec == 3 || ivdw == 3)
185 nnn = table_nelements*n0;
188 /* Coulomb interaction. ielec==0 means no interaction */
196 /* Vanilla cutoff coulomb */
198 fscal = vcoul*rinvsq;
204 vcoul = qq*(rinv+krsq-fr->c_rf);
205 fscal = qq*(rinv-2.0*krsq)*rinvsq;
209 /* Tabulated coulomb */
212 Geps = eps*VFtab[nnn+2];
213 Heps2 = eps2*VFtab[nnn+3];
217 FF = Fp+Geps+2.0*Heps2;
219 fscal = -qq*FF*tabscale*rinv;
224 gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
228 gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
232 } /* End of coulomb interactions */
235 /* VdW interaction. ivdw==0 means no interaction */
238 tj = nti+nvdwparam*type[aj];
243 /* Vanilla Lennard-Jones cutoff */
245 c12 = vdwparam[tj+1];
247 rinvsix = rinvsq*rinvsq*rinvsq;
248 Vvdw_disp = c6*rinvsix;
249 Vvdw_rep = c12*rinvsix*rinvsix;
250 fscal += (12.0*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
251 Vvdwtot = Vvdwtot+Vvdw_rep-Vvdw_disp;
257 cexp1 = vdwparam[tj+1];
258 cexp2 = vdwparam[tj+2];
260 rinvsix = rinvsq*rinvsq*rinvsq;
261 Vvdw_disp = c6*rinvsix;
263 Vvdw_rep = cexp1*exp(-br);
264 fscal += (br*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
265 Vvdwtot = Vvdwtot+Vvdw_rep-Vvdw_disp;
271 c12 = vdwparam[tj+1];
275 Geps = eps*VFtab[nnn+2];
276 Heps2 = eps2*VFtab[nnn+3];
279 FF = Fp+Geps+2.0*Heps2;
285 Geps = eps*VFtab[nnn+2];
286 Heps2 = eps2*VFtab[nnn+3];
289 FF = Fp+Geps+2.0*Heps2;
292 fscal += -(fijD+fijR)*tabscale*rinv;
293 Vvdwtot = Vvdwtot + Vvdw_disp + Vvdw_rep;
297 gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
300 } /* end VdW interactions */
320 fshift[is3+1] += fiy;
321 fshift[is3+2] += fiz;
322 ggid = nlist->gid[n];
324 Vvdw[ggid] += Vvdwtot;
326 /* Estimate flops, average for generic cg kernel:
327 * 12 flops per outer iteration
328 * 100 flops per inner iteration
330 inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC_CG, nlist->nri*12 + nlist->jindex[n]*100);