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.
43 #include "gromacs/math/vec.h"
45 #include "nonbonded.h"
46 #include "nb_kernel.h"
49 #include "nb_free_energy.h"
51 #include "gromacs/utility/fatalerror.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;
179 sh_invrc6 = fr->ic->sh_invrc6;
181 if (fr->coulomb_modifier == eintmodPOTSWITCH)
183 d = fr->rcoulomb-fr->rcoulomb_switch;
184 elec_swV3 = -10.0/(d*d*d);
185 elec_swV4 = 15.0/(d*d*d*d);
186 elec_swV5 = -6.0/(d*d*d*d*d);
187 elec_swF2 = -30.0/(d*d*d);
188 elec_swF3 = 60.0/(d*d*d*d);
189 elec_swF4 = -30.0/(d*d*d*d*d);
193 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
194 elec_swV3 = elec_swV4 = elec_swV5 = elec_swF2 = elec_swF3 = elec_swF4 = 0.0;
197 if (fr->vdw_modifier == eintmodPOTSWITCH)
199 d = fr->rvdw-fr->rvdw_switch;
200 vdw_swV3 = -10.0/(d*d*d);
201 vdw_swV4 = 15.0/(d*d*d*d);
202 vdw_swV5 = -6.0/(d*d*d*d*d);
203 vdw_swF2 = -30.0/(d*d*d);
204 vdw_swF3 = 60.0/(d*d*d*d);
205 vdw_swF4 = -30.0/(d*d*d*d*d);
209 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
210 vdw_swV3 = vdw_swV4 = vdw_swV5 = vdw_swF2 = vdw_swF3 = vdw_swF4 = 0.0;
213 if (fr->cutoff_scheme == ecutsVERLET)
215 const interaction_const_t *ic;
218 if (EVDW_PME(ic->vdwtype))
220 ivdw = GMX_NBKERNEL_VDW_LJEWALD;
224 ivdw = GMX_NBKERNEL_VDW_LENNARDJONES;
227 if (ic->eeltype == eelCUT || EEL_RF(ic->eeltype))
229 icoul = GMX_NBKERNEL_ELEC_REACTIONFIELD;
231 else if (EEL_PME_EWALD(ic->eeltype))
233 icoul = GMX_NBKERNEL_ELEC_EWALD;
237 gmx_incons("Unsupported eeltype with Verlet and free-energy");
240 bExactElecCutoff = TRUE;
241 bExactVdwCutoff = TRUE;
245 bExactElecCutoff = (fr->coulomb_modifier != eintmodNONE) || fr->eeltype == eelRF_ZERO;
246 bExactVdwCutoff = (fr->vdw_modifier != eintmodNONE);
249 bExactCutoffAll = (bExactElecCutoff && bExactVdwCutoff);
250 rcutoff_max2 = max(fr->rcoulomb, fr->rvdw);
251 rcutoff_max2 = rcutoff_max2*rcutoff_max2;
253 bEwald = (icoul == GMX_NBKERNEL_ELEC_EWALD);
254 bEwaldLJ = (ivdw == GMX_NBKERNEL_VDW_LJEWALD);
256 /* For Ewald/PME interactions we cannot easily apply the soft-core component to
257 * reciprocal space. When we use vanilla (not switch/shift) Ewald interactions, we
258 * can apply the small trick of subtracting the _reciprocal_ space contribution
259 * in this kernel, and instead apply the free energy interaction to the 1/r
260 * (standard coulomb) interaction.
262 * However, we cannot use this approach for switch-modified since we would then
263 * effectively end up evaluating a significantly different interaction here compared to the
264 * normal (non-free-energy) kernels, either by applying a cutoff at a different
265 * position than what the user requested, or by switching different
266 * things (1/r rather than short-range Ewald). For these settings, we just
267 * use the traditional short-range Ewald interaction in that case.
269 bConvertEwaldToCoulomb = (bEwald && (fr->coulomb_modifier != eintmodPOTSWITCH));
270 /* For now the below will always be true (since LJ-PME only works with Shift in Gromacs-5.0),
271 * but writing it this way means we stay in sync with coulomb, and it avoids future bugs.
273 bConvertLJEwaldToLJ6 = (bEwaldLJ && (fr->vdw_modifier != eintmodPOTSWITCH));
275 /* fix compiler warnings */
284 /* Lambda factor for state A, 1-lambda*/
285 LFC[STATE_A] = 1.0 - lambda_coul;
286 LFV[STATE_A] = 1.0 - lambda_vdw;
288 /* Lambda factor for state B, lambda*/
289 LFC[STATE_B] = lambda_coul;
290 LFV[STATE_B] = lambda_vdw;
292 /*derivative of the lambda factor for state A and B */
296 for (i = 0; i < NSTATES; i++)
298 lfac_coul[i] = (lam_power == 2 ? (1-LFC[i])*(1-LFC[i]) : (1-LFC[i]));
299 dlfac_coul[i] = DLF[i]*lam_power/sc_r_power*(lam_power == 2 ? (1-LFC[i]) : 1);
300 lfac_vdw[i] = (lam_power == 2 ? (1-LFV[i])*(1-LFV[i]) : (1-LFV[i]));
301 dlfac_vdw[i] = DLF[i]*lam_power/sc_r_power*(lam_power == 2 ? (1-LFV[i]) : 1);
304 sigma2_def = pow(sigma6_def, 1.0/3.0);
305 sigma2_min = pow(sigma6_min, 1.0/3.0);
307 /* Ewald (not PME) table is special (icoul==enbcoulFEWALD) */
309 do_tab = (icoul == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE ||
310 ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE);
313 tabscale = kernel_data->table_elec_vdw->scale;
314 VFtab = kernel_data->table_elec_vdw->data;
315 /* we always use the combined table here */
319 for (n = 0; (n < nri); n++)
321 int npair_within_cutoff;
323 npair_within_cutoff = 0;
327 shY = shiftvec[is3+1];
328 shZ = shiftvec[is3+2];
336 iqA = facel*chargeA[ii];
337 iqB = facel*chargeB[ii];
338 ntiA = 2*ntype*typeA[ii];
339 ntiB = 2*ntype*typeB[ii];
346 for (k = nj0; (k < nj1); k++)
353 rsq = dx*dx + dy*dy + dz*dz;
355 if (bExactCutoffAll && rsq >= rcutoff_max2)
357 /* We save significant time by skipping all code below.
358 * Note that with soft-core interactions, the actual cut-off
359 * check might be different. But since the soft-core distance
360 * is always larger than r, checking on r here is safe.
364 npair_within_cutoff++;
368 rinv = gmx_invsqrt(rsq);
373 /* The force at r=0 is zero, because of symmetry.
374 * But note that the potential is in general non-zero,
375 * since the soft-cored r will be non-zero.
381 if (sc_r_power == 6.0)
383 rpm2 = rsq*rsq; /* r4 */
384 rp = rpm2*rsq; /* r6 */
386 else if (sc_r_power == 48.0)
388 rp = rsq*rsq*rsq; /* r6 */
389 rp = rp*rp; /* r12 */
390 rp = rp*rp; /* r24 */
391 rp = rp*rp; /* r48 */
392 rpm2 = rp/rsq; /* r46 */
396 rp = pow(r, sc_r_power); /* not currently supported as input, but can handle it */
402 qq[STATE_A] = iqA*chargeA[jnr];
403 qq[STATE_B] = iqB*chargeB[jnr];
405 tj[STATE_A] = ntiA+2*typeA[jnr];
406 tj[STATE_B] = ntiB+2*typeB[jnr];
408 if (nlist->excl_fep == NULL || nlist->excl_fep[k])
410 c6[STATE_A] = nbfp[tj[STATE_A]];
411 c6[STATE_B] = nbfp[tj[STATE_B]];
413 for (i = 0; i < NSTATES; i++)
415 c12[i] = nbfp[tj[i]+1];
416 if ((c6[i] > 0) && (c12[i] > 0))
418 /* c12 is stored scaled with 12.0 and c6 is scaled with 6.0 - correct for this */
419 sigma6[i] = 0.5*c12[i]/c6[i];
420 sigma2[i] = pow(sigma6[i], 1.0/3.0);
421 /* should be able to get rid of this ^^^ internal pow call eventually. Will require agreement on
422 what data to store externally. Can't be fixed without larger scale changes, so not 4.6 */
423 if (sigma6[i] < sigma6_min) /* for disappearing coul and vdw with soft core at the same time */
425 sigma6[i] = sigma6_min;
426 sigma2[i] = sigma2_min;
431 sigma6[i] = sigma6_def;
432 sigma2[i] = sigma2_def;
434 if (sc_r_power == 6.0)
436 sigma_pow[i] = sigma6[i];
437 sigma_powm2[i] = sigma6[i]/sigma2[i];
439 else if (sc_r_power == 48.0)
441 sigma_pow[i] = sigma6[i]*sigma6[i]; /* sigma^12 */
442 sigma_pow[i] = sigma_pow[i]*sigma_pow[i]; /* sigma^24 */
443 sigma_pow[i] = sigma_pow[i]*sigma_pow[i]; /* sigma^48 */
444 sigma_powm2[i] = sigma_pow[i]/sigma2[i];
447 { /* not really supported as input, but in here for testing the general case*/
448 sigma_pow[i] = pow(sigma2[i], sc_r_power/2);
449 sigma_powm2[i] = sigma_pow[i]/(sigma2[i]);
453 /* only use softcore if one of the states has a zero endstate - softcore is for avoiding infinities!*/
454 if ((c12[STATE_A] > 0) && (c12[STATE_B] > 0))
461 alpha_vdw_eff = alpha_vdw;
462 alpha_coul_eff = alpha_coul;
465 for (i = 0; i < NSTATES; i++)
472 /* Only spend time on A or B state if it is non-zero */
473 if ( (qq[i] != 0) || (c6[i] != 0) || (c12[i] != 0) )
475 /* this section has to be inside the loop because of the dependence on sigma_pow */
476 rpinvC = 1.0/(alpha_coul_eff*lfac_coul[i]*sigma_pow[i]+rp);
477 rinvC = pow(rpinvC, 1.0/sc_r_power);
480 rpinvV = 1.0/(alpha_vdw_eff*lfac_vdw[i]*sigma_pow[i]+rp);
481 rinvV = pow(rpinvV, 1.0/sc_r_power);
490 n1C = tab_elemsize*n0;
496 n1V = tab_elemsize*n0;
499 /* Only process the coulomb interactions if we have charges,
500 * and if we either include all entries in the list (no cutoff
501 * used in the kernel), or if we are within the cutoff.
503 bComputeElecInteraction = !bExactElecCutoff ||
504 ( bConvertEwaldToCoulomb && r < rcoulomb) ||
505 (!bConvertEwaldToCoulomb && rC < rcoulomb);
507 if ( (qq[i] != 0) && bComputeElecInteraction)
511 case GMX_NBKERNEL_ELEC_COULOMB:
513 Vcoul[i] = qq[i]*rinvC;
514 FscalC[i] = Vcoul[i];
515 /* The shift for the Coulomb potential is stored in
516 * the RF parameter c_rf, which is 0 without shift
518 Vcoul[i] -= qq[i]*fr->ic->c_rf;
521 case GMX_NBKERNEL_ELEC_REACTIONFIELD:
523 Vcoul[i] = qq[i]*(rinvC + krf*rC*rC-crf);
524 FscalC[i] = qq[i]*(rinvC - 2.0*krf*rC*rC);
527 case GMX_NBKERNEL_ELEC_CUBICSPLINETABLE:
528 /* non-Ewald tabulated coulomb */
532 Geps = epsC*VFtab[nnn+2];
533 Heps2 = eps2C*VFtab[nnn+3];
536 FF = Fp+Geps+2.0*Heps2;
538 FscalC[i] = -qq[i]*tabscale*FF*rC;
541 case GMX_NBKERNEL_ELEC_GENERALIZEDBORN:
542 gmx_fatal(FARGS, "Free energy and GB not implemented.\n");
545 case GMX_NBKERNEL_ELEC_EWALD:
546 if (bConvertEwaldToCoulomb)
548 /* Ewald FEP is done only on the 1/r part */
549 Vcoul[i] = qq[i]*(rinvC-sh_ewald);
550 FscalC[i] = qq[i]*rinvC;
554 ewrt = rC*ewtabscale;
558 FscalC[i] = ewtab[ewitab]+eweps*ewtab[ewitab+1];
559 rinvcorr = rinvC-sh_ewald;
560 Vcoul[i] = qq[i]*(rinvcorr-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+FscalC[i])));
561 FscalC[i] = qq[i]*(rinvC-rC*FscalC[i]);
565 case GMX_NBKERNEL_ELEC_NONE:
571 gmx_incons("Invalid icoul in free energy kernel");
575 if (fr->coulomb_modifier == eintmodPOTSWITCH)
577 d = rC-fr->rcoulomb_switch;
578 d = (d > 0.0) ? d : 0.0;
580 sw = 1.0+d2*d*(elec_swV3+d*(elec_swV4+d*elec_swV5));
581 dsw = d2*(elec_swF2+d*(elec_swF3+d*elec_swF4));
583 FscalC[i] = FscalC[i]*sw - rC*Vcoul[i]*dsw;
586 FscalC[i] = (rC < rcoulomb) ? FscalC[i] : 0.0;
587 Vcoul[i] = (rC < rcoulomb) ? Vcoul[i] : 0.0;
591 /* Only process the VDW interactions if we have
592 * some non-zero parameters, and if we either
593 * include all entries in the list (no cutoff used
594 * in the kernel), or if we are within the cutoff.
596 bComputeVdwInteraction = !bExactVdwCutoff ||
597 ( bConvertLJEwaldToLJ6 && r < rvdw) ||
598 (!bConvertLJEwaldToLJ6 && rV < rvdw);
599 if ((c6[i] != 0 || c12[i] != 0) && bComputeVdwInteraction)
603 case GMX_NBKERNEL_VDW_LENNARDJONES:
604 case GMX_NBKERNEL_VDW_LJEWALD:
606 if (sc_r_power == 6.0)
612 rinv6 = pow(rinvV, 6.0);
615 Vvdw12 = c12[i]*rinv6*rinv6;
616 if (fr->vdw_modifier == eintmodPOTSHIFT)
618 Vvdw[i] = ( (Vvdw12-c12[i]*sh_invrc6*sh_invrc6)*(1.0/12.0)
619 -(Vvdw6-c6[i]*sh_invrc6)*(1.0/6.0));
623 Vvdw[i] = Vvdw12*(1.0/12.0) - Vvdw6*(1.0/6.0);
625 FscalV[i] = Vvdw12 - Vvdw6;
628 case GMX_NBKERNEL_VDW_BUCKINGHAM:
629 gmx_fatal(FARGS, "Buckingham free energy not supported.");
632 case GMX_NBKERNEL_VDW_CUBICSPLINETABLE:
638 Geps = epsV*VFtab[nnn+2];
639 Heps2 = eps2V*VFtab[nnn+3];
642 FF = Fp+Geps+2.0*Heps2;
644 FscalV[i] -= c6[i]*tabscale*FF*rV;
649 Geps = epsV*VFtab[nnn+6];
650 Heps2 = eps2V*VFtab[nnn+7];
653 FF = Fp+Geps+2.0*Heps2;
654 Vvdw[i] += c12[i]*VV;
655 FscalV[i] -= c12[i]*tabscale*FF*rV;
658 case GMX_NBKERNEL_VDW_NONE:
664 gmx_incons("Invalid ivdw in free energy kernel");
668 if (fr->vdw_modifier == eintmodPOTSWITCH)
670 d = rV-fr->rvdw_switch;
671 d = (d > 0.0) ? d : 0.0;
673 sw = 1.0+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
674 dsw = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
676 FscalV[i] = FscalV[i]*sw - rV*Vvdw[i]*dsw;
679 FscalV[i] = (rV < rvdw) ? FscalV[i] : 0.0;
680 Vvdw[i] = (rV < rvdw) ? Vvdw[i] : 0.0;
684 /* FscalC (and FscalV) now contain: dV/drC * rC
685 * Now we multiply by rC^-p, so it will be: dV/drC * rC^1-p
686 * Further down we first multiply by r^p-2 and then by
687 * the vector r, which in total gives: dV/drC * (r/rC)^1-p
694 /* Assemble A and B states */
695 for (i = 0; i < NSTATES; i++)
697 vctot += LFC[i]*Vcoul[i];
698 vvtot += LFV[i]*Vvdw[i];
700 Fscal += LFC[i]*FscalC[i]*rpm2;
701 Fscal += LFV[i]*FscalV[i]*rpm2;
703 dvdl_coul += Vcoul[i]*DLF[i] + LFC[i]*alpha_coul_eff*dlfac_coul[i]*FscalC[i]*sigma_pow[i];
704 dvdl_vdw += Vvdw[i]*DLF[i] + LFV[i]*alpha_vdw_eff*dlfac_vdw[i]*FscalV[i]*sigma_pow[i];
707 else if (icoul == GMX_NBKERNEL_ELEC_REACTIONFIELD)
709 /* For excluded pairs, which are only in this pair list when
710 * using the Verlet scheme, we don't use soft-core.
711 * The group scheme also doesn't soft-core for these.
712 * As there is no singularity, there is no need for soft-core.
722 for (i = 0; i < NSTATES; i++)
724 vctot += LFC[i]*qq[i]*VV;
725 Fscal += LFC[i]*qq[i]*FF;
726 dvdl_coul += DLF[i]*qq[i]*VV;
730 if (bConvertEwaldToCoulomb && ( !bExactElecCutoff || r < rcoulomb ) )
732 /* See comment in the preamble. When using Ewald interactions
733 * (unless we use a switch modifier) we subtract the reciprocal-space
734 * Ewald component here which made it possible to apply the free
735 * energy interaction to 1/r (vanilla coulomb short-range part)
736 * above. This gets us closer to the ideal case of applying
737 * the softcore to the entire electrostatic interaction,
738 * including the reciprocal-space component.
746 f_lr = ewtab[ewitab]+eweps*ewtab[ewitab+1];
747 v_lr = (ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+f_lr));
752 /* If we get here, the i particle (ii) has itself (jnr)
753 * in its neighborlist. This can only happen with the Verlet
754 * scheme, and corresponds to a self-interaction that will
755 * occur twice. Scale it down by 50% to only include it once.
760 for (i = 0; i < NSTATES; i++)
762 vctot -= LFC[i]*qq[i]*v_lr;
763 Fscal -= LFC[i]*qq[i]*f_lr;
764 dvdl_coul -= (DLF[i]*qq[i])*v_lr;
768 if (bConvertLJEwaldToLJ6 && (!bExactVdwCutoff || r < rvdw))
770 /* See comment in the preamble. When using LJ-Ewald interactions
771 * (unless we use a switch modifier) we subtract the reciprocal-space
772 * Ewald component here which made it possible to apply the free
773 * energy interaction to r^-6 (vanilla LJ6 short-range part)
774 * above. This gets us closer to the ideal case of applying
775 * the softcore to the entire VdW interaction,
776 * including the reciprocal-space component.
781 rs = rsq*rinv*ewtabscale;
784 f_lr = (1 - frac)*tab_ewald_F_lj[ri] + frac*tab_ewald_F_lj[ri+1];
786 VV = tab_ewald_V_lj[ri] - ewtabhalfspace*frac*(tab_ewald_F_lj[ri] + f_lr);
790 /* If we get here, the i particle (ii) has itself (jnr)
791 * in its neighborlist. This can only happen with the Verlet
792 * scheme, and corresponds to a self-interaction that will
793 * occur twice. Scale it down by 50% to only include it once.
798 for (i = 0; i < NSTATES; i++)
800 c6grid = nbfp_grid[tj[i]];
801 vvtot += LFV[i]*c6grid*VV*(1.0/6.0);
802 Fscal += LFV[i]*c6grid*FF*(1.0/6.0);
803 dvdl_vdw += (DLF[i]*c6grid)*VV*(1.0/6.0);
816 /* OpenMP atomics are expensive, but this kernels is also
817 * expensive, so we can take this hit, instead of using
818 * thread-local output buffers and extra reduction.
829 /* The atomics below are expensive with many OpenMP threads.
830 * Here unperturbed i-particles will usually only have a few
831 * (perturbed) j-particles in the list. Thus with a buffered list
832 * we can skip a significant number of i-reductions with a check.
834 if (npair_within_cutoff > 0)
850 fshift[is3+1] += fiy;
852 fshift[is3+2] += fiz;
866 dvdl[efptCOUL] += dvdl_coul;
868 dvdl[efptVDW] += dvdl_vdw;
870 /* Estimate flops, average for free energy stuff:
871 * 12 flops per outer iteration
872 * 150 flops per inner iteration
875 inc_nrnb(nrnb, eNR_NBKERNEL_FREE_ENERGY, nlist->nri*12 + nlist->jindex[n]*150);
879 nb_free_energy_evaluate_single(real r2, real sc_r_power, real alpha_coul, real alpha_vdw,
880 real tabscale, real *vftab,
881 real qqA, real c6A, real c12A, real qqB, real c6B, real c12B,
882 real LFC[2], real LFV[2], real DLF[2],
883 real lfac_coul[2], real lfac_vdw[2], real dlfac_coul[2], real dlfac_vdw[2],
884 real sigma6_def, real sigma6_min, real sigma2_def, real sigma2_min,
885 real *velectot, real *vvdwtot, real *dvdl)
887 real r, rp, rpm2, rtab, eps, eps2, Y, F, Geps, Heps2, Fp, VV, FF, fscal;
888 real qq[2], c6[2], c12[2], sigma6[2], sigma2[2], sigma_pow[2], sigma_powm2[2];
889 real alpha_coul_eff, alpha_vdw_eff, dvdl_coul, dvdl_vdw;
890 real rpinv, r_coul, r_vdw, velecsum, vvdwsum;
891 real fscal_vdw[2], fscal_elec[2];
892 real velec[2], vvdw[2];
902 if (sc_r_power == 6.0)
904 rpm2 = r2*r2; /* r4 */
905 rp = rpm2*r2; /* r6 */
907 else if (sc_r_power == 48.0)
909 rp = r2*r2*r2; /* r6 */
910 rp = rp*rp; /* r12 */
911 rp = rp*rp; /* r24 */
912 rp = rp*rp; /* r48 */
913 rpm2 = rp/r2; /* r46 */
917 rp = pow(r2, 0.5*sc_r_power); /* not currently supported as input, but can handle it */
921 /* Loop over state A(0) and B(1) */
922 for (i = 0; i < 2; i++)
924 if ((c6[i] > 0) && (c12[i] > 0))
926 /* The c6 & c12 coefficients now contain the constants 6.0 and 12.0, respectively.
927 * Correct for this by multiplying with (1/12.0)/(1/6.0)=6.0/12.0=0.5.
929 sigma6[i] = 0.5*c12[i]/c6[i];
930 sigma2[i] = pow(0.5*c12[i]/c6[i], 1.0/3.0);
931 /* should be able to get rid of this ^^^ internal pow call eventually. Will require agreement on
932 what data to store externally. Can't be fixed without larger scale changes, so not 5.0 */
933 if (sigma6[i] < sigma6_min) /* for disappearing coul and vdw with soft core at the same time */
935 sigma6[i] = sigma6_min;
936 sigma2[i] = sigma2_min;
941 sigma6[i] = sigma6_def;
942 sigma2[i] = sigma2_def;
944 if (sc_r_power == 6.0)
946 sigma_pow[i] = sigma6[i];
947 sigma_powm2[i] = sigma6[i]/sigma2[i];
949 else if (sc_r_power == 48.0)
951 sigma_pow[i] = sigma6[i]*sigma6[i]; /* sigma^12 */
952 sigma_pow[i] = sigma_pow[i]*sigma_pow[i]; /* sigma^24 */
953 sigma_pow[i] = sigma_pow[i]*sigma_pow[i]; /* sigma^48 */
954 sigma_powm2[i] = sigma_pow[i]/sigma2[i];
957 { /* not really supported as input, but in here for testing the general case*/
958 sigma_pow[i] = pow(sigma2[i], sc_r_power/2);
959 sigma_powm2[i] = sigma_pow[i]/(sigma2[i]);
963 /* only use softcore if one of the states has a zero endstate - softcore is for avoiding infinities!*/
964 if ((c12[0] > 0) && (c12[1] > 0))
971 alpha_vdw_eff = alpha_vdw;
972 alpha_coul_eff = alpha_coul;
975 /* Loop over A and B states again */
976 for (i = 0; i < 2; i++)
983 /* Only spend time on A or B state if it is non-zero */
984 if ( (qq[i] != 0) || (c6[i] != 0) || (c12[i] != 0) )
987 rpinv = 1.0/(alpha_coul_eff*lfac_coul[i]*sigma_pow[i]+rp);
988 r_coul = pow(rpinv, -1.0/sc_r_power);
990 /* Electrostatics table lookup data */
991 rtab = r_coul*tabscale;
999 Geps = eps*vftab[ntab+2];
1000 Heps2 = eps2*vftab[ntab+3];
1003 FF = Fp+Geps+2.0*Heps2;
1004 velec[i] = qq[i]*VV;
1005 fscal_elec[i] = -qq[i]*FF*r_coul*rpinv*tabscale;
1008 rpinv = 1.0/(alpha_vdw_eff*lfac_vdw[i]*sigma_pow[i]+rp);
1009 r_vdw = pow(rpinv, -1.0/sc_r_power);
1010 /* Vdw table lookup data */
1011 rtab = r_vdw*tabscale;
1019 Geps = eps*vftab[ntab+6];
1020 Heps2 = eps2*vftab[ntab+7];
1023 FF = Fp+Geps+2.0*Heps2;
1025 fscal_vdw[i] = -c6[i]*FF;
1030 Geps = eps*vftab[ntab+10];
1031 Heps2 = eps2*vftab[ntab+11];
1034 FF = Fp+Geps+2.0*Heps2;
1035 vvdw[i] += c12[i]*VV;
1036 fscal_vdw[i] -= c12[i]*FF;
1037 fscal_vdw[i] *= r_vdw*rpinv*tabscale;
1040 /* Now we have velec[i], vvdw[i], and fscal[i] for both states */
1041 /* Assemble A and B states */
1047 for (i = 0; i < 2; i++)
1049 velecsum += LFC[i]*velec[i];
1050 vvdwsum += LFV[i]*vvdw[i];
1052 fscal += (LFC[i]*fscal_elec[i]+LFV[i]*fscal_vdw[i])*rpm2;
1054 dvdl_coul += velec[i]*DLF[i] + LFC[i]*alpha_coul_eff*dlfac_coul[i]*fscal_elec[i]*sigma_pow[i];
1055 dvdl_vdw += vvdw[i]*DLF[i] + LFV[i]*alpha_vdw_eff*dlfac_vdw[i]*fscal_vdw[i]*sigma_pow[i];
1058 dvdl[efptCOUL] += dvdl_coul;
1059 dvdl[efptVDW] += dvdl_vdw;
1061 *velectot = velecsum;