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,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.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/mdtypes/md_enums.h"
49 #include "gromacs/utility/fatalerror.h"
52 gmx_nb_generic_kernel(t_nblist * nlist,
57 nb_kernel_data_t * kernel_data,
60 int ntype, table_nelements, ielec, ivdw;
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;
90 real ewtabscale, eweps, ewrt, ewtabhalfspace;
92 real rcoulomb2, rvdw, rvdw2, sh_dispersion, sh_repulsion;
93 real rcutoff, rcutoff2;
94 real d, d2, sw, dsw, rinvcorr;
95 real elec_swV3, elec_swV4, elec_swV5, elec_swF2, elec_swF3, elec_swF4;
96 real vdw_swV3, vdw_swV4, vdw_swV5, vdw_swF2, vdw_swF3, vdw_swF4;
97 real ewclj, ewclj2, ewclj6, ewcljrsq, poly, exponent, sh_lj_ewald;
98 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;
110 do_tab = (ielec == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE ||
111 ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE);
114 tabscale = kernel_data->table_elec_vdw->scale;
115 VFtab = kernel_data->table_elec_vdw->data;
123 const interaction_const_t *ic = fr->ic;
125 ewtab = ic->tabq_coul_FDV0;
126 ewtabscale = ic->tabq_scale;
127 ewtabhalfspace = 0.5/ewtabscale;
129 rcoulomb2 = ic->rcoulomb*ic->rcoulomb;
132 sh_dispersion = ic->dispersion_shift.cpot;
133 sh_repulsion = ic->repulsion_shift.cpot;
134 sh_lj_ewald = ic->sh_lj_ewald;
136 ewclj = ic->ewaldcoeff_lj;
137 ewclj2 = ewclj*ewclj;
138 ewclj6 = ewclj2*ewclj2*ewclj2;
140 if (ic->coulomb_modifier == eintmodPOTSWITCH)
142 d = ic->rcoulomb - ic->rcoulomb_switch;
143 elec_swV3 = -10.0/(d*d*d);
144 elec_swV4 = 15.0/(d*d*d*d);
145 elec_swV5 = -6.0/(d*d*d*d*d);
146 elec_swF2 = -30.0/(d*d*d);
147 elec_swF3 = 60.0/(d*d*d*d);
148 elec_swF4 = -30.0/(d*d*d*d*d);
152 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
153 elec_swV3 = elec_swV4 = elec_swV5 = elec_swF2 = elec_swF3 = elec_swF4 = 0.0;
155 if (ic->vdw_modifier == eintmodPOTSWITCH)
157 d = ic->rvdw - ic->rvdw_switch;
158 vdw_swV3 = -10.0/(d*d*d);
159 vdw_swV4 = 15.0/(d*d*d*d);
160 vdw_swV5 = -6.0/(d*d*d*d*d);
161 vdw_swF2 = -30.0/(d*d*d);
162 vdw_swF3 = 60.0/(d*d*d*d);
163 vdw_swF4 = -30.0/(d*d*d*d*d);
167 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
168 vdw_swV3 = vdw_swV4 = vdw_swV5 = vdw_swF2 = vdw_swF3 = vdw_swF4 = 0.0;
171 bExactElecCutoff = (ic->coulomb_modifier != eintmodNONE) || ic->eeltype == eelRF_ZERO;
172 bExactVdwCutoff = (ic->vdw_modifier != eintmodNONE);
173 bExactCutoff = bExactElecCutoff && bExactVdwCutoff;
177 rcutoff = ( ic->rcoulomb > ic->rvdw ) ? ic->rcoulomb : ic->rvdw;
178 rcutoff2 = rcutoff*rcutoff;
182 /* Fix warnings for stupid compilers */
186 /* avoid compiler warnings for cases that cannot happen */
191 /* 3 VdW parameters for Buckingham, otherwise 2 */
192 nvdwparam = (ivdw == GMX_NBKERNEL_VDW_BUCKINGHAM) ? 3 : 2;
193 table_nelements = 12;
195 charge = mdatoms->chargeA;
196 type = mdatoms->typeA;
197 facel = fr->ic->epsfac;
198 shiftvec = fr->shift_vec[0];
201 vdwgridparam = fr->ljpme_c6grid;
203 for (n = 0; (n < nlist->nri); n++)
205 is3 = 3*nlist->shift[n];
207 shY = shiftvec[is3+1];
208 shZ = shiftvec[is3+2];
209 nj0 = nlist->jindex[n];
210 nj1 = nlist->jindex[n+1];
216 iq = facel*charge[ii];
217 nti = nvdwparam*ntype*type[ii];
224 for (k = nj0; (k < nj1); k++)
226 jnr = nlist->jjnr[k];
234 rsq = dx*dx+dy*dy+dz*dz;
235 rinv = gmx::invsqrt(rsq);
242 if (bExactCutoff && rsq >= rcutoff2)
247 if (ielec == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE || ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE)
254 nnn = table_nelements*n0;
257 /* Coulomb interaction. ielec==0 means no interaction */
258 if (ielec != GMX_NBKERNEL_ELEC_NONE)
264 case GMX_NBKERNEL_ELEC_NONE:
267 case GMX_NBKERNEL_ELEC_COULOMB:
268 /* Vanilla cutoff coulomb */
270 felec = velec*rinvsq;
271 /* The shift for the Coulomb potential is stored in
272 * the RF parameter c_rf, which is 0 without shift
274 velec -= qq*ic->c_rf;
277 case GMX_NBKERNEL_ELEC_REACTIONFIELD:
279 velec = qq*(rinv + ic->k_rf*rsq-ic->c_rf);
280 felec = qq*(rinv*rinvsq - 2.0*ic->k_rf);
283 case GMX_NBKERNEL_ELEC_CUBICSPLINETABLE:
284 /* Tabulated coulomb */
287 Geps = eps*VFtab[nnn+2];
288 Heps2 = eps2*VFtab[nnn+3];
291 FF = Fp+Geps+2.0*Heps2;
293 felec = -qq*FF*tabscale*rinv;
296 case GMX_NBKERNEL_ELEC_EWALD:
297 ewrt = rsq*rinv*ewtabscale;
301 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
302 rinvcorr = (ic->coulomb_modifier == eintmodPOTSHIFT) ? rinv - ic->sh_ewald : rinv;
303 velec = qq*(rinvcorr-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
304 felec = qq*rinv*(rinvsq-felec);
308 gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
311 if (ic->coulomb_modifier == eintmodPOTSWITCH)
313 d = rsq*rinv - ic->rcoulomb_switch;
314 d = (d > 0.0) ? d : 0.0;
316 sw = 1.0+d2*d*(elec_swV3+d*(elec_swV4+d*elec_swV5));
317 dsw = d2*(elec_swF2+d*(elec_swF3+d*elec_swF4));
318 /* Apply switch function. Note that felec=f/r since it will be multiplied
319 * by the i-j displacement vector. This means felec'=f'/r=-(v*sw)'/r=
320 * -(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=felec*sw-v*dsw/r
322 felec = felec*sw - rinv*velec*dsw;
323 /* Once we have used velec to update felec we can modify velec too */
326 if (bExactElecCutoff)
328 felec = (rsq < rcoulomb2) ? felec : 0.0;
329 velec = (rsq < rcoulomb2) ? velec : 0.0;
332 } /* End of coulomb interactions */
335 /* VdW interaction. ivdw==0 means no interaction */
336 if (ivdw != GMX_NBKERNEL_VDW_NONE)
338 tj = nti+nvdwparam*type[jnr];
342 case GMX_NBKERNEL_VDW_NONE:
345 case GMX_NBKERNEL_VDW_LENNARDJONES:
346 /* Vanilla Lennard-Jones cutoff */
348 c12 = vdwparam[tj+1];
349 rinvsix = rinvsq*rinvsq*rinvsq;
350 vvdw_disp = c6*rinvsix;
351 vvdw_rep = c12*rinvsix*rinvsix;
352 fvdw = (vvdw_rep-vvdw_disp)*rinvsq;
353 if (ic->vdw_modifier == eintmodPOTSHIFT)
355 vvdw = (vvdw_rep + c12*sh_repulsion)/12.0 - (vvdw_disp + c6*sh_dispersion)/6.0;
359 vvdw = vvdw_rep/12.0-vvdw_disp/6.0;
363 case GMX_NBKERNEL_VDW_BUCKINGHAM:
366 cexp1 = vdwparam[tj+1];
367 cexp2 = vdwparam[tj+2];
369 rinvsix = rinvsq*rinvsq*rinvsq;
370 vvdw_disp = c6*rinvsix;
372 vvdw_rep = cexp1*std::exp(-br);
373 fvdw = (br*vvdw_rep-vvdw_disp)*rinvsq;
374 if (ic->vdw_modifier == eintmodPOTSHIFT)
376 vvdw = (vvdw_rep-cexp1*std::exp(-cexp2*rvdw))-(vvdw_disp + c6*sh_dispersion)/6.0;
380 vvdw = vvdw_rep-vvdw_disp/6.0;
384 case GMX_NBKERNEL_VDW_CUBICSPLINETABLE:
387 c12 = vdwparam[tj+1];
390 Geps = eps*VFtab[nnn+6];
391 Heps2 = eps2*VFtab[nnn+7];
394 FF = Fp+Geps+2.0*Heps2;
399 Geps = eps*VFtab[nnn+10];
400 Heps2 = eps2*VFtab[nnn+11];
403 FF = Fp+Geps+2.0*Heps2;
406 fvdw = -(fijD+fijR)*tabscale*rinv;
407 vvdw = vvdw_disp + vvdw_rep;
411 case GMX_NBKERNEL_VDW_LJEWALD:
413 rinvsix = rinvsq*rinvsq*rinvsq;
414 ewcljrsq = ewclj2*rsq;
415 exponent = std::exp(-ewcljrsq);
416 poly = exponent*(1.0 + ewcljrsq + ewcljrsq*ewcljrsq*0.5);
418 c12 = vdwparam[tj+1];
419 c6grid = vdwgridparam[tj];
420 vvdw_disp = (c6-c6grid*(1.0-poly))*rinvsix;
421 vvdw_rep = c12*rinvsix*rinvsix;
422 fvdw = (vvdw_rep - vvdw_disp - c6grid*(1.0/6.0)*exponent*ewclj6)*rinvsq;
423 if (ic->vdw_modifier == eintmodPOTSHIFT)
425 vvdw = (vvdw_rep + c12*sh_repulsion)/12.0 - (vvdw_disp + c6*sh_dispersion - c6grid*sh_lj_ewald)/6.0;
429 vvdw = vvdw_rep/12.0-vvdw_disp/6.0;
434 gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
437 if (ic->vdw_modifier == eintmodPOTSWITCH)
439 d = rsq*rinv - ic->rvdw_switch;
440 d = (d > 0.0) ? d : 0.0;
442 sw = 1.0+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
443 dsw = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
444 /* See coulomb interaction for the force-switch formula */
445 fvdw = fvdw*sw - rinv*vvdw*dsw;
450 fvdw = (rsq < rvdw2) ? fvdw : 0.0;
451 vvdw = (rsq < rvdw2) ? vvdw : 0.0;
454 } /* end VdW interactions */
464 f[j3+0] = f[j3+0] - tx;
465 f[j3+1] = f[j3+1] - ty;
466 f[j3+2] = f[j3+2] - tz;
469 f[ii3+0] = f[ii3+0] + fix;
470 f[ii3+1] = f[ii3+1] + fiy;
471 f[ii3+2] = f[ii3+2] + fiz;
472 fshift[is3] = fshift[is3]+fix;
473 fshift[is3+1] = fshift[is3+1]+fiy;
474 fshift[is3+2] = fshift[is3+2]+fiz;
475 ggid = nlist->gid[n];
476 velecgrp[ggid] += vctot;
477 vvdwgrp[ggid] += vvdwtot;
479 /* Estimate flops, average for generic kernel:
480 * 12 flops per outer iteration
481 * 50 flops per inner iteration
483 inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC, nlist->nri*12 + nlist->jindex[n]*50);