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) 2012, 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.h"
49 #include "nonbonded.h"
50 #include "nb_kernel.h"
54 gmx_nb_generic_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, ii, is3, ii3, k, nj0, nj1, jnr, j3, ggid, nnn, n0;
66 real fscal, felec, fvdw, velec, vvdw, tx, ty, tz;
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;
93 real ewtabscale, eweps, sh_ewald, ewrt, ewtabhalfspace;
95 real rcoulomb2, rvdw, rvdw2, sh_invrc6;
96 real rcutoff, rcutoff2;
97 real rswitch_elec, rswitch_vdw, d, d2, sw, dsw, rinvcorr;
98 real elec_swV3, elec_swV4, elec_swV5, elec_swF2, elec_swF3, elec_swF4;
99 real vdw_swV3, vdw_swV4, vdw_swV5, vdw_swF2, vdw_swF3, vdw_swF4;
100 gmx_bool bExactElecCutoff, bExactVdwCutoff, bExactCutoff;
104 ielec = nlist->ielec;
107 fshift = fr->fshift[0];
108 velecgrp = kernel_data->energygrp_elec;
109 vvdwgrp = kernel_data->energygrp_vdw;
110 tabscale = kernel_data->table_elec_vdw->scale;
111 VFtab = kernel_data->table_elec_vdw->data;
113 sh_ewald = fr->ic->sh_ewald;
114 ewtab = fr->ic->tabq_coul_FDV0;
115 ewtabscale = fr->ic->tabq_scale;
116 ewtabhalfspace = 0.5/ewtabscale;
118 rcoulomb2 = fr->rcoulomb*fr->rcoulomb;
121 sh_invrc6 = fr->ic->sh_invrc6;
123 if (fr->coulomb_modifier == eintmodPOTSWITCH)
125 d = fr->rcoulomb-fr->rcoulomb_switch;
126 elec_swV3 = -10.0/(d*d*d);
127 elec_swV4 = 15.0/(d*d*d*d);
128 elec_swV5 = -6.0/(d*d*d*d*d);
129 elec_swF2 = -30.0/(d*d*d);
130 elec_swF3 = 60.0/(d*d*d*d);
131 elec_swF4 = -30.0/(d*d*d*d*d);
135 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
136 elec_swV3 = elec_swV4 = elec_swV5 = elec_swF2 = elec_swF3 = elec_swF4 = 0.0;
138 if (fr->vdw_modifier == eintmodPOTSWITCH)
140 d = fr->rvdw-fr->rvdw_switch;
141 vdw_swV3 = -10.0/(d*d*d);
142 vdw_swV4 = 15.0/(d*d*d*d);
143 vdw_swV5 = -6.0/(d*d*d*d*d);
144 vdw_swF2 = -30.0/(d*d*d);
145 vdw_swF3 = 60.0/(d*d*d*d);
146 vdw_swF4 = -30.0/(d*d*d*d*d);
150 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
151 vdw_swV3 = vdw_swV4 = vdw_swV5 = vdw_swF2 = vdw_swF3 = vdw_swF4 = 0.0;
154 bExactElecCutoff = (fr->coulomb_modifier != eintmodNONE) || fr->eeltype == eelRF_ZERO;
155 bExactVdwCutoff = (fr->vdw_modifier != eintmodNONE);
156 bExactCutoff = bExactElecCutoff || bExactVdwCutoff;
160 rcutoff = ( fr->rcoulomb > fr->rvdw ) ? fr->rcoulomb : fr->rvdw;
161 rcutoff2 = rcutoff*rcutoff;
165 /* Fix warnings for stupid compilers */
166 rcutoff = rcutoff2 = 1e30;
169 /* avoid compiler warnings for cases that cannot happen */
174 /* 3 VdW parameters for buckingham, otherwise 2 */
175 nvdwparam = (ivdw == GMX_NBKERNEL_VDW_BUCKINGHAM) ? 3 : 2;
176 table_nelements = 12;
178 charge = mdatoms->chargeA;
179 type = mdatoms->typeA;
181 shiftvec = fr->shift_vec[0];
185 for (n = 0; (n < nlist->nri); n++)
187 is3 = 3*nlist->shift[n];
189 shY = shiftvec[is3+1];
190 shZ = shiftvec[is3+2];
191 nj0 = nlist->jindex[n];
192 nj1 = nlist->jindex[n+1];
198 iq = facel*charge[ii];
199 nti = nvdwparam*ntype*type[ii];
206 for (k = nj0; (k < nj1); k++)
208 jnr = nlist->jjnr[k];
216 rsq = dx*dx+dy*dy+dz*dz;
217 rinv = gmx_invsqrt(rsq);
224 if (bExactCutoff && rsq > rcutoff2)
229 if (ielec == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE || ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE)
236 nnn = table_nelements*n0;
239 /* Coulomb interaction. ielec==0 means no interaction */
240 if (ielec != GMX_NBKERNEL_ELEC_NONE)
246 case GMX_NBKERNEL_ELEC_NONE:
249 case GMX_NBKERNEL_ELEC_COULOMB:
250 /* Vanilla cutoff coulomb */
252 felec = velec*rinvsq;
255 case GMX_NBKERNEL_ELEC_REACTIONFIELD:
257 velec = qq*(rinv+fr->k_rf*rsq-fr->c_rf);
258 felec = qq*(rinv*rinvsq-2.0*fr->k_rf);
261 case GMX_NBKERNEL_ELEC_CUBICSPLINETABLE:
262 /* Tabulated coulomb */
265 Geps = eps*VFtab[nnn+2];
266 Heps2 = eps2*VFtab[nnn+3];
269 FF = Fp+Geps+2.0*Heps2;
271 felec = -qq*FF*tabscale*rinv;
274 case GMX_NBKERNEL_ELEC_GENERALIZEDBORN:
276 gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
279 case GMX_NBKERNEL_ELEC_EWALD:
280 ewrt = rsq*rinv*ewtabscale;
284 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
285 rinvcorr = (fr->coulomb_modifier == eintmodPOTSHIFT) ? rinv-fr->ic->sh_ewald : rinv;
286 velec = qq*(rinvcorr-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
287 felec = qq*rinv*(rinvsq-felec);
291 gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
294 if (fr->coulomb_modifier == eintmodPOTSWITCH)
296 d = rsq*rinv-fr->rcoulomb_switch;
297 d = (d > 0.0) ? d : 0.0;
299 sw = 1.0+d2*d*(elec_swV3+d*(elec_swV4+d*elec_swV5));
300 dsw = d2*(elec_swF2+d*(elec_swF3+d*elec_swF4));
301 /* Apply switch function. Note that felec=f/r since it will be multiplied
302 * by the i-j displacement vector. This means felec'=f'/r=-(v*sw)'/r=
303 * -(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=felec*sw-v*dsw/r
305 felec = felec*sw - rinv*velec*dsw;
306 /* Once we have used velec to update felec we can modify velec too */
309 if (bExactElecCutoff)
311 felec = (rsq <= rcoulomb2) ? felec : 0.0;
312 velec = (rsq <= rcoulomb2) ? velec : 0.0;
315 } /* End of coulomb interactions */
318 /* VdW interaction. ivdw==0 means no interaction */
319 if (ivdw != GMX_NBKERNEL_VDW_NONE)
321 tj = nti+nvdwparam*type[jnr];
325 case GMX_NBKERNEL_VDW_NONE:
328 case GMX_NBKERNEL_VDW_LENNARDJONES:
329 /* Vanilla Lennard-Jones cutoff */
331 c12 = vdwparam[tj+1];
332 rinvsix = rinvsq*rinvsq*rinvsq;
333 vvdw_disp = c6*rinvsix;
334 vvdw_rep = c12*rinvsix*rinvsix;
335 fvdw = (vvdw_rep-vvdw_disp)*rinvsq;
336 if (fr->vdw_modifier == eintmodPOTSHIFT)
338 vvdw = (vvdw_rep-c12*sh_invrc6*sh_invrc6)*(1.0/12.0)-(vvdw_disp-c6*sh_invrc6)*(1.0/6.0);
342 vvdw = vvdw_rep/12.0-vvdw_disp/6.0;
346 case GMX_NBKERNEL_VDW_BUCKINGHAM:
349 cexp1 = vdwparam[tj+1];
350 cexp2 = vdwparam[tj+2];
352 rinvsix = rinvsq*rinvsq*rinvsq;
353 vvdw_disp = c6*rinvsix;
355 vvdw_rep = cexp1*exp(-br);
356 fvdw = (br*vvdw_rep-vvdw_disp)*rinvsq;
357 if (fr->vdw_modifier == eintmodPOTSHIFT)
359 vvdw = (vvdw_rep-cexp1*exp(-cexp2*rvdw))-(vvdw_disp-c6*sh_invrc6)/6.0;
363 vvdw = vvdw_rep-vvdw_disp/6.0;
367 case GMX_NBKERNEL_VDW_CUBICSPLINETABLE:
370 c12 = vdwparam[tj+1];
373 Geps = eps*VFtab[nnn+6];
374 Heps2 = eps2*VFtab[nnn+7];
377 FF = Fp+Geps+2.0*Heps2;
382 Geps = eps*VFtab[nnn+10];
383 Heps2 = eps2*VFtab[nnn+11];
386 FF = Fp+Geps+2.0*Heps2;
389 fvdw = -(fijD+fijR)*tabscale*rinv;
390 vvdw = vvdw_disp + vvdw_rep;
394 gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
397 if (fr->vdw_modifier == eintmodPOTSWITCH)
399 d = rsq*rinv-fr->rvdw_switch;
400 d = (d > 0.0) ? d : 0.0;
402 sw = 1.0+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
403 dsw = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
404 /* See coulomb interaction for the force-switch formula */
405 fvdw = fvdw*sw - rinv*vvdw*dsw;
410 fvdw = (rsq <= rvdw2) ? fvdw : 0.0;
411 vvdw = (rsq <= rvdw2) ? vvdw : 0.0;
414 } /* end VdW interactions */
424 f[j3+0] = f[j3+0] - tx;
425 f[j3+1] = f[j3+1] - ty;
426 f[j3+2] = f[j3+2] - tz;
429 f[ii3+0] = f[ii3+0] + fix;
430 f[ii3+1] = f[ii3+1] + fiy;
431 f[ii3+2] = f[ii3+2] + fiz;
432 fshift[is3] = fshift[is3]+fix;
433 fshift[is3+1] = fshift[is3+1]+fiy;
434 fshift[is3+2] = fshift[is3+2]+fiz;
435 ggid = nlist->gid[n];
436 velecgrp[ggid] += vctot;
437 vvdwgrp[ggid] += vvdwtot;
439 /* Estimate flops, average for generic kernel:
440 * 12 flops per outer iteration
441 * 50 flops per inner iteration
443 inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC, nlist->nri*12 + nlist->jindex[n]*50);