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) 2012, The GROMACS development team,
7 * check out http://www.gromacs.org for more information.
8 * Copyright (c) 2012, by the GROMACS development team, led by
9 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
10 * others, as listed in the AUTHORS file in the top-level source
11 * directory and at http://www.gromacs.org.
13 * GROMACS is free software; you can redistribute it and/or
14 * modify it under the terms of the GNU Lesser General Public License
15 * as published by the Free Software Foundation; either version 2.1
16 * of the License, or (at your option) any later version.
18 * GROMACS is distributed in the hope that it will be useful,
19 * but WITHOUT ANY WARRANTY; without even the implied warranty of
20 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
21 * Lesser General Public License for more details.
23 * You should have received a copy of the GNU Lesser General Public
24 * License along with GROMACS; if not, see
25 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
26 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
28 * If you want to redistribute modifications to GROMACS, please
29 * consider that scientific software is very special. Version
30 * control is crucial - bugs must be traceable. We will be happy to
31 * consider code for inclusion in the official distribution, but
32 * derived work must not be called official GROMACS. Details are found
33 * in the README & COPYING files - if they are missing, get the
34 * official version at http://www.gromacs.org.
36 * To help us fund GROMACS development, we humbly ask that you cite
37 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "types/simple.h"
48 #include "nb_generic_adress.h"
51 #include "nonbonded.h"
52 #include "nb_kernel.h"
54 #define ALMOST_ZERO 1e-30
55 #define ALMOST_ONE 1-(1e-30)
57 gmx_nb_generic_adress_kernel(t_nblist * nlist,
62 nb_kernel_data_t * kernel_data,
65 int nri, ntype, table_nelements, ielec, ivdw;
66 real facel, gbtabscale;
67 int n, ii, is3, ii3, k, nj0, nj1, jnr, j3, ggid, nnn, n0;
69 real fscal, felec, fvdw, velec, vvdw, tx, ty, tz;
75 real rt, r, eps, eps2, Y, F, Geps, Heps2, VV, FF, Fp, fijD, fijR;
78 real vvdw_rep, vvdw_disp;
79 real ix, iy, iz, fix, fiy, fiz;
81 real dx, dy, dz, rsq, rinv;
82 real c6, c12, cexp1, cexp2, br;
96 real ewtabscale, eweps, sh_ewald, ewrt, ewtabhalfspace;
98 real rcoulomb2, rvdw, rvdw2, sh_invrc6;
99 real rcutoff, rcutoff2;
100 real rswitch_elec, rswitch_vdw, d, d2, sw, dsw, rinvcorr;
101 real elec_swV3, elec_swV4, elec_swV5, elec_swF2, elec_swF3, elec_swF4;
102 real vdw_swV3, vdw_swV4, vdw_swV5, vdw_swF2, vdw_swF3, vdw_swF4;
103 gmx_bool bExactElecCutoff, bExactVdwCutoff, bExactCutoff;
109 real hybscal; /* the multiplicator to the force for hybrid interactions*/
116 force_cap = fr->adress_ex_forcecap;
120 ielec = nlist->ielec;
123 fshift = fr->fshift[0];
124 velecgrp = kernel_data->energygrp_elec;
125 vvdwgrp = kernel_data->energygrp_vdw;
126 tabscale = kernel_data->table_elec_vdw->scale;
127 VFtab = kernel_data->table_elec_vdw->data;
129 sh_ewald = fr->ic->sh_ewald;
130 ewtab = fr->ic->tabq_coul_FDV0;
131 ewtabscale = fr->ic->tabq_scale;
132 ewtabhalfspace = 0.5/ewtabscale;
134 rcoulomb2 = fr->rcoulomb*fr->rcoulomb;
137 sh_invrc6 = fr->ic->sh_invrc6;
139 if (fr->coulomb_modifier == eintmodPOTSWITCH)
141 d = fr->rcoulomb-fr->rcoulomb_switch;
142 elec_swV3 = -10.0/(d*d*d);
143 elec_swV4 = 15.0/(d*d*d*d);
144 elec_swV5 = -6.0/(d*d*d*d*d);
145 elec_swF2 = -30.0/(d*d*d);
146 elec_swF3 = 60.0/(d*d*d*d);
147 elec_swF4 = -30.0/(d*d*d*d*d);
151 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
152 elec_swV3 = elec_swV4 = elec_swV5 = elec_swF2 = elec_swF3 = elec_swF4 = 0.0;
154 if (fr->vdw_modifier == eintmodPOTSWITCH)
156 d = fr->rvdw-fr->rvdw_switch;
157 vdw_swV3 = -10.0/(d*d*d);
158 vdw_swV4 = 15.0/(d*d*d*d);
159 vdw_swV5 = -6.0/(d*d*d*d*d);
160 vdw_swF2 = -30.0/(d*d*d);
161 vdw_swF3 = 60.0/(d*d*d*d);
162 vdw_swF4 = -30.0/(d*d*d*d*d);
166 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
167 vdw_swV3 = vdw_swV4 = vdw_swV5 = vdw_swF2 = vdw_swF3 = vdw_swF4 = 0.0;
170 bExactElecCutoff = (fr->coulomb_modifier != eintmodNONE) || fr->eeltype == eelRF_ZERO;
171 bExactVdwCutoff = (fr->vdw_modifier != eintmodNONE);
172 bExactCutoff = bExactElecCutoff || bExactVdwCutoff;
176 rcutoff = ( fr->rcoulomb > fr->rvdw ) ? fr->rcoulomb : fr->rvdw;
177 rcutoff2 = rcutoff*rcutoff;
181 /* Fix warnings for stupid compilers */
182 rcutoff = rcutoff2 = 1e30;
185 /* avoid compiler warnings for cases that cannot happen */
190 /* 3 VdW parameters for buckingham, otherwise 2 */
191 nvdwparam = (ivdw == GMX_NBKERNEL_VDW_BUCKINGHAM) ? 3 : 2;
192 table_nelements = 12;
194 charge = mdatoms->chargeA;
195 type = mdatoms->typeA;
197 shiftvec = fr->shift_vec[0];
201 for (n = 0; (n < nlist->nri); n++)
203 is3 = 3*nlist->shift[n];
205 shY = shiftvec[is3+1];
206 shZ = shiftvec[is3+2];
207 nj0 = nlist->jindex[n];
208 nj1 = nlist->jindex[n+1];
214 iq = facel*charge[ii];
215 nti = nvdwparam*ntype*type[ii];
222 /* We need to find out if this i atom is part of an
223 all-atom or CG energy group */
224 egp_nr = mdatoms->cENER[ii];
225 bCG = !fr->adress_group_explicit[egp_nr];
229 if ((!bCG) && weight_cg1 < ALMOST_ZERO)
234 for (k = nj0; (k < nj1); k++)
236 jnr = nlist->jjnr[k];
237 weight_cg2 = wf[jnr];
238 weight_product = weight_cg1*weight_cg2;
240 if (weight_product < ALMOST_ZERO)
242 /* if it's a explicit loop, skip this atom */
247 else /* if it's a coarse grained loop, include this atom */
252 else if (weight_product >= ALMOST_ONE)
255 /* if it's a explicit loop, include this atom */
260 else /* if it's a coarse grained loop, skip this atom */
265 /* both have double identity, get hybrid scaling factor */
268 hybscal = weight_product;
272 hybscal = 1.0 - hybscal;
283 rsq = dx*dx+dy*dy+dz*dz;
284 rinv = gmx_invsqrt(rsq);
291 if (bExactCutoff && rsq > rcutoff2)
296 if (ielec == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE || ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE)
303 nnn = table_nelements*n0;
306 /* Coulomb interaction. ielec==0 means no interaction */
307 if (ielec != GMX_NBKERNEL_ELEC_NONE)
313 case GMX_NBKERNEL_ELEC_NONE:
316 case GMX_NBKERNEL_ELEC_COULOMB:
317 /* Vanilla cutoff coulomb */
319 felec = velec*rinvsq;
322 case GMX_NBKERNEL_ELEC_REACTIONFIELD:
324 velec = qq*(rinv+fr->k_rf*rsq-fr->c_rf);
325 felec = qq*(rinv*rinvsq-2.0*fr->k_rf);
328 case GMX_NBKERNEL_ELEC_CUBICSPLINETABLE:
329 /* Tabulated coulomb */
332 Geps = eps*VFtab[nnn+2];
333 Heps2 = eps2*VFtab[nnn+3];
336 FF = Fp+Geps+2.0*Heps2;
338 felec = -qq*FF*tabscale*rinv;
341 case GMX_NBKERNEL_ELEC_GENERALIZEDBORN:
343 gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
346 case GMX_NBKERNEL_ELEC_EWALD:
347 ewrt = rsq*rinv*ewtabscale;
351 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
352 rinvcorr = (fr->coulomb_modifier == eintmodPOTSHIFT) ? rinv-fr->ic->sh_ewald : rinv;
353 velec = qq*(rinvcorr-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
354 felec = qq*rinv*(rinvsq-felec);
358 gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
361 if (fr->coulomb_modifier == eintmodPOTSWITCH)
363 d = rsq*rinv-fr->rcoulomb_switch;
364 d = (d > 0.0) ? d : 0.0;
366 sw = 1.0+d2*d*(elec_swV3+d*(elec_swV4+d*elec_swV5));
367 dsw = d2*(elec_swF2+d*(elec_swF3+d*elec_swF4));
368 /* Apply switch function. Note that felec=f/r since it will be multiplied
369 * by the i-j displacement vector. This means felec'=f'/r=-(v*sw)'/r=
370 * -(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=felec*sw-v*dsw/r
372 felec = felec*sw - rinv*velec*dsw;
373 /* Once we have used velec to update felec we can modify velec too */
376 if (bExactElecCutoff)
378 felec = (rsq <= rcoulomb2) ? felec : 0.0;
379 velec = (rsq <= rcoulomb2) ? velec : 0.0;
382 } /* End of coulomb interactions */
385 /* VdW interaction. ivdw==0 means no interaction */
386 if (ivdw != GMX_NBKERNEL_VDW_NONE)
388 tj = nti+nvdwparam*type[jnr];
392 case GMX_NBKERNEL_VDW_NONE:
395 case GMX_NBKERNEL_VDW_LENNARDJONES:
396 /* Vanilla Lennard-Jones cutoff */
398 c12 = vdwparam[tj+1];
399 rinvsix = rinvsq*rinvsq*rinvsq;
400 vvdw_disp = c6*rinvsix;
401 vvdw_rep = c12*rinvsix*rinvsix;
402 fvdw = (vvdw_rep-vvdw_disp)*rinvsq;
403 if (fr->vdw_modifier == eintmodPOTSHIFT)
405 vvdw = (vvdw_rep-c12*sh_invrc6*sh_invrc6)*(1.0/12.0)-(vvdw_disp-c6*sh_invrc6)*(1.0/6.0);
409 vvdw = vvdw_rep/12.0-vvdw_disp/6.0;
413 case GMX_NBKERNEL_VDW_BUCKINGHAM:
416 cexp1 = vdwparam[tj+1];
417 cexp2 = vdwparam[tj+2];
419 rinvsix = rinvsq*rinvsq*rinvsq;
420 vvdw_disp = c6*rinvsix;
422 vvdw_rep = cexp1*exp(-br);
423 fvdw = (br*vvdw_rep-vvdw_disp)*rinvsq;
424 if (fr->vdw_modifier == eintmodPOTSHIFT)
426 vvdw = (vvdw_rep-cexp1*exp(-cexp2*rvdw))-(vvdw_disp-c6*sh_invrc6)/6.0;
430 vvdw = vvdw_rep-vvdw_disp/6.0;
434 case GMX_NBKERNEL_VDW_CUBICSPLINETABLE:
437 c12 = vdwparam[tj+1];
440 Geps = eps*VFtab[nnn+6];
441 Heps2 = eps2*VFtab[nnn+7];
444 FF = Fp+Geps+2.0*Heps2;
449 Geps = eps*VFtab[nnn+10];
450 Heps2 = eps2*VFtab[nnn+11];
453 FF = Fp+Geps+2.0*Heps2;
456 fvdw = -(fijD+fijR)*tabscale*rinv;
457 vvdw = vvdw_disp + vvdw_rep;
461 gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
464 if (fr->vdw_modifier == eintmodPOTSWITCH)
466 d = rsq*rinv-fr->rvdw_switch;
467 d = (d > 0.0) ? d : 0.0;
469 sw = 1.0+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
470 dsw = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
471 /* See coulomb interaction for the force-switch formula */
472 fvdw = fvdw*sw - rinv*vvdw*dsw;
477 fvdw = (rsq <= rvdw2) ? fvdw : 0.0;
478 vvdw = (rsq <= rvdw2) ? vvdw : 0.0;
481 } /* end VdW interactions */
485 if (!bCG && force_cap > 0 && (fabs(fscal) > force_cap))
487 fscal = force_cap*fscal/fabs(fscal);
498 f[j3+0] = f[j3+0] - tx;
499 f[j3+1] = f[j3+1] - ty;
500 f[j3+2] = f[j3+2] - tz;
503 f[ii3+0] = f[ii3+0] + fix;
504 f[ii3+1] = f[ii3+1] + fiy;
505 f[ii3+2] = f[ii3+2] + fiz;
506 fshift[is3] = fshift[is3]+fix;
507 fshift[is3+1] = fshift[is3+1]+fiy;
508 fshift[is3+2] = fshift[is3+2]+fiz;
509 ggid = nlist->gid[n];
510 velecgrp[ggid] += vctot;
511 vvdwgrp[ggid] += vvdwtot;
513 /* Estimate flops, average for generic kernel:
514 * 12 flops per outer iteration
515 * 50 flops per inner iteration
517 inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC, nlist->nri*12 + nlist->jindex[n]*50);