2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009 Christoph Junghans, Brad Lambeth.
5 * Copyright (c) 2011 Christoph Junghans, Sebastian Fritsch.
6 * Copyright (c) 2011,2012,2013,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.
41 #include "gromacs/legacyheaders/types/simple.h"
42 #include "gromacs/math/vec.h"
43 #include "gromacs/legacyheaders/typedefs.h"
44 #include "nb_generic_adress.h"
45 #include "gromacs/legacyheaders/nrnb.h"
47 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/legacyheaders/nonbonded.h"
50 #include "nb_kernel.h"
52 #define ALMOST_ZERO 1e-30
53 #define ALMOST_ONE 1-(1e-30)
55 gmx_nb_generic_adress_kernel(t_nblist * nlist,
60 nb_kernel_data_t * kernel_data,
63 int nri, ntype, table_nelements, ielec, ivdw;
64 real facel, gbtabscale;
65 int n, ii, is3, ii3, k, nj0, nj1, jnr, j3, ggid, nnn, n0;
67 real fscal, felec, fvdw, velec, vvdw, tx, ty, tz;
73 real rt, r, eps, eps2, Y, F, Geps, Heps2, VV, FF, Fp, fijD, fijR;
76 real vvdw_rep, vvdw_disp;
77 real ix, iy, iz, fix, fiy, fiz;
79 real dx, dy, dz, rsq, rinv;
80 real c6, c12, cexp1, cexp2, br;
94 real ewtabscale, eweps, sh_ewald, ewrt, ewtabhalfspace;
96 real rcoulomb2, rvdw, rvdw2, sh_dispersion, sh_repulsion;
97 real rcutoff, rcutoff2;
98 real rswitch_elec, rswitch_vdw, d, d2, sw, dsw, rinvcorr;
99 real elec_swV3, elec_swV4, elec_swV5, elec_swF2, elec_swF3, elec_swF4;
100 real vdw_swV3, vdw_swV4, vdw_swV5, vdw_swF2, vdw_swF3, vdw_swF4;
101 gmx_bool bExactElecCutoff, bExactVdwCutoff, bExactCutoff;
107 real hybscal; /* the multiplicator to the force for hybrid interactions*/
114 force_cap = fr->adress_ex_forcecap;
118 ielec = nlist->ielec;
121 fshift = fr->fshift[0];
122 velecgrp = kernel_data->energygrp_elec;
123 vvdwgrp = kernel_data->energygrp_vdw;
124 tabscale = kernel_data->table_elec_vdw->scale;
125 VFtab = kernel_data->table_elec_vdw->data;
127 sh_ewald = fr->ic->sh_ewald;
128 ewtab = fr->ic->tabq_coul_FDV0;
129 ewtabscale = fr->ic->tabq_scale;
130 ewtabhalfspace = 0.5/ewtabscale;
132 rcoulomb2 = fr->rcoulomb*fr->rcoulomb;
135 sh_dispersion = fr->ic->dispersion_shift.cpot;
136 sh_repulsion = fr->ic->repulsion_shift.cpot;
138 if (fr->coulomb_modifier == eintmodPOTSWITCH)
140 d = fr->rcoulomb-fr->rcoulomb_switch;
141 elec_swV3 = -10.0/(d*d*d);
142 elec_swV4 = 15.0/(d*d*d*d);
143 elec_swV5 = -6.0/(d*d*d*d*d);
144 elec_swF2 = -30.0/(d*d*d);
145 elec_swF3 = 60.0/(d*d*d*d);
146 elec_swF4 = -30.0/(d*d*d*d*d);
150 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
151 elec_swV3 = elec_swV4 = elec_swV5 = elec_swF2 = elec_swF3 = elec_swF4 = 0.0;
153 if (fr->vdw_modifier == eintmodPOTSWITCH)
155 d = fr->rvdw-fr->rvdw_switch;
156 vdw_swV3 = -10.0/(d*d*d);
157 vdw_swV4 = 15.0/(d*d*d*d);
158 vdw_swV5 = -6.0/(d*d*d*d*d);
159 vdw_swF2 = -30.0/(d*d*d);
160 vdw_swF3 = 60.0/(d*d*d*d);
161 vdw_swF4 = -30.0/(d*d*d*d*d);
165 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
166 vdw_swV3 = vdw_swV4 = vdw_swV5 = vdw_swF2 = vdw_swF3 = vdw_swF4 = 0.0;
169 bExactElecCutoff = (fr->coulomb_modifier != eintmodNONE) || fr->eeltype == eelRF_ZERO;
170 bExactVdwCutoff = (fr->vdw_modifier != eintmodNONE);
171 bExactCutoff = bExactElecCutoff || bExactVdwCutoff;
175 rcutoff = ( fr->rcoulomb > fr->rvdw ) ? fr->rcoulomb : fr->rvdw;
176 rcutoff2 = rcutoff*rcutoff;
180 /* Fix warnings for stupid compilers */
181 rcutoff = rcutoff2 = 1e30;
184 /* avoid compiler warnings for cases that cannot happen */
189 /* 3 VdW parameters for buckingham, otherwise 2 */
190 nvdwparam = (ivdw == GMX_NBKERNEL_VDW_BUCKINGHAM) ? 3 : 2;
191 table_nelements = 12;
193 charge = mdatoms->chargeA;
194 type = mdatoms->typeA;
196 shiftvec = fr->shift_vec[0];
200 for (n = 0; (n < nlist->nri); n++)
202 is3 = 3*nlist->shift[n];
204 shY = shiftvec[is3+1];
205 shZ = shiftvec[is3+2];
206 nj0 = nlist->jindex[n];
207 nj1 = nlist->jindex[n+1];
213 iq = facel*charge[ii];
214 nti = nvdwparam*ntype*type[ii];
221 /* We need to find out if this i atom is part of an
222 all-atom or CG energy group */
223 egp_nr = mdatoms->cENER[ii];
224 bCG = !fr->adress_group_explicit[egp_nr];
228 if ((!bCG) && weight_cg1 < ALMOST_ZERO)
233 for (k = nj0; (k < nj1); k++)
235 jnr = nlist->jjnr[k];
236 weight_cg2 = wf[jnr];
237 weight_product = weight_cg1*weight_cg2;
239 if (weight_product < ALMOST_ZERO)
241 /* if it's a explicit loop, skip this atom */
246 else /* if it's a coarse grained loop, include this atom */
251 else if (weight_product >= ALMOST_ONE)
254 /* if it's a explicit loop, include this atom */
259 else /* if it's a coarse grained loop, skip this atom */
264 /* both have double identity, get hybrid scaling factor */
267 hybscal = weight_product;
271 hybscal = 1.0 - hybscal;
282 rsq = dx*dx+dy*dy+dz*dz;
283 rinv = gmx_invsqrt(rsq);
290 if (bExactCutoff && rsq > rcutoff2)
295 if (ielec == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE || ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE)
302 nnn = table_nelements*n0;
305 /* Coulomb interaction. ielec==0 means no interaction */
306 if (ielec != GMX_NBKERNEL_ELEC_NONE)
312 case GMX_NBKERNEL_ELEC_NONE:
315 case GMX_NBKERNEL_ELEC_COULOMB:
316 /* Vanilla cutoff coulomb */
318 felec = velec*rinvsq;
321 case GMX_NBKERNEL_ELEC_REACTIONFIELD:
323 velec = qq*(rinv+fr->k_rf*rsq-fr->c_rf);
324 felec = qq*(rinv*rinvsq-2.0*fr->k_rf);
327 case GMX_NBKERNEL_ELEC_CUBICSPLINETABLE:
328 /* Tabulated coulomb */
331 Geps = eps*VFtab[nnn+2];
332 Heps2 = eps2*VFtab[nnn+3];
335 FF = Fp+Geps+2.0*Heps2;
337 felec = -qq*FF*tabscale*rinv;
340 case GMX_NBKERNEL_ELEC_GENERALIZEDBORN:
342 gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
345 case GMX_NBKERNEL_ELEC_EWALD:
346 ewrt = rsq*rinv*ewtabscale;
350 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
351 rinvcorr = (fr->coulomb_modifier == eintmodPOTSHIFT) ? rinv-fr->ic->sh_ewald : rinv;
352 velec = qq*(rinvcorr-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
353 felec = qq*rinv*(rinvsq-felec);
357 gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
360 if (fr->coulomb_modifier == eintmodPOTSWITCH)
362 d = rsq*rinv-fr->rcoulomb_switch;
363 d = (d > 0.0) ? d : 0.0;
365 sw = 1.0+d2*d*(elec_swV3+d*(elec_swV4+d*elec_swV5));
366 dsw = d2*(elec_swF2+d*(elec_swF3+d*elec_swF4));
367 /* Apply switch function. Note that felec=f/r since it will be multiplied
368 * by the i-j displacement vector. This means felec'=f'/r=-(v*sw)'/r=
369 * -(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=felec*sw-v*dsw/r
371 felec = felec*sw - rinv*velec*dsw;
372 /* Once we have used velec to update felec we can modify velec too */
375 if (bExactElecCutoff)
377 felec = (rsq <= rcoulomb2) ? felec : 0.0;
378 velec = (rsq <= rcoulomb2) ? velec : 0.0;
381 } /* End of coulomb interactions */
384 /* VdW interaction. ivdw==0 means no interaction */
385 if (ivdw != GMX_NBKERNEL_VDW_NONE)
387 tj = nti+nvdwparam*type[jnr];
391 case GMX_NBKERNEL_VDW_NONE:
394 case GMX_NBKERNEL_VDW_LENNARDJONES:
395 /* Vanilla Lennard-Jones cutoff */
397 c12 = vdwparam[tj+1];
398 rinvsix = rinvsq*rinvsq*rinvsq;
399 vvdw_disp = c6*rinvsix;
400 vvdw_rep = c12*rinvsix*rinvsix;
401 fvdw = (vvdw_rep-vvdw_disp)*rinvsq;
402 if (fr->vdw_modifier == eintmodPOTSHIFT)
404 vvdw = (vvdw_rep + c12*sh_repulsion)/12.0 - (vvdw_disp + c6*sh_dispersion)/6.0;
408 vvdw = vvdw_rep/12.0-vvdw_disp/6.0;
412 case GMX_NBKERNEL_VDW_BUCKINGHAM:
415 cexp1 = vdwparam[tj+1];
416 cexp2 = vdwparam[tj+2];
418 rinvsix = rinvsq*rinvsq*rinvsq;
419 vvdw_disp = c6*rinvsix;
421 vvdw_rep = cexp1*exp(-br);
422 fvdw = (br*vvdw_rep-vvdw_disp)*rinvsq;
423 if (fr->vdw_modifier == eintmodPOTSHIFT)
425 vvdw = (vvdw_rep-cexp1*exp(-cexp2*rvdw)) - (vvdw_disp + c6*sh_dispersion)/6.0;
429 vvdw = vvdw_rep-vvdw_disp/6.0;
433 case GMX_NBKERNEL_VDW_CUBICSPLINETABLE:
436 c12 = vdwparam[tj+1];
439 Geps = eps*VFtab[nnn+6];
440 Heps2 = eps2*VFtab[nnn+7];
443 FF = Fp+Geps+2.0*Heps2;
448 Geps = eps*VFtab[nnn+10];
449 Heps2 = eps2*VFtab[nnn+11];
452 FF = Fp+Geps+2.0*Heps2;
455 fvdw = -(fijD+fijR)*tabscale*rinv;
456 vvdw = vvdw_disp + vvdw_rep;
460 gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
463 if (fr->vdw_modifier == eintmodPOTSWITCH)
465 d = rsq*rinv-fr->rvdw_switch;
466 d = (d > 0.0) ? d : 0.0;
468 sw = 1.0+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
469 dsw = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
470 /* See coulomb interaction for the force-switch formula */
471 fvdw = fvdw*sw - rinv*vvdw*dsw;
476 fvdw = (rsq <= rvdw2) ? fvdw : 0.0;
477 vvdw = (rsq <= rvdw2) ? vvdw : 0.0;
480 } /* end VdW interactions */
484 if (!bCG && force_cap > 0 && (fabs(fscal) > force_cap))
486 fscal = force_cap*fscal/fabs(fscal);
497 f[j3+0] = f[j3+0] - tx;
498 f[j3+1] = f[j3+1] - ty;
499 f[j3+2] = f[j3+2] - tz;
502 f[ii3+0] = f[ii3+0] + fix;
503 f[ii3+1] = f[ii3+1] + fiy;
504 f[ii3+2] = f[ii3+2] + fiz;
505 fshift[is3] = fshift[is3]+fix;
506 fshift[is3+1] = fshift[is3+1]+fiy;
507 fshift[is3+2] = fshift[is3+2]+fiz;
508 ggid = nlist->gid[n];
509 velecgrp[ggid] += vctot;
510 vvdwgrp[ggid] += vvdwtot;
512 /* Estimate flops, average for generic adress kernel:
513 * 14 flops per outer iteration
514 * 54 flops per inner iteration
516 inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC_ADRESS, nlist->nri*14 + nlist->jindex[n]*54);