Separate CPU NB kernel and buffer clearing subcounters
[alexxy/gromacs.git] / src / gromacs / nbnxm / kerneldispatch.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
36 #include "gmxpre.h"
37
38 #include "gromacs/gmxlib/nrnb.h"
39 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
40 #include "gromacs/gmxlib/nonbonded/nb_kernel.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/force_flags.h"
46 #include "gromacs/mdlib/gmx_omp_nthreads.h"
47 #include "gromacs/mdtypes/enerdata.h"
48 #include "gromacs/mdtypes/inputrec.h"
49 #include "gromacs/mdtypes/interaction_const.h"
50 #include "gromacs/mdtypes/md_enums.h"
51 #include "gromacs/mdtypes/mdatom.h"
52 #include "gromacs/nbnxm/gpu_data_mgmt.h"
53 #include "gromacs/nbnxm/nbnxm.h"
54 #include "gromacs/nbnxm/nbnxm_simd.h"
55 #include "gromacs/nbnxm/kernels_reference/kernel_gpu_ref.h"
56 #include "gromacs/simd/simd.h"
57 #include "gromacs/timing/wallcycle.h"
58 #include "gromacs/utility/gmxassert.h"
59 #include "gromacs/utility/real.h"
60
61 #include "kernel_common.h"
62 #include "pairlistset.h"
63 #include "pairlistsets.h"
64 #define INCLUDE_KERNELFUNCTION_TABLES
65 #include "gromacs/nbnxm/kernels_reference/kernel_ref.h"
66 #ifdef GMX_NBNXN_SIMD_2XNN
67 #include "gromacs/nbnxm/kernels_simd_2xmm/kernels.h"
68 #endif
69 #ifdef GMX_NBNXN_SIMD_4XN
70 #include "gromacs/nbnxm/kernels_simd_4xm/kernels.h"
71 #endif
72 #undef INCLUDE_FUNCTION_TABLES
73
74 /*! \brief Clears the energy group output buffers
75  *
76  * \param[in,out] out  nbnxn kernel output struct
77  */
78 static void clearGroupEnergies(nbnxn_atomdata_output_t *out)
79 {
80     std::fill(out->Vvdw.begin(), out->Vvdw.end(), 0.0_real);
81     std::fill(out->Vc.begin(), out->Vc.end(), 0.0_real);
82     std::fill(out->VSvdw.begin(), out->VSvdw.end(), 0.0_real);
83     std::fill(out->VSc.begin(), out->VSc.end(), 0.0_real);
84 }
85
86 /*! \brief Reduce the group-pair energy buffers produced by a SIMD kernel
87  * to single terms in the output buffers.
88  *
89  * The SIMD kernels produce a large number of energy buffer in SIMD registers
90  * to avoid scattered reads and writes.
91  *
92  * \tparam        unrollj         The unroll size for j-particles in the SIMD kernel
93  * \param[in]     numGroups       The number of energy groups
94  * \param[in]     numGroups_2log  Log2 of numGroups, rounded up
95  * \param[in,out] out             Struct with energy buffers
96  */
97 template <int unrollj> static void
98 reduceGroupEnergySimdBuffers(int                       numGroups,
99                              int                       numGroups_2log,
100                              nbnxn_atomdata_output_t  *out)
101 {
102     const int                 unrollj_half     = unrollj/2;
103     /* Energies are stored in SIMD registers with size 2^numGroups_2log */
104     const int                 numGroupsStorage = (1 << numGroups_2log);
105
106     const real * gmx_restrict vVdwSimd     = out->VSvdw.data();
107     const real * gmx_restrict vCoulombSimd = out->VSc.data();
108     real * gmx_restrict       vVdw         = out->Vvdw.data();
109     real * gmx_restrict       vCoulomb     = out->Vc.data();
110
111     /* The size of the SIMD energy group buffer array is:
112      * numGroups*numGroups*numGroupsStorage*unrollj_half*simd_width
113      */
114     for (int i = 0; i < numGroups; i++)
115     {
116         for (int j1 = 0; j1 < numGroups; j1++)
117         {
118             for (int j0 = 0; j0 < numGroups; j0++)
119             {
120                 int c = ((i*numGroups + j1)*numGroupsStorage + j0)*unrollj_half*unrollj;
121                 for (int s = 0; s < unrollj_half; s++)
122                 {
123                     vVdw    [i*numGroups + j0] += vVdwSimd    [c + 0];
124                     vVdw    [i*numGroups + j1] += vVdwSimd    [c + 1];
125                     vCoulomb[i*numGroups + j0] += vCoulombSimd[c + 0];
126                     vCoulomb[i*numGroups + j1] += vCoulombSimd[c + 1];
127                     c                          += unrollj + 2;
128                 }
129             }
130         }
131     }
132 }
133
134 /*! \brief Dispatches the non-bonded N versus M atom cluster CPU kernels.
135  *
136  * OpenMP parallelization is performed within this function.
137  * Energy reduction, but not force and shift force reduction, is performed
138  * within this function.
139  *
140  * \param[in]     pairlistSet   Pairlists with local or non-local interactions to compute
141  * \param[in]     kernelSetup   The non-bonded kernel setup
142  * \param[in,out] nbat          The atomdata for the interactions
143  * \param[in]     ic            Non-bonded interaction constants
144  * \param[in]     shiftVectors  The PBC shift vectors
145  * \param[in]     forceFlags    Flags that tell what to compute
146  * \param[in]     clearF        Enum that tells if to clear the force output buffer
147  * \param[out]    vCoulomb      Output buffer for Coulomb energies
148  * \param[out]    vVdw          Output buffer for Van der Waals energies
149  * \param[in]     wcycle        Pointer to cycle counting data structure.
150  */
151 static void
152 nbnxn_kernel_cpu(const PairlistSet              &pairlistSet,
153                  const Nbnxm::KernelSetup       &kernelSetup,
154                  nbnxn_atomdata_t               *nbat,
155                  const interaction_const_t      &ic,
156                  rvec                           *shiftVectors,
157                  int                             forceFlags,
158                  int                             clearF,
159                  real                           *vCoulomb,
160                  real                           *vVdw,
161                  gmx_wallcycle                  *wcycle)
162 {
163
164     int                      coulkt;
165     if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
166     {
167         coulkt = coulktRF;
168     }
169     else
170     {
171         if (kernelSetup.ewaldExclusionType == Nbnxm::EwaldExclusionType::Table)
172         {
173             if (ic.rcoulomb == ic.rvdw)
174             {
175                 coulkt = coulktTAB;
176             }
177             else
178             {
179                 coulkt = coulktTAB_TWIN;
180             }
181         }
182         else
183         {
184             if (ic.rcoulomb == ic.rvdw)
185             {
186                 coulkt = coulktEWALD;
187             }
188             else
189             {
190                 coulkt = coulktEWALD_TWIN;
191             }
192         }
193     }
194
195     const nbnxn_atomdata_t::Params &nbatParams = nbat->params();
196
197     int vdwkt = 0;
198     if (ic.vdwtype == evdwCUT)
199     {
200         switch (ic.vdw_modifier)
201         {
202             case eintmodNONE:
203             case eintmodPOTSHIFT:
204                 switch (nbatParams.comb_rule)
205                 {
206                     case ljcrGEOM: vdwkt = vdwktLJCUT_COMBGEOM; break;
207                     case ljcrLB:   vdwkt = vdwktLJCUT_COMBLB;   break;
208                     case ljcrNONE: vdwkt = vdwktLJCUT_COMBNONE; break;
209                     default:
210                         GMX_RELEASE_ASSERT(false, "Unknown combination rule");
211                 }
212                 break;
213             case eintmodFORCESWITCH:
214                 vdwkt = vdwktLJFORCESWITCH;
215                 break;
216             case eintmodPOTSWITCH:
217                 vdwkt = vdwktLJPOTSWITCH;
218                 break;
219             default:
220                 GMX_RELEASE_ASSERT(false, "Unsupported VdW interaction modifier");
221         }
222     }
223     else if (ic.vdwtype == evdwPME)
224     {
225         if (ic.ljpme_comb_rule == eljpmeGEOM)
226         {
227             vdwkt = vdwktLJEWALDCOMBGEOM;
228         }
229         else
230         {
231             vdwkt = vdwktLJEWALDCOMBLB;
232             /* At setup we (should have) selected the C reference kernel */
233             GMX_RELEASE_ASSERT(kernelSetup.kernelType == Nbnxm::KernelType::Cpu4x4_PlainC, "Only the C reference nbnxn SIMD kernel supports LJ-PME with LB combination rules");
234         }
235     }
236     else
237     {
238         GMX_RELEASE_ASSERT(false, "Unsupported VdW interaction type");
239     }
240
241     gmx::ArrayRef<const NbnxnPairlistCpu> pairlists = pairlistSet.cpuLists();
242
243     int gmx_unused                        nthreads = gmx_omp_nthreads_get(emntNonbonded);
244     wallcycle_sub_start(wcycle, ewcsNBFCLEARBUF);
245 #pragma omp parallel for schedule(static) num_threads(nthreads)
246     for (int nb = 0; nb < pairlists.ssize(); nb++)
247     {
248         // Presently, the kernels do not call C++ code that can throw,
249         // so no need for a try/catch pair in this OpenMP region.
250         nbnxn_atomdata_output_t *out = &nbat->out[nb];
251
252         if (clearF == enbvClearFYes)
253         {
254             clearForceBuffer(nbat, nb);
255
256             clear_fshift(out->fshift.data());
257         }
258
259         if (nb == 0)
260         {
261             wallcycle_sub_stop(wcycle, ewcsNBFCLEARBUF);
262             wallcycle_sub_start(wcycle, ewcsNBFKERNEL);
263         }
264
265         // TODO: Change to reference
266         const NbnxnPairlistCpu *pairlist = &pairlists[nb];
267
268         if (!(forceFlags & GMX_FORCE_ENERGY))
269         {
270             /* Don't calculate energies */
271             switch (kernelSetup.kernelType)
272             {
273                 case Nbnxm::KernelType::Cpu4x4_PlainC:
274                     nbnxn_kernel_noener_ref[coulkt][vdwkt](pairlist, nbat,
275                                                            &ic,
276                                                            shiftVectors,
277                                                            out);
278                     break;
279 #ifdef GMX_NBNXN_SIMD_2XNN
280                 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
281                     nbnxm_kernel_noener_simd_2xmm[coulkt][vdwkt](pairlist, nbat,
282                                                                  &ic,
283                                                                  shiftVectors,
284                                                                  out);
285                     break;
286 #endif
287 #ifdef GMX_NBNXN_SIMD_4XN
288                 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
289                     nbnxm_kernel_noener_simd_4xm[coulkt][vdwkt](pairlist, nbat,
290                                                                 &ic,
291                                                                 shiftVectors,
292                                                                 out);
293                     break;
294 #endif
295                 default:
296                     GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
297             }
298         }
299         else if (out->Vvdw.size() == 1)
300         {
301             /* A single energy group (pair) */
302             out->Vvdw[0] = 0;
303             out->Vc[0]   = 0;
304
305             switch (kernelSetup.kernelType)
306             {
307                 case Nbnxm::KernelType::Cpu4x4_PlainC:
308                     nbnxn_kernel_ener_ref[coulkt][vdwkt](pairlist, nbat,
309                                                          &ic,
310                                                          shiftVectors,
311                                                          out);
312                     break;
313 #ifdef GMX_NBNXN_SIMD_2XNN
314                 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
315                     nbnxm_kernel_ener_simd_2xmm[coulkt][vdwkt](pairlist, nbat,
316                                                                &ic,
317                                                                shiftVectors,
318                                                                out);
319                     break;
320 #endif
321 #ifdef GMX_NBNXN_SIMD_4XN
322                 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
323                     nbnxm_kernel_ener_simd_4xm[coulkt][vdwkt](pairlist, nbat,
324                                                               &ic,
325                                                               shiftVectors,
326                                                               out);
327                     break;
328 #endif
329                 default:
330                     GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
331             }
332         }
333         else
334         {
335             /* Calculate energy group contributions */
336             clearGroupEnergies(out);
337
338             int unrollj = 0;
339
340             switch (kernelSetup.kernelType)
341             {
342                 case Nbnxm::KernelType::Cpu4x4_PlainC:
343                     unrollj = c_nbnxnCpuIClusterSize;
344                     nbnxn_kernel_energrp_ref[coulkt][vdwkt](pairlist, nbat,
345                                                             &ic,
346                                                             shiftVectors,
347                                                             out);
348                     break;
349 #ifdef GMX_NBNXN_SIMD_2XNN
350                 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
351                     unrollj = GMX_SIMD_REAL_WIDTH/2;
352                     nbnxm_kernel_energrp_simd_2xmm[coulkt][vdwkt](pairlist, nbat,
353                                                                   &ic,
354                                                                   shiftVectors,
355                                                                   out);
356                     break;
357 #endif
358 #ifdef GMX_NBNXN_SIMD_4XN
359                 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
360                     unrollj = GMX_SIMD_REAL_WIDTH;
361                     nbnxm_kernel_energrp_simd_4xm[coulkt][vdwkt](pairlist, nbat,
362                                                                  &ic,
363                                                                  shiftVectors,
364                                                                  out);
365                     break;
366 #endif
367                 default:
368                     GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
369             }
370
371             if (kernelSetup.kernelType != Nbnxm::KernelType::Cpu4x4_PlainC)
372             {
373                 switch (unrollj)
374                 {
375                     case 2:
376                         reduceGroupEnergySimdBuffers<2>(nbatParams.nenergrp,
377                                                         nbatParams.neg_2log,
378                                                         out);
379                         break;
380                     case 4:
381                         reduceGroupEnergySimdBuffers<4>(nbatParams.nenergrp,
382                                                         nbatParams.neg_2log,
383                                                         out);
384                         break;
385                     case 8:
386                         reduceGroupEnergySimdBuffers<8>(nbatParams.nenergrp,
387                                                         nbatParams.neg_2log,
388                                                         out);
389                         break;
390                     default:
391                         GMX_RELEASE_ASSERT(false, "Unsupported j-unroll size");
392                 }
393             }
394         }
395     }
396     wallcycle_sub_stop(wcycle, ewcsNBFKERNEL);
397
398     if (forceFlags & GMX_FORCE_ENERGY)
399     {
400         reduce_energies_over_lists(nbat, pairlists.ssize(), vVdw, vCoulomb);
401     }
402 }
403
404 static void accountFlops(t_nrnb                           *nrnb,
405                          const PairlistSet                &pairlistSet,
406                          const nonbonded_verlet_t         &nbv,
407                          const interaction_const_t        &ic,
408                          const int                         forceFlags)
409 {
410     const bool usingGpuKernels = nbv.useGpu();
411
412     int        enr_nbnxn_kernel_ljc;
413     if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
414     {
415         enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
416     }
417     else if ((!usingGpuKernels && nbv.kernelSetup().ewaldExclusionType == Nbnxm::EwaldExclusionType::Analytical) ||
418              (usingGpuKernels && Nbnxm::gpu_is_kernel_ewald_analytical(nbv.gpu_nbv)))
419     {
420         enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
421     }
422     else
423     {
424         enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
425     }
426     int enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
427     if (forceFlags & GMX_FORCE_ENERGY)
428     {
429         /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
430         enr_nbnxn_kernel_ljc += 1;
431         enr_nbnxn_kernel_lj  += 1;
432     }
433
434     inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
435              pairlistSet.natpair_ljq_);
436     inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
437              pairlistSet.natpair_lj_);
438     /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
439     inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
440              pairlistSet.natpair_q_);
441
442     const bool calcEnergy = ((forceFlags & GMX_FORCE_ENERGY) != 0);
443     if (ic.vdw_modifier == eintmodFORCESWITCH)
444     {
445         /* We add up the switch cost separately */
446         inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW + (calcEnergy ? 1 : 0),
447                  pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
448     }
449     if (ic.vdw_modifier == eintmodPOTSWITCH)
450     {
451         /* We add up the switch cost separately */
452         inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW + (calcEnergy ? 1 : 0),
453                  pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
454     }
455     if (ic.vdwtype == evdwPME)
456     {
457         /* We add up the LJ Ewald cost separately */
458         inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD + (calcEnergy ? 1 : 0),
459                  pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
460     }
461 }
462
463 void
464 nonbonded_verlet_t::dispatchNonbondedKernel(Nbnxm::InteractionLocality iLocality,
465                                             const interaction_const_t &ic,
466                                             int                        forceFlags,
467                                             int                        clearF,
468                                             const t_forcerec          &fr,
469                                             gmx_enerdata_t            *enerd,
470                                             t_nrnb                    *nrnb,
471                                             gmx_wallcycle             *wcycle)
472 {
473     const PairlistSet &pairlistSet = pairlistSets().pairlistSet(iLocality);
474
475     switch (kernelSetup().kernelType)
476     {
477         case Nbnxm::KernelType::Cpu4x4_PlainC:
478         case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
479         case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
480             nbnxn_kernel_cpu(pairlistSet,
481                              kernelSetup(),
482                              nbat.get(),
483                              ic,
484                              fr.shift_vec,
485                              forceFlags,
486                              clearF,
487                              enerd->grpp.ener[egCOULSR].data(),
488                              fr.bBHAM ?
489                              enerd->grpp.ener[egBHAMSR].data() :
490                              enerd->grpp.ener[egLJSR].data(),
491                              wcycle);
492             break;
493
494         case Nbnxm::KernelType::Gpu8x8x8:
495             Nbnxm::gpu_launch_kernel(gpu_nbv, forceFlags, iLocality);
496             break;
497
498         case Nbnxm::KernelType::Cpu8x8x8_PlainC:
499             nbnxn_kernel_gpu_ref(pairlistSet.gpuList(),
500                                  nbat.get(), &ic,
501                                  fr.shift_vec,
502                                  forceFlags,
503                                  clearF,
504                                  nbat->out[0].f,
505                                  nbat->out[0].fshift.data(),
506                                  enerd->grpp.ener[egCOULSR].data(),
507                                  fr.bBHAM ?
508                                  enerd->grpp.ener[egBHAMSR].data() :
509                                  enerd->grpp.ener[egLJSR].data());
510             break;
511
512         default:
513             GMX_RELEASE_ASSERT(false, "Invalid nonbonded kernel type passed!");
514
515     }
516
517     accountFlops(nrnb, pairlistSet, *this, ic, forceFlags);
518 }
519
520 void
521 nonbonded_verlet_t::dispatchFreeEnergyKernel(Nbnxm::InteractionLocality  iLocality,
522                                              t_forcerec                 *fr,
523                                              rvec                        x[],
524                                              rvec                        f[],
525                                              const t_mdatoms            &mdatoms,
526                                              t_lambda                   *fepvals,
527                                              real                       *lambda,
528                                              gmx_enerdata_t             *enerd,
529                                              const int                   forceFlags,
530                                              t_nrnb                     *nrnb)
531 {
532     const auto nbl_fep = pairlistSets().pairlistSet(iLocality).fepLists();
533
534     /* When the first list is empty, all are empty and there is nothing to do */
535     if (!pairlistSets().params().haveFep || nbl_fep[0]->nrj == 0)
536     {
537         return;
538     }
539
540     int donb_flags = 0;
541     /* Add short-range interactions */
542     donb_flags |= GMX_NONBONDED_DO_SR;
543
544     /* Currently all group scheme kernels always calculate (shift-)forces */
545     if (forceFlags & GMX_FORCE_FORCES)
546     {
547         donb_flags |= GMX_NONBONDED_DO_FORCE;
548     }
549     if (forceFlags & GMX_FORCE_VIRIAL)
550     {
551         donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
552     }
553     if (forceFlags & GMX_FORCE_ENERGY)
554     {
555         donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
556     }
557
558     nb_kernel_data_t kernel_data;
559     real             dvdl_nb[efptNR] = { 0 };
560     kernel_data.flags  = donb_flags;
561     kernel_data.lambda = lambda;
562     kernel_data.dvdl   = dvdl_nb;
563
564     kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR].data();
565     kernel_data.energygrp_vdw  = enerd->grpp.ener[egLJSR].data();
566
567     GMX_ASSERT(gmx_omp_nthreads_get(emntNonbonded) == nbl_fep.ssize(), "Number of lists should be same as number of NB threads");
568
569 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
570     for (int th = 0; th < nbl_fep.ssize(); th++)
571     {
572         try
573         {
574             gmx_nb_free_energy_kernel(nbl_fep[th].get(),
575                                       x, f, fr, &mdatoms, &kernel_data, nrnb);
576         }
577         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
578     }
579
580     if (fepvals->sc_alpha != 0)
581     {
582         enerd->dvdl_nonlin[efptVDW]  += dvdl_nb[efptVDW];
583         enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
584     }
585     else
586     {
587         enerd->dvdl_lin[efptVDW]  += dvdl_nb[efptVDW];
588         enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
589     }
590
591     /* If we do foreign lambda and we have soft-core interactions
592      * we have to recalculate the (non-linear) energies contributions.
593      */
594     if (fepvals->n_lambda > 0 && (forceFlags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
595     {
596         real lam_i[efptNR];
597         kernel_data.flags          = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
598         kernel_data.lambda         = lam_i;
599         kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR].data();
600         kernel_data.energygrp_vdw  = enerd->foreign_grpp.ener[egLJSR].data();
601         /* Note that we add to kernel_data.dvdl, but ignore the result */
602
603         for (size_t i = 0; i < enerd->enerpart_lambda.size(); i++)
604         {
605             for (int j = 0; j < efptNR; j++)
606             {
607                 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
608             }
609             reset_foreign_enerdata(enerd);
610 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
611             for (int th = 0; th < nbl_fep.ssize(); th++)
612             {
613                 try
614                 {
615                     gmx_nb_free_energy_kernel(nbl_fep[th].get(),
616                                               x, f, fr, &mdatoms, &kernel_data, nrnb);
617                 }
618                 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
619             }
620
621             sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
622             enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
623         }
624     }
625 }