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,2014, 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"
51 #include "gmx_fatal.h"
54 gmx_nb_generic_cg_kernel(t_nblist * nlist,
59 nb_kernel_data_t * kernel_data,
62 int nri, ntype, table_nelements, ielec, ivdw;
63 real facel, gbtabscale;
64 int n, is3, i3, k, nj0, nj1, j3, ggid, nnn, n0;
65 int ai0, ai1, ai, aj0, aj1, aj;
67 real fscal, tx, ty, tz;
70 real qq, vcoul, krsq, vctot;
73 real rt, r, eps, eps2, Y, F, Geps, Heps2, VV, FF, Fp, fijD, fijR;
76 real Vvdw_rep, Vvdw_disp;
77 real ix, iy, iz, fix, fiy, fiz;
79 real dx, dy, dz, rsq, rinv;
80 real c6, c12, cexp1, cexp2, br;
100 fshift = fr->fshift[0];
101 Vc = kernel_data->energygrp_elec;
102 Vvdw = kernel_data->energygrp_vdw;
103 tabscale = kernel_data->table_elec_vdw->scale;
104 VFtab = kernel_data->table_elec_vdw->data;
106 /* avoid compiler warnings for cases that cannot happen */
112 /* 3 VdW parameters for buckingham, otherwise 2 */
113 nvdwparam = (nlist->ivdw == 2) ? 3 : 2;
114 table_nelements = (ielec == 3) ? 4 : 0;
115 table_nelements += (ivdw == 3) ? 8 : 0;
117 charge = mdatoms->chargeA;
118 type = mdatoms->typeA;
120 shiftvec = fr->shift_vec[0];
124 for (n = 0; (n < nlist->nri); n++)
126 is3 = 3*nlist->shift[n];
128 shY = shiftvec[is3+1];
129 shZ = shiftvec[is3+2];
130 nj0 = nlist->jindex[n];
131 nj1 = nlist->jindex[n+1];
132 ai0 = nlist->iinr[n];
133 ai1 = nlist->iinr_end[n];
140 for (k = nj0; (k < nj1); k++)
142 aj0 = nlist->jjnr[k];
143 aj1 = nlist->jjnr_end[k];
144 excl = &nlist->excl[k*MAX_CGCGSIZE];
146 for (ai = ai0; (ai < ai1); ai++)
152 iq = facel*charge[ai];
153 nti = nvdwparam*ntype*type[ai];
155 /* Note that this code currently calculates
156 * all LJ and Coulomb interactions,
157 * even if the LJ parameters or charges are zero.
158 * If required, this can be optimized.
161 for (aj = aj0; (aj < aj1); aj++)
163 /* Check if this interaction is excluded */
164 if (excl[aj-aj0] & (1<<(ai-ai0)))
176 rsq = dx*dx+dy*dy+dz*dz;
177 rinv = gmx_invsqrt(rsq);
181 if (ielec == 3 || ivdw == 3)
188 nnn = table_nelements*n0;
191 /* Coulomb interaction. ielec==0 means no interaction */
199 /* Vanilla cutoff coulomb */
201 fscal = vcoul*rinvsq;
207 vcoul = qq*(rinv+krsq-fr->c_rf);
208 fscal = qq*(rinv-2.0*krsq)*rinvsq;
212 /* Tabulated coulomb */
215 Geps = eps*VFtab[nnn+2];
216 Heps2 = eps2*VFtab[nnn+3];
220 FF = Fp+Geps+2.0*Heps2;
222 fscal = -qq*FF*tabscale*rinv;
227 gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
231 gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
235 } /* End of coulomb interactions */
238 /* VdW interaction. ivdw==0 means no interaction */
241 tj = nti+nvdwparam*type[aj];
246 /* Vanilla Lennard-Jones cutoff */
248 c12 = vdwparam[tj+1];
250 rinvsix = rinvsq*rinvsq*rinvsq;
251 Vvdw_disp = c6*rinvsix;
252 Vvdw_rep = c12*rinvsix*rinvsix;
253 fscal += (12.0*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
254 Vvdwtot = Vvdwtot+Vvdw_rep-Vvdw_disp;
260 cexp1 = vdwparam[tj+1];
261 cexp2 = vdwparam[tj+2];
263 rinvsix = rinvsq*rinvsq*rinvsq;
264 Vvdw_disp = c6*rinvsix;
266 Vvdw_rep = cexp1*exp(-br);
267 fscal += (br*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
268 Vvdwtot = Vvdwtot+Vvdw_rep-Vvdw_disp;
274 c12 = vdwparam[tj+1];
278 Geps = eps*VFtab[nnn+2];
279 Heps2 = eps2*VFtab[nnn+3];
282 FF = Fp+Geps+2.0*Heps2;
288 Geps = eps*VFtab[nnn+2];
289 Heps2 = eps2*VFtab[nnn+3];
292 FF = Fp+Geps+2.0*Heps2;
295 fscal += -(fijD+fijR)*tabscale*rinv;
296 Vvdwtot = Vvdwtot + Vvdw_disp + Vvdw_rep;
300 gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
303 } /* end VdW interactions */
323 fshift[is3+1] += fiy;
324 fshift[is3+2] += fiz;
325 ggid = nlist->gid[n];
327 Vvdw[ggid] += Vvdwtot;
329 /* Estimate flops, average for generic cg kernel:
330 * 12 flops per outer iteration
331 * 100 flops per inner iteration
333 inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC_CG, nlist->nri*12 + nlist->jindex[n]*100);