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_gpu_data_mgmt.h"
67 #include "nbnxm_simd.h"
68 #include "pairlistset.h"
69 #include "pairlistsets.h"
70 #include "kernels_reference/kernel_gpu_ref.h"
71 #define INCLUDE_KERNELFUNCTION_TABLES
72 #include "kernels_reference/kernel_ref.h"
73 #ifdef GMX_NBNXN_SIMD_2XNN
74 # include "kernels_simd_2xmm/kernels.h"
76 #ifdef GMX_NBNXN_SIMD_4XN
77 # include "kernels_simd_4xm/kernels.h"
79 #undef INCLUDE_FUNCTION_TABLES
81 /*! \brief Clears the energy group output buffers
83 * \param[in,out] out nbnxn kernel output struct
85 static void clearGroupEnergies(nbnxn_atomdata_output_t* out)
87 std::fill(out->Vvdw.begin(), out->Vvdw.end(), 0.0_real);
88 std::fill(out->Vc.begin(), out->Vc.end(), 0.0_real);
89 std::fill(out->VSvdw.begin(), out->VSvdw.end(), 0.0_real);
90 std::fill(out->VSc.begin(), out->VSc.end(), 0.0_real);
93 /*! \brief Reduce the group-pair energy buffers produced by a SIMD kernel
94 * to single terms in the output buffers.
96 * The SIMD kernels produce a large number of energy buffer in SIMD registers
97 * to avoid scattered reads and writes.
99 * \tparam unrollj The unroll size for j-particles in the SIMD kernel
100 * \param[in] numGroups The number of energy groups
101 * \param[in] numGroups_2log Log2 of numGroups, rounded up
102 * \param[in,out] out Struct with energy buffers
104 template<int unrollj>
105 static void reduceGroupEnergySimdBuffers(int numGroups, int numGroups_2log, nbnxn_atomdata_output_t* out)
107 const int unrollj_half = unrollj / 2;
108 /* Energies are stored in SIMD registers with size 2^numGroups_2log */
109 const int numGroupsStorage = (1 << numGroups_2log);
111 const real* gmx_restrict vVdwSimd = out->VSvdw.data();
112 const real* gmx_restrict vCoulombSimd = out->VSc.data();
113 real* gmx_restrict vVdw = out->Vvdw.data();
114 real* gmx_restrict vCoulomb = out->Vc.data();
116 /* The size of the SIMD energy group buffer array is:
117 * numGroups*numGroups*numGroupsStorage*unrollj_half*simd_width
119 for (int i = 0; i < numGroups; i++)
121 for (int j1 = 0; j1 < numGroups; j1++)
123 for (int j0 = 0; j0 < numGroups; j0++)
125 int c = ((i * numGroups + j1) * numGroupsStorage + j0) * unrollj_half * unrollj;
126 for (int s = 0; s < unrollj_half; s++)
128 vVdw[i * numGroups + j0] += vVdwSimd[c + 0];
129 vVdw[i * numGroups + j1] += vVdwSimd[c + 1];
130 vCoulomb[i * numGroups + j0] += vCoulombSimd[c + 0];
131 vCoulomb[i * numGroups + j1] += vCoulombSimd[c + 1];
139 static int getCoulombKernelType(const Nbnxm::KernelSetup& kernelSetup, const interaction_const_t& ic)
142 if (EEL_RF(ic.eeltype) || ic.eeltype == CoulombInteractionType::Cut)
148 if (kernelSetup.ewaldExclusionType == Nbnxm::EwaldExclusionType::Table)
150 if (ic.rcoulomb == ic.rvdw)
156 return coulktTAB_TWIN;
161 if (ic.rcoulomb == ic.rvdw)
167 return coulktEWALD_TWIN;
173 static int getVdwKernelType(const Nbnxm::KernelSetup& kernelSetup,
174 const nbnxn_atomdata_t::Params& nbatParams,
175 const interaction_const_t& ic)
177 if (ic.vdwtype == VanDerWaalsType::Cut)
179 switch (ic.vdw_modifier)
181 case InteractionModifiers::None:
182 case InteractionModifiers::PotShift:
183 switch (nbatParams.ljCombinationRule)
185 case LJCombinationRule::Geometric: return vdwktLJCUT_COMBGEOM;
186 case LJCombinationRule::LorentzBerthelot: return vdwktLJCUT_COMBLB;
187 case LJCombinationRule::None: return vdwktLJCUT_COMBNONE;
188 default: gmx_incons("Unknown combination rule");
190 case InteractionModifiers::ForceSwitch: return vdwktLJFORCESWITCH;
191 case InteractionModifiers::PotSwitch: return vdwktLJPOTSWITCH;
193 std::string errorMsg =
194 gmx::formatString("Unsupported VdW interaction modifier %s (%d)",
195 enumValueToString(ic.vdw_modifier),
196 static_cast<int>(ic.vdw_modifier));
197 gmx_incons(errorMsg);
200 else if (ic.vdwtype == VanDerWaalsType::Pme)
202 if (ic.ljpme_comb_rule == LongRangeVdW::Geom)
204 return vdwktLJEWALDCOMBGEOM;
208 /* At setup we (should have) selected the C reference kernel */
209 GMX_RELEASE_ASSERT(kernelSetup.kernelType == Nbnxm::KernelType::Cpu4x4_PlainC,
210 "Only the C reference nbnxn SIMD kernel supports LJ-PME with LB "
211 "combination rules");
212 return vdwktLJEWALDCOMBLB;
217 std::string errorMsg = gmx::formatString("Unsupported VdW interaction type %s (%d)",
218 enumValueToString(ic.vdwtype),
219 static_cast<int>(ic.vdwtype));
220 gmx_incons(errorMsg);
224 /*! \brief Dispatches the non-bonded N versus M atom cluster CPU kernels.
226 * OpenMP parallelization is performed within this function.
227 * Energy reduction, but not force and shift force reduction, is performed
228 * within this function.
230 * \param[in] pairlistSet Pairlists with local or non-local interactions to compute
231 * \param[in] kernelSetup The non-bonded kernel setup
232 * \param[in,out] nbat The atomdata for the interactions
233 * \param[in] ic Non-bonded interaction constants
234 * \param[in] shiftVectors The PBC shift vectors
235 * \param[in] stepWork Flags that tell what to compute
236 * \param[in] clearF Enum that tells if to clear the force output buffer
237 * \param[out] vCoulomb Output buffer for Coulomb energies
238 * \param[out] vVdw Output buffer for Van der Waals energies
239 * \param[in] wcycle Pointer to cycle counting data structure.
241 static void nbnxn_kernel_cpu(const PairlistSet& pairlistSet,
242 const Nbnxm::KernelSetup& kernelSetup,
243 nbnxn_atomdata_t* nbat,
244 const interaction_const_t& ic,
246 const gmx::StepWorkload& stepWork,
250 gmx_wallcycle* wcycle)
253 const nbnxn_atomdata_t::Params& nbatParams = nbat->params();
255 const int coulkt = getCoulombKernelType(kernelSetup, ic);
256 const int vdwkt = getVdwKernelType(kernelSetup, nbatParams, ic);
258 gmx::ArrayRef<const NbnxnPairlistCpu> pairlists = pairlistSet.cpuLists();
260 int gmx_unused nthreads = gmx_omp_nthreads_get(emntNonbonded);
261 wallcycle_sub_start(wcycle, ewcsNONBONDED_CLEAR);
262 #pragma omp parallel for schedule(static) num_threads(nthreads)
263 for (gmx::index nb = 0; nb < pairlists.ssize(); nb++)
265 // Presently, the kernels do not call C++ code that can throw,
266 // so no need for a try/catch pair in this OpenMP region.
267 nbnxn_atomdata_output_t* out = &nbat->out[nb];
269 if (clearF == enbvClearFYes)
271 clearForceBuffer(nbat, nb);
273 clear_fshift(out->fshift.data());
278 wallcycle_sub_stop(wcycle, ewcsNONBONDED_CLEAR);
279 wallcycle_sub_start(wcycle, ewcsNONBONDED_KERNEL);
282 // TODO: Change to reference
283 const NbnxnPairlistCpu* pairlist = &pairlists[nb];
285 if (!stepWork.computeEnergy)
287 /* Don't calculate energies */
288 switch (kernelSetup.kernelType)
290 case Nbnxm::KernelType::Cpu4x4_PlainC:
291 nbnxn_kernel_noener_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
293 #ifdef GMX_NBNXN_SIMD_2XNN
294 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
295 nbnxm_kernel_noener_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
298 #ifdef GMX_NBNXN_SIMD_4XN
299 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
300 nbnxm_kernel_noener_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
303 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
306 else if (out->Vvdw.size() == 1)
308 /* A single energy group (pair) */
312 switch (kernelSetup.kernelType)
314 case Nbnxm::KernelType::Cpu4x4_PlainC:
315 nbnxn_kernel_ener_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
317 #ifdef GMX_NBNXN_SIMD_2XNN
318 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
319 nbnxm_kernel_ener_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
322 #ifdef GMX_NBNXN_SIMD_4XN
323 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
324 nbnxm_kernel_ener_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
327 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
332 /* Calculate energy group contributions */
333 clearGroupEnergies(out);
337 switch (kernelSetup.kernelType)
339 case Nbnxm::KernelType::Cpu4x4_PlainC:
340 unrollj = c_nbnxnCpuIClusterSize;
341 nbnxn_kernel_energrp_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
343 #ifdef GMX_NBNXN_SIMD_2XNN
344 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
345 unrollj = GMX_SIMD_REAL_WIDTH / 2;
346 nbnxm_kernel_energrp_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
349 #ifdef GMX_NBNXN_SIMD_4XN
350 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
351 unrollj = GMX_SIMD_REAL_WIDTH;
352 nbnxm_kernel_energrp_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
355 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
358 if (kernelSetup.kernelType != Nbnxm::KernelType::Cpu4x4_PlainC)
363 reduceGroupEnergySimdBuffers<2>(nbatParams.nenergrp, nbatParams.neg_2log, out);
366 reduceGroupEnergySimdBuffers<4>(nbatParams.nenergrp, nbatParams.neg_2log, out);
369 reduceGroupEnergySimdBuffers<8>(nbatParams.nenergrp, nbatParams.neg_2log, out);
371 default: GMX_RELEASE_ASSERT(false, "Unsupported j-unroll size");
376 wallcycle_sub_stop(wcycle, ewcsNONBONDED_KERNEL);
378 if (stepWork.computeEnergy)
380 reduce_energies_over_lists(nbat, pairlists.ssize(), vVdw, vCoulomb);
384 static void accountFlops(t_nrnb* nrnb,
385 const PairlistSet& pairlistSet,
386 const nonbonded_verlet_t& nbv,
387 const interaction_const_t& ic,
388 const gmx::StepWorkload& stepWork)
390 const bool usingGpuKernels = nbv.useGpu();
392 int enr_nbnxn_kernel_ljc = eNRNB;
393 if (EEL_RF(ic.eeltype) || ic.eeltype == CoulombInteractionType::Cut)
395 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
397 else if ((!usingGpuKernels && nbv.kernelSetup().ewaldExclusionType == Nbnxm::EwaldExclusionType::Analytical)
398 || (usingGpuKernels && Nbnxm::gpu_is_kernel_ewald_analytical(nbv.gpu_nbv)))
400 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
404 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
406 int enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
407 if (stepWork.computeEnergy)
409 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
410 enr_nbnxn_kernel_ljc += 1;
411 enr_nbnxn_kernel_lj += 1;
414 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc, pairlistSet.natpair_ljq_);
415 inc_nrnb(nrnb, enr_nbnxn_kernel_lj, pairlistSet.natpair_lj_);
416 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
417 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc - eNR_NBNXN_LJ_RF + eNR_NBNXN_RF, pairlistSet.natpair_q_);
419 if (ic.vdw_modifier == InteractionModifiers::ForceSwitch)
421 /* We add up the switch cost separately */
423 eNR_NBNXN_ADD_LJ_FSW + (stepWork.computeEnergy ? 1 : 0),
424 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
426 if (ic.vdw_modifier == InteractionModifiers::PotSwitch)
428 /* We add up the switch cost separately */
430 eNR_NBNXN_ADD_LJ_PSW + (stepWork.computeEnergy ? 1 : 0),
431 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
433 if (ic.vdwtype == VanDerWaalsType::Pme)
435 /* We add up the LJ Ewald cost separately */
437 eNR_NBNXN_ADD_LJ_EWALD + (stepWork.computeEnergy ? 1 : 0),
438 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
442 void nonbonded_verlet_t::dispatchNonbondedKernel(gmx::InteractionLocality iLocality,
443 const interaction_const_t& ic,
444 const gmx::StepWorkload& stepWork,
446 const t_forcerec& fr,
447 gmx_enerdata_t* enerd,
450 const PairlistSet& pairlistSet = pairlistSets().pairlistSet(iLocality);
452 switch (kernelSetup().kernelType)
454 case Nbnxm::KernelType::Cpu4x4_PlainC:
455 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
456 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
465 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
466 fr.bBHAM ? enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::BuckinghamSR].data()
467 : enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
471 case Nbnxm::KernelType::Gpu8x8x8:
472 Nbnxm::gpu_launch_kernel(gpu_nbv, stepWork, iLocality);
475 case Nbnxm::KernelType::Cpu8x8x8_PlainC:
476 nbnxn_kernel_gpu_ref(
477 pairlistSet.gpuList(),
484 nbat->out[0].fshift.data(),
485 enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
486 fr.bBHAM ? enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::BuckinghamSR].data()
487 : enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data());
490 default: GMX_RELEASE_ASSERT(false, "Invalid nonbonded kernel type passed!");
493 accountFlops(nrnb, pairlistSet, *this, ic, stepWork);
496 void nonbonded_verlet_t::dispatchFreeEnergyKernel(gmx::InteractionLocality iLocality,
497 const t_forcerec& fr,
498 gmx::ArrayRef<const gmx::RVec> coords,
499 gmx::ForceWithShiftForces* forceWithShiftForces,
500 gmx::ArrayRef<const real> chargeA,
501 gmx::ArrayRef<const real> chargeB,
502 gmx::ArrayRef<const int> typeA,
503 gmx::ArrayRef<const int> typeB,
505 gmx::ArrayRef<const real> lambda,
506 gmx_enerdata_t* enerd,
507 const gmx::StepWorkload& stepWork,
510 const auto nbl_fep = pairlistSets().pairlistSet(iLocality).fepLists();
512 /* When the first list is empty, all are empty and there is nothing to do */
513 if (!pairlistSets().params().haveFep || nbl_fep[0]->nrj == 0)
519 /* Add short-range interactions */
520 donb_flags |= GMX_NONBONDED_DO_SR;
522 if (stepWork.computeForces)
524 donb_flags |= GMX_NONBONDED_DO_FORCE;
526 if (stepWork.computeVirial)
528 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
530 if (stepWork.computeEnergy)
532 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
535 gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, real> dvdl_nb = { 0 };
536 int kernelFlags = donb_flags;
537 gmx::ArrayRef<const real> kernelLambda = lambda;
538 gmx::ArrayRef<real> kernelDvdl = dvdl_nb;
540 gmx::ArrayRef<real> energygrp_elec = enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR];
541 gmx::ArrayRef<real> energygrp_vdw = enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR];
543 GMX_ASSERT(gmx_omp_nthreads_get(emntNonbonded) == nbl_fep.ssize(),
544 "Number of lists should be same as number of NB threads");
546 wallcycle_sub_start(wcycle_, ewcsNONBONDED_FEP);
547 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
548 for (gmx::index th = 0; th < nbl_fep.ssize(); th++)
552 gmx_nb_free_energy_kernel(*nbl_fep[th],
554 forceWithShiftForces,
567 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
570 if (fepvals->sc_alpha != 0)
572 enerd->dvdl_nonlin[FreeEnergyPerturbationCouplingType::Vdw] +=
573 dvdl_nb[FreeEnergyPerturbationCouplingType::Vdw];
574 enerd->dvdl_nonlin[FreeEnergyPerturbationCouplingType::Coul] +=
575 dvdl_nb[FreeEnergyPerturbationCouplingType::Coul];
579 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] +=
580 dvdl_nb[FreeEnergyPerturbationCouplingType::Vdw];
581 enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Coul] +=
582 dvdl_nb[FreeEnergyPerturbationCouplingType::Coul];
585 /* If we do foreign lambda and we have soft-core interactions
586 * we have to recalculate the (non-linear) energies contributions.
588 if (fepvals->n_lambda > 0 && stepWork.computeDhdl && fepvals->sc_alpha != 0)
590 gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, real> lam_i;
591 kernelFlags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE))
592 | GMX_NONBONDED_DO_FOREIGNLAMBDA;
593 kernelLambda = lam_i;
594 kernelDvdl = dvdl_nb;
595 gmx::ArrayRef<real> energygrp_elec =
596 enerd->foreign_grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR];
597 gmx::ArrayRef<real> energygrp_vdw =
598 enerd->foreign_grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR];
600 for (gmx::index i = 0; i < 1 + enerd->foreignLambdaTerms.numLambdas(); i++)
602 std::fill(std::begin(dvdl_nb), std::end(dvdl_nb), 0);
603 for (int j = 0; j < static_cast<int>(FreeEnergyPerturbationCouplingType::Count); j++)
605 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i - 1]);
607 reset_foreign_enerdata(enerd);
608 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
609 for (gmx::index th = 0; th < nbl_fep.ssize(); th++)
613 gmx_nb_free_energy_kernel(*nbl_fep[th],
615 forceWithShiftForces,
628 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
631 sum_epot(enerd->foreign_grpp, enerd->foreign_term);
632 enerd->foreignLambdaTerms.accumulate(
634 enerd->foreign_term[F_EPOT],
635 dvdl_nb[FreeEnergyPerturbationCouplingType::Vdw]
636 + dvdl_nb[FreeEnergyPerturbationCouplingType::Coul]);
639 wallcycle_sub_stop(wcycle_, ewcsNONBONDED_FEP);