3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * GROningen Mixture of Alchemy and Childrens' Stories
41 #include "types/simple.h"
44 #include "nb_generic.h"
47 #include "nonbonded.h"
48 #include "nb_kernel.h"
52 gmx_nb_generic_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, ii, is3, ii3, k, nj0, nj1, jnr, j3, ggid, nnn, n0;
64 real fscal, felec, fvdw, velec, vvdw, tx, ty, tz;
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;
91 real ewtabscale, eweps, sh_ewald, ewrt, ewtabhalfspace;
93 real rcoulomb2, rvdw, rvdw2, sh_invrc6;
94 real rcutoff, rcutoff2;
95 real rswitch_elec, rswitch_vdw, d, d2, sw, dsw, rinvcorr;
96 real elec_swV3, elec_swV4, elec_swV5, elec_swF2, elec_swF3, elec_swF4;
97 real vdw_swV3, vdw_swV4, vdw_swV5, vdw_swF2, vdw_swF3, vdw_swF4;
98 gmx_bool bExactElecCutoff, bExactVdwCutoff, bExactCutoff;
102 ielec = nlist->ielec;
105 fshift = fr->fshift[0];
106 velecgrp = kernel_data->energygrp_elec;
107 vvdwgrp = kernel_data->energygrp_vdw;
108 tabscale = kernel_data->table_elec_vdw->scale;
109 VFtab = kernel_data->table_elec_vdw->data;
111 sh_ewald = fr->ic->sh_ewald;
112 ewtab = fr->ic->tabq_coul_FDV0;
113 ewtabscale = fr->ic->tabq_scale;
114 ewtabhalfspace = 0.5/ewtabscale;
116 rcoulomb2 = fr->rcoulomb*fr->rcoulomb;
119 sh_invrc6 = fr->ic->sh_invrc6;
121 if (fr->coulomb_modifier == eintmodPOTSWITCH)
123 d = fr->rcoulomb-fr->rcoulomb_switch;
124 elec_swV3 = -10.0/(d*d*d);
125 elec_swV4 = 15.0/(d*d*d*d);
126 elec_swV5 = -6.0/(d*d*d*d*d);
127 elec_swF2 = -30.0/(d*d*d);
128 elec_swF3 = 60.0/(d*d*d*d);
129 elec_swF4 = -30.0/(d*d*d*d*d);
133 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
134 elec_swV3 = elec_swV4 = elec_swV5 = elec_swF2 = elec_swF3 = elec_swF4 = 0.0;
136 if (fr->vdw_modifier == eintmodPOTSWITCH)
138 d = fr->rvdw-fr->rvdw_switch;
139 vdw_swV3 = -10.0/(d*d*d);
140 vdw_swV4 = 15.0/(d*d*d*d);
141 vdw_swV5 = -6.0/(d*d*d*d*d);
142 vdw_swF2 = -30.0/(d*d*d);
143 vdw_swF3 = 60.0/(d*d*d*d);
144 vdw_swF4 = -30.0/(d*d*d*d*d);
148 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
149 vdw_swV3 = vdw_swV4 = vdw_swV5 = vdw_swF2 = vdw_swF3 = vdw_swF4 = 0.0;
152 bExactElecCutoff = (fr->coulomb_modifier != eintmodNONE) || fr->eeltype == eelRF_ZERO;
153 bExactVdwCutoff = (fr->vdw_modifier != eintmodNONE);
154 bExactCutoff = bExactElecCutoff || bExactVdwCutoff;
158 rcutoff = ( fr->rcoulomb > fr->rvdw ) ? fr->rcoulomb : fr->rvdw;
159 rcutoff2 = rcutoff*rcutoff;
163 /* Fix warnings for stupid compilers */
164 rcutoff = rcutoff2 = 1e30;
167 /* avoid compiler warnings for cases that cannot happen */
172 /* 3 VdW parameters for buckingham, otherwise 2 */
173 nvdwparam = (ivdw == GMX_NBKERNEL_VDW_BUCKINGHAM) ? 3 : 2;
174 table_nelements = 12;
176 charge = mdatoms->chargeA;
177 type = mdatoms->typeA;
179 shiftvec = fr->shift_vec[0];
183 for (n = 0; (n < nlist->nri); n++)
185 is3 = 3*nlist->shift[n];
187 shY = shiftvec[is3+1];
188 shZ = shiftvec[is3+2];
189 nj0 = nlist->jindex[n];
190 nj1 = nlist->jindex[n+1];
196 iq = facel*charge[ii];
197 nti = nvdwparam*ntype*type[ii];
204 for (k = nj0; (k < nj1); k++)
206 jnr = nlist->jjnr[k];
214 rsq = dx*dx+dy*dy+dz*dz;
215 rinv = gmx_invsqrt(rsq);
222 if (bExactCutoff && rsq > rcutoff2)
227 if (ielec == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE || ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE)
234 nnn = table_nelements*n0;
237 /* Coulomb interaction. ielec==0 means no interaction */
238 if (ielec != GMX_NBKERNEL_ELEC_NONE)
244 case GMX_NBKERNEL_ELEC_NONE:
247 case GMX_NBKERNEL_ELEC_COULOMB:
248 /* Vanilla cutoff coulomb */
250 felec = velec*rinvsq;
253 case GMX_NBKERNEL_ELEC_REACTIONFIELD:
255 velec = qq*(rinv+fr->k_rf*rsq-fr->c_rf);
256 felec = qq*(rinv*rinvsq-2.0*fr->k_rf);
259 case GMX_NBKERNEL_ELEC_CUBICSPLINETABLE:
260 /* Tabulated coulomb */
263 Geps = eps*VFtab[nnn+2];
264 Heps2 = eps2*VFtab[nnn+3];
267 FF = Fp+Geps+2.0*Heps2;
269 felec = -qq*FF*tabscale*rinv;
272 case GMX_NBKERNEL_ELEC_GENERALIZEDBORN:
274 gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
277 case GMX_NBKERNEL_ELEC_EWALD:
278 ewrt = rsq*rinv*ewtabscale;
282 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
283 rinvcorr = (fr->coulomb_modifier == eintmodPOTSHIFT) ? rinv-fr->ic->sh_ewald : rinv;
284 velec = qq*(rinvcorr-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
285 felec = qq*rinv*(rinvsq-felec);
289 gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
292 if (fr->coulomb_modifier == eintmodPOTSWITCH)
294 d = rsq*rinv-fr->rcoulomb_switch;
295 d = (d > 0.0) ? d : 0.0;
297 sw = 1.0+d2*d*(elec_swV3+d*(elec_swV4+d*elec_swV5));
298 dsw = d2*(elec_swF2+d*(elec_swF3+d*elec_swF4));
299 /* Apply switch function. Note that felec=f/r since it will be multiplied
300 * by the i-j displacement vector. This means felec'=f'/r=-(v*sw)'/r=
301 * -(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=felec*sw-v*dsw/r
303 felec = felec*sw - rinv*velec*dsw;
304 /* Once we have used velec to update felec we can modify velec too */
307 if (bExactElecCutoff)
309 felec = (rsq <= rcoulomb2) ? felec : 0.0;
310 velec = (rsq <= rcoulomb2) ? velec : 0.0;
313 } /* End of coulomb interactions */
316 /* VdW interaction. ivdw==0 means no interaction */
317 if (ivdw != GMX_NBKERNEL_VDW_NONE)
319 tj = nti+nvdwparam*type[jnr];
323 case GMX_NBKERNEL_VDW_NONE:
326 case GMX_NBKERNEL_VDW_LENNARDJONES:
327 /* Vanilla Lennard-Jones cutoff */
329 c12 = vdwparam[tj+1];
330 rinvsix = rinvsq*rinvsq*rinvsq;
331 vvdw_disp = c6*rinvsix;
332 vvdw_rep = c12*rinvsix*rinvsix;
333 fvdw = (vvdw_rep-vvdw_disp)*rinvsq;
334 if (fr->vdw_modifier == eintmodPOTSHIFT)
336 vvdw = (vvdw_rep-c12*sh_invrc6*sh_invrc6)*(1.0/12.0)-(vvdw_disp-c6*sh_invrc6)*(1.0/6.0);
340 vvdw = vvdw_rep/12.0-vvdw_disp/6.0;
344 case GMX_NBKERNEL_VDW_BUCKINGHAM:
347 cexp1 = vdwparam[tj+1];
348 cexp2 = vdwparam[tj+2];
350 rinvsix = rinvsq*rinvsq*rinvsq;
351 vvdw_disp = c6*rinvsix;
353 vvdw_rep = cexp1*exp(-br);
354 fvdw = (br*vvdw_rep-vvdw_disp)*rinvsq;
355 if (fr->vdw_modifier == eintmodPOTSHIFT)
357 vvdw = (vvdw_rep-cexp1*exp(-cexp2*rvdw))-(vvdw_disp-c6*sh_invrc6)/6.0;
361 vvdw = vvdw_rep-vvdw_disp/6.0;
365 case GMX_NBKERNEL_VDW_CUBICSPLINETABLE:
368 c12 = vdwparam[tj+1];
371 Geps = eps*VFtab[nnn+6];
372 Heps2 = eps2*VFtab[nnn+7];
375 FF = Fp+Geps+2.0*Heps2;
380 Geps = eps*VFtab[nnn+10];
381 Heps2 = eps2*VFtab[nnn+11];
384 FF = Fp+Geps+2.0*Heps2;
387 fvdw = -(fijD+fijR)*tabscale*rinv;
388 vvdw = vvdw_disp + vvdw_rep;
392 gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
395 if (fr->vdw_modifier == eintmodPOTSWITCH)
397 d = rsq*rinv-fr->rvdw_switch;
398 d = (d > 0.0) ? d : 0.0;
400 sw = 1.0+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
401 dsw = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
402 /* See coulomb interaction for the force-switch formula */
403 fvdw = fvdw*sw - rinv*vvdw*dsw;
408 fvdw = (rsq <= rvdw2) ? fvdw : 0.0;
409 vvdw = (rsq <= rvdw2) ? vvdw : 0.0;
412 } /* end VdW interactions */
422 f[j3+0] = f[j3+0] - tx;
423 f[j3+1] = f[j3+1] - ty;
424 f[j3+2] = f[j3+2] - tz;
427 f[ii3+0] = f[ii3+0] + fix;
428 f[ii3+1] = f[ii3+1] + fiy;
429 f[ii3+2] = f[ii3+2] + fiz;
430 fshift[is3] = fshift[is3]+fix;
431 fshift[is3+1] = fshift[is3+1]+fiy;
432 fshift[is3+2] = fshift[is3+2]+fiz;
433 ggid = nlist->gid[n];
434 velecgrp[ggid] += vctot;
435 vvdwgrp[ggid] += vvdwtot;
437 /* Estimate flops, average for generic kernel:
438 * 12 flops per outer iteration
439 * 50 flops per inner iteration
441 inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC, nlist->nri*12 + nlist->jindex[n]*50);