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,2015,2017,2018, 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.
39 #include "nb_generic_cg.h"
43 #include "gromacs/gmxlib/nrnb.h"
44 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
45 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
46 #include "gromacs/math/functions.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/utility/fatalerror.h"
51 gmx_nb_generic_cg_kernel(t_nblist * nlist,
56 nb_kernel_data_t * kernel_data,
59 int ntype, table_nelements, ielec, ivdw;
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;
101 do_tab = (ielec == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE ||
102 ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE);
105 tabscale = kernel_data->table_elec_vdw->scale;
106 VFtab = kernel_data->table_elec_vdw->data;
114 /* avoid compiler warnings for cases that cannot happen */
119 /* 3 VdW parameters for buckingham, otherwise 2 */
120 nvdwparam = (nlist->ivdw == 2) ? 3 : 2;
121 table_nelements = (ielec == 3) ? 4 : 0;
122 table_nelements += (ivdw == 3) ? 8 : 0;
124 charge = mdatoms->chargeA;
125 type = mdatoms->typeA;
126 facel = fr->ic->epsfac;
127 shiftvec = fr->shift_vec[0];
131 for (n = 0; (n < nlist->nri); n++)
133 is3 = 3*nlist->shift[n];
135 shY = shiftvec[is3+1];
136 shZ = shiftvec[is3+2];
137 nj0 = nlist->jindex[n];
138 nj1 = nlist->jindex[n+1];
139 ai0 = nlist->iinr[n];
140 ai1 = nlist->iinr_end[n];
147 for (k = nj0; (k < nj1); k++)
149 aj0 = nlist->jjnr[k];
150 aj1 = nlist->jjnr_end[k];
151 excl = &nlist->excl[k*MAX_CGCGSIZE];
153 for (ai = ai0; (ai < ai1); ai++)
159 iq = facel*charge[ai];
160 nti = nvdwparam*ntype*type[ai];
162 /* Note that this code currently calculates
163 * all LJ and Coulomb interactions,
164 * even if the LJ parameters or charges are zero.
165 * If required, this can be optimized.
168 for (aj = aj0; (aj < aj1); aj++)
170 /* Check if this interaction is excluded */
171 if (excl[aj-aj0] & (1<<(ai-ai0)))
183 rsq = dx*dx+dy*dy+dz*dz;
184 rinv = gmx::invsqrt(rsq);
188 if (ielec == 3 || ivdw == 3)
195 nnn = table_nelements*n0;
198 /* Coulomb interaction. ielec==0 means no interaction */
206 /* Vanilla cutoff coulomb */
208 fscal = vcoul*rinvsq;
213 krsq = fr->ic->k_rf*rsq;
214 vcoul = qq*(rinv + krsq - fr->ic->c_rf);
215 fscal = qq*(rinv-2.0*krsq)*rinvsq;
219 /* Tabulated coulomb */
222 Geps = eps*VFtab[nnn+2];
223 Heps2 = eps2*VFtab[nnn+3];
227 FF = Fp+Geps+2.0*Heps2;
229 fscal = -qq*FF*tabscale*rinv;
234 gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
236 gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
239 } /* End of coulomb interactions */
242 /* VdW interaction. ivdw==0 means no interaction */
245 tj = nti+nvdwparam*type[aj];
250 /* Vanilla Lennard-Jones cutoff */
252 c12 = vdwparam[tj+1];
254 rinvsix = rinvsq*rinvsq*rinvsq;
255 Vvdw_disp = c6*rinvsix;
256 Vvdw_rep = c12*rinvsix*rinvsix;
257 fscal += (12.0*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
258 Vvdwtot = Vvdwtot+Vvdw_rep-Vvdw_disp;
264 cexp1 = vdwparam[tj+1];
265 cexp2 = vdwparam[tj+2];
267 rinvsix = rinvsq*rinvsq*rinvsq;
268 Vvdw_disp = c6*rinvsix;
270 Vvdw_rep = cexp1*std::exp(-br);
271 fscal += (br*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
272 Vvdwtot = Vvdwtot+Vvdw_rep-Vvdw_disp;
278 c12 = vdwparam[tj+1];
282 Geps = eps*VFtab[nnn+2];
283 Heps2 = eps2*VFtab[nnn+3];
286 FF = Fp+Geps+2.0*Heps2;
292 Geps = eps*VFtab[nnn+2];
293 Heps2 = eps2*VFtab[nnn+3];
296 FF = Fp+Geps+2.0*Heps2;
299 fscal += -(fijD+fijR)*tabscale*rinv;
300 Vvdwtot = Vvdwtot + Vvdw_disp + Vvdw_rep;
304 gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
306 } /* end VdW interactions */
326 fshift[is3+1] += fiy;
327 fshift[is3+2] += fiz;
328 ggid = nlist->gid[n];
330 Vvdw[ggid] += Vvdwtot;
332 /* Estimate flops, average for generic cg kernel:
333 * 12 flops per outer iteration
334 * 100 flops per inner iteration
336 inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC_CG, nlist->nri*12 + nlist->jindex[n]*100);