2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5 * Copyright (c) 2017,2018,2019,2020,2021, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
39 #include "gromacs/gmxlib/nrnb.h"
40 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
41 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/mdlib/enerdata_utils.h"
44 #include "gromacs/mdlib/force.h"
45 #include "gromacs/mdlib/gmx_omp_nthreads.h"
46 #include "gromacs/mdtypes/enerdata.h"
47 #include "gromacs/mdtypes/forceoutput.h"
48 #include "gromacs/mdtypes/forcerec.h"
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/mdtypes/interaction_const.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/mdtypes/mdatom.h"
53 #include "gromacs/mdtypes/nblist.h"
54 #include "gromacs/mdtypes/simulation_workload.h"
55 #include "gromacs/nbnxm/gpu_data_mgmt.h"
56 #include "gromacs/nbnxm/nbnxm.h"
57 #include "gromacs/simd/simd.h"
58 #include "gromacs/timing/wallcycle.h"
59 #include "gromacs/utility/enumerationhelpers.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/gmxassert.h"
62 #include "gromacs/utility/real.h"
64 #include "kernel_common.h"
65 #include "nbnxm_gpu.h"
66 #include "nbnxm_simd.h"
67 #include "pairlistset.h"
68 #include "pairlistsets.h"
69 #include "kernels_reference/kernel_gpu_ref.h"
70 #define INCLUDE_KERNELFUNCTION_TABLES
71 #include "kernels_reference/kernel_ref.h"
72 #ifdef GMX_NBNXN_SIMD_2XNN
73 # include "kernels_simd_2xmm/kernels.h"
75 #ifdef GMX_NBNXN_SIMD_4XN
76 # include "kernels_simd_4xm/kernels.h"
78 #undef INCLUDE_FUNCTION_TABLES
80 /*! \brief Clears the energy group output buffers
82 * \param[in,out] out nbnxn kernel output struct
84 static void clearGroupEnergies(nbnxn_atomdata_output_t* out)
86 std::fill(out->Vvdw.begin(), out->Vvdw.end(), 0.0_real);
87 std::fill(out->Vc.begin(), out->Vc.end(), 0.0_real);
88 std::fill(out->VSvdw.begin(), out->VSvdw.end(), 0.0_real);
89 std::fill(out->VSc.begin(), out->VSc.end(), 0.0_real);
92 /*! \brief Reduce the group-pair energy buffers produced by a SIMD kernel
93 * to single terms in the output buffers.
95 * The SIMD kernels produce a large number of energy buffer in SIMD registers
96 * to avoid scattered reads and writes.
98 * \tparam unrollj The unroll size for j-particles in the SIMD kernel
99 * \param[in] numGroups The number of energy groups
100 * \param[in] numGroups_2log Log2 of numGroups, rounded up
101 * \param[in,out] out Struct with energy buffers
103 template<int unrollj>
104 static void reduceGroupEnergySimdBuffers(int numGroups, int numGroups_2log, nbnxn_atomdata_output_t* out)
106 const int unrollj_half = unrollj / 2;
107 /* Energies are stored in SIMD registers with size 2^numGroups_2log */
108 const int numGroupsStorage = (1 << numGroups_2log);
110 const real* gmx_restrict vVdwSimd = out->VSvdw.data();
111 const real* gmx_restrict vCoulombSimd = out->VSc.data();
112 real* gmx_restrict vVdw = out->Vvdw.data();
113 real* gmx_restrict vCoulomb = out->Vc.data();
115 /* The size of the SIMD energy group buffer array is:
116 * numGroups*numGroups*numGroupsStorage*unrollj_half*simd_width
118 for (int i = 0; i < numGroups; i++)
120 for (int j1 = 0; j1 < numGroups; j1++)
122 for (int j0 = 0; j0 < numGroups; j0++)
124 int c = ((i * numGroups + j1) * numGroupsStorage + j0) * unrollj_half * unrollj;
125 for (int s = 0; s < unrollj_half; s++)
127 vVdw[i * numGroups + j0] += vVdwSimd[c + 0];
128 vVdw[i * numGroups + j1] += vVdwSimd[c + 1];
129 vCoulomb[i * numGroups + j0] += vCoulombSimd[c + 0];
130 vCoulomb[i * numGroups + j1] += vCoulombSimd[c + 1];
138 static int getCoulombKernelType(const Nbnxm::KernelSetup& kernelSetup, const interaction_const_t& ic)
141 if (EEL_RF(ic.eeltype) || ic.eeltype == CoulombInteractionType::Cut)
147 if (kernelSetup.ewaldExclusionType == Nbnxm::EwaldExclusionType::Table)
149 if (ic.rcoulomb == ic.rvdw)
155 return coulktTAB_TWIN;
160 if (ic.rcoulomb == ic.rvdw)
166 return coulktEWALD_TWIN;
172 static int getVdwKernelType(const Nbnxm::KernelSetup& kernelSetup,
173 const nbnxn_atomdata_t::Params& nbatParams,
174 const interaction_const_t& ic)
176 if (ic.vdwtype == VanDerWaalsType::Cut)
178 switch (ic.vdw_modifier)
180 case InteractionModifiers::None:
181 case InteractionModifiers::PotShift:
182 switch (nbatParams.ljCombinationRule)
184 case LJCombinationRule::Geometric: return vdwktLJCUT_COMBGEOM;
185 case LJCombinationRule::LorentzBerthelot: return vdwktLJCUT_COMBLB;
186 case LJCombinationRule::None: return vdwktLJCUT_COMBNONE;
187 default: gmx_incons("Unknown combination rule");
189 case InteractionModifiers::ForceSwitch: return vdwktLJFORCESWITCH;
190 case InteractionModifiers::PotSwitch: return vdwktLJPOTSWITCH;
192 std::string errorMsg =
193 gmx::formatString("Unsupported VdW interaction modifier %s (%d)",
194 enumValueToString(ic.vdw_modifier),
195 static_cast<int>(ic.vdw_modifier));
196 gmx_incons(errorMsg);
199 else if (ic.vdwtype == VanDerWaalsType::Pme)
201 if (ic.ljpme_comb_rule == LongRangeVdW::Geom)
203 return vdwktLJEWALDCOMBGEOM;
207 /* At setup we (should have) selected the C reference kernel */
208 GMX_RELEASE_ASSERT(kernelSetup.kernelType == Nbnxm::KernelType::Cpu4x4_PlainC,
209 "Only the C reference nbnxn SIMD kernel supports LJ-PME with LB "
210 "combination rules");
211 return vdwktLJEWALDCOMBLB;
216 std::string errorMsg = gmx::formatString("Unsupported VdW interaction type %s (%d)",
217 enumValueToString(ic.vdwtype),
218 static_cast<int>(ic.vdwtype));
219 gmx_incons(errorMsg);
223 /*! \brief Dispatches the non-bonded N versus M atom cluster CPU kernels.
225 * OpenMP parallelization is performed within this function.
226 * Energy reduction, but not force and shift force reduction, is performed
227 * within this function.
229 * \param[in] pairlistSet Pairlists with local or non-local interactions to compute
230 * \param[in] kernelSetup The non-bonded kernel setup
231 * \param[in,out] nbat The atomdata for the interactions
232 * \param[in] ic Non-bonded interaction constants
233 * \param[in] shiftVectors The PBC shift vectors
234 * \param[in] stepWork Flags that tell what to compute
235 * \param[in] clearF Enum that tells if to clear the force output buffer
236 * \param[out] vCoulomb Output buffer for Coulomb energies
237 * \param[out] vVdw Output buffer for Van der Waals energies
238 * \param[in] wcycle Pointer to cycle counting data structure.
240 static void nbnxn_kernel_cpu(const PairlistSet& pairlistSet,
241 const Nbnxm::KernelSetup& kernelSetup,
242 nbnxn_atomdata_t* nbat,
243 const interaction_const_t& ic,
244 gmx::ArrayRef<const gmx::RVec> shiftVectors,
245 const gmx::StepWorkload& stepWork,
249 gmx_wallcycle* wcycle)
252 const nbnxn_atomdata_t::Params& nbatParams = nbat->params();
254 const int coulkt = getCoulombKernelType(kernelSetup, ic);
255 const int vdwkt = getVdwKernelType(kernelSetup, nbatParams, ic);
257 gmx::ArrayRef<const NbnxnPairlistCpu> pairlists = pairlistSet.cpuLists();
259 auto shiftVecPointer = as_rvec_array(shiftVectors.data());
261 int gmx_unused nthreads = gmx_omp_nthreads_get(emntNonbonded);
262 wallcycle_sub_start(wcycle, ewcsNONBONDED_CLEAR);
263 #pragma omp parallel for schedule(static) num_threads(nthreads)
264 for (gmx::index nb = 0; nb < pairlists.ssize(); nb++)
266 // Presently, the kernels do not call C++ code that can throw,
267 // so no need for a try/catch pair in this OpenMP region.
268 nbnxn_atomdata_output_t* out = &nbat->out[nb];
270 if (clearF == enbvClearFYes)
272 clearForceBuffer(nbat, nb);
274 clear_fshift(out->fshift.data());
279 wallcycle_sub_stop(wcycle, ewcsNONBONDED_CLEAR);
280 wallcycle_sub_start(wcycle, ewcsNONBONDED_KERNEL);
283 // TODO: Change to reference
284 const NbnxnPairlistCpu* pairlist = &pairlists[nb];
286 if (!stepWork.computeEnergy)
288 /* Don't calculate energies */
289 switch (kernelSetup.kernelType)
291 case Nbnxm::KernelType::Cpu4x4_PlainC:
292 nbnxn_kernel_noener_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
294 #ifdef GMX_NBNXN_SIMD_2XNN
295 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
296 nbnxm_kernel_noener_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
299 #ifdef GMX_NBNXN_SIMD_4XN
300 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
301 nbnxm_kernel_noener_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
304 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
307 else if (out->Vvdw.size() == 1)
309 /* A single energy group (pair) */
313 switch (kernelSetup.kernelType)
315 case Nbnxm::KernelType::Cpu4x4_PlainC:
316 nbnxn_kernel_ener_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
318 #ifdef GMX_NBNXN_SIMD_2XNN
319 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
320 nbnxm_kernel_ener_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
323 #ifdef GMX_NBNXN_SIMD_4XN
324 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
325 nbnxm_kernel_ener_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
328 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
333 /* Calculate energy group contributions */
334 clearGroupEnergies(out);
338 switch (kernelSetup.kernelType)
340 case Nbnxm::KernelType::Cpu4x4_PlainC:
341 unrollj = c_nbnxnCpuIClusterSize;
342 nbnxn_kernel_energrp_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
344 #ifdef GMX_NBNXN_SIMD_2XNN
345 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
346 unrollj = GMX_SIMD_REAL_WIDTH / 2;
347 nbnxm_kernel_energrp_simd_2xmm[coulkt][vdwkt](
348 pairlist, nbat, &ic, shiftVecPointer, out);
351 #ifdef GMX_NBNXN_SIMD_4XN
352 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
353 unrollj = GMX_SIMD_REAL_WIDTH;
354 nbnxm_kernel_energrp_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
357 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
360 if (kernelSetup.kernelType != Nbnxm::KernelType::Cpu4x4_PlainC)
365 reduceGroupEnergySimdBuffers<2>(nbatParams.nenergrp, nbatParams.neg_2log, out);
368 reduceGroupEnergySimdBuffers<4>(nbatParams.nenergrp, nbatParams.neg_2log, out);
371 reduceGroupEnergySimdBuffers<8>(nbatParams.nenergrp, nbatParams.neg_2log, out);
373 default: GMX_RELEASE_ASSERT(false, "Unsupported j-unroll size");
378 wallcycle_sub_stop(wcycle, ewcsNONBONDED_KERNEL);
380 if (stepWork.computeEnergy)
382 reduce_energies_over_lists(nbat, pairlists.ssize(), vVdw, vCoulomb);
386 static void accountFlops(t_nrnb* nrnb,
387 const PairlistSet& pairlistSet,
388 const nonbonded_verlet_t& nbv,
389 const interaction_const_t& ic,
390 const gmx::StepWorkload& stepWork)
392 const bool usingGpuKernels = nbv.useGpu();
394 int enr_nbnxn_kernel_ljc = eNRNB;
395 if (EEL_RF(ic.eeltype) || ic.eeltype == CoulombInteractionType::Cut)
397 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
399 else if ((!usingGpuKernels && nbv.kernelSetup().ewaldExclusionType == Nbnxm::EwaldExclusionType::Analytical)
400 || (usingGpuKernels && Nbnxm::gpu_is_kernel_ewald_analytical(nbv.gpu_nbv)))
402 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
406 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
408 int enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
409 if (stepWork.computeEnergy)
411 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
412 enr_nbnxn_kernel_ljc += 1;
413 enr_nbnxn_kernel_lj += 1;
416 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc, pairlistSet.natpair_ljq_);
417 inc_nrnb(nrnb, enr_nbnxn_kernel_lj, pairlistSet.natpair_lj_);
418 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
419 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc - eNR_NBNXN_LJ_RF + eNR_NBNXN_RF, pairlistSet.natpair_q_);
421 if (ic.vdw_modifier == InteractionModifiers::ForceSwitch)
423 /* We add up the switch cost separately */
425 eNR_NBNXN_ADD_LJ_FSW + (stepWork.computeEnergy ? 1 : 0),
426 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
428 if (ic.vdw_modifier == InteractionModifiers::PotSwitch)
430 /* We add up the switch cost separately */
432 eNR_NBNXN_ADD_LJ_PSW + (stepWork.computeEnergy ? 1 : 0),
433 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
435 if (ic.vdwtype == VanDerWaalsType::Pme)
437 /* We add up the LJ Ewald cost separately */
439 eNR_NBNXN_ADD_LJ_EWALD + (stepWork.computeEnergy ? 1 : 0),
440 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
444 void nonbonded_verlet_t::dispatchNonbondedKernel(gmx::InteractionLocality iLocality,
445 const interaction_const_t& ic,
446 const gmx::StepWorkload& stepWork,
448 const t_forcerec& fr,
449 gmx_enerdata_t* enerd,
452 const PairlistSet& pairlistSet = pairlistSets().pairlistSet(iLocality);
454 switch (kernelSetup().kernelType)
456 case Nbnxm::KernelType::Cpu4x4_PlainC:
457 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
458 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
467 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
469 ? enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::BuckinghamSR].data()
470 : enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
474 case Nbnxm::KernelType::Gpu8x8x8:
475 Nbnxm::gpu_launch_kernel(gpu_nbv, stepWork, iLocality);
478 case Nbnxm::KernelType::Cpu8x8x8_PlainC:
479 nbnxn_kernel_gpu_ref(
480 pairlistSet.gpuList(),
487 nbat->out[0].fshift.data(),
488 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
490 ? enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::BuckinghamSR].data()
491 : enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data());
494 default: GMX_RELEASE_ASSERT(false, "Invalid nonbonded kernel type passed!");
497 accountFlops(nrnb, pairlistSet, *this, ic, stepWork);
500 void nonbonded_verlet_t::dispatchFreeEnergyKernel(gmx::InteractionLocality iLocality,
501 const t_forcerec& fr,
502 gmx::ArrayRef<const gmx::RVec> coords,
503 gmx::ForceWithShiftForces* forceWithShiftForces,
504 gmx::ArrayRef<const real> chargeA,
505 gmx::ArrayRef<const real> chargeB,
506 gmx::ArrayRef<const int> typeA,
507 gmx::ArrayRef<const int> typeB,
509 gmx::ArrayRef<const real> lambda,
510 gmx_enerdata_t* enerd,
511 const gmx::StepWorkload& stepWork,
514 const auto nbl_fep = pairlistSets().pairlistSet(iLocality).fepLists();
516 /* When the first list is empty, all are empty and there is nothing to do */
517 if (!pairlistSets().params().haveFep || nbl_fep[0]->nrj == 0)
523 /* Add short-range interactions */
524 donb_flags |= GMX_NONBONDED_DO_SR;
526 if (stepWork.computeForces)
528 donb_flags |= GMX_NONBONDED_DO_FORCE;
530 if (stepWork.computeVirial)
532 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
534 if (stepWork.computeEnergy)
536 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
539 gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, real> dvdl_nb = { 0 };
540 int kernelFlags = donb_flags;
541 gmx::ArrayRef<const real> kernelLambda = lambda;
542 gmx::ArrayRef<real> kernelDvdl = dvdl_nb;
544 gmx::ArrayRef<real> energygrp_elec = enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR];
545 gmx::ArrayRef<real> energygrp_vdw = enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR];
547 GMX_ASSERT(gmx_omp_nthreads_get(emntNonbonded) == nbl_fep.ssize(),
548 "Number of lists should be same as number of NB threads");
550 wallcycle_sub_start(wcycle_, ewcsNONBONDED_FEP);
551 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
552 for (gmx::index th = 0; th < nbl_fep.ssize(); th++)
556 gmx_nb_free_energy_kernel(*nbl_fep[th],
558 forceWithShiftForces,
571 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
574 if (fepvals->sc_alpha != 0)
576 enerd->dvdl_nonlin[FreeEnergyPerturbationCouplingType::Vdw] +=
577 dvdl_nb[FreeEnergyPerturbationCouplingType::Vdw];
578 enerd->dvdl_nonlin[FreeEnergyPerturbationCouplingType::Coul] +=
579 dvdl_nb[FreeEnergyPerturbationCouplingType::Coul];
583 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] +=
584 dvdl_nb[FreeEnergyPerturbationCouplingType::Vdw];
585 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Coul] +=
586 dvdl_nb[FreeEnergyPerturbationCouplingType::Coul];
589 /* If we do foreign lambda and we have soft-core interactions
590 * we have to recalculate the (non-linear) energies contributions.
592 if (fepvals->n_lambda > 0 && stepWork.computeDhdl && fepvals->sc_alpha != 0)
594 gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, real> lam_i;
595 kernelFlags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE))
596 | GMX_NONBONDED_DO_FOREIGNLAMBDA;
597 kernelLambda = lam_i;
598 kernelDvdl = dvdl_nb;
599 gmx::ArrayRef<real> energygrp_elec =
600 enerd->foreign_grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR];
601 gmx::ArrayRef<real> energygrp_vdw =
602 enerd->foreign_grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR];
604 for (gmx::index i = 0; i < 1 + enerd->foreignLambdaTerms.numLambdas(); i++)
606 std::fill(std::begin(dvdl_nb), std::end(dvdl_nb), 0);
607 for (int j = 0; j < static_cast<int>(FreeEnergyPerturbationCouplingType::Count); j++)
609 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i - 1]);
611 reset_foreign_enerdata(enerd);
612 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
613 for (gmx::index th = 0; th < nbl_fep.ssize(); th++)
617 gmx_nb_free_energy_kernel(*nbl_fep[th],
619 forceWithShiftForces,
632 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
635 sum_epot(enerd->foreign_grpp, enerd->foreign_term.data());
636 enerd->foreignLambdaTerms.accumulate(
638 enerd->foreign_term[F_EPOT],
639 dvdl_nb[FreeEnergyPerturbationCouplingType::Vdw]
640 + dvdl_nb[FreeEnergyPerturbationCouplingType::Coul]);
643 wallcycle_sub_stop(wcycle_, ewcsNONBONDED_FEP);