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,2015,2016,2017,2018,2019, 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_free_energy.h"
45 #include "gromacs/gmxlib/nrnb.h"
46 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
47 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
48 #include "gromacs/math/functions.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/mdtypes/forceoutput.h"
51 #include "gromacs/mdtypes/forcerec.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/utility/fatalerror.h"
56 //! Enum for templating the soft-core treatment in the kernel
57 enum class SoftCoreTreatment
59 None, //!< No soft-core
60 RPower6, //!< Soft-core with r-power = 6
61 RPower48 //!< Soft-core with r-power = 48
64 //! Most treatments are fine with float in mixed-precision mode.
65 template <SoftCoreTreatment softCoreTreatment>
68 //! Real type for soft-core calculations
72 //! This treatment requires double precision for some computations.
74 struct SoftCoreReal<SoftCoreTreatment::RPower48>
76 //! Real type for soft-core calculations
80 //! Computes r^(1/p) and 1/r^(1/p) for the standard p=6
81 template <SoftCoreTreatment softCoreTreatment>
82 static inline void pthRoot(const real r,
86 *invPthRoot = gmx::invsqrt(std::cbrt(r));
87 *pthRoot = 1/(*invPthRoot);
90 // We need a double version to make the specialization below work
92 //! Computes r^(1/p) and 1/r^(1/p) for the standard p=6
93 template <SoftCoreTreatment softCoreTreatment>
94 static inline void pthRoot(const double r,
98 *invPthRoot = gmx::invsqrt(std::cbrt(r));
99 *pthRoot = 1/(*invPthRoot);
103 //! Computes r^(1/p) and 1/r^(1/p) for p=48
105 inline void pthRoot<SoftCoreTreatment::RPower48>(const double r,
109 *pthRoot = std::pow(r, 1.0/48.0);
110 *invPthRoot = 1/(*pthRoot);
113 //! Templated free-energy non-bonded kernel
114 template<SoftCoreTreatment softCoreTreatment, bool scLambdasOrAlphasDiffer>
116 nb_free_energy_kernel(const t_nblist * gmx_restrict nlist,
117 rvec * gmx_restrict xx,
118 gmx::ForceWithShiftForces * forceWithShiftForces,
119 const t_forcerec * gmx_restrict fr,
120 const t_mdatoms * gmx_restrict mdatoms,
121 nb_kernel_data_t * gmx_restrict kernel_data,
122 t_nrnb * gmx_restrict nrnb)
124 using SCReal = typename SoftCoreReal<softCoreTreatment>::Real;
126 constexpr bool useSoftCore = (softCoreTreatment != SoftCoreTreatment::None);
131 int i, n, ii, is3, ii3, k, nj0, nj1, jnr, j3, ggid;
133 real tx, ty, tz, Fscal;
134 SCReal FscalC[NSTATES], FscalV[NSTATES]; /* Needs double for sc_power==48 */
135 real Vcoul[NSTATES], Vvdw[NSTATES];
138 real qq[NSTATES], vctot;
139 int ntiA, ntiB, tj[NSTATES];
140 real Vvdw6, Vvdw12, vvtot;
141 real ix, iy, iz, fix, fiy, fiz;
142 real dx, dy, dz, rsq, rinv;
143 real c6[NSTATES], c12[NSTATES], c6grid;
144 real LFC[NSTATES], LFV[NSTATES], DLF[NSTATES];
145 SCReal dvdl_coul, dvdl_vdw;
146 real lfac_coul[NSTATES], dlfac_coul[NSTATES], lfac_vdw[NSTATES], dlfac_vdw[NSTATES];
147 real sigma6[NSTATES], alpha_vdw_eff, alpha_coul_eff;
149 SCReal rp, rpm2, rC, rV, rpinvC, rpinvV; /* Needs double for sc_power==48 */
150 real sigma_pow[NSTATES];
162 const real * shiftvec;
163 const real * chargeA;
164 const real * chargeB;
165 real sigma6_min, sigma6_def, lam_power;
166 real alpha_coul, alpha_vdw, lambda_coul, lambda_vdw;
167 const real * nbfp, *nbfp_grid;
171 gmx_bool bDoForces, bDoShiftForces, bDoPotential;
172 gmx_bool bEwald, bEwaldLJ;
174 const real * tab_ewald_F_lj = nullptr;
175 const real * tab_ewald_V_lj = nullptr;
177 real vdw_swV3, vdw_swV4, vdw_swV5, vdw_swF2, vdw_swF3, vdw_swF4;
178 gmx_bool bComputeVdwInteraction, bComputeElecInteraction;
179 const real * ewtab = nullptr;
181 real ewrt, eweps, ewtabscale = 0, ewtabhalfspace = 0, sh_ewald = 0;
183 const real onetwelfth = 1.0/12.0;
184 const real onesixth = 1.0/6.0;
185 const real zero = 0.0;
186 const real half = 0.5;
187 const real one = 1.0;
188 const real two = 2.0;
189 const real six = 6.0;
191 /* Extract pointer to non-bonded interaction constants */
192 const interaction_const_t *ic = fr->ic;
194 // TODO: We should get rid of using pointers to real
195 const real *x = xx[0];
196 real * gmx_restrict f = &(forceWithShiftForces->force()[0][0]);
198 real * gmx_restrict fshift = &(forceWithShiftForces->shiftForces()[0][0]);
200 // Extract pair list data
203 jindex = nlist->jindex;
205 shift = nlist->shift;
208 shiftvec = fr->shift_vec[0];
209 chargeA = mdatoms->chargeA;
210 chargeB = mdatoms->chargeB;
211 Vc = kernel_data->energygrp_elec;
212 typeA = mdatoms->typeA;
213 typeB = mdatoms->typeB;
216 nbfp_grid = fr->ljpme_c6grid;
217 Vv = kernel_data->energygrp_vdw;
218 lambda_coul = kernel_data->lambda[efptCOUL];
219 lambda_vdw = kernel_data->lambda[efptVDW];
220 dvdl = kernel_data->dvdl;
221 alpha_coul = fr->sc_alphacoul;
222 alpha_vdw = fr->sc_alphavdw;
223 lam_power = fr->sc_power;
224 sigma6_def = fr->sc_sigma6_def;
225 sigma6_min = fr->sc_sigma6_min;
226 bDoForces = ((kernel_data->flags & GMX_NONBONDED_DO_FORCE) != 0);
227 bDoShiftForces = ((kernel_data->flags & GMX_NONBONDED_DO_SHIFTFORCE) != 0);
228 bDoPotential = ((kernel_data->flags & GMX_NONBONDED_DO_POTENTIAL) != 0);
230 // Extract data from interaction_const_t
231 const real facel = ic->epsfac;
232 const real rcoulomb = ic->rcoulomb;
233 const real krf = ic->k_rf;
234 const real crf = ic->c_rf;
235 const real sh_lj_ewald = ic->sh_lj_ewald;
236 const real rvdw = ic->rvdw;
237 const real dispersionShift = ic->dispersion_shift.cpot;
238 const real repulsionShift = ic->repulsion_shift.cpot;
240 // Note that the nbnxm kernels do not support Coulomb potential switching at all
241 GMX_ASSERT(ic->coulomb_modifier != eintmodPOTSWITCH,
242 "Potential switching is not supported for Coulomb with FEP");
244 if (ic->vdw_modifier == eintmodPOTSWITCH)
246 d = ic->rvdw - ic->rvdw_switch;
247 vdw_swV3 = -10.0/(d*d*d);
248 vdw_swV4 = 15.0/(d*d*d*d);
249 vdw_swV5 = -6.0/(d*d*d*d*d);
250 vdw_swF2 = -30.0/(d*d*d);
251 vdw_swF3 = 60.0/(d*d*d*d);
252 vdw_swF4 = -30.0/(d*d*d*d*d);
256 /* Avoid warnings from stupid compilers (looking at you, Clang!) */
257 vdw_swV3 = vdw_swV4 = vdw_swV5 = vdw_swF2 = vdw_swF3 = vdw_swF4 = 0.0;
260 if (EVDW_PME(ic->vdwtype))
262 ivdw = GMX_NBKERNEL_VDW_LJEWALD;
266 ivdw = GMX_NBKERNEL_VDW_LENNARDJONES;
269 if (ic->eeltype == eelCUT || EEL_RF(ic->eeltype))
271 icoul = GMX_NBKERNEL_ELEC_REACTIONFIELD;
273 else if (EEL_PME_EWALD(ic->eeltype))
275 icoul = GMX_NBKERNEL_ELEC_EWALD;
279 gmx_incons("Unsupported eeltype with Verlet and free-energy");
282 rcutoff_max2 = std::max(ic->rcoulomb, ic->rvdw);
283 rcutoff_max2 = rcutoff_max2*rcutoff_max2;
285 bEwald = (icoul == GMX_NBKERNEL_ELEC_EWALD);
286 bEwaldLJ = (ivdw == GMX_NBKERNEL_VDW_LJEWALD);
288 if (bEwald || bEwaldLJ)
290 const auto &tables = *ic->coulombEwaldTables;
291 sh_ewald = ic->sh_ewald;
292 ewtab = tables.tableFDV0.data();
293 ewtabscale = tables.scale;
294 ewtabhalfspace = half/ewtabscale;
295 tab_ewald_F_lj = tables.tableF.data();
296 tab_ewald_V_lj = tables.tableV.data();
299 /* For Ewald/PME interactions we cannot easily apply the soft-core component to
300 * reciprocal space. When we use non-switched Ewald interactions, we
301 * assume the soft-coring does not significantly affect the grid contribution
302 * and apply the soft-core only to the full 1/r (- shift) pair contribution.
304 * However, we cannot use this approach for switch-modified since we would then
305 * effectively end up evaluating a significantly different interaction here compared to the
306 * normal (non-free-energy) kernels, either by applying a cutoff at a different
307 * position than what the user requested, or by switching different
308 * things (1/r rather than short-range Ewald). For these settings, we just
309 * use the traditional short-range Ewald interaction in that case.
311 GMX_RELEASE_ASSERT(!(bEwald && ic->coulomb_modifier == eintmodPOTSWITCH) &&
312 !(bEwaldLJ && ic->vdw_modifier == eintmodPOTSWITCH),
313 "Can not apply soft-core to switched Ewald potentials");
318 /* Lambda factor for state A, 1-lambda*/
319 LFC[STATE_A] = one - lambda_coul;
320 LFV[STATE_A] = one - lambda_vdw;
322 /* Lambda factor for state B, lambda*/
323 LFC[STATE_B] = lambda_coul;
324 LFV[STATE_B] = lambda_vdw;
326 /*derivative of the lambda factor for state A and B */
330 constexpr real sc_r_power = (softCoreTreatment == SoftCoreTreatment::RPower48 ? 48.0_real : 6.0_real);
331 for (i = 0; i < NSTATES; i++)
333 lfac_coul[i] = (lam_power == 2 ? (1-LFC[i])*(1-LFC[i]) : (1-LFC[i]));
334 dlfac_coul[i] = DLF[i]*lam_power/sc_r_power*(lam_power == 2 ? (1-LFC[i]) : 1);
335 lfac_vdw[i] = (lam_power == 2 ? (1-LFV[i])*(1-LFV[i]) : (1-LFV[i]));
336 dlfac_vdw[i] = DLF[i]*lam_power/sc_r_power*(lam_power == 2 ? (1-LFV[i]) : 1);
339 for (n = 0; (n < nri); n++)
341 int npair_within_cutoff;
343 npair_within_cutoff = 0;
347 shY = shiftvec[is3+1];
348 shZ = shiftvec[is3+2];
356 iqA = facel*chargeA[ii];
357 iqB = facel*chargeB[ii];
358 ntiA = 2*ntype*typeA[ii];
359 ntiB = 2*ntype*typeB[ii];
366 for (k = nj0; (k < nj1); k++)
373 rsq = dx*dx + dy*dy + dz*dz;
375 if (rsq >= rcutoff_max2)
377 /* We save significant time by skipping all code below.
378 * Note that with soft-core interactions, the actual cut-off
379 * check might be different. But since the soft-core distance
380 * is always larger than r, checking on r here is safe.
384 npair_within_cutoff++;
388 /* Note that unlike in the nbnxn kernels, we do not need
389 * to clamp the value of rsq before taking the invsqrt
390 * to avoid NaN in the LJ calculation, since here we do
391 * not calculate LJ interactions when C6 and C12 are zero.
394 rinv = gmx::invsqrt(rsq);
399 /* The force at r=0 is zero, because of symmetry.
400 * But note that the potential is in general non-zero,
401 * since the soft-cored r will be non-zero.
407 if (softCoreTreatment == SoftCoreTreatment::None)
409 /* The soft-core power p will not affect the results
410 * with not using soft-core, so we use power of 0 which gives
411 * the simplest math and cheapest code.
416 if (softCoreTreatment == SoftCoreTreatment::RPower6)
418 rpm2 = rsq*rsq; /* r4 */
419 rp = rpm2*rsq; /* r6 */
421 if (softCoreTreatment == SoftCoreTreatment::RPower48)
423 rp = rsq*rsq*rsq; /* r6 */
424 rp = rp*rp; /* r12 */
425 rp = rp*rp; /* r24 */
426 rp = rp*rp; /* r48 */
427 rpm2 = rp/rsq; /* r46 */
432 qq[STATE_A] = iqA*chargeA[jnr];
433 qq[STATE_B] = iqB*chargeB[jnr];
435 tj[STATE_A] = ntiA+2*typeA[jnr];
436 tj[STATE_B] = ntiB+2*typeB[jnr];
438 if (nlist->excl_fep == nullptr || nlist->excl_fep[k])
440 c6[STATE_A] = nbfp[tj[STATE_A]];
441 c6[STATE_B] = nbfp[tj[STATE_B]];
443 for (i = 0; i < NSTATES; i++)
445 c12[i] = nbfp[tj[i]+1];
448 if ((c6[i] > 0) && (c12[i] > 0))
450 /* c12 is stored scaled with 12.0 and c6 is scaled with 6.0 - correct for this */
451 sigma6[i] = half*c12[i]/c6[i];
452 if (sigma6[i] < sigma6_min) /* for disappearing coul and vdw with soft core at the same time */
454 sigma6[i] = sigma6_min;
459 sigma6[i] = sigma6_def;
461 if (softCoreTreatment == SoftCoreTreatment::RPower6)
463 sigma_pow[i] = sigma6[i];
467 sigma_pow[i] = sigma6[i]*sigma6[i]; /* sigma^12 */
468 sigma_pow[i] = sigma_pow[i]*sigma_pow[i]; /* sigma^24 */
469 sigma_pow[i] = sigma_pow[i]*sigma_pow[i]; /* sigma^48 */
476 /* only use softcore if one of the states has a zero endstate - softcore is for avoiding infinities!*/
477 if ((c12[STATE_A] > 0) && (c12[STATE_B] > 0))
484 alpha_vdw_eff = alpha_vdw;
485 alpha_coul_eff = alpha_coul;
489 for (i = 0; i < NSTATES; i++)
496 /* Only spend time on A or B state if it is non-zero */
497 if ( (qq[i] != 0) || (c6[i] != 0) || (c12[i] != 0) )
499 /* this section has to be inside the loop because of the dependence on sigma_pow */
502 rpinvC = one/(alpha_coul_eff*lfac_coul[i]*sigma_pow[i]+rp);
503 pthRoot<softCoreTreatment>(rpinvC, &rinvC, &rC);
505 if (scLambdasOrAlphasDiffer)
507 rpinvV = one/(alpha_vdw_eff*lfac_vdw[i]*sigma_pow[i]+rp);
508 pthRoot<softCoreTreatment>(rpinvV, &rinvV, &rV);
512 /* We can avoid one expensive pow and one / operation */
529 /* Only process the coulomb interactions if we have charges,
530 * and if we either include all entries in the list (no cutoff
531 * used in the kernel), or if we are within the cutoff.
533 bComputeElecInteraction =
534 ( bEwald && r < rcoulomb) ||
535 (!bEwald && rC < rcoulomb);
537 if ( (qq[i] != 0) && bComputeElecInteraction)
541 /* Ewald FEP is done only on the 1/r part */
542 Vcoul[i] = qq[i]*(rinvC-sh_ewald);
543 FscalC[i] = qq[i]*rinvC;
548 Vcoul[i] = qq[i]*(rinvC + krf*rC*rC-crf);
549 FscalC[i] = qq[i]*(rinvC - two*krf*rC*rC);
553 /* Only process the VDW interactions if we have
554 * some non-zero parameters, and if we either
555 * include all entries in the list (no cutoff used
556 * in the kernel), or if we are within the cutoff.
558 bComputeVdwInteraction =
559 ( bEwaldLJ && r < rvdw) ||
560 (!bEwaldLJ && rV < rvdw);
561 if ((c6[i] != 0 || c12[i] != 0) && bComputeVdwInteraction)
563 /* cutoff LJ, also handles part of Ewald LJ */
564 if (softCoreTreatment == SoftCoreTreatment::RPower6)
571 rinv6 = rinv6*rinv6*rinv6;
574 Vvdw12 = c12[i]*rinv6*rinv6;
576 Vvdw[i] = ( (Vvdw12 + c12[i]*repulsionShift)*onetwelfth
577 - (Vvdw6 + c6[i]*dispersionShift)*onesixth);
578 FscalV[i] = Vvdw12 - Vvdw6;
582 /* Subtract the grid potential at the cut-off */
583 c6grid = nbfp_grid[tj[i]];
584 Vvdw[i] += c6grid*sh_lj_ewald*onesixth;
587 if (ic->vdw_modifier == eintmodPOTSWITCH)
589 d = rV - ic->rvdw_switch;
590 d = (d > zero) ? d : zero;
592 sw = one+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
593 dsw = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
595 FscalV[i] = FscalV[i]*sw - rV*Vvdw[i]*dsw;
598 FscalV[i] = (rV < rvdw) ? FscalV[i] : zero;
599 Vvdw[i] = (rV < rvdw) ? Vvdw[i] : zero;
603 /* FscalC (and FscalV) now contain: dV/drC * rC
604 * Now we multiply by rC^-p, so it will be: dV/drC * rC^1-p
605 * Further down we first multiply by r^p-2 and then by
606 * the vector r, which in total gives: dV/drC * (r/rC)^1-p
613 /* Assemble A and B states */
614 for (i = 0; i < NSTATES; i++)
616 vctot += LFC[i]*Vcoul[i];
617 vvtot += LFV[i]*Vvdw[i];
619 Fscal += LFC[i]*FscalC[i]*rpm2;
620 Fscal += LFV[i]*FscalV[i]*rpm2;
624 dvdl_coul += Vcoul[i]*DLF[i] + LFC[i]*alpha_coul_eff*dlfac_coul[i]*FscalC[i]*sigma_pow[i];
625 dvdl_vdw += Vvdw[i]*DLF[i] + LFV[i]*alpha_vdw_eff*dlfac_vdw[i]*FscalV[i]*sigma_pow[i];
629 dvdl_coul += Vcoul[i]*DLF[i];
630 dvdl_vdw += Vvdw[i]*DLF[i];
634 else if (icoul == GMX_NBKERNEL_ELEC_REACTIONFIELD)
636 /* For excluded pairs, which are only in this pair list when
637 * using the Verlet scheme, we don't use soft-core.
638 * The group scheme also doesn't soft-core for these.
639 * As there is no singularity, there is no need for soft-core.
649 for (i = 0; i < NSTATES; i++)
651 vctot += LFC[i]*qq[i]*VV;
652 Fscal += LFC[i]*qq[i]*FF;
653 dvdl_coul += DLF[i]*qq[i]*VV;
657 if (bEwald && r < rcoulomb)
659 /* See comment in the preamble. When using Ewald interactions
660 * (unless we use a switch modifier) we subtract the reciprocal-space
661 * Ewald component here which made it possible to apply the free
662 * energy interaction to 1/r (vanilla coulomb short-range part)
663 * above. This gets us closer to the ideal case of applying
664 * the softcore to the entire electrostatic interaction,
665 * including the reciprocal-space component.
670 ewitab = static_cast<int>(ewrt);
673 f_lr = ewtab[ewitab]+eweps*ewtab[ewitab+1];
674 v_lr = (ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+f_lr));
677 /* Note that any possible Ewald shift has already been applied in
678 * the normal interaction part above.
683 /* If we get here, the i particle (ii) has itself (jnr)
684 * in its neighborlist. This can only happen with the Verlet
685 * scheme, and corresponds to a self-interaction that will
686 * occur twice. Scale it down by 50% to only include it once.
691 for (i = 0; i < NSTATES; i++)
693 vctot -= LFC[i]*qq[i]*v_lr;
694 Fscal -= LFC[i]*qq[i]*f_lr;
695 dvdl_coul -= (DLF[i]*qq[i])*v_lr;
699 if (bEwaldLJ && r < rvdw)
701 /* See comment in the preamble. When using LJ-Ewald interactions
702 * (unless we use a switch modifier) we subtract the reciprocal-space
703 * Ewald component here which made it possible to apply the free
704 * energy interaction to r^-6 (vanilla LJ6 short-range part)
705 * above. This gets us closer to the ideal case of applying
706 * the softcore to the entire VdW interaction,
707 * including the reciprocal-space component.
709 /* We could also use the analytical form here
710 * iso a table, but that can cause issues for
711 * r close to 0 for non-interacting pairs.
716 rs = rsq*rinv*ewtabscale;
717 ri = static_cast<int>(rs);
719 f_lr = (1 - frac)*tab_ewald_F_lj[ri] + frac*tab_ewald_F_lj[ri+1];
720 /* TODO: Currently the Ewald LJ table does not contain
721 * the factor 1/6, we should add this.
724 VV = (tab_ewald_V_lj[ri] - ewtabhalfspace*frac*(tab_ewald_F_lj[ri] + f_lr))/six;
728 /* If we get here, the i particle (ii) has itself (jnr)
729 * in its neighborlist. This can only happen with the Verlet
730 * scheme, and corresponds to a self-interaction that will
731 * occur twice. Scale it down by 50% to only include it once.
736 for (i = 0; i < NSTATES; i++)
738 c6grid = nbfp_grid[tj[i]];
739 vvtot += LFV[i]*c6grid*VV;
740 Fscal += LFV[i]*c6grid*FF;
741 dvdl_vdw += (DLF[i]*c6grid)*VV;
753 /* OpenMP atomics are expensive, but this kernels is also
754 * expensive, so we can take this hit, instead of using
755 * thread-local output buffers and extra reduction.
757 * All the OpenMP regions in this file are trivial and should
758 * not throw, so no need for try/catch.
769 /* The atomics below are expensive with many OpenMP threads.
770 * Here unperturbed i-particles will usually only have a few
771 * (perturbed) j-particles in the list. Thus with a buffered list
772 * we can skip a significant number of i-reductions with a check.
774 if (npair_within_cutoff > 0)
790 fshift[is3+1] += fiy;
792 fshift[is3+2] += fiz;
806 dvdl[efptCOUL] += dvdl_coul;
808 dvdl[efptVDW] += dvdl_vdw;
810 /* Estimate flops, average for free energy stuff:
811 * 12 flops per outer iteration
812 * 150 flops per inner iteration
815 inc_nrnb(nrnb, eNR_NBKERNEL_FREE_ENERGY, nlist->nri*12 + nlist->jindex[n]*150);
818 void gmx_nb_free_energy_kernel(const t_nblist *nlist,
820 gmx::ForceWithShiftForces *ff,
821 const t_forcerec *fr,
822 const t_mdatoms *mdatoms,
823 nb_kernel_data_t *kernel_data,
826 if (fr->sc_alphacoul == 0 && fr->sc_alphavdw == 0)
828 nb_free_energy_kernel<SoftCoreTreatment::None, false>(nlist, xx, ff, fr, mdatoms, kernel_data, nrnb);
830 else if (fr->sc_r_power == 6.0_real)
832 if (kernel_data->lambda[efptCOUL] == kernel_data->lambda[efptVDW] &&
833 fr->sc_alphacoul == fr->sc_alphavdw)
835 nb_free_energy_kernel<SoftCoreTreatment::RPower6, false>(nlist, xx, ff, fr, mdatoms, kernel_data, nrnb);
839 nb_free_energy_kernel<SoftCoreTreatment::RPower6, true>(nlist, xx, ff, fr, mdatoms, kernel_data, nrnb);
842 else if (fr->sc_r_power == 48.0_real)
844 nb_free_energy_kernel<SoftCoreTreatment::RPower48, true>(nlist, xx, ff, fr, mdatoms, kernel_data, nrnb);
848 GMX_RELEASE_ASSERT(false, "Unsupported soft-core r-power");