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,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.
39 #include "nb_generic.h"
43 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
44 #include "gromacs/legacyheaders/nonbonded.h"
45 #include "gromacs/legacyheaders/nrnb.h"
46 #include "gromacs/legacyheaders/typedefs.h"
47 #include "gromacs/legacyheaders/types/simple.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/utility/fatalerror.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, c6grid, cexp1, cexp2, br;
80 real * vdwparam, *vdwgridparam;
91 real ewtabscale, eweps, sh_ewald, ewrt, ewtabhalfspace;
93 real rcoulomb2, rvdw, rvdw2, sh_dispersion, sh_repulsion;
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 real ewclj, ewclj2, ewclj6, ewcljrsq, poly, exponent, sh_lj_ewald;
99 gmx_bool bExactElecCutoff, bExactVdwCutoff, bExactCutoff;
103 ielec = nlist->ielec;
106 fshift = fr->fshift[0];
107 velecgrp = kernel_data->energygrp_elec;
108 vvdwgrp = kernel_data->energygrp_vdw;
109 tabscale = kernel_data->table_elec_vdw->scale;
110 VFtab = kernel_data->table_elec_vdw->data;
112 sh_ewald = fr->ic->sh_ewald;
113 ewtab = fr->ic->tabq_coul_FDV0;
114 ewtabscale = fr->ic->tabq_scale;
115 ewtabhalfspace = 0.5/ewtabscale;
117 rcoulomb2 = fr->rcoulomb*fr->rcoulomb;
120 sh_dispersion = fr->ic->dispersion_shift.cpot;
121 sh_repulsion = fr->ic->repulsion_shift.cpot;
122 sh_lj_ewald = fr->ic->sh_lj_ewald;
124 ewclj = fr->ewaldcoeff_lj;
125 ewclj2 = ewclj*ewclj;
126 ewclj6 = ewclj2*ewclj2*ewclj2;
128 if (fr->coulomb_modifier == eintmodPOTSWITCH)
130 d = fr->rcoulomb-fr->rcoulomb_switch;
131 elec_swV3 = -10.0/(d*d*d);
132 elec_swV4 = 15.0/(d*d*d*d);
133 elec_swV5 = -6.0/(d*d*d*d*d);
134 elec_swF2 = -30.0/(d*d*d);
135 elec_swF3 = 60.0/(d*d*d*d);
136 elec_swF4 = -30.0/(d*d*d*d*d);
140 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
141 elec_swV3 = elec_swV4 = elec_swV5 = elec_swF2 = elec_swF3 = elec_swF4 = 0.0;
143 if (fr->vdw_modifier == eintmodPOTSWITCH)
145 d = fr->rvdw-fr->rvdw_switch;
146 vdw_swV3 = -10.0/(d*d*d);
147 vdw_swV4 = 15.0/(d*d*d*d);
148 vdw_swV5 = -6.0/(d*d*d*d*d);
149 vdw_swF2 = -30.0/(d*d*d);
150 vdw_swF3 = 60.0/(d*d*d*d);
151 vdw_swF4 = -30.0/(d*d*d*d*d);
155 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
156 vdw_swV3 = vdw_swV4 = vdw_swV5 = vdw_swF2 = vdw_swF3 = vdw_swF4 = 0.0;
159 bExactElecCutoff = (fr->coulomb_modifier != eintmodNONE) || fr->eeltype == eelRF_ZERO;
160 bExactVdwCutoff = (fr->vdw_modifier != eintmodNONE);
161 bExactCutoff = bExactElecCutoff && bExactVdwCutoff;
165 rcutoff = ( fr->rcoulomb > fr->rvdw ) ? fr->rcoulomb : fr->rvdw;
166 rcutoff2 = rcutoff*rcutoff;
170 /* Fix warnings for stupid compilers */
171 rcutoff = rcutoff2 = 1e30;
174 /* avoid compiler warnings for cases that cannot happen */
179 /* 3 VdW parameters for Buckingham, otherwise 2 */
180 nvdwparam = (ivdw == GMX_NBKERNEL_VDW_BUCKINGHAM) ? 3 : 2;
181 table_nelements = 12;
183 charge = mdatoms->chargeA;
184 type = mdatoms->typeA;
186 shiftvec = fr->shift_vec[0];
189 vdwgridparam = fr->ljpme_c6grid;
191 for (n = 0; (n < nlist->nri); n++)
193 is3 = 3*nlist->shift[n];
195 shY = shiftvec[is3+1];
196 shZ = shiftvec[is3+2];
197 nj0 = nlist->jindex[n];
198 nj1 = nlist->jindex[n+1];
204 iq = facel*charge[ii];
205 nti = nvdwparam*ntype*type[ii];
212 for (k = nj0; (k < nj1); k++)
214 jnr = nlist->jjnr[k];
222 rsq = dx*dx+dy*dy+dz*dz;
223 rinv = gmx_invsqrt(rsq);
230 if (bExactCutoff && rsq >= rcutoff2)
235 if (ielec == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE || ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE)
242 nnn = table_nelements*n0;
245 /* Coulomb interaction. ielec==0 means no interaction */
246 if (ielec != GMX_NBKERNEL_ELEC_NONE)
252 case GMX_NBKERNEL_ELEC_NONE:
255 case GMX_NBKERNEL_ELEC_COULOMB:
256 /* Vanilla cutoff coulomb */
258 felec = velec*rinvsq;
259 /* The shift for the Coulomb potential is stored in
260 * the RF parameter c_rf, which is 0 without shift
262 velec -= qq*fr->ic->c_rf;
265 case GMX_NBKERNEL_ELEC_REACTIONFIELD:
267 velec = qq*(rinv+fr->k_rf*rsq-fr->c_rf);
268 felec = qq*(rinv*rinvsq-2.0*fr->k_rf);
271 case GMX_NBKERNEL_ELEC_CUBICSPLINETABLE:
272 /* Tabulated coulomb */
275 Geps = eps*VFtab[nnn+2];
276 Heps2 = eps2*VFtab[nnn+3];
279 FF = Fp+Geps+2.0*Heps2;
281 felec = -qq*FF*tabscale*rinv;
284 case GMX_NBKERNEL_ELEC_GENERALIZEDBORN:
286 gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
289 case GMX_NBKERNEL_ELEC_EWALD:
290 ewrt = rsq*rinv*ewtabscale;
294 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
295 rinvcorr = (fr->coulomb_modifier == eintmodPOTSHIFT) ? rinv-fr->ic->sh_ewald : rinv;
296 velec = qq*(rinvcorr-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
297 felec = qq*rinv*(rinvsq-felec);
301 gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
304 if (fr->coulomb_modifier == eintmodPOTSWITCH)
306 d = rsq*rinv-fr->rcoulomb_switch;
307 d = (d > 0.0) ? d : 0.0;
309 sw = 1.0+d2*d*(elec_swV3+d*(elec_swV4+d*elec_swV5));
310 dsw = d2*(elec_swF2+d*(elec_swF3+d*elec_swF4));
311 /* Apply switch function. Note that felec=f/r since it will be multiplied
312 * by the i-j displacement vector. This means felec'=f'/r=-(v*sw)'/r=
313 * -(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=felec*sw-v*dsw/r
315 felec = felec*sw - rinv*velec*dsw;
316 /* Once we have used velec to update felec we can modify velec too */
319 if (bExactElecCutoff)
321 felec = (rsq < rcoulomb2) ? felec : 0.0;
322 velec = (rsq < rcoulomb2) ? velec : 0.0;
325 } /* End of coulomb interactions */
328 /* VdW interaction. ivdw==0 means no interaction */
329 if (ivdw != GMX_NBKERNEL_VDW_NONE)
331 tj = nti+nvdwparam*type[jnr];
335 case GMX_NBKERNEL_VDW_NONE:
338 case GMX_NBKERNEL_VDW_LENNARDJONES:
339 /* Vanilla Lennard-Jones cutoff */
341 c12 = vdwparam[tj+1];
342 rinvsix = rinvsq*rinvsq*rinvsq;
343 vvdw_disp = c6*rinvsix;
344 vvdw_rep = c12*rinvsix*rinvsix;
345 fvdw = (vvdw_rep-vvdw_disp)*rinvsq;
346 if (fr->vdw_modifier == eintmodPOTSHIFT)
348 vvdw = (vvdw_rep + c12*sh_repulsion)/12.0 - (vvdw_disp + c6*sh_dispersion)/6.0;
352 vvdw = vvdw_rep/12.0-vvdw_disp/6.0;
356 case GMX_NBKERNEL_VDW_BUCKINGHAM:
359 cexp1 = vdwparam[tj+1];
360 cexp2 = vdwparam[tj+2];
362 rinvsix = rinvsq*rinvsq*rinvsq;
363 vvdw_disp = c6*rinvsix;
365 vvdw_rep = cexp1*exp(-br);
366 fvdw = (br*vvdw_rep-vvdw_disp)*rinvsq;
367 if (fr->vdw_modifier == eintmodPOTSHIFT)
369 vvdw = (vvdw_rep-cexp1*exp(-cexp2*rvdw))-(vvdw_disp + c6*sh_dispersion)/6.0;
373 vvdw = vvdw_rep-vvdw_disp/6.0;
377 case GMX_NBKERNEL_VDW_CUBICSPLINETABLE:
380 c12 = vdwparam[tj+1];
383 Geps = eps*VFtab[nnn+6];
384 Heps2 = eps2*VFtab[nnn+7];
387 FF = Fp+Geps+2.0*Heps2;
392 Geps = eps*VFtab[nnn+10];
393 Heps2 = eps2*VFtab[nnn+11];
396 FF = Fp+Geps+2.0*Heps2;
399 fvdw = -(fijD+fijR)*tabscale*rinv;
400 vvdw = vvdw_disp + vvdw_rep;
404 case GMX_NBKERNEL_VDW_LJEWALD:
406 rinvsix = rinvsq*rinvsq*rinvsq;
407 ewcljrsq = ewclj2*rsq;
408 exponent = exp(-ewcljrsq);
409 poly = exponent*(1.0 + ewcljrsq + ewcljrsq*ewcljrsq*0.5);
411 c12 = vdwparam[tj+1];
412 c6grid = vdwgridparam[tj];
413 vvdw_disp = (c6-c6grid*(1.0-poly))*rinvsix;
414 vvdw_rep = c12*rinvsix*rinvsix;
415 fvdw = (vvdw_rep - vvdw_disp - c6grid*(1.0/6.0)*exponent*ewclj6)*rinvsq;
416 if (fr->vdw_modifier == eintmodPOTSHIFT)
418 vvdw = (vvdw_rep + c12*sh_repulsion)/12.0 - (vvdw_disp + c6*sh_dispersion - c6grid*sh_lj_ewald)/6.0;
422 vvdw = vvdw_rep/12.0-vvdw_disp/6.0;
427 gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
430 if (fr->vdw_modifier == eintmodPOTSWITCH)
432 d = rsq*rinv-fr->rvdw_switch;
433 d = (d > 0.0) ? d : 0.0;
435 sw = 1.0+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
436 dsw = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
437 /* See coulomb interaction for the force-switch formula */
438 fvdw = fvdw*sw - rinv*vvdw*dsw;
443 fvdw = (rsq < rvdw2) ? fvdw : 0.0;
444 vvdw = (rsq < rvdw2) ? vvdw : 0.0;
447 } /* end VdW interactions */
457 f[j3+0] = f[j3+0] - tx;
458 f[j3+1] = f[j3+1] - ty;
459 f[j3+2] = f[j3+2] - tz;
462 f[ii3+0] = f[ii3+0] + fix;
463 f[ii3+1] = f[ii3+1] + fiy;
464 f[ii3+2] = f[ii3+2] + fiz;
465 fshift[is3] = fshift[is3]+fix;
466 fshift[is3+1] = fshift[is3+1]+fiy;
467 fshift[is3+2] = fshift[is3+2]+fiz;
468 ggid = nlist->gid[n];
469 velecgrp[ggid] += vctot;
470 vvdwgrp[ggid] += vvdwtot;
472 /* Estimate flops, average for generic kernel:
473 * 12 flops per outer iteration
474 * 50 flops per inner iteration
476 inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC, nlist->nri*12 + nlist->jindex[n]*50);