2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
43 #include "types/simple.h"
46 #include "nb_generic_cg.h"
47 #include "nonbonded.h"
48 #include "nb_kernel.h"
52 gmx_nb_generic_cg_kernel(t_nblist * nlist,
57 nb_kernel_data_t * kernel_data,
60 int nri, ntype, table_nelements, ielec, ivdw;
61 real facel, gbtabscale;
62 int n, is3, i3, k, nj0, nj1, j3, ggid, nnn, n0;
63 int ai0, ai1, ai, aj0, aj1, aj;
65 real fscal, tx, ty, tz;
68 real qq, vcoul, krsq, vctot;
71 real rt, r, eps, eps2, Y, F, Geps, Heps2, VV, FF, Fp, fijD, fijR;
74 real Vvdw_rep, Vvdw_disp;
75 real ix, iy, iz, fix, fiy, fiz;
77 real dx, dy, dz, rsq, rinv;
78 real c6, c12, cexp1, cexp2, br;
98 fshift = fr->fshift[0];
99 Vc = kernel_data->energygrp_elec;
100 Vvdw = kernel_data->energygrp_vdw;
101 tabscale = kernel_data->table_elec_vdw->scale;
102 VFtab = kernel_data->table_elec_vdw->data;
104 /* avoid compiler warnings for cases that cannot happen */
110 /* 3 VdW parameters for buckingham, otherwise 2 */
111 nvdwparam = (nlist->ivdw == 2) ? 3 : 2;
112 table_nelements = (ielec == 3) ? 4 : 0;
113 table_nelements += (ivdw == 3) ? 8 : 0;
115 charge = mdatoms->chargeA;
116 type = mdatoms->typeA;
118 shiftvec = fr->shift_vec[0];
122 for (n = 0; (n < nlist->nri); n++)
124 is3 = 3*nlist->shift[n];
126 shY = shiftvec[is3+1];
127 shZ = shiftvec[is3+2];
128 nj0 = nlist->jindex[n];
129 nj1 = nlist->jindex[n+1];
130 ai0 = nlist->iinr[n];
131 ai1 = nlist->iinr_end[n];
138 for (k = nj0; (k < nj1); k++)
140 aj0 = nlist->jjnr[k];
141 aj1 = nlist->jjnr_end[k];
142 excl = &nlist->excl[k*MAX_CGCGSIZE];
144 for (ai = ai0; (ai < ai1); ai++)
150 iq = facel*charge[ai];
151 nti = nvdwparam*ntype*type[ai];
153 /* Note that this code currently calculates
154 * all LJ and Coulomb interactions,
155 * even if the LJ parameters or charges are zero.
156 * If required, this can be optimized.
159 for (aj = aj0; (aj < aj1); aj++)
161 /* Check if this interaction is excluded */
162 if (excl[aj-aj0] & (1<<(ai-ai0)))
174 rsq = dx*dx+dy*dy+dz*dz;
175 rinv = gmx_invsqrt(rsq);
179 if (ielec == 3 || ivdw == 3)
186 nnn = table_nelements*n0;
189 /* Coulomb interaction. ielec==0 means no interaction */
197 /* Vanilla cutoff coulomb */
199 fscal = vcoul*rinvsq;
205 vcoul = qq*(rinv+krsq-fr->c_rf);
206 fscal = qq*(rinv-2.0*krsq)*rinvsq;
210 /* Tabulated coulomb */
213 Geps = eps*VFtab[nnn+2];
214 Heps2 = eps2*VFtab[nnn+3];
218 FF = Fp+Geps+2.0*Heps2;
220 fscal = -qq*FF*tabscale*rinv;
225 gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
229 gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
233 } /* End of coulomb interactions */
236 /* VdW interaction. ivdw==0 means no interaction */
239 tj = nti+nvdwparam*type[aj];
244 /* Vanilla Lennard-Jones cutoff */
246 c12 = vdwparam[tj+1];
248 rinvsix = rinvsq*rinvsq*rinvsq;
249 Vvdw_disp = c6*rinvsix;
250 Vvdw_rep = c12*rinvsix*rinvsix;
251 fscal += (12.0*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
252 Vvdwtot = Vvdwtot+Vvdw_rep-Vvdw_disp;
258 cexp1 = vdwparam[tj+1];
259 cexp2 = vdwparam[tj+2];
261 rinvsix = rinvsq*rinvsq*rinvsq;
262 Vvdw_disp = c6*rinvsix;
264 Vvdw_rep = cexp1*exp(-br);
265 fscal += (br*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
266 Vvdwtot = Vvdwtot+Vvdw_rep-Vvdw_disp;
272 c12 = vdwparam[tj+1];
276 Geps = eps*VFtab[nnn+2];
277 Heps2 = eps2*VFtab[nnn+3];
280 FF = Fp+Geps+2.0*Heps2;
286 Geps = eps*VFtab[nnn+2];
287 Heps2 = eps2*VFtab[nnn+3];
290 FF = Fp+Geps+2.0*Heps2;
293 fscal += -(fijD+fijR)*tabscale*rinv;
294 Vvdwtot = Vvdwtot + Vvdw_disp + Vvdw_rep;
298 gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
301 } /* end VdW interactions */
321 fshift[is3+1] += fiy;
322 fshift[is3+2] += fiz;
323 ggid = nlist->gid[n];
325 Vvdw[ggid] += Vvdwtot;
327 /* Estimate flops, average for generic cg kernel:
328 * 12 flops per outer iteration
329 * 100 flops per inner iteration
331 inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC_CG, nlist->nri*12 + nlist->jindex[n]*100);