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) 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.
45 #include "nonbonded.h"
46 #include "nb_kernel.h"
49 #include "nb_free_energy.h"
51 #include "gmx_fatal.h"
54 gmx_nb_free_energy_kernel(const t_nblist * gmx_restrict nlist,
55 rvec * gmx_restrict xx,
56 rvec * gmx_restrict ff,
57 t_forcerec * gmx_restrict fr,
58 const t_mdatoms * gmx_restrict mdatoms,
59 nb_kernel_data_t * gmx_restrict kernel_data,
60 t_nrnb * gmx_restrict nrnb)
66 int i, j, n, ii, is3, ii3, k, nj0, nj1, jnr, j3, ggid;
68 real Fscal, FscalC[NSTATES], FscalV[NSTATES], tx, ty, tz;
69 real Vcoul[NSTATES], Vvdw[NSTATES];
70 real rinv6, r, rt, rtC, rtV;
72 real qq[NSTATES], vctot, krsq;
73 int ntiA, ntiB, tj[NSTATES];
74 real Vvdw6, Vvdw12, vvtot;
75 real ix, iy, iz, fix, fiy, fiz;
76 real dx, dy, dz, rsq, rinv;
77 real c6[NSTATES], c12[NSTATES], c6grid;
78 real LFC[NSTATES], LFV[NSTATES], DLF[NSTATES];
79 double dvdl_coul, dvdl_vdw;
80 real lfac_coul[NSTATES], dlfac_coul[NSTATES], lfac_vdw[NSTATES], dlfac_vdw[NSTATES];
81 real sigma6[NSTATES], alpha_vdw_eff, alpha_coul_eff, sigma2_def, sigma2_min;
82 real rp, rpm2, rC, rV, rinvC, rpinvC, rinvV, rpinvV;
83 real sigma2[NSTATES], sigma_pow[NSTATES], sigma_powm2[NSTATES], rs, rs2;
84 int do_tab, tab_elemsize;
85 int n0, n1C, n1V, nnn;
86 real Y, F, G, H, Fp, Geps, Heps2, epsC, eps2C, epsV, eps2V, VV, FF;
97 const real * shiftvec;
101 const real * VFtab = NULL;
104 real facel, krf, crf;
105 const real * chargeA;
106 const real * chargeB;
107 real sigma6_min, sigma6_def, lam_power, sc_power, sc_r_power;
108 real alpha_coul, alpha_vdw, lambda_coul, lambda_vdw, ewc_lj;
109 const real * nbfp, *nbfp_grid;
113 gmx_bool bDoForces, bDoShiftForces, bDoPotential;
114 real rcoulomb, rvdw, sh_invrc6;
115 gmx_bool bExactElecCutoff, bExactVdwCutoff, bExactCutoffAll;
116 gmx_bool bEwald, bEwaldLJ;
118 const real * tab_ewald_F_lj;
119 const real * tab_ewald_V_lj;
120 real d, d2, sw, dsw, rinvcorr;
121 real elec_swV3, elec_swV4, elec_swV5, elec_swF2, elec_swF3, elec_swF4;
122 real vdw_swV3, vdw_swV4, vdw_swV5, vdw_swF2, vdw_swF3, vdw_swF4;
123 gmx_bool bConvertEwaldToCoulomb, bConvertLJEwaldToLJ6;
124 gmx_bool bComputeVdwInteraction, bComputeElecInteraction;
127 real ewrt, eweps, ewtabscale, ewtabhalfspace, sh_ewald;
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;
133 tab_ewald_F_lj = fr->ic->tabq_vdw_F;
134 tab_ewald_V_lj = fr->ic->tabq_vdw_V;
139 fshift = fr->fshift[0];
143 jindex = nlist->jindex;
145 icoul = nlist->ielec;
147 shift = nlist->shift;
150 shiftvec = fr->shift_vec[0];
151 chargeA = mdatoms->chargeA;
152 chargeB = mdatoms->chargeB;
156 ewc_lj = fr->ewaldcoeff_lj;
157 Vc = kernel_data->energygrp_elec;
158 typeA = mdatoms->typeA;
159 typeB = mdatoms->typeB;
162 nbfp_grid = fr->ljpme_c6grid;
163 Vv = kernel_data->energygrp_vdw;
164 lambda_coul = kernel_data->lambda[efptCOUL];
165 lambda_vdw = kernel_data->lambda[efptVDW];
166 dvdl = kernel_data->dvdl;
167 alpha_coul = fr->sc_alphacoul;
168 alpha_vdw = fr->sc_alphavdw;
169 lam_power = fr->sc_power;
170 sc_r_power = fr->sc_r_power;
171 sigma6_def = fr->sc_sigma6_def;
172 sigma6_min = fr->sc_sigma6_min;
173 bDoForces = kernel_data->flags & GMX_NONBONDED_DO_FORCE;
174 bDoShiftForces = kernel_data->flags & GMX_NONBONDED_DO_SHIFTFORCE;
175 bDoPotential = kernel_data->flags & GMX_NONBONDED_DO_POTENTIAL;
177 rcoulomb = fr->rcoulomb;
178 sh_ewald = fr->ic->sh_ewald;
180 sh_invrc6 = fr->ic->sh_invrc6;
182 if (fr->coulomb_modifier == eintmodPOTSWITCH)
184 d = fr->rcoulomb-fr->rcoulomb_switch;
185 elec_swV3 = -10.0/(d*d*d);
186 elec_swV4 = 15.0/(d*d*d*d);
187 elec_swV5 = -6.0/(d*d*d*d*d);
188 elec_swF2 = -30.0/(d*d*d);
189 elec_swF3 = 60.0/(d*d*d*d);
190 elec_swF4 = -30.0/(d*d*d*d*d);
194 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
195 elec_swV3 = elec_swV4 = elec_swV5 = elec_swF2 = elec_swF3 = elec_swF4 = 0.0;
198 if (fr->vdw_modifier == eintmodPOTSWITCH)
200 d = fr->rvdw-fr->rvdw_switch;
201 vdw_swV3 = -10.0/(d*d*d);
202 vdw_swV4 = 15.0/(d*d*d*d);
203 vdw_swV5 = -6.0/(d*d*d*d*d);
204 vdw_swF2 = -30.0/(d*d*d);
205 vdw_swF3 = 60.0/(d*d*d*d);
206 vdw_swF4 = -30.0/(d*d*d*d*d);
210 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
211 vdw_swV3 = vdw_swV4 = vdw_swV5 = vdw_swF2 = vdw_swF3 = vdw_swF4 = 0.0;
214 if (fr->cutoff_scheme == ecutsVERLET)
216 const interaction_const_t *ic;
219 if (EVDW_PME(ic->vdwtype))
221 ivdw = GMX_NBKERNEL_VDW_LJEWALD;
225 ivdw = GMX_NBKERNEL_VDW_LENNARDJONES;
228 if (ic->eeltype == eelCUT || EEL_RF(ic->eeltype))
230 icoul = GMX_NBKERNEL_ELEC_REACTIONFIELD;
232 else if (EEL_PME_EWALD(ic->eeltype))
234 icoul = GMX_NBKERNEL_ELEC_EWALD;
238 gmx_incons("Unsupported eeltype with Verlet and free-energy");
241 bExactElecCutoff = TRUE;
242 bExactVdwCutoff = TRUE;
246 bExactElecCutoff = (fr->coulomb_modifier != eintmodNONE) || fr->eeltype == eelRF_ZERO;
247 bExactVdwCutoff = (fr->vdw_modifier != eintmodNONE);
250 bExactCutoffAll = (bExactElecCutoff && bExactVdwCutoff);
251 rcutoff_max2 = max(fr->rcoulomb, fr->rvdw);
252 rcutoff_max2 = rcutoff_max2*rcutoff_max2;
254 bEwald = (icoul == GMX_NBKERNEL_ELEC_EWALD);
255 bEwaldLJ = (ivdw == GMX_NBKERNEL_VDW_LJEWALD);
257 /* For Ewald/PME interactions we cannot easily apply the soft-core component to
258 * reciprocal space. When we use vanilla (not switch/shift) Ewald interactions, we
259 * can apply the small trick of subtracting the _reciprocal_ space contribution
260 * in this kernel, and instead apply the free energy interaction to the 1/r
261 * (standard coulomb) interaction.
263 * However, we cannot use this approach for switch-modified since we would then
264 * effectively end up evaluating a significantly different interaction here compared to the
265 * normal (non-free-energy) kernels, either by applying a cutoff at a different
266 * position than what the user requested, or by switching different
267 * things (1/r rather than short-range Ewald). For these settings, we just
268 * use the traditional short-range Ewald interaction in that case.
270 bConvertEwaldToCoulomb = (bEwald && (fr->coulomb_modifier != eintmodPOTSWITCH));
271 /* For now the below will always be true (since LJ-PME only works with Shift in Gromacs-5.0),
272 * but writing it this way means we stay in sync with coulomb, and it avoids future bugs.
274 bConvertLJEwaldToLJ6 = (bEwaldLJ && (fr->vdw_modifier != eintmodPOTSWITCH));
276 /* fix compiler warnings */
285 /* Lambda factor for state A, 1-lambda*/
286 LFC[STATE_A] = 1.0 - lambda_coul;
287 LFV[STATE_A] = 1.0 - lambda_vdw;
289 /* Lambda factor for state B, lambda*/
290 LFC[STATE_B] = lambda_coul;
291 LFV[STATE_B] = lambda_vdw;
293 /*derivative of the lambda factor for state A and B */
297 for (i = 0; i < NSTATES; i++)
299 lfac_coul[i] = (lam_power == 2 ? (1-LFC[i])*(1-LFC[i]) : (1-LFC[i]));
300 dlfac_coul[i] = DLF[i]*lam_power/sc_r_power*(lam_power == 2 ? (1-LFC[i]) : 1);
301 lfac_vdw[i] = (lam_power == 2 ? (1-LFV[i])*(1-LFV[i]) : (1-LFV[i]));
302 dlfac_vdw[i] = DLF[i]*lam_power/sc_r_power*(lam_power == 2 ? (1-LFV[i]) : 1);
305 sigma2_def = pow(sigma6_def, 1.0/3.0);
306 sigma2_min = pow(sigma6_min, 1.0/3.0);
308 /* Ewald (not PME) table is special (icoul==enbcoulFEWALD) */
310 do_tab = (icoul == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE ||
311 ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE);
314 tabscale = kernel_data->table_elec_vdw->scale;
315 VFtab = kernel_data->table_elec_vdw->data;
316 /* we always use the combined table here */
320 for (n = 0; (n < nri); n++)
322 int npair_within_cutoff;
324 npair_within_cutoff = 0;
328 shY = shiftvec[is3+1];
329 shZ = shiftvec[is3+2];
337 iqA = facel*chargeA[ii];
338 iqB = facel*chargeB[ii];
339 ntiA = 2*ntype*typeA[ii];
340 ntiB = 2*ntype*typeB[ii];
347 for (k = nj0; (k < nj1); k++)
354 rsq = dx*dx + dy*dy + dz*dz;
356 if (bExactCutoffAll && rsq >= rcutoff_max2)
358 /* We save significant time by skipping all code below.
359 * Note that with soft-core interactions, the actual cut-off
360 * check might be different. But since the soft-core distance
361 * is always larger than r, checking on r here is safe.
365 npair_within_cutoff++;
369 rinv = gmx_invsqrt(rsq);
374 /* The force at r=0 is zero, because of symmetry.
375 * But note that the potential is in general non-zero,
376 * since the soft-cored r will be non-zero.
382 if (sc_r_power == 6.0)
384 rpm2 = rsq*rsq; /* r4 */
385 rp = rpm2*rsq; /* r6 */
387 else if (sc_r_power == 48.0)
389 rp = rsq*rsq*rsq; /* r6 */
390 rp = rp*rp; /* r12 */
391 rp = rp*rp; /* r24 */
392 rp = rp*rp; /* r48 */
393 rpm2 = rp/rsq; /* r46 */
397 rp = pow(r, sc_r_power); /* not currently supported as input, but can handle it */
403 qq[STATE_A] = iqA*chargeA[jnr];
404 qq[STATE_B] = iqB*chargeB[jnr];
406 tj[STATE_A] = ntiA+2*typeA[jnr];
407 tj[STATE_B] = ntiB+2*typeB[jnr];
409 if (nlist->excl_fep == NULL || nlist->excl_fep[k])
411 c6[STATE_A] = nbfp[tj[STATE_A]];
412 c6[STATE_B] = nbfp[tj[STATE_B]];
414 for (i = 0; i < NSTATES; i++)
416 c12[i] = nbfp[tj[i]+1];
417 if ((c6[i] > 0) && (c12[i] > 0))
419 /* c12 is stored scaled with 12.0 and c6 is scaled with 6.0 - correct for this */
420 sigma6[i] = 0.5*c12[i]/c6[i];
421 sigma2[i] = pow(sigma6[i], 1.0/3.0);
422 /* should be able to get rid of this ^^^ internal pow call eventually. Will require agreement on
423 what data to store externally. Can't be fixed without larger scale changes, so not 4.6 */
424 if (sigma6[i] < sigma6_min) /* for disappearing coul and vdw with soft core at the same time */
426 sigma6[i] = sigma6_min;
427 sigma2[i] = sigma2_min;
432 sigma6[i] = sigma6_def;
433 sigma2[i] = sigma2_def;
435 if (sc_r_power == 6.0)
437 sigma_pow[i] = sigma6[i];
438 sigma_powm2[i] = sigma6[i]/sigma2[i];
440 else if (sc_r_power == 48.0)
442 sigma_pow[i] = sigma6[i]*sigma6[i]; /* sigma^12 */
443 sigma_pow[i] = sigma_pow[i]*sigma_pow[i]; /* sigma^24 */
444 sigma_pow[i] = sigma_pow[i]*sigma_pow[i]; /* sigma^48 */
445 sigma_powm2[i] = sigma_pow[i]/sigma2[i];
448 { /* not really supported as input, but in here for testing the general case*/
449 sigma_pow[i] = pow(sigma2[i], sc_r_power/2);
450 sigma_powm2[i] = sigma_pow[i]/(sigma2[i]);
454 /* only use softcore if one of the states has a zero endstate - softcore is for avoiding infinities!*/
455 if ((c12[STATE_A] > 0) && (c12[STATE_B] > 0))
462 alpha_vdw_eff = alpha_vdw;
463 alpha_coul_eff = alpha_coul;
466 for (i = 0; i < NSTATES; i++)
473 /* Only spend time on A or B state if it is non-zero */
474 if ( (qq[i] != 0) || (c6[i] != 0) || (c12[i] != 0) )
476 /* this section has to be inside the loop because of the dependence on sigma_pow */
477 rpinvC = 1.0/(alpha_coul_eff*lfac_coul[i]*sigma_pow[i]+rp);
478 rinvC = pow(rpinvC, 1.0/sc_r_power);
481 rpinvV = 1.0/(alpha_vdw_eff*lfac_vdw[i]*sigma_pow[i]+rp);
482 rinvV = pow(rpinvV, 1.0/sc_r_power);
491 n1C = tab_elemsize*n0;
497 n1V = tab_elemsize*n0;
500 /* Only process the coulomb interactions if we have charges,
501 * and if we either include all entries in the list (no cutoff
502 * used in the kernel), or if we are within the cutoff.
504 bComputeElecInteraction = !bExactElecCutoff ||
505 ( bConvertEwaldToCoulomb && r < rcoulomb) ||
506 (!bConvertEwaldToCoulomb && rC < rcoulomb);
508 if ( (qq[i] != 0) && bComputeElecInteraction)
512 case GMX_NBKERNEL_ELEC_COULOMB:
514 Vcoul[i] = qq[i]*rinvC;
515 FscalC[i] = Vcoul[i];
516 /* The shift for the Coulomb potential is stored in
517 * the RF parameter c_rf, which is 0 without shift
519 Vcoul[i] -= qq[i]*fr->ic->c_rf;
522 case GMX_NBKERNEL_ELEC_REACTIONFIELD:
524 Vcoul[i] = qq[i]*(rinvC + krf*rC*rC-crf);
525 FscalC[i] = qq[i]*(rinvC - 2.0*krf*rC*rC);
528 case GMX_NBKERNEL_ELEC_CUBICSPLINETABLE:
529 /* non-Ewald tabulated coulomb */
533 Geps = epsC*VFtab[nnn+2];
534 Heps2 = eps2C*VFtab[nnn+3];
537 FF = Fp+Geps+2.0*Heps2;
539 FscalC[i] = -qq[i]*tabscale*FF*rC;
542 case GMX_NBKERNEL_ELEC_GENERALIZEDBORN:
543 gmx_fatal(FARGS, "Free energy and GB not implemented.\n");
546 case GMX_NBKERNEL_ELEC_EWALD:
547 if (bConvertEwaldToCoulomb)
549 /* Ewald FEP is done only on the 1/r part */
550 Vcoul[i] = qq[i]*(rinvC-sh_ewald);
551 FscalC[i] = qq[i]*rinvC;
555 ewrt = rC*ewtabscale;
559 FscalC[i] = ewtab[ewitab]+eweps*ewtab[ewitab+1];
560 rinvcorr = rinvC-sh_ewald;
561 Vcoul[i] = qq[i]*(rinvcorr-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+FscalC[i])));
562 FscalC[i] = qq[i]*(rinvC-rC*FscalC[i]);
566 case GMX_NBKERNEL_ELEC_NONE:
572 gmx_incons("Invalid icoul in free energy kernel");
576 if (fr->coulomb_modifier == eintmodPOTSWITCH)
578 d = rC-fr->rcoulomb_switch;
579 d = (d > 0.0) ? d : 0.0;
581 sw = 1.0+d2*d*(elec_swV3+d*(elec_swV4+d*elec_swV5));
582 dsw = d2*(elec_swF2+d*(elec_swF3+d*elec_swF4));
584 FscalC[i] = FscalC[i]*sw - rC*Vcoul[i]*dsw;
587 FscalC[i] = (rC < rcoulomb) ? FscalC[i] : 0.0;
588 Vcoul[i] = (rC < rcoulomb) ? Vcoul[i] : 0.0;
592 /* Only process the VDW interactions if we have
593 * some non-zero parameters, and if we either
594 * include all entries in the list (no cutoff used
595 * in the kernel), or if we are within the cutoff.
597 bComputeVdwInteraction = !bExactVdwCutoff ||
598 ( bConvertLJEwaldToLJ6 && r < rvdw) ||
599 (!bConvertLJEwaldToLJ6 && rV < rvdw);
600 if ((c6[i] != 0 || c12[i] != 0) && bComputeVdwInteraction)
604 case GMX_NBKERNEL_VDW_LENNARDJONES:
605 case GMX_NBKERNEL_VDW_LJEWALD:
607 if (sc_r_power == 6.0)
613 rinv6 = pow(rinvV, 6.0);
616 Vvdw12 = c12[i]*rinv6*rinv6;
617 if (fr->vdw_modifier == eintmodPOTSHIFT)
619 Vvdw[i] = ( (Vvdw12-c12[i]*sh_invrc6*sh_invrc6)*(1.0/12.0)
620 -(Vvdw6-c6[i]*sh_invrc6)*(1.0/6.0));
624 Vvdw[i] = Vvdw12*(1.0/12.0) - Vvdw6*(1.0/6.0);
626 FscalV[i] = Vvdw12 - Vvdw6;
629 case GMX_NBKERNEL_VDW_BUCKINGHAM:
630 gmx_fatal(FARGS, "Buckingham free energy not supported.");
633 case GMX_NBKERNEL_VDW_CUBICSPLINETABLE:
639 Geps = epsV*VFtab[nnn+2];
640 Heps2 = eps2V*VFtab[nnn+3];
643 FF = Fp+Geps+2.0*Heps2;
645 FscalV[i] -= c6[i]*tabscale*FF*rV;
650 Geps = epsV*VFtab[nnn+6];
651 Heps2 = eps2V*VFtab[nnn+7];
654 FF = Fp+Geps+2.0*Heps2;
655 Vvdw[i] += c12[i]*VV;
656 FscalV[i] -= c12[i]*tabscale*FF*rV;
659 case GMX_NBKERNEL_VDW_NONE:
665 gmx_incons("Invalid ivdw in free energy kernel");
669 if (fr->vdw_modifier == eintmodPOTSWITCH)
671 d = rV-fr->rvdw_switch;
672 d = (d > 0.0) ? d : 0.0;
674 sw = 1.0+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
675 dsw = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
677 FscalV[i] = FscalV[i]*sw - rV*Vvdw[i]*dsw;
680 FscalV[i] = (rV < rvdw) ? FscalV[i] : 0.0;
681 Vvdw[i] = (rV < rvdw) ? Vvdw[i] : 0.0;
685 /* FscalC (and FscalV) now contain: dV/drC * rC
686 * Now we multiply by rC^-p, so it will be: dV/drC * rC^1-p
687 * Further down we first multiply by r^p-2 and then by
688 * the vector r, which in total gives: dV/drC * (r/rC)^1-p
695 /* Assemble A and B states */
696 for (i = 0; i < NSTATES; i++)
698 vctot += LFC[i]*Vcoul[i];
699 vvtot += LFV[i]*Vvdw[i];
701 Fscal += LFC[i]*FscalC[i]*rpm2;
702 Fscal += LFV[i]*FscalV[i]*rpm2;
704 dvdl_coul += Vcoul[i]*DLF[i] + LFC[i]*alpha_coul_eff*dlfac_coul[i]*FscalC[i]*sigma_pow[i];
705 dvdl_vdw += Vvdw[i]*DLF[i] + LFV[i]*alpha_vdw_eff*dlfac_vdw[i]*FscalV[i]*sigma_pow[i];
708 else if (icoul == GMX_NBKERNEL_ELEC_REACTIONFIELD)
710 /* For excluded pairs, which are only in this pair list when
711 * using the Verlet scheme, we don't use soft-core.
712 * The group scheme also doesn't soft-core for these.
713 * As there is no singularity, there is no need for soft-core.
723 for (i = 0; i < NSTATES; i++)
725 vctot += LFC[i]*qq[i]*VV;
726 Fscal += LFC[i]*qq[i]*FF;
727 dvdl_coul += DLF[i]*qq[i]*VV;
731 if (bConvertEwaldToCoulomb && ( !bExactElecCutoff || r < rcoulomb ) )
733 /* See comment in the preamble. When using Ewald interactions
734 * (unless we use a switch modifier) we subtract the reciprocal-space
735 * Ewald component here which made it possible to apply the free
736 * energy interaction to 1/r (vanilla coulomb short-range part)
737 * above. This gets us closer to the ideal case of applying
738 * the softcore to the entire electrostatic interaction,
739 * including the reciprocal-space component.
747 f_lr = ewtab[ewitab]+eweps*ewtab[ewitab+1];
748 v_lr = (ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+f_lr));
753 /* If we get here, the i particle (ii) has itself (jnr)
754 * in its neighborlist. This can only happen with the Verlet
755 * scheme, and corresponds to a self-interaction that will
756 * occur twice. Scale it down by 50% to only include it once.
761 for (i = 0; i < NSTATES; i++)
763 vctot -= LFC[i]*qq[i]*v_lr;
764 Fscal -= LFC[i]*qq[i]*f_lr;
765 dvdl_coul -= (DLF[i]*qq[i])*v_lr;
769 if (bConvertLJEwaldToLJ6 && (!bExactVdwCutoff || r < rvdw))
771 /* See comment in the preamble. When using LJ-Ewald interactions
772 * (unless we use a switch modifier) we subtract the reciprocal-space
773 * Ewald component here which made it possible to apply the free
774 * energy interaction to r^-6 (vanilla LJ6 short-range part)
775 * above. This gets us closer to the ideal case of applying
776 * the softcore to the entire VdW interaction,
777 * including the reciprocal-space component.
782 rs = rsq*rinv*ewtabscale;
785 f_lr = (1 - frac)*tab_ewald_F_lj[ri] + frac*tab_ewald_F_lj[ri+1];
787 VV = tab_ewald_V_lj[ri] - ewtabhalfspace*frac*(tab_ewald_F_lj[ri] + f_lr);
791 /* If we get here, the i particle (ii) has itself (jnr)
792 * in its neighborlist. This can only happen with the Verlet
793 * scheme, and corresponds to a self-interaction that will
794 * occur twice. Scale it down by 50% to only include it once.
799 for (i = 0; i < NSTATES; i++)
801 c6grid = nbfp_grid[tj[i]];
802 vvtot += LFV[i]*c6grid*VV*(1.0/6.0);
803 Fscal += LFV[i]*c6grid*FF*(1.0/6.0);
804 dvdl_vdw += (DLF[i]*c6grid)*VV*(1.0/6.0);
817 /* OpenMP atomics are expensive, but this kernels is also
818 * expensive, so we can take this hit, instead of using
819 * thread-local output buffers and extra reduction.
830 /* The atomics below are expensive with many OpenMP threads.
831 * Here unperturbed i-particles will usually only have a few
832 * (perturbed) j-particles in the list. Thus with a buffered list
833 * we can skip a significant number of i-reductions with a check.
835 if (npair_within_cutoff > 0)
851 fshift[is3+1] += fiy;
853 fshift[is3+2] += fiz;
867 dvdl[efptCOUL] += dvdl_coul;
869 dvdl[efptVDW] += dvdl_vdw;
871 /* Estimate flops, average for free energy stuff:
872 * 12 flops per outer iteration
873 * 150 flops per inner iteration
876 inc_nrnb(nrnb, eNR_NBKERNEL_FREE_ENERGY, nlist->nri*12 + nlist->jindex[n]*150);
880 nb_free_energy_evaluate_single(real r2, real sc_r_power, real alpha_coul, real alpha_vdw,
881 real tabscale, real *vftab,
882 real qqA, real c6A, real c12A, real qqB, real c6B, real c12B,
883 real LFC[2], real LFV[2], real DLF[2],
884 real lfac_coul[2], real lfac_vdw[2], real dlfac_coul[2], real dlfac_vdw[2],
885 real sigma6_def, real sigma6_min, real sigma2_def, real sigma2_min,
886 real *velectot, real *vvdwtot, real *dvdl)
888 real r, rp, rpm2, rtab, eps, eps2, Y, F, Geps, Heps2, Fp, VV, FF, fscal;
889 real qq[2], c6[2], c12[2], sigma6[2], sigma2[2], sigma_pow[2], sigma_powm2[2];
890 real alpha_coul_eff, alpha_vdw_eff, dvdl_coul, dvdl_vdw;
891 real rpinv, r_coul, r_vdw, velecsum, vvdwsum;
892 real fscal_vdw[2], fscal_elec[2];
893 real velec[2], vvdw[2];
903 if (sc_r_power == 6.0)
905 rpm2 = r2*r2; /* r4 */
906 rp = rpm2*r2; /* r6 */
908 else if (sc_r_power == 48.0)
910 rp = r2*r2*r2; /* r6 */
911 rp = rp*rp; /* r12 */
912 rp = rp*rp; /* r24 */
913 rp = rp*rp; /* r48 */
914 rpm2 = rp/r2; /* r46 */
918 rp = pow(r2, 0.5*sc_r_power); /* not currently supported as input, but can handle it */
922 /* Loop over state A(0) and B(1) */
923 for (i = 0; i < 2; i++)
925 if ((c6[i] > 0) && (c12[i] > 0))
927 /* The c6 & c12 coefficients now contain the constants 6.0 and 12.0, respectively.
928 * Correct for this by multiplying with (1/12.0)/(1/6.0)=6.0/12.0=0.5.
930 sigma6[i] = 0.5*c12[i]/c6[i];
931 sigma2[i] = pow(0.5*c12[i]/c6[i], 1.0/3.0);
932 /* should be able to get rid of this ^^^ internal pow call eventually. Will require agreement on
933 what data to store externally. Can't be fixed without larger scale changes, so not 5.0 */
934 if (sigma6[i] < sigma6_min) /* for disappearing coul and vdw with soft core at the same time */
936 sigma6[i] = sigma6_min;
937 sigma2[i] = sigma2_min;
942 sigma6[i] = sigma6_def;
943 sigma2[i] = sigma2_def;
945 if (sc_r_power == 6.0)
947 sigma_pow[i] = sigma6[i];
948 sigma_powm2[i] = sigma6[i]/sigma2[i];
950 else if (sc_r_power == 48.0)
952 sigma_pow[i] = sigma6[i]*sigma6[i]; /* sigma^12 */
953 sigma_pow[i] = sigma_pow[i]*sigma_pow[i]; /* sigma^24 */
954 sigma_pow[i] = sigma_pow[i]*sigma_pow[i]; /* sigma^48 */
955 sigma_powm2[i] = sigma_pow[i]/sigma2[i];
958 { /* not really supported as input, but in here for testing the general case*/
959 sigma_pow[i] = pow(sigma2[i], sc_r_power/2);
960 sigma_powm2[i] = sigma_pow[i]/(sigma2[i]);
964 /* only use softcore if one of the states has a zero endstate - softcore is for avoiding infinities!*/
965 if ((c12[0] > 0) && (c12[1] > 0))
972 alpha_vdw_eff = alpha_vdw;
973 alpha_coul_eff = alpha_coul;
976 /* Loop over A and B states again */
977 for (i = 0; i < 2; i++)
984 /* Only spend time on A or B state if it is non-zero */
985 if ( (qq[i] != 0) || (c6[i] != 0) || (c12[i] != 0) )
988 rpinv = 1.0/(alpha_coul_eff*lfac_coul[i]*sigma_pow[i]+rp);
989 r_coul = pow(rpinv, -1.0/sc_r_power);
991 /* Electrostatics table lookup data */
992 rtab = r_coul*tabscale;
1000 Geps = eps*vftab[ntab+2];
1001 Heps2 = eps2*vftab[ntab+3];
1004 FF = Fp+Geps+2.0*Heps2;
1005 velec[i] = qq[i]*VV;
1006 fscal_elec[i] = -qq[i]*FF*r_coul*rpinv*tabscale;
1009 rpinv = 1.0/(alpha_vdw_eff*lfac_vdw[i]*sigma_pow[i]+rp);
1010 r_vdw = pow(rpinv, -1.0/sc_r_power);
1011 /* Vdw table lookup data */
1012 rtab = r_vdw*tabscale;
1020 Geps = eps*vftab[ntab+6];
1021 Heps2 = eps2*vftab[ntab+7];
1024 FF = Fp+Geps+2.0*Heps2;
1026 fscal_vdw[i] = -c6[i]*FF;
1031 Geps = eps*vftab[ntab+10];
1032 Heps2 = eps2*vftab[ntab+11];
1035 FF = Fp+Geps+2.0*Heps2;
1036 vvdw[i] += c12[i]*VV;
1037 fscal_vdw[i] -= c12[i]*FF;
1038 fscal_vdw[i] *= r_vdw*rpinv*tabscale;
1041 /* Now we have velec[i], vvdw[i], and fscal[i] for both states */
1042 /* Assemble A and B states */
1048 for (i = 0; i < 2; i++)
1050 velecsum += LFC[i]*velec[i];
1051 vvdwsum += LFV[i]*vvdw[i];
1053 fscal += (LFC[i]*fscal_elec[i]+LFV[i]*fscal_vdw[i])*rpm2;
1055 dvdl_coul += velec[i]*DLF[i] + LFC[i]*alpha_coul_eff*dlfac_coul[i]*fscal_elec[i]*sigma_pow[i];
1056 dvdl_vdw += vvdw[i]*DLF[i] + LFV[i]*alpha_vdw_eff*dlfac_vdw[i]*fscal_vdw[i]*sigma_pow[i];
1059 dvdl[efptCOUL] += dvdl_coul;
1060 dvdl[efptVDW] += dvdl_vdw;
1062 *velectot = velecsum;